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

Return the natural logarithm of the square-root of the determinant of the input positive-definite square matrix. More...

Detailed Description

Return the natural logarithm of 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.
Computing the logarithm of the square-root of the determinant is preferrable for high rank matrices where there is a possibility of overflow, although it comes with the additional computational cost of computing size(mat, 1) extra logarithms.

Note
This generic interface is meant to facilitate simultaneous computation of the Cholesky factorization and determinant of a input positive definite matrix.
If the Cholesky factorization is already computed, simply compute the determinant as the square of the multiplicative trace of the Cholesky factor.
Parameters
[in]mat: The input contiguous matrix of shape (1:ndim, 1:ndim) of
  1. type complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128),
  2. type real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing a subset of the positive-definite matrix whose log(sqrt(determinant)) is to be computed.
[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.
[out]detSqrtLog: The output scalar of type real of the same kind as the inputmatrepresenting the natural logarithm of the square root of its determinant.<br> \param[out] info : The output scalarintegerof default kind [IK](@ref pm_kind::IK) that is0<b>if and only if</b> the Cholesky factorization succeeds.<br> Otherwise, it is set to the order of the leading minor of the specified input subset ofmat` that is not positive definite, indicating the occurrence of an error and that the factorization could not be completed.
A non-zero value implies failure of the Cholesky factorization.
A non-zero info implied the input mat is not positive definite.
[in,out]chol: The input/output matrix of the same type, kind, and shape as the input mat containing the Cholesky factorization in its triangular subset dictated by the input arguments subset and operation.
See the documentation of the input argument mat above for more information.
[in]operation: The input scalar that can be either,
  1. the constant nothing implying that the same subset of the output chol as the input subset must be populated with the Cholesky factorization.
  2. the constant transHerm implying that the complex transpose of the computed Cholesky factorization must be written to the output chol.


Possible calling interfaces

call setMatDetSqrtLog(mat, subset, detSqrtLog, info, chol, operation)
Return the natural logarithm of the square-root of the determinant of the input positive-definite squ...
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, in which case, the procedure will return the control with info /= 0.
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.
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
5 use pm_matrixCopy, only: transHerm
10 use pm_matrixChol, only: uppDia, lowDia
11 use pm_matrixChol, only: getMatChol
12 use pm_distCov, only: getCovRand
13 use pm_arrayResize, only: setResized
14 use pm_distUnif, only: getUnifRand
15 use pm_matrixInv, only: getMatInv
16 use pm_io, only: display_type
17 use pm_io, only: getFormat
18
19 implicit none
20
21 integer(IK) :: info, ndim, itry, ntry = 10
22 character(:, SK), allocatable :: format
23 type(display_type) :: disp
24
25 disp = display_type(file = "main.out.F90")
26
27 call disp%skip()
28 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 call disp%show("! Compute the sqrt of the determinant of the positive definite matrix.")
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%skip()
32
33 block
34 use pm_kind, only: TKG => RKS
35 real(TKG), allocatable :: mat(:,:)
36 real(TKG) :: detSqrtLog
37 format = getFormat(mold = [0._TKG], ed = SK_"e")
38 do itry = 1, ntry
39 call disp%skip()
40 call disp%show("ndim = getUnifRand(1, 5)")
41 ndim = getUnifRand(1, 5)
42 call disp%show("ndim")
43 call disp%show( ndim )
44 call disp%show("call setResized(mat, [ndim, ndim + 1])")
45 call setResized(mat, [ndim, ndim + 1])
46 call disp%show("mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)")
47 mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
48 call disp%show("mat(:,1:ndim)")
49 call disp%show( mat(:,1:ndim) , format = format )
50 call disp%show("mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.")
51 mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG)
52 call disp%show("mat(:,1:ndim)")
53 call disp%show( mat(:,1:ndim) , format = format )
54 call disp%show("call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)")
55 call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
56 call disp%show("if (info /= 0) error stop 'Cholesky factorization failed.'")
57 if (info /= 0) error stop 'Cholesky factorization failed.'
58 call disp%show("detSqrtLog")
59 call disp%show( detSqrtLog )
60 call disp%show("mat")
61 call disp%show( mat , format = format )
62 call disp%show("getMatChol(mat(:,1:ndim), lowDia)")
63 call disp%show( getMatChol(mat(:,1:ndim), lowDia) , format = format )
64 call disp%show("getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.")
65 call disp%show( getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) )
66 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.")
67 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) )
68 call disp%skip()
69 call disp%show("mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)")
70 mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
71 call disp%show("mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.")
72 mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG)
73 call disp%show("mat(:,2:ndim+1)")
74 call disp%show( mat(:,2:ndim+1) , format = format )
75 call disp%show("call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)")
76 call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
77 call disp%show("if (info /= 0) error stop 'Cholesky factorization failed.'")
78 if (info /= 0) error stop 'Cholesky factorization failed.'
79 call disp%show("detSqrtLog")
80 call disp%show( detSqrtLog )
81 call disp%show("mat")
82 call disp%show( mat , format = format )
83 call disp%show("getMatChol(mat(:,2:ndim+1), uppDia)")
84 call disp%show( getMatChol(mat(:,2:ndim+1), uppDia) , format = format )
85 call disp%show("getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.")
86 call disp%show( getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) )
87 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.")
88 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) )
89 call disp%skip()
90 end do
91 end block
92
93 block
94 use pm_kind, only: TKG => CKS
95 integer(IK), parameter :: ndim = 3
96 complex(TKG), allocatable :: tmp(:,:)
97 complex(TKG), parameter :: mat(*,*) = reshape( [ (9.0, 0.0), (3.0, 3.0), (3.0, -3.0) &
98 , (3.0, -3.0),(18.0, 0.0), (8.0, -6.0) &
99 , (3.0, 3.0), (8.0, 6.0),(43.0, 0.0) &
100 ], shape = [ndim, ndim], order = [2, 1])
101 real(TKG) :: detSqrtLog
102 format = getFormat(mold = [(0._TKG, 0._TKG)], ed = SK_"e")
103 call disp%skip()
104 call disp%show("mat")
105 call disp%show( mat , format = format )
106 call disp%show("call setResized(tmp, [ndim, ndim + 1])")
107 call setResized(tmp, [ndim, ndim + 1])
108 call disp%show("tmp(:,1:ndim) = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
109 tmp(:,1:ndim) = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG))
110 call disp%show("tmp(:,1:ndim)")
111 call disp%show( tmp(:,1:ndim) , format = format )
112 call disp%show("call setMatDetSqrtLog(tmp(:,1:ndim), lowDia, detSqrtLog, info, tmp(:,2:ndim+1), transHerm)")
113 call setMatDetSqrtLog(tmp(:,1:ndim), lowDia, detSqrtLog, info, tmp(:,2:ndim+1), transHerm)
114 call disp%show("if (info /= 0) error stop 'Cholesky factorization failed.'")
115 if (info /= 0) error stop 'Cholesky factorization failed.'
116 call disp%show("detSqrtLog")
117 call disp%show( detSqrtLog )
118 call disp%show("tmp")
119 call disp%show( tmp , format = format )
120 call disp%show("getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! for comparison.")
121 call disp%show( getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) )
122 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! must be one")
123 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) )
124 call disp%skip()
125 call disp%show("tmp(:,2:ndim+1) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
126 tmp(:,2:ndim+1) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, uppDia, init = (0._TKG, 0._TKG))
127 call disp%show("tmp(:,2:ndim+1)")
128 call disp%show( tmp(:,2:ndim+1) , format = format )
129 call disp%show("call setMatDetSqrtLog(tmp(:,2:ndim+1), uppDia, detSqrtLog, info, tmp(:,1:ndim), transHerm)")
130 call setMatDetSqrtLog(tmp(:,2:ndim+1), uppDia, detSqrtLog, info, tmp(:,1:ndim), transHerm)
131 call disp%show("if (info /= 0) error stop 'Cholesky factorization failed.'")
132 if (info /= 0) error stop 'Cholesky factorization failed.'
133 call disp%show("detSqrtLog")
134 call disp%show( detSqrtLog )
135 call disp%show("tmp")
136 call disp%show( tmp , format = format )
137 call disp%show("getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! for comparison.")
138 call disp%show( getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) )
139 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! must be one.")
140 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) )
141 call disp%skip()
142 end block
143
144 block
145 use pm_kind, only: TKG => CKS
146 integer(IK), parameter :: ndim = 3
147 complex(TKG), allocatable :: tmp(:,:)
148 complex(TKG), parameter :: mat(*,*) = reshape( [ (25.0, 0.0), (-5.0, -5.0), (10.0, 5.0) &
149 , (-5.0, 5.0), (51.0, 0.0), (4.0, -6.0) &
150 , (10.0, -5.0), (4.0, 6.0), (71.0, 0.0) &
151 ], shape = [ndim, ndim], order = [2, 1])
152 real(TKG) :: detSqrtLog
153 format = getFormat(mold = [(0._TKG, 0._TKG)], ed = SK_"e")
154 call disp%skip()
155 call disp%show("mat")
156 call disp%show( mat , format = format )
157 call disp%show("call setResized(tmp, [ndim, ndim + 1])")
158 call setResized(tmp, [ndim, ndim + 1])
159 call disp%show("tmp(:,1:ndim) = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
160 tmp(:,1:ndim) = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG))
161 call disp%show("tmp(:,1:ndim)")
162 call disp%show( tmp(:,1:ndim) , format = format )
163 call disp%show("call setMatDetSqrtLog(tmp(:,1:ndim), lowDia, detSqrtLog, info, tmp(:,2:ndim+1), transHerm)")
164 call setMatDetSqrtLog(tmp(:,1:ndim), lowDia, detSqrtLog, info, tmp(:,2:ndim+1), transHerm)
165 call disp%show("if (info /= 0) error stop 'Cholesky factorization failed.'")
166 if (info /= 0) error stop 'Cholesky factorization failed.'
167 call disp%show("detSqrtLog")
168 call disp%show( detSqrtLog )
169 call disp%show("tmp")
170 call disp%show( tmp , format = format )
171 call disp%show("getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! for comparison.")
172 call disp%show( getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) )
173 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! must be one")
174 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) )
175 call disp%skip()
176 call disp%show("tmp(:,2:ndim+1) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.")
177 tmp(:,2:ndim+1) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, uppDia, init = (0._TKG, 0._TKG))
178 call disp%show("tmp(:,2:ndim+1)")
179 call disp%show( tmp(:,2:ndim+1) , format = format )
180 call disp%show("call setMatDetSqrtLog(tmp(:,2:ndim+1), uppDia, detSqrtLog, info, tmp(:,1:ndim), transHerm)")
181 call setMatDetSqrtLog(tmp(:,2:ndim+1), uppDia, detSqrtLog, info, tmp(:,1:ndim), transHerm)
182 call disp%show("if (info /= 0) error stop 'Cholesky factorization failed.'")
183 if (info /= 0) error stop 'Cholesky factorization failed.'
184 call disp%show("detSqrtLog")
185 call disp%show( detSqrtLog )
186 call disp%show("tmp")
187 call disp%show( tmp , format = format )
188 call disp%show("getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! for comparison.")
189 call disp%show( getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) )
190 call disp%show("detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! must be one.")
191 call disp%show( detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) )
192 call disp%skip()
193 end block
194
195end 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...
Generate and return a generic or type/kind-specific IO format with the requested specifications that ...
Definition: pm_io.F90:18485
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 natural logarithm of the square-root of the determinant of the input positive...
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
10call setResized(mat, [ndim, ndim + 1])
11mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
12mat(:,1:ndim)
13 0.100000E+01, -0.476027E+00
14-0.476027E+00, 0.100000E+01
15mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
16mat(:,1:ndim)
17 0.100000E+01, 0.000000E+00
18-0.476027E+00, 0.100000E+01
19call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
20if (info /= 0) error stop 'Cholesky factorization failed.'
21detSqrtLog
22-0.128480628
23mat
24 0.100000E+01, 0.100000E+01, -0.476027E+00
25-0.476027E+00, 0.100000E+01, 0.879431E+00
26getMatChol(mat(:,1:ndim), lowDia)
27 0.100000E+01, 0.000000E+00
28-0.476027E+00, 0.879431E+00
29getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
30-0.128480628
31detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
32+0.00000000
33
34mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
35mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
36mat(:,2:ndim+1)
37 0.100000E+01, 0.811936E+00
38 0.000000E+00, 0.100000E+01
39call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
40if (info /= 0) error stop 'Cholesky factorization failed.'
41detSqrtLog
42-0.538289726
43mat
44 0.100000E+01, 0.100000E+01, 0.811936E+00
45 0.811936E+00, 0.583746E+00, 0.100000E+01
46getMatChol(mat(:,2:ndim+1), uppDia)
47 0.100000E+01, 0.811936E+00
48 0.000000E+00, 0.583746E+00
49getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
50-0.538289726
51detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
52+0.00000000
53
54
55ndim = getUnifRand(1, 5)
56ndim
57+3
58call setResized(mat, [ndim, ndim + 1])
59mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
60mat(:,1:ndim)
61 0.100000E+01, -0.647685E-01, -0.286596E+00
62-0.647685E-01, 0.100000E+01, -0.317350E+00
63-0.286596E+00, -0.317350E+00, 0.100000E+01
64mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
65mat(:,1:ndim)
66 0.100000E+01, 0.000000E+00, 0.000000E+00
67-0.647685E-01, 0.100000E+01, 0.000000E+00
68-0.286596E+00, -0.317350E+00, 0.100000E+01
69call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
70if (info /= 0) error stop 'Cholesky factorization failed.'
71detSqrtLog
72-0.110837713
73mat
74 0.100000E+01, 0.100000E+01, -0.647685E-01, -0.286596E+00
75-0.647685E-01, 0.100000E+01, 0.997900E+00, -0.336619E+00
76-0.286596E+00, -0.317350E+00, 0.100000E+01, 0.896967E+00
77getMatChol(mat(:,1:ndim), lowDia)
78 0.100000E+01, 0.000000E+00, 0.000000E+00
79-0.647685E-01, 0.997900E+00, 0.000000E+00
80-0.286596E+00, -0.336619E+00, 0.896967E+00
81getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
82-0.110837713
83detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
84+0.00000000
85
86mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
87mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
88mat(:,2:ndim+1)
89 0.100000E+01, 0.278812E+00, 0.185030E+00
90 0.000000E+00, 0.100000E+01, -0.621793E+00
91 0.000000E+00, 0.000000E+00, 0.100000E+01
92call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
93if (info /= 0) error stop 'Cholesky factorization failed.'
94detSqrtLog
95-0.413629055
96mat
97 0.100000E+01, 0.100000E+01, 0.278812E+00, 0.185030E+00
98 0.278812E+00, 0.960346E+00, 0.100000E+01, -0.621793E+00
99 0.185030E+00, -0.701186E+00, 0.688550E+00, 0.100000E+01
100getMatChol(mat(:,2:ndim+1), uppDia)
101 0.100000E+01, 0.278812E+00, 0.185030E+00
102 0.000000E+00, 0.960346E+00, -0.701186E+00
103 0.000000E+00, 0.000000E+00, 0.688550E+00
104getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
105-0.413629055
106detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
107+0.00000000
108
109
110ndim = getUnifRand(1, 5)
111ndim
112+1
113call setResized(mat, [ndim, ndim + 1])
114mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
115mat(:,1:ndim)
116 0.100000E+01
117mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
118mat(:,1:ndim)
119 0.100000E+01
120call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
121if (info /= 0) error stop 'Cholesky factorization failed.'
122detSqrtLog
123+0.00000000
124mat
125 0.100000E+01, 0.100000E+01
126getMatChol(mat(:,1:ndim), lowDia)
127 0.100000E+01
128getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
129+0.00000000
130detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
131+0.00000000
132
133mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
134mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
135mat(:,2:ndim+1)
136 0.100000E+01
137call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
138if (info /= 0) error stop 'Cholesky factorization failed.'
139detSqrtLog
140+0.00000000
141mat
142 0.100000E+01, 0.100000E+01
143getMatChol(mat(:,2:ndim+1), uppDia)
144 0.100000E+01
145getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
146+0.00000000
147detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
148+0.00000000
149
150
151ndim = getUnifRand(1, 5)
152ndim
153+1
154call setResized(mat, [ndim, ndim + 1])
155mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
156mat(:,1:ndim)
157 0.100000E+01
158mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
159mat(:,1:ndim)
160 0.100000E+01
161call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
162if (info /= 0) error stop 'Cholesky factorization failed.'
163detSqrtLog
164+0.00000000
165mat
166 0.100000E+01, 0.100000E+01
167getMatChol(mat(:,1:ndim), lowDia)
168 0.100000E+01
169getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
170+0.00000000
171detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
172+0.00000000
173
174mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
175mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
176mat(:,2:ndim+1)
177 0.100000E+01
178call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
179if (info /= 0) error stop 'Cholesky factorization failed.'
180detSqrtLog
181+0.00000000
182mat
183 0.100000E+01, 0.100000E+01
184getMatChol(mat(:,2:ndim+1), uppDia)
185 0.100000E+01
186getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
187+0.00000000
188detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
189+0.00000000
190
191
192ndim = getUnifRand(1, 5)
193ndim
194+1
195call setResized(mat, [ndim, ndim + 1])
196mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
197mat(:,1:ndim)
198 0.100000E+01
199mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
200mat(:,1:ndim)
201 0.100000E+01
202call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
203if (info /= 0) error stop 'Cholesky factorization failed.'
204detSqrtLog
205+0.00000000
206mat
207 0.100000E+01, 0.100000E+01
208getMatChol(mat(:,1:ndim), lowDia)
209 0.100000E+01
210getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
211+0.00000000
212detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
213+0.00000000
214
215mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
216mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
217mat(:,2:ndim+1)
218 0.100000E+01
219call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
220if (info /= 0) error stop 'Cholesky factorization failed.'
221detSqrtLog
222+0.00000000
223mat
224 0.100000E+01, 0.100000E+01
225getMatChol(mat(:,2:ndim+1), uppDia)
226 0.100000E+01
227getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
228+0.00000000
229detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
230+0.00000000
231
232
233ndim = getUnifRand(1, 5)
234ndim
235+1
236call setResized(mat, [ndim, ndim + 1])
237mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
238mat(:,1:ndim)
239 0.100000E+01
240mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
241mat(:,1:ndim)
242 0.100000E+01
243call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
244if (info /= 0) error stop 'Cholesky factorization failed.'
245detSqrtLog
246+0.00000000
247mat
248 0.100000E+01, 0.100000E+01
249getMatChol(mat(:,1:ndim), lowDia)
250 0.100000E+01
251getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
252+0.00000000
253detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
254+0.00000000
255
256mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
257mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
258mat(:,2:ndim+1)
259 0.100000E+01
260call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
261if (info /= 0) error stop 'Cholesky factorization failed.'
262detSqrtLog
263+0.00000000
264mat
265 0.100000E+01, 0.100000E+01
266getMatChol(mat(:,2:ndim+1), uppDia)
267 0.100000E+01
268getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
269+0.00000000
270detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
271+0.00000000
272
273
274ndim = getUnifRand(1, 5)
275ndim
276+2
277call setResized(mat, [ndim, ndim + 1])
278mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
279mat(:,1:ndim)
280 0.100000E+01, 0.988230E+00
281 0.988230E+00, 0.100000E+01
282mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
283mat(:,1:ndim)
284 0.100000E+01, 0.000000E+00
285 0.988230E+00, 0.100000E+01
286call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
287if (info /= 0) error stop 'Cholesky factorization failed.'
288detSqrtLog
289-1.87749791
290mat
291 0.100000E+01, 0.100000E+01, 0.988230E+00
292 0.988230E+00, 0.100000E+01, 0.152972E+00
293getMatChol(mat(:,1:ndim), lowDia)
294 0.100000E+01, 0.000000E+00
295 0.988230E+00, 0.152972E+00
296getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
297-1.87749791
298detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
299+0.00000000
300
301mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
302mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
303mat(:,2:ndim+1)
304 0.100000E+01, 0.989782E+00
305 0.000000E+00, 0.100000E+01
306call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
307if (info /= 0) error stop 'Cholesky factorization failed.'
308detSqrtLog
309-1.94780946
310mat
311 0.100000E+01, 0.100000E+01, 0.989782E+00
312 0.989782E+00, 0.142586E+00, 0.100000E+01
313getMatChol(mat(:,2:ndim+1), uppDia)
314 0.100000E+01, 0.989782E+00
315 0.000000E+00, 0.142586E+00
316getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
317-1.94780946
318detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
319+0.00000000
320
321
322ndim = getUnifRand(1, 5)
323ndim
324+1
325call setResized(mat, [ndim, ndim + 1])
326mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
327mat(:,1:ndim)
328 0.100000E+01
329mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
330mat(:,1:ndim)
331 0.100000E+01
332call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
333if (info /= 0) error stop 'Cholesky factorization failed.'
334detSqrtLog
335+0.00000000
336mat
337 0.100000E+01, 0.100000E+01
338getMatChol(mat(:,1:ndim), lowDia)
339 0.100000E+01
340getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
341+0.00000000
342detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
343+0.00000000
344
345mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
346mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
347mat(:,2:ndim+1)
348 0.100000E+01
349call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
350if (info /= 0) error stop 'Cholesky factorization failed.'
351detSqrtLog
352+0.00000000
353mat
354 0.100000E+01, 0.100000E+01
355getMatChol(mat(:,2:ndim+1), uppDia)
356 0.100000E+01
357getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
358+0.00000000
359detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
360+0.00000000
361
362
363ndim = getUnifRand(1, 5)
364ndim
365+3
366call setResized(mat, [ndim, ndim + 1])
367mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
368mat(:,1:ndim)
369 0.100000E+01, 0.997602E+00, -0.160301E+00
370 0.997602E+00, 0.100000E+01, -0.202446E+00
371-0.160301E+00, -0.202446E+00, 0.100000E+01
372mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
373mat(:,1:ndim)
374 0.100000E+01, 0.000000E+00, 0.000000E+00
375 0.997602E+00, 0.100000E+01, 0.000000E+00
376-0.160301E+00, -0.202446E+00, 0.100000E+01
377call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
378if (info /= 0) error stop 'Cholesky factorization failed.'
379detSqrtLog
380-2.92877221
381mat
382 0.100000E+01, 0.100000E+01, 0.997602E+00, -0.160301E+00
383 0.997602E+00, 0.100000E+01, 0.692103E-01, -0.614491E+00
384-0.160301E+00, -0.202446E+00, 0.100000E+01, 0.772466E+00
385getMatChol(mat(:,1:ndim), lowDia)
386 0.100000E+01, 0.000000E+00, 0.000000E+00
387 0.997602E+00, 0.692103E-01, 0.000000E+00
388-0.160301E+00, -0.614491E+00, 0.772466E+00
389getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
390-2.92877221
391detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
392+0.00000000
393
394mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
395mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
396mat(:,2:ndim+1)
397 0.100000E+01, 0.434083E-01, 0.354631E+00
398 0.000000E+00, 0.100000E+01, 0.877353E+00
399 0.000000E+00, 0.000000E+00, 0.100000E+01
400call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
401if (info /= 0) error stop 'Cholesky factorization failed.'
402detSqrtLog
403-1.02159023
404mat
405 0.100000E+01, 0.100000E+01, 0.434083E-01, 0.354631E+00
406 0.434083E-01, 0.999057E+00, 0.100000E+01, 0.877353E+00
407 0.354631E+00, 0.862772E+00, 0.360362E+00, 0.100000E+01
408getMatChol(mat(:,2:ndim+1), uppDia)
409 0.100000E+01, 0.434083E-01, 0.354631E+00
410 0.000000E+00, 0.999057E+00, 0.862772E+00
411 0.000000E+00, 0.000000E+00, 0.360362E+00
412getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
413-1.02159023
414detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
415+0.00000000
416
417
418ndim = getUnifRand(1, 5)
419ndim
420+3
421call setResized(mat, [ndim, ndim + 1])
422mat(:,1:ndim) = getCovRand(mold = 1._TKG, ndim = ndim)
423mat(:,1:ndim)
424 0.100000E+01, -0.784509E+00, 0.397969E+00
425-0.784509E+00, 0.100000E+01, 0.142067E+00
426 0.397969E+00, 0.142067E+00, 0.100000E+01
427mat(:,1:ndim) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, lowDia, init = 0._TKG) ! reset the upper.
428mat(:,1:ndim)
429 0.100000E+01, 0.000000E+00, 0.000000E+00
430-0.784509E+00, 0.100000E+01, 0.000000E+00
431 0.397969E+00, 0.142067E+00, 0.100000E+01
432call setMatDetSqrtLog(mat(:,1:ndim), lowDia, detSqrtLog, info, mat(:,2:ndim+1), transHerm)
433if (info /= 0) error stop 'Cholesky factorization failed.'
434detSqrtLog
435-1.07162166
436mat
437 0.100000E+01, 0.100000E+01, -0.784509E+00, 0.397969E+00
438-0.784509E+00, 0.100000E+01, 0.620118E+00, 0.732567E+00
439 0.397969E+00, 0.142067E+00, 0.100000E+01, 0.552238E+00
440getMatChol(mat(:,1:ndim), lowDia)
441 0.100000E+01, 0.000000E+00, 0.000000E+00
442-0.784509E+00, 0.620118E+00, 0.000000E+00
443 0.397969E+00, 0.732567E+00, 0.552238E+00
444getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! for comparison.
445-1.07162166
446detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,1:ndim), lowDia)) ! must be one.
447+0.00000000
448
449mat(:,2:ndim+1) = getCovRand(mold = 1._TKG, ndim = ndim)
450mat(:,2:ndim+1) = getMatCopy(rdpack, mat(:,2:ndim+1), rdpack, uppDia, init = 0._TKG) ! reset the lower.
451mat(:,2:ndim+1)
452 0.100000E+01, 0.987103E+00, -0.774583E+00
453 0.000000E+00, 0.100000E+01, -0.838769E+00
454 0.000000E+00, 0.000000E+00, 0.100000E+01
455call setMatDetSqrtLog(mat(:,2:ndim+1), uppDia, detSqrtLog, info, mat(:,1:ndim), transHerm)
456if (info /= 0) error stop 'Cholesky factorization failed.'
457detSqrtLog
458-2.67486763
459mat
460 0.100000E+01, 0.100000E+01, 0.987103E+00, -0.774583E+00
461 0.987103E+00, 0.160085E+00, 0.100000E+01, -0.838769E+00
462-0.774583E+00, -0.463351E+00, 0.430495E+00, 0.100000E+01
463getMatChol(mat(:,2:ndim+1), uppDia)
464 0.100000E+01, 0.987103E+00, -0.774583E+00
465 0.000000E+00, 0.160085E+00, -0.463351E+00
466 0.000000E+00, 0.000000E+00, 0.430495E+00
467getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! for comparison.
468-2.67486763
469detSqrtLog - getMatMulTraceLog(getMatChol(mat(:,2:ndim+1), uppDia)) ! must be one.
470+0.00000000
471
472
473mat
474( 0.900000E+01, 0.000000E+00), ( 0.300000E+01, 0.300000E+01), ( 0.300000E+01, -0.300000E+01)
475( 0.300000E+01, -0.300000E+01), ( 0.180000E+02, 0.000000E+00), ( 0.800000E+01, -0.600000E+01)
476( 0.300000E+01, 0.300000E+01), ( 0.800000E+01, 0.600000E+01), ( 0.430000E+02, 0.000000E+00)
477call setResized(tmp, [ndim, ndim + 1])
478tmp(:,1:ndim) = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
479tmp(:,1:ndim)
480( 0.900000E+01, 0.000000E+00), ( 0.000000E+00, 0.000000E+00), ( 0.000000E+00, 0.000000E+00)
481( 0.300000E+01, -0.300000E+01), ( 0.180000E+02, 0.000000E+00), ( 0.000000E+00, 0.000000E+00)
482( 0.300000E+01, 0.300000E+01), ( 0.800000E+01, 0.600000E+01), ( 0.430000E+02, 0.000000E+00)
483call setMatDetSqrtLog(tmp(:,1:ndim), lowDia, detSqrtLog, info, tmp(:,2:ndim+1), transHerm)
484if (info /= 0) error stop 'Cholesky factorization failed.'
485detSqrtLog
486+4.27666616
487tmp
488( 0.900000E+01, 0.000000E+00), ( 0.300000E+01, 0.000000E+00), ( 0.100000E+01, 0.100000E+01), ( 0.100000E+01, -0.100000E+01)
489( 0.300000E+01, -0.300000E+01), ( 0.180000E+02, 0.000000E+00), ( 0.400000E+01, 0.000000E+00), ( 0.200000E+01, -0.100000E+01)
490( 0.300000E+01, 0.300000E+01), ( 0.800000E+01, 0.600000E+01), ( 0.430000E+02, 0.000000E+00), ( 0.600000E+01, 0.000000E+00)
491getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! for comparison.
492(+4.27666616, +0.00000000)
493detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! must be one
494(+0.00000000, +0.00000000)
495
496tmp(:,2:ndim+1) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
497tmp(:,2:ndim+1)
498( 0.900000E+01, 0.000000E+00), ( 0.300000E+01, 0.300000E+01), ( 0.300000E+01, -0.300000E+01)
499( 0.000000E+00, 0.000000E+00), ( 0.180000E+02, 0.000000E+00), ( 0.800000E+01, -0.600000E+01)
500( 0.000000E+00, 0.000000E+00), ( 0.000000E+00, 0.000000E+00), ( 0.430000E+02, 0.000000E+00)
501call setMatDetSqrtLog(tmp(:,2:ndim+1), uppDia, detSqrtLog, info, tmp(:,1:ndim), transHerm)
502if (info /= 0) error stop 'Cholesky factorization failed.'
503detSqrtLog
504+4.27666616
505tmp
506( 0.300000E+01, 0.000000E+00), ( 0.900000E+01, 0.000000E+00), ( 0.300000E+01, 0.300000E+01), ( 0.300000E+01, -0.300000E+01)
507( 0.100000E+01, -0.100000E+01), ( 0.400000E+01, 0.000000E+00), ( 0.180000E+02, 0.000000E+00), ( 0.800000E+01, -0.600000E+01)
508( 0.100000E+01, 0.100000E+01), ( 0.200000E+01, 0.100000E+01), ( 0.600000E+01, 0.000000E+00), ( 0.430000E+02, 0.000000E+00)
509getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! for comparison.
510(+4.27666616, +0.00000000)
511detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! must be one.
512(+0.00000000, +0.00000000)
513
514
515mat
516( 0.250000E+02, 0.000000E+00), (-0.500000E+01, -0.500000E+01), ( 0.100000E+02, 0.500000E+01)
517(-0.500000E+01, 0.500000E+01), ( 0.510000E+02, 0.000000E+00), ( 0.400000E+01, -0.600000E+01)
518( 0.100000E+02, -0.500000E+01), ( 0.400000E+01, 0.600000E+01), ( 0.710000E+02, 0.000000E+00)
519call setResized(tmp, [ndim, ndim + 1])
520tmp(:,1:ndim) = getMatCopy(rdpack, mat, rdpack, lowDia, init = (0._TKG, 0._TKG)) ! reset the upper.
521tmp(:,1:ndim)
522( 0.250000E+02, 0.000000E+00), ( 0.000000E+00, 0.000000E+00), ( 0.000000E+00, 0.000000E+00)
523(-0.500000E+01, 0.500000E+01), ( 0.510000E+02, 0.000000E+00), ( 0.000000E+00, 0.000000E+00)
524( 0.100000E+02, -0.500000E+01), ( 0.400000E+01, 0.600000E+01), ( 0.710000E+02, 0.000000E+00)
525call setMatDetSqrtLog(tmp(:,1:ndim), lowDia, detSqrtLog, info, tmp(:,2:ndim+1), transHerm)
526if (info /= 0) error stop 'Cholesky factorization failed.'
527detSqrtLog
528+5.63478947
529tmp
530( 0.250000E+02, 0.000000E+00), ( 0.500000E+01, 0.000000E+00), (-0.100000E+01, -0.100000E+01), ( 0.200000E+01, 0.100000E+01)
531(-0.500000E+01, 0.500000E+01), ( 0.510000E+02, 0.000000E+00), ( 0.700000E+01, 0.000000E+00), ( 0.100000E+01, -0.100000E+01)
532( 0.100000E+02, -0.500000E+01), ( 0.400000E+01, 0.600000E+01), ( 0.710000E+02, 0.000000E+00), ( 0.800000E+01, 0.000000E+00)
533getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! for comparison.
534(+5.63478947, +0.00000000)
535detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,1:ndim), lowDia)) ! must be one
536(+0.00000000, +0.00000000)
537
538tmp(:,2:ndim+1) = getMatCopy(rdpack, mat(:,1:ndim), rdpack, uppDia, init = (0._TKG, 0._TKG)) ! reset the upper.
539tmp(:,2:ndim+1)
540( 0.250000E+02, 0.000000E+00), (-0.500000E+01, -0.500000E+01), ( 0.100000E+02, 0.500000E+01)
541( 0.000000E+00, 0.000000E+00), ( 0.510000E+02, 0.000000E+00), ( 0.400000E+01, -0.600000E+01)
542( 0.000000E+00, 0.000000E+00), ( 0.000000E+00, 0.000000E+00), ( 0.710000E+02, 0.000000E+00)
543call setMatDetSqrtLog(tmp(:,2:ndim+1), uppDia, detSqrtLog, info, tmp(:,1:ndim), transHerm)
544if (info /= 0) error stop 'Cholesky factorization failed.'
545detSqrtLog
546+5.63478947
547tmp
548( 0.500000E+01, 0.000000E+00), ( 0.250000E+02, 0.000000E+00), (-0.500000E+01, -0.500000E+01), ( 0.100000E+02, 0.500000E+01)
549(-0.100000E+01, 0.100000E+01), ( 0.700000E+01, 0.000000E+00), ( 0.510000E+02, 0.000000E+00), ( 0.400000E+01, -0.600000E+01)
550( 0.200000E+01, -0.100000E+01), ( 0.100000E+01, 0.100000E+01), ( 0.800000E+01, 0.000000E+00), ( 0.710000E+02, 0.000000E+00)
551getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! for comparison.
552(+5.63478947, +0.00000000)
553detSqrtLog - getMatMulTraceLog(getMatChol(tmp(:,2:ndim+1), uppDia)) ! must be one.
554(+0.00000000, +0.00000000)
555
556
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 1747 of file pm_matrixDet.F90.


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