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+2
10mat = getCovRand(mold = 1._TKG, ndim = ndim)
11mat
12+1.00000000, -0.769710481
13-0.769710481, +1.00000000
14detSqrtLog = getMatDetSqrtLog(mat)
15detSqrtLog
16-0.448800951
17getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
18-0.448800951
19detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
20+0.00000000
21
22mat = getCovRand(mold = 1._TKG, ndim = ndim)
23mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
24mat
25+1.00000000, +0.00000000
26-0.805478334, +1.00000000
27detSqrtLog = getMatDetSqrtLog(mat, lowDia)
28detSqrtLog
29-0.523193121
30getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
31-0.523193121
32detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
33+0.00000000
34
35mat = getCovRand(mold = 1._TKG, ndim = ndim)
36mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
37mat
38+1.00000000, +0.780387580
39+0.00000000, +1.00000000
40detSqrtLog = getMatDetSqrtLog(mat, uppDia)
41detSqrtLog
42-0.469529957
43getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
44-0.469529957
45detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
46+0.00000000
47
48
49ndim = getUnifRand(1, 5)
50ndim
51+3
52mat = getCovRand(mold = 1._TKG, ndim = ndim)
53mat
54+1.00000000, -0.657553613, +0.756684393E-1
55-0.657553613, +1.00000000, -0.712460756
56+0.756684393E-1, -0.712460756, +1.00000000
57detSqrtLog = getMatDetSqrtLog(mat)
58detSqrtLog
59-1.03893852
60getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
61-1.03893852
62detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
63+0.357627869E-6
64
65mat = getCovRand(mold = 1._TKG, ndim = ndim)
66mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
67mat
68+1.00000000, +0.00000000, +0.00000000
69-0.170223251, +1.00000000, +0.00000000
70+0.764763653, -0.103227898, +1.00000000
71detSqrtLog = getMatDetSqrtLog(mat, lowDia)
72detSqrtLog
73-0.455177784
74getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
75-0.455177784
76detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
77+0.00000000
78
79mat = getCovRand(mold = 1._TKG, ndim = ndim)
80mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
81mat
82+1.00000000, +0.389872521, +0.513094306
83+0.00000000, +1.00000000, -0.162081793
84+0.00000000, +0.00000000, +1.00000000
85detSqrtLog = getMatDetSqrtLog(mat, uppDia)
86detSqrtLog
87-0.352997601
88getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
89-0.352997601
90detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
91+0.00000000
92
93
94ndim = getUnifRand(1, 5)
95ndim
96+5
97mat = getCovRand(mold = 1._TKG, ndim = ndim)
98mat
99+1.00000000, +0.169951934E-1, -0.652485669, -0.326669872, +0.265041351
100+0.169951934E-1, +1.00000000, +0.715217888, +0.583260298, +0.473803759
101-0.652485669, +0.715217888, +1.00000000, +0.546844840, +0.188181654
102-0.326669872, +0.583260298, +0.546844840, +1.00000000, -0.150400907
103+0.265041351, +0.473803759, +0.188181654, -0.150400907, +1.00000000
104detSqrtLog = getMatDetSqrtLog(mat)
105detSqrtLog
106-2.44936061
107getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
108-2.44936061
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, +0.00000000, +0.00000000, +0.00000000
116+0.389324874, +1.00000000, +0.00000000, +0.00000000, +0.00000000
117-0.945386231, -0.709943175E-1, +1.00000000, +0.00000000, +0.00000000
118-0.788472772, -0.583545208, +0.623233497, +1.00000000, +0.00000000
119-0.619876146, -0.204184517, +0.632815540, +0.953209102E-1, +1.00000000
120detSqrtLog = getMatDetSqrtLog(mat, lowDia)
121detSqrtLog
122-8.90280724
123getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
124-8.90280724
125detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
126+0.00000000
127
128mat = getCovRand(mold = 1._TKG, ndim = ndim)
129mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
130mat
131+1.00000000, -0.773320675, +0.491137594, +0.563879728, +0.479513854
132+0.00000000, +1.00000000, -0.565559864, -0.144912094, -0.666484237
133+0.00000000, +0.00000000, +1.00000000, +0.475127637, +0.667961717
134+0.00000000, +0.00000000, +0.00000000, +1.00000000, -0.203014016E-1
135+0.00000000, +0.00000000, +0.00000000, +0.00000000, +1.00000000
136detSqrtLog = getMatDetSqrtLog(mat, uppDia)
137detSqrtLog
138-1.91700339
139getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
140-1.91700339
141detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
142+0.00000000
143
144
145ndim = getUnifRand(1, 5)
146ndim
147+2
148mat = getCovRand(mold = 1._TKG, ndim = ndim)
149mat
150+1.00000000, -0.945053518
151-0.945053518, +1.00000000
152detSqrtLog = getMatDetSqrtLog(mat)
153detSqrtLog
154-1.11805296
155getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
156-1.11805296
157detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
158+0.00000000
159
160mat = getCovRand(mold = 1._TKG, ndim = ndim)
161mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
162mat
163+1.00000000, +0.00000000
164-0.679010868, +1.00000000
165detSqrtLog = getMatDetSqrtLog(mat, lowDia)
166detSqrtLog
167-0.309071571
168getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
169-0.309071571
170detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
171+0.00000000
172
173mat = getCovRand(mold = 1._TKG, ndim = ndim)
174mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
175mat
176+1.00000000, +0.335241526
177+0.00000000, +1.00000000
178detSqrtLog = getMatDetSqrtLog(mat, uppDia)
179detSqrtLog
180-0.596096367E-1
181getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
182-0.596096367E-1
183detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
184+0.00000000
185
186
187ndim = getUnifRand(1, 5)
188ndim
189+5
190mat = getCovRand(mold = 1._TKG, ndim = ndim)
191mat
192+1.00000000, +0.985954642, -0.155992419, -0.960156500, +0.104709864
193+0.985954642, +1.00000000, -0.844101608E-2, -0.906984329, +0.768529624E-1
194-0.155992419, -0.844101608E-2, +1.00000000, +0.310564876, -0.528498709
195-0.960156500, -0.906984329, +0.310564876, +1.00000000, -0.118123330
196+0.104709864, +0.768529624E-1, -0.528498709, -0.118123330, +1.00000000
197detSqrtLog = getMatDetSqrtLog(mat)
198detSqrtLog
199-6.55869579
200getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
201-6.55869579
202detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
203+0.203084946E-2
204
205mat = getCovRand(mold = 1._TKG, ndim = ndim)
206mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
207mat
208+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
209-0.999829948, +1.00000000, +0.00000000, +0.00000000, +0.00000000
210+0.324458808, -0.308008075, +1.00000000, +0.00000000, +0.00000000
211-0.250718981, +0.262643188, +0.391957313, +1.00000000, +0.00000000
212+0.499530770E-1, -0.393406563E-1, +0.265812188, +0.620156586, +1.00000000
213detSqrtLog = getMatDetSqrtLog(mat, lowDia)
214detSqrtLog
215-7.81679249
216getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
217-7.81679249
218detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
219+0.00000000
220
221mat = getCovRand(mold = 1._TKG, ndim = ndim)
222mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
223mat
224+1.00000000, +0.316341460, +0.436131597, -0.646341324, +0.582020283
225+0.00000000, +1.00000000, -0.441888034, -0.748386204, +0.689044893
226+0.00000000, +0.00000000, +1.00000000, +0.259849489, -0.289439380
227+0.00000000, +0.00000000, +0.00000000, +1.00000000, -0.727643311
228+0.00000000, +0.00000000, +0.00000000, +0.00000000, +1.00000000
229detSqrtLog = getMatDetSqrtLog(mat, uppDia)
230detSqrtLog
231-2.03846312
232getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
233-2.03846312
234detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
235+0.00000000
236
237
238ndim = getUnifRand(1, 5)
239ndim
240+3
241mat = getCovRand(mold = 1._TKG, ndim = ndim)
242mat
243+1.00000000, +0.759554446, -0.321252376
244+0.759554446, +1.00000000, -0.651270628
245-0.321252376, -0.651270628, +1.00000000
246detSqrtLog = getMatDetSqrtLog(mat)
247detSqrtLog
248-0.771937728
249getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
250-0.771937728
251detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
252+0.119209290E-6
253
254mat = getCovRand(mold = 1._TKG, ndim = ndim)
255mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
256mat
257+1.00000000, +0.00000000, +0.00000000
258-0.818485796, +1.00000000, +0.00000000
259+0.103240877, +0.126732647, +1.00000000
260detSqrtLog = getMatDetSqrtLog(mat, lowDia)
261detSqrtLog
262-0.633025169
263getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
264-0.633025169
265detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
266+0.00000000
267
268mat = getCovRand(mold = 1._TKG, ndim = ndim)
269mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
270mat
271+1.00000000, -0.339436799, -0.874218404
272+0.00000000, +1.00000000, +0.589697182
273+0.00000000, +0.00000000, +1.00000000
274detSqrtLog = getMatDetSqrtLog(mat, uppDia)
275detSqrtLog
276-1.04877090
277getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
278-1.04877090
279detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
280+0.00000000
281
282
283ndim = getUnifRand(1, 5)
284ndim
285+2
286mat = getCovRand(mold = 1._TKG, ndim = ndim)
287mat
288+1.00000000, -0.801043093
289-0.801043093, +1.00000000
290detSqrtLog = getMatDetSqrtLog(mat)
291detSqrtLog
292-0.513150573
293getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
294-0.513150573
295detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
296-0.238418579E-6
297
298mat = getCovRand(mold = 1._TKG, ndim = ndim)
299mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
300mat
301+1.00000000, +0.00000000
302+0.840460896, +1.00000000
303detSqrtLog = getMatDetSqrtLog(mat, lowDia)
304detSqrtLog
305-0.612725079
306getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
307-0.612725079
308detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
309+0.00000000
310
311mat = getCovRand(mold = 1._TKG, ndim = ndim)
312mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
313mat
314+1.00000000, +0.999919057
315+0.00000000, +1.00000000
316detSqrtLog = getMatDetSqrtLog(mat, uppDia)
317detSqrtLog
318-4.36430836
319getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
320-4.36430836
321detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
322+0.00000000
323
324
325ndim = getUnifRand(1, 5)
326ndim
327+5
328mat = getCovRand(mold = 1._TKG, ndim = ndim)
329mat
330+1.00000000, -0.887362540, -0.839923441, +0.942252576E-1, -0.455317855
331-0.887362540, +1.00000000, +0.992714047, +0.151003286, +0.887530744E-1
332-0.839923441, +0.992714047, +1.00000000, +0.251991391, -0.142256729E-2
333+0.942252576E-1, +0.151003286, +0.251991391, +1.00000000, -0.569754124
334-0.455317855, +0.887530744E-1, -0.142256729E-2, -0.569754124, +1.00000000
335detSqrtLog = getMatDetSqrtLog(mat)
336detSqrtLog
337-4.67533445
338getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
339-4.67533445
340detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
341-0.886917114E-4
342
343mat = getCovRand(mold = 1._TKG, ndim = ndim)
344mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
345mat
346+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
347+0.709386289, +1.00000000, +0.00000000, +0.00000000, +0.00000000
348-0.984122574, -0.815454483, +1.00000000, +0.00000000, +0.00000000
349-0.198942199, +0.341859162E-1, +0.195321798, +1.00000000, +0.00000000
350-0.421742588, -0.762103319, +0.549212575, -0.593156368E-1, +1.00000000
351detSqrtLog = getMatDetSqrtLog(mat, lowDia)
352detSqrtLog
353-4.71597624
354getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
355-4.71597624
356detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
357+0.00000000
358
359mat = getCovRand(mold = 1._TKG, ndim = ndim)
360mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
361mat
362+1.00000000, -0.777927637, +0.596161306, +0.455557078, +0.522729397
363+0.00000000, +1.00000000, -0.734426975E-1, -0.745029569, -0.257589757
364+0.00000000, +0.00000000, +1.00000000, -0.207445621, +0.488719970
365+0.00000000, +0.00000000, +0.00000000, +1.00000000, -0.225611702
366+0.00000000, +0.00000000, +0.00000000, +0.00000000, +1.00000000
367detSqrtLog = getMatDetSqrtLog(mat, uppDia)
368detSqrtLog
369-2.07221103
370getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
371-2.07221103
372detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
373+0.00000000
374
375
376ndim = getUnifRand(1, 5)
377ndim
378+1
379mat = getCovRand(mold = 1._TKG, ndim = ndim)
380mat
381+1.00000000
382detSqrtLog = getMatDetSqrtLog(mat)
383detSqrtLog
384+0.00000000
385getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
386+0.00000000
387detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
388+0.00000000
389
390mat = getCovRand(mold = 1._TKG, ndim = ndim)
391mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
392mat
393+1.00000000
394detSqrtLog = getMatDetSqrtLog(mat, lowDia)
395detSqrtLog
396+0.00000000
397getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
398+0.00000000
399detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
400+0.00000000
401
402mat = getCovRand(mold = 1._TKG, ndim = ndim)
403mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
404mat
405+1.00000000
406detSqrtLog = getMatDetSqrtLog(mat, uppDia)
407detSqrtLog
408+0.00000000
409getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
410+0.00000000
411detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
412+0.00000000
413
414
415ndim = getUnifRand(1, 5)
416ndim
417+5
418mat = getCovRand(mold = 1._TKG, ndim = ndim)
419mat
420+1.00000000, +0.888132334, -0.480939627, -0.556586027, -0.458889484
421+0.888132334, +1.00000000, -0.220319897, -0.509577692, -0.179750189
422-0.480939627, -0.220319897, +1.00000000, -0.246100575, +0.113878071
423-0.556586027, -0.509577692, -0.246100575, +1.00000000, +0.816105723
424-0.458889484, -0.179750189, +0.113878071, +0.816105723, +1.00000000
425detSqrtLog = getMatDetSqrtLog(mat)
426detSqrtLog
427-3.64937711
428getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
429-3.64937711
430detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
431-0.643730164E-5
432
433mat = getCovRand(mold = 1._TKG, ndim = ndim)
434mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
435mat
436+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
437+0.597684562, +1.00000000, +0.00000000, +0.00000000, +0.00000000
438+0.621274710, +0.810606062, +1.00000000, +0.00000000, +0.00000000
439+0.691059828E-1, +0.668436289, +0.816653013, +1.00000000, +0.00000000
440+0.509898245, +0.649358869, +0.861838222, +0.682172537, +1.00000000
441detSqrtLog = getMatDetSqrtLog(mat, lowDia)
442detSqrtLog
443-5.65679836
444getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
445-5.65679836
446detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
447+0.00000000
448
449mat = getCovRand(mold = 1._TKG, ndim = ndim)
450mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
451mat
452+1.00000000, -0.329446524, +0.820383579E-1, +0.232297480, -0.694908261
453+0.00000000, +1.00000000, -0.552640676, -0.634252608, -0.110254869
454+0.00000000, +0.00000000, +1.00000000, +0.297738194, +0.580430105E-1
455+0.00000000, +0.00000000, +0.00000000, +1.00000000, +0.514770687
456+0.00000000, +0.00000000, +0.00000000, +0.00000000, +1.00000000
457detSqrtLog = getMatDetSqrtLog(mat, uppDia)
458detSqrtLog
459-2.36555791
460getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
461-2.36555791
462detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
463+0.00000000
464
465
466mat
467(+9.00000000, +0.00000000), (+3.00000000, +3.00000000), (+3.00000000, -3.00000000)
468(+3.00000000, -3.00000000), (+18.0000000, +0.00000000), (+8.00000000, -6.00000000)
469(+3.00000000, +3.00000000), (+8.00000000, +6.00000000), (+43.0000000, +0.00000000)
470detSqrtLog = getMatDetSqrtLog(mat)
471detSqrtLog
472+4.27666616
473getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
474(+4.27666616, +0.00000000)
475detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.
476-18.2898731
477
478tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
479tmp
480(+9.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
481(+3.00000000, -3.00000000), (+18.0000000, +0.00000000), (+0.00000000, +0.00000000)
482(+3.00000000, +3.00000000), (+8.00000000, +6.00000000), (+43.0000000, +0.00000000)
483detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
484detSqrtLog
485+4.27666616
486getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.
487(+4.27666616, +0.00000000)
488detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one
489(+0.00000000, +0.00000000)
490
491tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
492tmp
493(+9.00000000, +0.00000000), (+3.00000000, +3.00000000), (+3.00000000, -3.00000000)
494(+0.00000000, +0.00000000), (+18.0000000, +0.00000000), (+8.00000000, -6.00000000)
495(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+43.0000000, +0.00000000)
496detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
497detSqrtLog
498+4.27666616
499getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.
500(+4.27666616, +0.00000000)
501detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.
502(+0.00000000, +0.00000000)
503
504
505mat
506(+25.0000000, +0.00000000), (-5.00000000, -5.00000000), (+10.0000000, +5.00000000)
507(-5.00000000, +5.00000000), (+51.0000000, +0.00000000), (+4.00000000, -6.00000000)
508(+10.0000000, -5.00000000), (+4.00000000, +6.00000000), (+71.0000000, +0.00000000)
509detSqrtLog = getMatDetSqrtLog(mat)
510detSqrtLog
511+5.63478947
512getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
513(+5.63478947, +0.00000000)
514detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.
515-31.7508545
516
517tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
518tmp
519(+25.0000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
520(-5.00000000, +5.00000000), (+51.0000000, +0.00000000), (+0.00000000, +0.00000000)
521(+10.0000000, -5.00000000), (+4.00000000, +6.00000000), (+71.0000000, +0.00000000)
522detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
523detSqrtLog
524+5.63478947
525getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.
526(+5.63478947, +0.00000000)
527detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one.
528(+0.00000000, +0.00000000)
529
530tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
531tmp
532(+25.0000000, +0.00000000), (-5.00000000, -5.00000000), (+10.0000000, +5.00000000)
533(+0.00000000, +0.00000000), (+51.0000000, +0.00000000), (+4.00000000, -6.00000000)
534(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+71.0000000, +0.00000000)
535detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
536detSqrtLog
537+5.63478947
538getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.
539(+5.63478947, +0.00000000)
540detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.
541(+0.00000000, +0.00000000)
542
543
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: