ParaMonte Fortran 2.0.0
Parallel Monte Carlo and Machine Learning Library
See the latest version documentation.
pm_matrixDet::getMatDetSqrt Interface Reference

Generate and return the square-root of the determinant of the input positive-definite square matrix. More...

Detailed Description

Generate and return the square-root of the determinant of the input positive-definite square matrix.

The positive-definiteness guarantee by the user enables a fast method of computing the determinant of the input positive-definite square matrix via its Cholesky Factorization.
For high rank matrices, there is an overflow possibility for the output of this generic interface.
The overflow risk can be avoided by calling getMatDetSqrtLog instead.

Parameters
[in]mat: The input contiguous positive-definite square matrix of shape (ndim,ndim) of type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128). Only the upper triangle and diagonals of the matrix are needed and used to compute the Cholesky Factorization.
[in]subset: The input scalar that can be either,
  1. the constant uppDia implying that the upper-diagonal triangular block of mat should be used for building the lower-diagonal Cholesky factor in the output chol.
  2. the constant lowDia implying that the lower-diagonal triangular block of mat should be used for building the upper-diagonal Cholesky factor in the output chol.
This argument is merely a convenience to differentiate the different procedure functionalities within this generic interface.
(optional, default = uppDia)
Returns
detSqrtLog : The output scalar of type real of the same kind as the input mat containing the square-root of the determinant of the input positive-definite square matrix.
If the input mat is not positive-definite, the algorithm will halt the program by calling error stop.


Possible calling interfaces

detSqrtLog = getMatDetSqrtLog(mat, subset = subset)
Generate and return the natural logarithm of the square-root of the determinant of the input positive...
This module contains procedures and generic interfaces relevant to the computation of the determinant...
Warning
The condition size(mat, 1) == size(mat, 2) must hold for the corresponding input arguments.
This condition is verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The input mat is assumed to be positive-definite, otherwise the Cholesky Factorization within the procedure will fail, leading to the invocation of the error stop Fortran statement.
See setMatDetSqrtLog for the faster fault-tolerant subroutine versions of these procedures.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
Remarks
The procedures under this generic interface are potentially computationally demanding compared to setMatDetSqrtLog.
See setMatDetSqrtLog for a more performant subroutine version of this generic interface for repeated calls.
See also
getMatDet
setMatDet
getMatDetSqrt
setMatDetSqrt
getMatDetSqrtLog
setMatDetSqrtLog
getMatMulTrace


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_matrixCopy, only: rdpack
8 use pm_matrixChol, only: uppDia, lowDia
10 use pm_distCov, only: getCovRand
11 use pm_arrayResize, only: setResized
12 use pm_distUnif, only: getUnifRand
13 use pm_matrixInv, only: getMatInv
14 use pm_io, only: display_type
15
16 implicit none
17
18 integer(IK) :: ndim, itry, ntry = 10
19 type(display_type) :: disp
20
21 disp = display_type(file = "main.out.F90")
22
23 call disp%skip()
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%show("! Compute the sqrt of the determinant of the positive definite matrix.")
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%skip()
28
29 block
30 use pm_kind, only: TKG => RKS
31 real(TKG), allocatable :: mat(:,:)
32 real(TKG) :: detSqrtLog
33 do itry = 1, ntry
34 call disp%skip()
35 call disp%show("ndim = getUnifRand(1, 5)")
36 ndim = getUnifRand(1, 5)
37 call disp%show("ndim")
38 call disp%show( ndim )
39 call disp%show("mat = getCovRand(mold = 1._TKG, ndim = ndim)")
40 mat = getCovRand(mold = 1._TKG, ndim = ndim)
41 call disp%show("mat")
42 call disp%show( mat )
43 call disp%show("detSqrtLog = getMatDetSqrtLog(mat)")
44 detSqrtLog = getMatDetSqrtLog(mat)
45 call disp%show("detSqrtLog")
46 call disp%show( detSqrtLog )
47 call disp%show("getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.")
48 call disp%show( getMatMulTraceLog(getMatChol(mat, uppDia)) )
49 call disp%show("detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.")
50 call disp%show( detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) )
51 call disp%skip()
52 call disp%show("mat = getCovRand(mold = 1._TKG, ndim = ndim)")
53 mat = getCovRand(mold = 1._TKG, ndim = ndim)
54 call disp%show("mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.")
55 mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG)
56 call disp%show("mat")
57 call disp%show( mat )
58 call disp%show("detSqrtLog = getMatDetSqrtLog(mat, lowDia)")
59 detSqrtLog = getMatDetSqrtLog(mat, lowDia)
60 call disp%show("detSqrtLog")
61 call disp%show( detSqrtLog )
62 call disp%show("getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.")
63 call disp%show( getMatMulTraceLog(getMatChol(mat, lowDia)) )
64 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.")
65 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) )
66 call disp%skip()
67 call disp%show("mat = getCovRand(mold = 1._TKG, ndim = ndim)")
68 mat = getCovRand(mold = 1._TKG, ndim = ndim)
69 call disp%show("mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.")
70 mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG)
71 call disp%show("mat")
72 call disp%show( mat )
73 call disp%show("detSqrtLog = getMatDetSqrtLog(mat, uppDia)")
74 detSqrtLog = getMatDetSqrtLog(mat, uppDia)
75 call disp%show("detSqrtLog")
76 call disp%show( detSqrtLog )
77 call disp%show("getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.")
78 call disp%show( getMatMulTraceLog(getMatChol(mat, uppDia)) )
79 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.")
80 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) )
81 call disp%skip()
82 end do
83 end block
84
85 block
86 use pm_kind, only: TKG => CKS
87 complex(TKG), allocatable :: tmp(:,:)
88 complex(TKG), parameter :: mat(*,*) = reshape( [ (9.0, 0.0), (3.0, 3.0), (3.0, -3.0) &
89 , (3.0, -3.0),(18.0, 0.0), (8.0, -6.0) &
90 , (3.0, 3.0), (8.0, 6.0),(43.0, 0.0) &
91 ], shape = [3, 3], order = [2, 1])
92 real(TKG) :: detSqrtLog
93 call disp%skip()
94 call disp%show("mat")
95 call disp%show( mat )
96 call disp%show("detSqrtLog = getMatDetSqrtLog(mat)")
97 detSqrtLog = getMatDetSqrtLog(mat)
98 call disp%show("detSqrtLog")
99 call disp%show( detSqrtLog )
100 call disp%show("getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.")
101 call disp%show( getMatMulTraceLog(getMatChol(mat, uppDia)) )
102 call disp%show("detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.")
103 call disp%show( detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) )
104 call disp%skip()
105 call disp%show("tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
106 tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG))
107 call disp%show("tmp")
108 call disp%show( tmp )
109 call disp%show("detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)")
110 detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
111 call disp%show("detSqrtLog")
112 call disp%show( detSqrtLog )
113 call disp%show("getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.")
114 call disp%show( getMatMulTraceLog(getMatChol(tmp, lowDia)) )
115 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one")
116 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) )
117 call disp%skip()
118 call disp%show("tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
119 tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG))
120 call disp%show("tmp")
121 call disp%show( tmp )
122 call disp%show("detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)")
123 detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
124 call disp%show("detSqrtLog")
125 call disp%show( detSqrtLog )
126 call disp%show("getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.")
127 call disp%show( getMatMulTraceLog(getMatChol(tmp, uppDia)) )
128 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.")
129 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) )
130 call disp%skip()
131 end block
132
133 block
134 use pm_kind, only: TKG => CKS
135 complex(TKG), allocatable :: tmp(:,:)
136 complex(TKG), parameter :: mat(*,*) = reshape( [ (25.0, 0.0), (-5.0, -5.0), (10.0, 5.0) &
137 , (-5.0, 5.0), (51.0, 0.0), (4.0, -6.0) &
138 , (10.0, -5.0), (4.0, 6.0), (71.0, 0.0) &
139 ], shape = [3, 3], order = [2, 1])
140 real(TKG) :: detSqrtLog
141 call disp%skip()
142 call disp%show("mat")
143 call disp%show( mat )
144 call disp%show("detSqrtLog = getMatDetSqrtLog(mat)")
145 detSqrtLog = getMatDetSqrtLog(mat)
146 call disp%show("detSqrtLog")
147 call disp%show( detSqrtLog )
148 call disp%show("getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.")
149 call disp%show( getMatMulTraceLog(getMatChol(mat, uppDia)) )
150 call disp%show("detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.")
151 call disp%show( detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) )
152 call disp%skip()
153 call disp%show("tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
154 tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG))
155 call disp%show("tmp")
156 call disp%show( tmp )
157 call disp%show("detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)")
158 detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
159 call disp%show("detSqrtLog")
160 call disp%show( detSqrtLog )
161 call disp%show("getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.")
162 call disp%show( getMatMulTraceLog(getMatChol(tmp, lowDia)) )
163 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one.")
164 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) )
165 call disp%skip()
166 call disp%show("tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
167 tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG))
168 call disp%show("tmp")
169 call disp%show( tmp )
170 call disp%show("detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)")
171 detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
172 call disp%show("detSqrtLog")
173 call disp%show( detSqrtLog )
174 call disp%show("getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.")
175 call disp%show( getMatMulTraceLog(getMatChol(tmp, uppDia)) )
176 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.")
177 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) )
178 call disp%skip()
179 end block
180
181end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return a random positive-definite (correlation or covariance) matrix using the Gram meth...
Definition: pm_distCov.F90:394
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
Generate and return the upper or the lower Cholesky factorization of the input symmetric positive-def...
Generate and return a copy of a desired subset of the input source matrix of arbitrary shape (:) or (...
Generate and return the full inverse of an input matrix of general or triangular form directly or thr...
Generate and return the natural logarithm of the multiplicative trace of an input square matrix of ty...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains classes and procedures for generating random matrices distributed on the space o...
Definition: pm_distCov.F90:72
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for input/output (IO) or generic display operations on st...
Definition: pm_io.F90:252
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
Definition: pm_io.F90:11393
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
Definition: pm_kind.F90:268
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
Definition: pm_kind.F90:541
integer, parameter CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
Definition: pm_kind.F90:570
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
Definition: pm_kind.F90:540
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
Definition: pm_kind.F90:539
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
Definition: pm_kind.F90:567
This module contains procedures and generic interfaces for computing the Cholesky factorization of po...
This module contains procedures and generic interfaces relevant to copying (diagonal or upper/lower t...
This module contains abstract and concrete derived types and procedures related to the inversion of s...
This module contains procedures and generic interfaces for computing the additive or multiplicative t...
Generate and return an object of type display_type.
Definition: pm_io.F90:10282

Example Unix compile command via Intel ifort compiler
1#!/usr/bin/env sh
2rm main.exe
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example Windows Batch compile command via Intel ifort compiler
1del main.exe
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
4main.exe

Example Unix / MinGW compile command via GNU gfortran compiler
1#!/usr/bin/env sh
2rm main.exe
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
4./main.exe

Example output
1
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3! Compute the sqrt of the determinant of the positive definite matrix.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = getUnifRand(1, 5)
8ndim
9+5
10mat = getCovRand(mold = 1._TKG, ndim = ndim)
11mat
12+1.00000000, -0.923925519, +0.795801759, +0.556083620, -0.439362645
13-0.923925519, +1.00000000, -0.648698449, -0.243134469, +0.260473728
14+0.795801759, -0.648698449, +1.00000000, +0.672725320, -0.362898499
15+0.556083620, -0.243134469, +0.672725320, +1.00000000, -0.787845373
16-0.439362645, +0.260473728, -0.362898499, -0.787845373, +1.00000000
17detSqrtLog = getMatDetSqrtLog(mat)
18detSqrtLog
19-3.32033920
20getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
21-3.32033920
22detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
23+0.548362732E-5
24
25mat = getCovRand(mold = 1._TKG, ndim = ndim)
26mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
27mat
28+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
29+0.473446935, +1.00000000, +0.00000000, +0.00000000, +0.00000000
30-0.769212842, -0.875687480, +1.00000000, +0.00000000, +0.00000000
31-0.630960226, +0.210506231, +0.242536217E-1, +1.00000000, +0.00000000
32+0.357581615, -0.503823459E-1, -0.267758846, -0.109704062E-1, +1.00000000
33detSqrtLog = getMatDetSqrtLog(mat, lowDia)
34detSqrtLog
35-3.66848445
36getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
37-3.66848445
38detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
39+0.00000000
40
41mat = getCovRand(mold = 1._TKG, ndim = ndim)
42mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
43mat
44+1.00000000, -0.987090841E-1, -0.388906419, -0.744392455, -0.202922136
45+0.00000000, +1.00000000, -0.677588224, -0.128657043, +0.555987418
46+0.00000000, +0.00000000, +1.00000000, +0.530256510, -0.634381950
47+0.00000000, +0.00000000, +0.00000000, +1.00000000, +0.253579259
48+0.00000000, +0.00000000, +0.00000000, +0.00000000, +1.00000000
49detSqrtLog = getMatDetSqrtLog(mat, uppDia)
50detSqrtLog
51-2.18715143
52getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
53-2.18715143
54detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
55+0.00000000
56
57
58ndim = getUnifRand(1, 5)
59ndim
60+1
61mat = getCovRand(mold = 1._TKG, ndim = ndim)
62mat
63+1.00000000
64detSqrtLog = getMatDetSqrtLog(mat)
65detSqrtLog
66+0.00000000
67getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
68+0.00000000
69detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
70+0.00000000
71
72mat = getCovRand(mold = 1._TKG, ndim = ndim)
73mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
74mat
75+1.00000000
76detSqrtLog = getMatDetSqrtLog(mat, lowDia)
77detSqrtLog
78+0.00000000
79getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
80+0.00000000
81detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
82+0.00000000
83
84mat = getCovRand(mold = 1._TKG, ndim = ndim)
85mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
86mat
87+1.00000000
88detSqrtLog = getMatDetSqrtLog(mat, uppDia)
89detSqrtLog
90+0.00000000
91getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
92+0.00000000
93detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
94+0.00000000
95
96
97ndim = getUnifRand(1, 5)
98ndim
99+2
100mat = getCovRand(mold = 1._TKG, ndim = ndim)
101mat
102+1.00000000, +0.636380613
103+0.636380613, +1.00000000
104detSqrtLog = getMatDetSqrtLog(mat)
105detSqrtLog
106-0.259580404
107getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
108-0.259580404
109detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
110+0.00000000
111
112mat = getCovRand(mold = 1._TKG, ndim = ndim)
113mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
114mat
115+1.00000000, +0.00000000
116+0.708982050, +1.00000000
117detSqrtLog = getMatDetSqrtLog(mat, lowDia)
118detSqrtLog
119-0.349236190
120getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
121-0.349236190
122detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
123+0.00000000
124
125mat = getCovRand(mold = 1._TKG, ndim = ndim)
126mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
127mat
128+1.00000000, -0.155172020
129+0.00000000, +1.00000000
130detSqrtLog = getMatDetSqrtLog(mat, uppDia)
131detSqrtLog
132-0.121865133E-1
133getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
134-0.121865133E-1
135detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
136+0.00000000
137
138
139ndim = getUnifRand(1, 5)
140ndim
141+4
142mat = getCovRand(mold = 1._TKG, ndim = ndim)
143mat
144+1.00000000, -0.998882890, +0.340588838, -0.817403376
145-0.998882890, +1.00000000, -0.382320315, +0.831684768
146+0.340588838, -0.382320315, +1.00000000, -0.430089980
147-0.817403376, +0.831684768, -0.430089980, +1.00000000
148detSqrtLog = getMatDetSqrtLog(mat)
149detSqrtLog
150-6.08370113
151getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
152-6.08370113
153detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
154+0.305414200E-2
155
156mat = getCovRand(mold = 1._TKG, ndim = ndim)
157mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
158mat
159+1.00000000, +0.00000000, +0.00000000, +0.00000000
160-0.985799372, +1.00000000, +0.00000000, +0.00000000
161+0.833298504, -0.911155999, +1.00000000, +0.00000000
162+0.430329829, -0.340287954, +0.135605678, +1.00000000
163detSqrtLog = getMatDetSqrtLog(mat, lowDia)
164detSqrtLog
165-4.10915709
166getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
167-4.10915709
168detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
169+0.00000000
170
171mat = getCovRand(mold = 1._TKG, ndim = ndim)
172mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
173mat
174+1.00000000, -0.942418635, -0.931140065, +0.504103959
175+0.00000000, +1.00000000, +0.832991242, -0.273008823
176+0.00000000, +0.00000000, +1.00000000, -0.750759542
177+0.00000000, +0.00000000, +0.00000000, +1.00000000
178detSqrtLog = getMatDetSqrtLog(mat, uppDia)
179detSqrtLog
180-3.91975594
181getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
182-3.91975594
183detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
184+0.00000000
185
186
187ndim = getUnifRand(1, 5)
188ndim
189+3
190mat = getCovRand(mold = 1._TKG, ndim = ndim)
191mat
192+1.00000000, -0.749091804, -0.691371620
193-0.749091804, +1.00000000, +0.855242133
194-0.691371620, +0.855242133, +1.00000000
195detSqrtLog = getMatDetSqrtLog(mat)
196detSqrtLog
197-1.08015716
198getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
199-1.08015716
200detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
201+0.119209290E-6
202
203mat = getCovRand(mold = 1._TKG, ndim = ndim)
204mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
205mat
206+1.00000000, +0.00000000, +0.00000000
207+0.437291205, +1.00000000, +0.00000000
208-0.210297331, +0.480057687, +1.00000000
209detSqrtLog = getMatDetSqrtLog(mat, lowDia)
210detSqrtLog
211-0.403939456
212getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
213-0.403939456
214detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
215+0.00000000
216
217mat = getCovRand(mold = 1._TKG, ndim = ndim)
218mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
219mat
220+1.00000000, +0.550339460, +0.333082527
221+0.00000000, +1.00000000, +0.257115960
222+0.00000000, +0.00000000, +1.00000000
223detSqrtLog = getMatDetSqrtLog(mat, uppDia)
224detSqrtLog
225-0.243605807
226getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
227-0.243605807
228detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
229+0.00000000
230
231
232ndim = getUnifRand(1, 5)
233ndim
234+2
235mat = getCovRand(mold = 1._TKG, ndim = ndim)
236mat
237+1.00000000, -0.631532431
238-0.631532431, +1.00000000
239detSqrtLog = getMatDetSqrtLog(mat)
240detSqrtLog
241-0.254441470
242getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
243-0.254441470
244detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
245-0.298023224E-7
246
247mat = getCovRand(mold = 1._TKG, ndim = ndim)
248mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
249mat
250+1.00000000, +0.00000000
251+0.804985464, +1.00000000
252detSqrtLog = getMatDetSqrtLog(mat, lowDia)
253detSqrtLog
254-0.522064388
255getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
256-0.522064388
257detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
258+0.00000000
259
260mat = getCovRand(mold = 1._TKG, ndim = ndim)
261mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
262mat
263+1.00000000, +0.626405716
264+0.00000000, +1.00000000
265detSqrtLog = getMatDetSqrtLog(mat, uppDia)
266detSqrtLog
267-0.249106169
268getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
269-0.249106169
270detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
271+0.00000000
272
273
274ndim = getUnifRand(1, 5)
275ndim
276+1
277mat = getCovRand(mold = 1._TKG, ndim = ndim)
278mat
279+1.00000000
280detSqrtLog = getMatDetSqrtLog(mat)
281detSqrtLog
282+0.00000000
283getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
284+0.00000000
285detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
286+0.00000000
287
288mat = getCovRand(mold = 1._TKG, ndim = ndim)
289mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
290mat
291+1.00000000
292detSqrtLog = getMatDetSqrtLog(mat, lowDia)
293detSqrtLog
294+0.00000000
295getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
296+0.00000000
297detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
298+0.00000000
299
300mat = getCovRand(mold = 1._TKG, ndim = ndim)
301mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
302mat
303+1.00000000
304detSqrtLog = getMatDetSqrtLog(mat, uppDia)
305detSqrtLog
306+0.00000000
307getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
308+0.00000000
309detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
310+0.00000000
311
312
313ndim = getUnifRand(1, 5)
314ndim
315+4
316mat = getCovRand(mold = 1._TKG, ndim = ndim)
317mat
318+1.00000000, -0.605744481, -0.340187997, -0.519873261
319-0.605744481, +1.00000000, +0.512276530, +0.256492853
320-0.340187997, +0.512276530, +1.00000000, -0.208593518
321-0.519873261, +0.256492853, -0.208593518, +1.00000000
322detSqrtLog = getMatDetSqrtLog(mat)
323detSqrtLog
324-0.679720879
325getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
326-0.679720879
327detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
328+0.00000000
329
330mat = getCovRand(mold = 1._TKG, ndim = ndim)
331mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
332mat
333+1.00000000, +0.00000000, +0.00000000, +0.00000000
334-0.740314424, +1.00000000, +0.00000000, +0.00000000
335-0.775035262, +0.898291111, +1.00000000, +0.00000000
336-0.540926516, +0.859361172, +0.914804459, +1.00000000
337detSqrtLog = getMatDetSqrtLog(mat, lowDia)
338detSqrtLog
339-2.58510637
340getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
341-2.58510637
342detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
343+0.00000000
344
345mat = getCovRand(mold = 1._TKG, ndim = ndim)
346mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
347mat
348+1.00000000, -0.161007509, -0.247560740, -0.552706540
349+0.00000000, +1.00000000, +0.971962869, +0.697107553
350+0.00000000, +0.00000000, +1.00000000, +0.829490900
351+0.00000000, +0.00000000, +0.00000000, +1.00000000
352detSqrtLog = getMatDetSqrtLog(mat, uppDia)
353detSqrtLog
354-3.01644516
355getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
356-3.01644516
357detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
358+0.00000000
359
360
361ndim = getUnifRand(1, 5)
362ndim
363+1
364mat = getCovRand(mold = 1._TKG, ndim = ndim)
365mat
366+1.00000000
367detSqrtLog = getMatDetSqrtLog(mat)
368detSqrtLog
369+0.00000000
370getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
371+0.00000000
372detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
373+0.00000000
374
375mat = getCovRand(mold = 1._TKG, ndim = ndim)
376mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
377mat
378+1.00000000
379detSqrtLog = getMatDetSqrtLog(mat, lowDia)
380detSqrtLog
381+0.00000000
382getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
383+0.00000000
384detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
385+0.00000000
386
387mat = getCovRand(mold = 1._TKG, ndim = ndim)
388mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
389mat
390+1.00000000
391detSqrtLog = getMatDetSqrtLog(mat, uppDia)
392detSqrtLog
393+0.00000000
394getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
395+0.00000000
396detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
397+0.00000000
398
399
400ndim = getUnifRand(1, 5)
401ndim
402+3
403mat = getCovRand(mold = 1._TKG, ndim = ndim)
404mat
405+1.00000000, +0.801576301E-2, -0.388933152
406+0.801576301E-2, +1.00000000, +0.116096601
407-0.388933152, +0.116096601, +1.00000000
408detSqrtLog = getMatDetSqrtLog(mat)
409detSqrtLog
410-0.904825479E-1
411getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
412-0.904825479E-1
413detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
414+0.447034836E-7
415
416mat = getCovRand(mold = 1._TKG, ndim = ndim)
417mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
418mat
419+1.00000000, +0.00000000, +0.00000000
420-0.461657077, +1.00000000, +0.00000000
421+0.480568051, +0.289506048, +1.00000000
422detSqrtLog = getMatDetSqrtLog(mat, lowDia)
423detSqrtLog
424-0.534058452
425getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
426-0.534058452
427detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
428+0.00000000
429
430mat = getCovRand(mold = 1._TKG, ndim = ndim)
431mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
432mat
433+1.00000000, +0.709700584, +0.335972309
434+0.00000000, +1.00000000, -0.263459623
435+0.00000000, +0.00000000, +1.00000000
436detSqrtLog = getMatDetSqrtLog(mat, uppDia)
437detSqrtLog
438-0.834598660
439getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
440-0.834598660
441detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
442+0.00000000
443
444
445mat
446(+9.00000000, +0.00000000), (+3.00000000, +3.00000000), (+3.00000000, -3.00000000)
447(+3.00000000, -3.00000000), (+18.0000000, +0.00000000), (+8.00000000, -6.00000000)
448(+3.00000000, +3.00000000), (+8.00000000, +6.00000000), (+43.0000000, +0.00000000)
449detSqrtLog = getMatDetSqrtLog(mat)
450detSqrtLog
451+4.27666616
452getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
453(+4.27666616, +0.00000000)
454detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.
455-18.2898731
456
457tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
458tmp
459(+9.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
460(+3.00000000, -3.00000000), (+18.0000000, +0.00000000), (+0.00000000, +0.00000000)
461(+3.00000000, +3.00000000), (+8.00000000, +6.00000000), (+43.0000000, +0.00000000)
462detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
463detSqrtLog
464+4.27666616
465getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.
466(+4.27666616, +0.00000000)
467detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one
468(+0.00000000, +0.00000000)
469
470tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
471tmp
472(+9.00000000, +0.00000000), (+3.00000000, +3.00000000), (+3.00000000, -3.00000000)
473(+0.00000000, +0.00000000), (+18.0000000, +0.00000000), (+8.00000000, -6.00000000)
474(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+43.0000000, +0.00000000)
475detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
476detSqrtLog
477+4.27666616
478getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.
479(+4.27666616, +0.00000000)
480detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.
481(+0.00000000, +0.00000000)
482
483
484mat
485(+25.0000000, +0.00000000), (-5.00000000, -5.00000000), (+10.0000000, +5.00000000)
486(-5.00000000, +5.00000000), (+51.0000000, +0.00000000), (+4.00000000, -6.00000000)
487(+10.0000000, -5.00000000), (+4.00000000, +6.00000000), (+71.0000000, +0.00000000)
488detSqrtLog = getMatDetSqrtLog(mat)
489detSqrtLog
490+5.63478947
491getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
492(+5.63478947, +0.00000000)
493detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.
494-31.7508545
495
496tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
497tmp
498(+25.0000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
499(-5.00000000, +5.00000000), (+51.0000000, +0.00000000), (+0.00000000, +0.00000000)
500(+10.0000000, -5.00000000), (+4.00000000, +6.00000000), (+71.0000000, +0.00000000)
501detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
502detSqrtLog
503+5.63478947
504getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.
505(+5.63478947, +0.00000000)
506detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one.
507(+0.00000000, +0.00000000)
508
509tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
510tmp
511(+25.0000000, +0.00000000), (-5.00000000, -5.00000000), (+10.0000000, +5.00000000)
512(+0.00000000, +0.00000000), (+51.0000000, +0.00000000), (+4.00000000, -6.00000000)
513(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+71.0000000, +0.00000000)
514detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
515detSqrtLog
516+5.63478947
517getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.
518(+5.63478947, +0.00000000)
519detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.
520(+0.00000000, +0.00000000)
521
522
Test:
test_pm_matrixDet


Final Remarks


If you believe this algorithm or its documentation can be improved, we appreciate your contribution and help to edit this page's documentation and source file on GitHub.
For details on the naming abbreviations, see this page.
For details on the naming conventions, see this page.
This software is distributed under the MIT license with additional terms outlined below.

  1. If you use any parts or concepts from this library to any extent, please acknowledge the usage by citing the relevant publications of the ParaMonte library.
  2. If you regenerate any parts/ideas from this library in a programming environment other than those currently supported by this ParaMonte library (i.e., other than C, C++, Fortran, MATLAB, Python, R), please also ask the end users to cite this original ParaMonte library.

This software is available to the public under a highly permissive license.
Help us justify its continued development and maintenance by acknowledging its benefit to society, distributing it, and contributing to it.

Author:
Amir Shahmoradi, Apr 21, 2017, 1:43 PM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 600 of file pm_matrixDet.F90.


The documentation for this interface was generated from the following file: