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+4
10mat = getCovRand(mold = 1._TKG, ndim = ndim)
11mat
12+1.00000000, +0.194989108E-1, -0.992550135, +0.636948466
13+0.194989108E-1, +1.00000000, -0.968626887E-1, +0.226403967
14-0.992550135, -0.968626887E-1, +1.00000000, -0.610434651
15+0.636948466, +0.226403967, -0.610434651, +1.00000000
16detSqrtLog = getMatDetSqrtLog(mat)
17detSqrtLog
18-2.84602022
19getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
20-2.84602022
21detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
22+0.357627869E-5
23
24mat = getCovRand(mold = 1._TKG, ndim = ndim)
25mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
26mat
27+1.00000000, +0.00000000, +0.00000000, +0.00000000
28+0.942800283, +1.00000000, +0.00000000, +0.00000000
29-0.929292202, -0.935344517, +1.00000000, +0.00000000
30+0.581445456, +0.564269543, -0.331789404, +1.00000000
31detSqrtLog = getMatDetSqrtLog(mat, lowDia)
32detSqrtLog
33-3.00610924
34getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
35-3.00610924
36detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
37+0.00000000
38
39mat = getCovRand(mold = 1._TKG, ndim = ndim)
40mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
41mat
42+1.00000000, +0.655721366, -0.663232148, +0.250736684
43+0.00000000, +1.00000000, +0.106879979, -0.486804247
44+0.00000000, +0.00000000, +1.00000000, -0.814963639
45+0.00000000, +0.00000000, +0.00000000, +1.00000000
46detSqrtLog = getMatDetSqrtLog(mat, uppDia)
47detSqrtLog
48-2.70481777
49getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
50-2.70481777
51detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
52+0.00000000
53
54
55ndim = getUnifRand(1, 5)
56ndim
57+3
58mat = getCovRand(mold = 1._TKG, ndim = ndim)
59mat
60+1.00000000, +0.605357707, +0.885512412
61+0.605357707, +1.00000000, +0.250792056
62+0.885512412, +0.250792056, +1.00000000
63detSqrtLog = getMatDetSqrtLog(mat)
64detSqrtLog
65-1.44669425
66getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
67-1.44669425
68detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
69+0.143051147E-5
70
71mat = getCovRand(mold = 1._TKG, ndim = ndim)
72mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
73mat
74+1.00000000, +0.00000000, +0.00000000
75+0.319143176, +1.00000000, +0.00000000
76-0.372258425, -0.799992144, +1.00000000
77detSqrtLog = getMatDetSqrtLog(mat, lowDia)
78detSqrtLog
79-0.586127341
80getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
81-0.586127341
82detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
83+0.00000000
84
85mat = getCovRand(mold = 1._TKG, ndim = ndim)
86mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
87mat
88+1.00000000, +0.378937036, +0.765523672
89+0.00000000, +1.00000000, +0.791987896
90+0.00000000, +0.00000000, +1.00000000
91detSqrtLog = getMatDetSqrtLog(mat, uppDia)
92detSqrtLog
93-1.13834429
94getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
95-1.13834429
96detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
97+0.00000000
98
99
100ndim = getUnifRand(1, 5)
101ndim
102+1
103mat = getCovRand(mold = 1._TKG, ndim = ndim)
104mat
105+1.00000000
106detSqrtLog = getMatDetSqrtLog(mat)
107detSqrtLog
108+0.00000000
109getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
110+0.00000000
111detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
112+0.00000000
113
114mat = getCovRand(mold = 1._TKG, ndim = ndim)
115mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
116mat
117+1.00000000
118detSqrtLog = getMatDetSqrtLog(mat, lowDia)
119detSqrtLog
120+0.00000000
121getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
122+0.00000000
123detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
124+0.00000000
125
126mat = getCovRand(mold = 1._TKG, ndim = ndim)
127mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
128mat
129+1.00000000
130detSqrtLog = getMatDetSqrtLog(mat, uppDia)
131detSqrtLog
132+0.00000000
133getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
134+0.00000000
135detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
136+0.00000000
137
138
139ndim = getUnifRand(1, 5)
140ndim
141+1
142mat = getCovRand(mold = 1._TKG, ndim = ndim)
143mat
144+1.00000000
145detSqrtLog = getMatDetSqrtLog(mat)
146detSqrtLog
147+0.00000000
148getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
149+0.00000000
150detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
151+0.00000000
152
153mat = getCovRand(mold = 1._TKG, ndim = ndim)
154mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
155mat
156+1.00000000
157detSqrtLog = getMatDetSqrtLog(mat, lowDia)
158detSqrtLog
159+0.00000000
160getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
161+0.00000000
162detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
163+0.00000000
164
165mat = getCovRand(mold = 1._TKG, ndim = ndim)
166mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
167mat
168+1.00000000
169detSqrtLog = getMatDetSqrtLog(mat, uppDia)
170detSqrtLog
171+0.00000000
172getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
173+0.00000000
174detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
175+0.00000000
176
177
178ndim = getUnifRand(1, 5)
179ndim
180+5
181mat = getCovRand(mold = 1._TKG, ndim = ndim)
182mat
183+1.00000000, +0.829913080, +0.426425874, -0.261845648, +0.348446034E-1
184+0.829913080, +1.00000000, -0.102190375E-1, -0.181813389, +0.399501026
185+0.426425874, -0.102190375E-1, +1.00000000, -0.739061296, -0.177195236
186-0.261845648, -0.181813389, -0.739061296, +1.00000000, -0.318746984
187+0.348446034E-1, +0.399501026, -0.177195236, -0.318746984, +1.00000000
188detSqrtLog = getMatDetSqrtLog(mat)
189detSqrtLog
190-2.97719312
191getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
192-2.97719312
193detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
194+0.810623169E-5
195
196mat = getCovRand(mold = 1._TKG, ndim = ndim)
197mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
198mat
199+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
200-0.649118006, +1.00000000, +0.00000000, +0.00000000, +0.00000000
201+0.391001225, +0.360414147, +1.00000000, +0.00000000, +0.00000000
202-0.649768710, +0.904200375, +0.889530480E-1, +1.00000000, +0.00000000
203-0.348328501, -0.168172196, -0.471172243, -0.263199210, +1.00000000
204detSqrtLog = getMatDetSqrtLog(mat, lowDia)
205detSqrtLog
206-3.45151424
207getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
208-3.45151424
209detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
210+0.00000000
211
212mat = getCovRand(mold = 1._TKG, ndim = ndim)
213mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
214mat
215+1.00000000, -0.941285908, -0.410704136, -0.235568985, +0.383620411
216+0.00000000, +1.00000000, +0.593273878, +0.115311585, -0.374017388
217+0.00000000, +0.00000000, +1.00000000, +0.179179847, +0.145443484
218+0.00000000, +0.00000000, +0.00000000, +1.00000000, +0.599070311
219+0.00000000, +0.00000000, +0.00000000, +0.00000000, +1.00000000
220detSqrtLog = getMatDetSqrtLog(mat, uppDia)
221detSqrtLog
222-2.31154060
223getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
224-2.31154060
225detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
226+0.00000000
227
228
229ndim = getUnifRand(1, 5)
230ndim
231+4
232mat = getCovRand(mold = 1._TKG, ndim = ndim)
233mat
234+1.00000000, -0.240141615, +0.644364297, +0.111268386
235-0.240141615, +1.00000000, +0.795950741E-1, +0.126492307
236+0.644364297, +0.795950741E-1, +1.00000000, -0.394532293
237+0.111268386, +0.126492307, -0.394532293, +1.00000000
238detSqrtLog = getMatDetSqrtLog(mat)
239detSqrtLog
240-0.717794120
241getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
242-0.717794120
243detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
244-0.178813934E-6
245
246mat = getCovRand(mold = 1._TKG, ndim = ndim)
247mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
248mat
249+1.00000000, +0.00000000, +0.00000000, +0.00000000
250+0.527687669, +1.00000000, +0.00000000, +0.00000000
251+0.167883843, +0.885345101, +1.00000000, +0.00000000
252+0.296970874, -0.388576269, -0.438351154, +1.00000000
253detSqrtLog = getMatDetSqrtLog(mat, lowDia)
254detSqrtLog
255-1.86873221
256getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
257-1.86873221
258detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
259+0.00000000
260
261mat = getCovRand(mold = 1._TKG, ndim = ndim)
262mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
263mat
264+1.00000000, +0.599177122, -0.455994099, -0.279145539
265+0.00000000, +1.00000000, +0.115580082, -0.766944706
266+0.00000000, +0.00000000, +1.00000000, -0.680036366
267+0.00000000, +0.00000000, +0.00000000, +1.00000000
268detSqrtLog = getMatDetSqrtLog(mat, uppDia)
269detSqrtLog
270-2.97581863
271getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
272-2.97581863
273detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
274+0.00000000
275
276
277ndim = getUnifRand(1, 5)
278ndim
279+1
280mat = getCovRand(mold = 1._TKG, ndim = ndim)
281mat
282+1.00000000
283detSqrtLog = getMatDetSqrtLog(mat)
284detSqrtLog
285+0.00000000
286getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
287+0.00000000
288detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
289+0.00000000
290
291mat = getCovRand(mold = 1._TKG, ndim = ndim)
292mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
293mat
294+1.00000000
295detSqrtLog = getMatDetSqrtLog(mat, lowDia)
296detSqrtLog
297+0.00000000
298getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
299+0.00000000
300detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
301+0.00000000
302
303mat = getCovRand(mold = 1._TKG, ndim = ndim)
304mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
305mat
306+1.00000000
307detSqrtLog = getMatDetSqrtLog(mat, uppDia)
308detSqrtLog
309+0.00000000
310getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
311+0.00000000
312detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
313+0.00000000
314
315
316ndim = getUnifRand(1, 5)
317ndim
318+2
319mat = getCovRand(mold = 1._TKG, ndim = ndim)
320mat
321+1.00000000, +0.702977657
322+0.702977657, +1.00000000
323detSqrtLog = getMatDetSqrtLog(mat)
324detSqrtLog
325-0.340784848
326getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
327-0.340784848
328detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
329-0.596046448E-7
330
331mat = getCovRand(mold = 1._TKG, ndim = ndim)
332mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
333mat
334+1.00000000, +0.00000000
335-0.963681042, +1.00000000
336detSqrtLog = getMatDetSqrtLog(mat, lowDia)
337detSqrtLog
338-1.32029712
339getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
340-1.32029712
341detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
342+0.00000000
343
344mat = getCovRand(mold = 1._TKG, ndim = ndim)
345mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
346mat
347+1.00000000, +0.876026988
348+0.00000000, +1.00000000
349detSqrtLog = getMatDetSqrtLog(mat, uppDia)
350detSqrtLog
351-0.729267538
352getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
353-0.729267538
354detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
355+0.00000000
356
357
358ndim = getUnifRand(1, 5)
359ndim
360+2
361mat = getCovRand(mold = 1._TKG, ndim = ndim)
362mat
363+1.00000000, -0.606280863
364-0.606280863, +1.00000000
365detSqrtLog = getMatDetSqrtLog(mat)
366detSqrtLog
367-0.229098007
368getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
369-0.229098007
370detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
371-0.894069672E-7
372
373mat = getCovRand(mold = 1._TKG, ndim = ndim)
374mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
375mat
376+1.00000000, +0.00000000
377-0.613612771, +1.00000000
378detSqrtLog = getMatDetSqrtLog(mat, lowDia)
379detSqrtLog
380-0.236219794
381getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
382-0.236219794
383detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
384+0.00000000
385
386mat = getCovRand(mold = 1._TKG, ndim = ndim)
387mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
388mat
389+1.00000000, +0.799726546
390+0.00000000, +1.00000000
391detSqrtLog = getMatDetSqrtLog(mat, uppDia)
392detSqrtLog
393-0.510218382
394getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
395-0.510218382
396detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
397+0.00000000
398
399
400ndim = getUnifRand(1, 5)
401ndim
402+2
403mat = getCovRand(mold = 1._TKG, ndim = ndim)
404mat
405+1.00000000, -0.973335564
406-0.973335564, +1.00000000
407detSqrtLog = getMatDetSqrtLog(mat)
408detSqrtLog
409-1.47234941
410getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
411-1.47234941
412detSqrtLog + getMatDetSqrtLog(getMatInv(mat)) ! must be one.
413-0.953674316E-6
414
415mat = getCovRand(mold = 1._TKG, ndim = ndim)
416mat = getMatCopy(rdpack, mat, rdpack, lowDia, init = 0._TKG) ! reset the upper.
417mat
418+1.00000000, +0.00000000
419+0.666825354, +1.00000000
420detSqrtLog = getMatDetSqrtLog(mat, lowDia)
421detSqrtLog
422-0.294083804
423getMatMulTraceLog(getMatChol(mat, lowDia)) ! for comparison.
424-0.294083804
425detSqrtLog - getMatMulTraceLog(getMatChol(mat, lowDia)) ! must be one.
426+0.00000000
427
428mat = getCovRand(mold = 1._TKG, ndim = ndim)
429mat = getMatCopy(rdpack, mat, rdpack, uppDia, init = 0._TKG) ! reset the lower.
430mat
431+1.00000000, +0.834811479E-1
432+0.00000000, +1.00000000
433detSqrtLog = getMatDetSqrtLog(mat, uppDia)
434detSqrtLog
435-0.349673326E-2
436getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
437-0.349673326E-2
438detSqrtLog - getMatMulTraceLog(getMatChol(mat, uppDia)) ! must be one.
439+0.00000000
440
441
442mat
443(+9.00000000, +0.00000000), (+3.00000000, +3.00000000), (+3.00000000, -3.00000000)
444(+3.00000000, -3.00000000), (+18.0000000, +0.00000000), (+8.00000000, -6.00000000)
445(+3.00000000, +3.00000000), (+8.00000000, +6.00000000), (+43.0000000, +0.00000000)
446detSqrtLog = getMatDetSqrtLog(mat)
447detSqrtLog
448+4.27666616
449getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
450(+4.27666616, +0.00000000)
451detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.
452-18.2898731
453
454tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
455tmp
456(+9.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
457(+3.00000000, -3.00000000), (+18.0000000, +0.00000000), (+0.00000000, +0.00000000)
458(+3.00000000, +3.00000000), (+8.00000000, +6.00000000), (+43.0000000, +0.00000000)
459detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
460detSqrtLog
461+4.27666616
462getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.
463(+4.27666616, +0.00000000)
464detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one
465(+0.00000000, +0.00000000)
466
467tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
468tmp
469(+9.00000000, +0.00000000), (+3.00000000, +3.00000000), (+3.00000000, -3.00000000)
470(+0.00000000, +0.00000000), (+18.0000000, +0.00000000), (+8.00000000, -6.00000000)
471(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+43.0000000, +0.00000000)
472detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
473detSqrtLog
474+4.27666616
475getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.
476(+4.27666616, +0.00000000)
477detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.
478(+0.00000000, +0.00000000)
479
480
481mat
482(+25.0000000, +0.00000000), (-5.00000000, -5.00000000), (+10.0000000, +5.00000000)
483(-5.00000000, +5.00000000), (+51.0000000, +0.00000000), (+4.00000000, -6.00000000)
484(+10.0000000, -5.00000000), (+4.00000000, +6.00000000), (+71.0000000, +0.00000000)
485detSqrtLog = getMatDetSqrtLog(mat)
486detSqrtLog
487+5.63478947
488getMatMulTraceLog(getMatChol(mat, uppDia)) ! for comparison.
489(+5.63478947, +0.00000000)
490detSqrtLog * getMatDetSqrtLog(getMatInv(mat)) ! must be one.
491-31.7508545
492
493tmp = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
494tmp
495(+25.0000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
496(-5.00000000, +5.00000000), (+51.0000000, +0.00000000), (+0.00000000, +0.00000000)
497(+10.0000000, -5.00000000), (+4.00000000, +6.00000000), (+71.0000000, +0.00000000)
498detSqrtLog = getMatDetSqrtLog(tmp, subset = lowDia)
499detSqrtLog
500+5.63478947
501getMatMulTraceLog(getMatChol(tmp, lowDia)) ! for comparison.
502(+5.63478947, +0.00000000)
503detSqrtLog - getMatMulTraceLog(getMatChol(tmp, lowDia)) ! must be one.
504(+0.00000000, +0.00000000)
505
506tmp = getMatCopy(rdpack, mat, rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
507tmp
508(+25.0000000, +0.00000000), (-5.00000000, -5.00000000), (+10.0000000, +5.00000000)
509(+0.00000000, +0.00000000), (+51.0000000, +0.00000000), (+4.00000000, -6.00000000)
510(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+71.0000000, +0.00000000)
511detSqrtLog = getMatDetSqrtLog(tmp, subset = uppDia)
512detSqrtLog
513+5.63478947
514getMatMulTraceLog(getMatChol(tmp, uppDia)) ! for comparison.
515(+5.63478947, +0.00000000)
516detSqrtLog - getMatMulTraceLog(getMatChol(tmp, uppDia)) ! must be one.
517(+0.00000000, +0.00000000)
518
519
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 at Austin

Definition at line 600 of file pm_matrixDet.F90.


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