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

Return a random positive-definite power-law-distributed (correlation) matrix.
More...

Detailed Description

Return a random positive-definite power-law-distributed (correlation) matrix.

See the documentation of pm_distCov for details.

Parameters
[in,out]rng: The input/output scalar that can be an object of,
  1. type rngf_type, implying the use of intrinsic Fortran uniform RNG.
  2. type xoshiro256ssw_type, implying the use of xoshiro256** uniform RNG.
[out]rand: The output 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 random (optionally power-law-distributed determinant) positive-definite matrix.
The output rand can of complex type if and only if the optional input argument method is missing.
[in,out]method: The input/output scalar constant that can be one of the following:
  1. The scalar input constant dvine implying the use of the Dvine algorithm for generating random covariance matrices whose determinants are power-law distributed with exponent eta.
    In this case, the argument method has intent(in).
  2. A scalar output variable of type onion_type such as onion implying the use of the Onion algorithm for generating random covariance matrices whose determinants are power-law distributed with exponent eta.
    In this case, the argument method has intent(out).
    If the Cholesky factorization within the Onion algorithm fails, methodinfo will be set to the order of the leading minor of the specified input subset of mat that is not positive definite, indicating the occurrence of an error and that the factorization could not be completed.
    Otherwise, the info component of the onion method is set to 0.
The resulting matrix distribution from dvine and onion are identically distributed but onion method tends to have slightly faster runtime.
The larger eta is, the more the output random matrix looks like the Identity matrix.
Setting eta = 0. corresponds to a uniform distribution of the output matrix over the space of positive-definite correlation matrices.
See the description of the output argument rand for more information on the effects of eta on the off-diagonal elements of the output positive-definite matrix.
(optional. If missing the Gram method is used for random matrix generation. It must be missing for output rand of type complex.)
[in]eta: The input non-negative scalar of type real of the same kind as the output argument rand.
The larger eta is, the more the output random matrix looks like the Identity matrix.
Setting eta = 0. corresponds to a uniform distribution of the output matrix over the space of positive-definite correlation matrices.
See the description of the output argument rand for more information on the effects of eta on the off-diagonal elements of the output positive-definite matrix.
(optional. It must be present if and only if the input argument method is also present.)
[in]scale: The input scalar or contiguous vector of size ndim of type real of the same kind as the output argument rand, representing the scale of the matrix (e.g., the standard deviation of a covariance matrix) along each dimension.
(optional. default = 1.)


Possible calling interfaces

! Default (Gram) method.
call setCovRand(rng, rand(1:ndim, 1:ndim))
call setCovRand(rng, rand(1:ndim, 1:ndim), scale)
call setCovRand(rng, rand(1:ndim, 1:ndim), scale(1:ndim))
! Other methods.
call setCovRand(rng, rand(1:ndim, 1:ndim), method, eta)
call setCovRand(rng, rand(1:ndim, 1:ndim), method, eta, scale)
call setCovRand(rng, rand(1:ndim, 1:ndim), method, eta, scale(1:ndim))
Return a random positive-definite power-law-distributed (correlation) matrix.
Definition: pm_distCov.F90:787
This module contains classes and procedures for generating random matrices distributed on the space o...
Definition: pm_distCov.F90:72
Warning
The condition 0 <= eta must hold for the corresponding input arguments.
The condition all([0 < scale]) must hold for the corresponding input arguments.
The condition size(rand, 1) == size(rand, 2) must hold for the corresponding input arguments.
The condition rank(scale) == 0 .or. all(size(scale) == shape(rand)) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Beware that when the input argument scale is missing, the diagonal elements of the output correlation matrix are not enforced to match 1.
As such, numerical matrix multiplication errors may lead to diagonal matrix values slightly deviating from 1. If you need such a guarantee on the diagonal elements of the output random correlation matrix, use getCovRand.


Example usage

1program example
2
3 use pm_kind, only: SK
4 use pm_kind, only: IK, LK
5 use pm_io, only: field_type
6 use pm_io, only: display_type
7 use pm_matrixChol, only: setChoLow
10 use pm_distCov, only: setCovRand, dvine, onion
12 use pm_arrayResize, only: setResized
13 use pm_matrixDet, only: getMatDet
14
15 implicit none
16
17 integer(IK) :: itry, ndim
18
19 type(display_type) :: disp
20 disp = display_type(file = SK_"main.out.F90", format = field_type(complex = SK_"math"))
21
22 call disp%skip()
23 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
24 call disp%show("!Gram method for real covariance.")
25 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
26 call disp%skip()
27
28 block
29
30 use pm_kind, only: TKG => RKS ! all kinds are supported.
31 real(TKG), allocatable :: scale(:)
32 real(TKG), allocatable :: rand(:,:)
33 do itry = 1, 5
34
35 call disp%skip()
36 call disp%show("ndim = getUnifRand(3, 9)")
37 ndim = getUnifRand(3, 9)
38 call disp%show("ndim")
39 call disp%show( ndim )
40 call disp%show("call setResized(rand, [ndim, ndim])")
41 call setResized(rand, [ndim, ndim])
42 call disp%show("scale = getUnifRand(1, 10, ndim)")
43 scale = getUnifRand(1, 10, ndim)
44 call disp%show("scale")
45 call disp%show( scale )
46 call disp%skip()
47
48 call disp%show("call setCovRand(rngf, rand)")
49 call setCovRand(rngf, rand)
50 call disp%show("rand")
51 call disp%show( rand )
52 call disp%show("isMatClass(rand, posdefmat)")
53 call disp%show( isMatClass(rand, posdefmat) )
54 call disp%show("isMatClass(rand, hermitian)")
55 call disp%show( isMatClass(rand, hermitian) )
56 call disp%show("getMatDet(rand)")
57 call disp%show( getMatDet(rand) )
58
59 call disp%show("call setCovRand(rngf, rand, scale(1))")
60 call setCovRand(rngf, rand, scale(1))
61 call disp%show("rand")
62 call disp%show( rand )
63 call disp%show("isMatClass(rand, posdefmat)")
64 call disp%show( isMatClass(rand, posdefmat) )
65 call disp%show("isMatClass(rand, hermitian)")
66 call disp%show( isMatClass(rand, hermitian) )
67 call disp%show("getMatDet(rand)")
68 call disp%show( getMatDet(rand) )
69 call disp%skip()
70
71 call disp%show("call setCovRand(rngf, rand, scale)")
72 call setCovRand(rngf, rand, scale)
73 call disp%show("rand")
74 call disp%show( rand )
75 call disp%show("isMatClass(rand, posdefmat)")
76 call disp%show( isMatClass(rand, posdefmat) )
77 call disp%show("isMatClass(rand, hermitian)")
78 call disp%show( isMatClass(rand, hermitian) )
79 call disp%show("getMatDet(rand)")
80 call disp%show( getMatDet(rand) )
81 call disp%skip()
82
83 end do
84
85 end block
86
87 call disp%skip()
88 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
89 call disp%show("!Gram method for complex covariance.")
90 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
91 call disp%skip()
92
93 block
94
95 use pm_kind, only: TKG => RKD ! all kinds are supported.
96 real(TKG), allocatable :: scale(:)
97 complex(TKG), allocatable :: rand(:,:)
98 do itry = 1, 5
99
100 call disp%skip()
101 call disp%show("ndim = getUnifRand(3, 9)")
102 ndim = getUnifRand(3, 9)
103 call disp%show("ndim")
104 call disp%show( ndim )
105 call disp%show("call setResized(rand, [ndim, ndim])")
106 call setResized(rand, [ndim, ndim])
107 call disp%show("scale = getUnifRand(1, 10, ndim)")
108 scale = getUnifRand(1, 10, ndim)
109 call disp%show("scale")
110 call disp%show( scale )
111 call disp%skip()
112
113 call disp%show("call setCovRand(rngf, rand)")
114 call setCovRand(rngf, rand)
115 call disp%show("rand")
116 call disp%show( rand )
117 call disp%show("isMatClass(rand, posdefmat)")
118 call disp%show( isMatClass(rand, posdefmat) )
119 call disp%show("isMatClass(rand, hermitian)")
120 call disp%show( isMatClass(rand, hermitian) )
121 call disp%show("getMatDet(rand)")
122 call disp%show( getMatDet(rand) )
123
124 call disp%show("call setCovRand(rngf, rand, scale(1))")
125 call setCovRand(rngf, rand, scale(1))
126 call disp%show("rand")
127 call disp%show( rand )
128 call disp%show("isMatClass(rand, posdefmat)")
129 call disp%show( isMatClass(rand, posdefmat) )
130 call disp%show("isMatClass(rand, hermitian)")
131 call disp%show( isMatClass(rand, hermitian) )
132 call disp%show("getMatDet(rand)")
133 call disp%show( getMatDet(rand) )
134 call disp%skip()
135
136 call disp%show("call setCovRand(rngf, rand, scale)")
137 call setCovRand(rngf, rand, scale)
138 call disp%show("rand")
139 call disp%show( rand )
140 call disp%show("isMatClass(rand, posdefmat)")
141 call disp%show( isMatClass(rand, posdefmat) )
142 call disp%show("isMatClass(rand, hermitian)")
143 call disp%show( isMatClass(rand, hermitian) )
144 call disp%show("getMatDet(rand)")
145 call disp%show( getMatDet(rand) )
146 call disp%skip()
147
148 end do
149
150 end block
151
152 call disp%skip()
153 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
154 call disp%show("!Dvine and Onion methods.")
155 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%")
156 call disp%skip()
157
158 block
159
160 use pm_kind, only: TKG => RKS ! all kinds are supported.
161 real(TKG), allocatable :: rand(:,:)
162 real(TKG) :: eta, scale
163 do itry = 1, 10
164
165 call disp%skip()
166 call disp%show("eta = getUnifRand(1, 10)")
167 eta = getUnifRand(1, 10)
168 call disp%show("ndim = getUnifRand(2, 5)")
169 ndim = getUnifRand(2, 5)
170 call disp%show("call setResized(rand, [ndim, ndim])")
171 call setResized(rand, [ndim, ndim])
172 call disp%show("scale = getUnifRand(1, 10)")
173 scale = getUnifRand(1, 10)
174 call disp%skip()
175
176 call disp%show("call setCovRand(rngf, rand, dvine, eta)")
177 call setCovRand(rngf, rand, dvine, eta)
178 call disp%show("onion%info")
179 call disp%show( onion%info )
180 call disp%show("rand")
181 call disp%show( rand )
182 call disp%show("isMatClass(rand, posdefmat)")
183 call disp%show( isMatClass(rand, posdefmat) )
184 call disp%show("isMatClass(rand, hermitian)")
185 call disp%show( isMatClass(rand, hermitian) )
186 call disp%show("getMatDet(rand)")
187 call disp%show( getMatDet(rand) )
188
189 call disp%show("call setCovRand(rngf, rand, onion, eta)")
190 call setCovRand(rngf, rand, onion, eta)
191 call disp%show("onion%info")
192 call disp%show( onion%info )
193 call disp%show("rand")
194 call disp%show( rand )
195 call disp%show("isMatClass(rand, posdefmat)")
196 call disp%show( isMatClass(rand, posdefmat) )
197 call disp%show("isMatClass(rand, hermitian)")
198 call disp%show( isMatClass(rand, hermitian) )
199 call disp%show("getMatDet(rand)")
200 call disp%show( getMatDet(rand) )
201
202 call disp%show("call setCovRand(rngf, rand, dvine, eta, scale)")
203 call setCovRand(rngf, rand, dvine, eta, scale)
204 call disp%show("rand")
205 call disp%show( rand )
206 call disp%show("isMatClass(rand, posdefmat)")
207 call disp%show( isMatClass(rand, posdefmat) )
208 call disp%show("isMatClass(rand, hermitian)")
209 call disp%show( isMatClass(rand, hermitian) )
210 call disp%show("getMatDet(rand)")
211 call disp%show( getMatDet(rand) )
212 call disp%skip()
213
214 call disp%show("call setCovRand(rngf, rand, onion, eta, scale)")
215 call setCovRand(rngf, rand, onion, eta, scale)
216 call disp%show("onion%info")
217 call disp%show( onion%info )
218 call disp%show("rand")
219 call disp%show( rand )
220 call disp%show("isMatClass(rand, posdefmat)")
221 call disp%show( isMatClass(rand, posdefmat) )
222 call disp%show("isMatClass(rand, hermitian)")
223 call disp%show( isMatClass(rand, hermitian) )
224 call disp%show("getMatDet(rand)")
225 call disp%show( getMatDet(rand) )
226 call disp%skip()
227
228 end do
229
230 call disp%skip()
231 call disp%show("ndim = getUnifRand(2, 10)")
232 ndim = getUnifRand(2, 10)
233 call disp%show("call setResized(rand, [ndim, ndim])")
234 call setResized(rand, [ndim, ndim])
235 call disp%skip()
236
237 call disp%show("call setCovRand(rngf, rand, dvine, eta = 0._TKG, scale = [(real(itry, TKG), itry = 1, ndim)])")
238 call setCovRand(rngf, rand, dvine, eta = 0._TKG, scale = [(real(itry, TKG), itry = 1, ndim)])
239 call disp%show("rand")
240 call disp%show( rand )
241 call disp%show("isMatClass(rand, posdefmat)")
242 call disp%show( isMatClass(rand, posdefmat) )
243 call disp%skip()
244
245 call disp%show("call setCovRand(rngf, rand, onion, eta = 0._TKG, scale = [(real(itry, TKG), itry = 1, ndim)])")
246 call setCovRand(rngf, rand, onion, eta = 0._TKG, scale = [(real(itry, TKG), itry = 1, ndim)])
247 call disp%show("onion%info")
248 call disp%show( onion%info )
249 call disp%show("rand")
250 call disp%show( rand )
251 call disp%show("isMatClass(rand, posdefmat)")
252 call disp%show( isMatClass(rand, posdefmat) )
253 call disp%skip()
254
255 end block
256
257end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate and return a (collection) of random vector(s) of size ndim from the ndim-dimensional MultiVa...
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
[LEGACY code] Return the lower-triangle of the Cholesky factorization of the symmetric positive-def...
Generate and return .true. if and only if the input matrix is of the specified input class.
Generate and return the determinant of the input general square matrix.
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
type(onion_type) onion
The scalar module variable object of type onion_type implying the use of the Onion algorithm for gene...
Definition: pm_distCov.F90:331
type(dvine_type), parameter dvine
The scalar constant of type dvine_type implying the use of the Dvine algorithm for generating random ...
Definition: pm_distCov.F90:238
This module contains classes and procedures for computing various statistical quantities related to t...
This module contains classes and procedures for computing various statistical quantities related to t...
type(rngf_type) rngf
The scalar constant object of type rngf_type whose presence signified the use of the Fortran intrinsi...
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 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 RKD
The double precision real kind in Fortran mode. On most platforms, this is an 64-bit real kind.
Definition: pm_kind.F90:568
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 abstract and concrete derived types that are required for compile-time resolutio...
type(posdefmat_type), parameter posdefmat
This is a scalar parameter object of type hermitian_type that is exclusively used to signify the Herm...
type(hermitian_type), parameter hermitian
This is a scalar parameter object of type hermitian_type that is exclusively used to signify the Herm...
This module contains procedures and generic interfaces relevant to the computation of the determinant...
Generate and return an object of type display_type.
Definition: pm_io.F90:10282
The derived type that can be used for constructing containers of format or left and right delimiters ...
Definition: pm_io.F90:482

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!Gram method for real covariance.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = getUnifRand(3, 9)
8ndim
9+9
10call setResized(rand, [ndim, ndim])
11scale = getUnifRand(1, 10, ndim)
12scale
13+8.00000000, +3.00000000, +4.00000000, +10.0000000, +6.00000000, +4.00000000, +10.0000000, +10.0000000, +3.00000000
14
15call setCovRand(rngf, rand)
16rand
17+1.00000000, +0.957517326, +0.568198740, +0.542673826, +0.489191800, -0.356931627, -0.962954666E-2, +0.654745340, -0.282974422
18+0.957517326, +1.00000012, +0.755919039, +0.337433100, +0.334102571, -0.458285511, -0.631121844E-1, +0.538206458, -0.418572485
19+0.568198740, +0.755919039, +1.00000012, -0.280721724, -0.250233740, -0.308606625, -0.270803452, +0.208552212, -0.540053070
20+0.542673826, +0.337433100, -0.280721724, +0.999999881, +0.564207017, -0.326139718, +0.437263787, +0.371268362, -0.449737161E-1
21+0.489191800, +0.334102571, -0.250233740, +0.564207017, +1.00000000, -0.233650595, -0.141660780, +0.511791646, +0.288811713
22-0.356931627, -0.458285511, -0.308606625, -0.326139718, -0.233650595, +1.00000000, -0.117863551, +0.504920855E-1, +0.496963412
23-0.962954666E-2, -0.631121844E-1, -0.270803452, +0.437263787, -0.141660780, -0.117863551, +1.00000000, +0.739401132E-1, +0.155677482
24+0.654745340, +0.538206458, +0.208552212, +0.371268362, +0.511791646, +0.504920855E-1, +0.739401132E-1, +1.00000000, +0.362090945
25-0.282974422, -0.418572485, -0.540053070, -0.449737161E-1, +0.288811713, +0.496963412, +0.155677482, +0.362090945, +1.00000024
27T
29T
30getMatDet(rand)
31+0.128518263E-7
32call setCovRand(rngf, rand, scale(1))
33rand
34+64.0000000, +50.4022102, +32.4713364, -21.6564598, +40.4578247, +27.3008213, -5.24287510, -15.2494507, +15.4966755
35+50.4022102, +64.0000076, +4.51029778, +5.01695251, +6.92782593, +29.0213833, -22.3060226, -17.9950905, +15.9546604
36+32.4713364, +4.51029778, +64.0000000, -4.66325188, +42.9756355, +19.6075058, +35.9222946, -10.8078575, +4.54169369
37-21.6564598, +5.01695251, -4.66325188, +63.9999962, -27.7601547, +27.1944046, -3.22864532, -22.4624348, +12.6389408
38+40.4578247, +6.92782593, +42.9756355, -27.7601547, +63.9999924, +32.8835220, -0.750291824, -15.6556320, +21.1342888
39+27.3008213, +29.0213833, +19.6075058, +27.1944046, +32.8835220, +63.9999886, -24.8746185, -32.8423042, +37.0633850
40-5.24287510, -22.3060226, +35.9222946, -3.22864532, -0.750291824, -24.8746185, +63.9999886, +16.4446545, -28.3457317
41-15.2494507, -17.9950905, -10.8078575, -22.4624348, -15.6556320, -32.8423042, +16.4446545, +64.0000000, -0.239372253E-1
42+15.4966755, +15.9546604, +4.54169369, +12.6389408, +21.1342888, +37.0633850, -28.3457317, -0.239372253E-1, +64.0000000
44T
46T
47getMatDet(rand)
48+0.855332506E+10
49
50call setCovRand(rngf, rand, scale)
51rand
52+64.0000000, +5.92940044, +14.2702913, +67.8983994, +6.26989841, -9.10739803, -30.5183334, +42.8082962, -4.16588497
53+5.92940044, +8.99999809, -1.28502917, +3.38465500, -1.56195664, -7.32579088, -5.83927631, -10.7693233, -0.629105568
54+14.2702913, -1.28502917, +15.9999981, +2.35771537, -12.0460091, -3.79404116, -11.1901169, +17.9625244, +3.52691841
55+67.8983994, +3.38465500, +2.35771537, +99.9999924, +10.3877478, +6.49343443, -22.6729660, +50.6163979, -12.6810360
56+6.26989841, -1.56195664, -12.0460091, +10.3877478, +36.0000000, -5.13571167, -13.6847868, -7.86860085, -2.44178915
57-9.10739803, -7.32579088, -3.79404116, +6.49343443, -5.13571167, +16.0000019, +9.04471588, +8.30318069, -4.32467318
58-30.5183334, -5.83927631, -11.1901169, -22.6729660, -13.6847868, +9.04471588, +100.000000, +28.7007694, -3.43619967
59+42.8082962, -10.7693233, +17.9625244, +50.6163979, -7.86860085, +8.30318069, +28.7007694, +99.9999847, -10.3450241
60-4.16588497, -0.629105568, +3.52691841, -12.6810360, -2.44178915, -4.32467318, -3.43619967, -10.3450241, +9.00000000
62T
64T
65getMatDet(rand)
66+1882902.62
67
68
69ndim = getUnifRand(3, 9)
70ndim
71+7
72call setResized(rand, [ndim, ndim])
73scale = getUnifRand(1, 10, ndim)
74scale
75+5.00000000, +8.00000000, +5.00000000, +2.00000000, +5.00000000, +1.00000000, +3.00000000
76
77call setCovRand(rngf, rand)
78rand
79+1.00000000, -0.968285501, -0.649314642, -0.912330821E-1, -0.155782253, -0.288999975, +0.498902321
80-0.968285501, +1.00000000, +0.440165997, +0.235607222, +0.324412793, +0.168063670, -0.375256956
81-0.649314642, +0.440165997, +0.999999881, -0.451365829, -0.476461202, +0.544231296, -0.678153098
82-0.912330821E-1, +0.235607222, -0.451365829, +0.999999881, +0.682478368, -0.156975895, +0.272120059
83-0.155782253, +0.324412793, -0.476461202, +0.682478368, +1.00000012, -0.472894877, +0.637292922
84-0.288999975, +0.168063670, +0.544231296, -0.156975895, -0.472894877, +1.00000000, -0.464530885
85+0.498902321, -0.375256956, -0.678153098, +0.272120059, +0.637292922, -0.464530885, +1.00000000
87T
89T
90getMatDet(rand)
91+0.588229021E-8
92call setCovRand(rngf, rand, scale(1))
93rand
94+25.0000000, -19.9442482, +24.4878025, -16.2620773, -17.8176136, -12.7186184, +12.6411219
95-19.9442482, +25.0000000, -19.2004433, +23.3937302, +21.6040516, +7.32900810, -1.02460456
96+24.4878025, -19.2004433, +24.9999981, -16.1335182, -18.4336433, -12.1916275, +11.3499441
97-16.2620773, +23.3937302, -16.1335182, +25.0000000, +23.6711597, +7.10887766, +2.22400022
98-17.8176136, +21.6040516, -18.4336433, +23.6711597, +25.0000000, +6.44416952, +0.361280888
99-12.7186184, +7.32900810, -12.1916275, +7.10887766, +6.44416952, +25.0000038, -19.8210487
100+12.6411219, -1.02460456, +11.3499441, +2.22400022, +0.361280888, -19.8210487, +24.9999943
102T
104T
105getMatDet(rand)
106+1528.03247
107
108call setCovRand(rngf, rand, scale)
109rand
110+25.0000000, +4.43194056, +7.39055586, -2.01091743, -14.0758848, +1.15525091, +2.65164709
111+4.43194056, +63.9999847, -35.8767853, +6.31929636, +19.9159508, +5.34445381, +14.0770884
112+7.39055586, -35.8767853, +25.0000000, -5.95658064, -20.0584373, -3.00715637, -6.84581470
113-2.01091743, +6.31929636, -5.95658064, +3.99999976, +8.52338696, +1.40324390, -0.838441014
114-14.0758848, +19.9159508, -20.0584373, +8.52338696, +25.0000000, +2.24430227, +1.21620166
115+1.15525091, +5.34445381, -3.00715637, +1.40324390, +2.24430227, +1.00000012, +0.244440839
116+2.65164709, +14.0770884, -6.84581470, -0.838441014, +1.21620166, +0.244440839, +9.00000000
118T
120T
121getMatDet(rand)
122+16.7363205
123
124
125ndim = getUnifRand(3, 9)
126ndim
127+4
128call setResized(rand, [ndim, ndim])
129scale = getUnifRand(1, 10, ndim)
130scale
131+1.00000000, +10.0000000, +6.00000000, +10.0000000
132
133call setCovRand(rngf, rand)
134rand
135+0.999999881, -0.983515084, +0.607585371, -0.476067483
136-0.983515084, +1.00000000, -0.696168244, +0.466376722
137+0.607585371, -0.696168244, +1.00000000, -0.625689209
138-0.476067483, +0.466376722, -0.625689209, +1.00000000
140T
142T
143getMatDet(rand)
144+0.460834661E-2
145call setCovRand(rngf, rand, scale(1))
146rand
147+1.00000000, +0.859935403, +0.795688748, -0.468120366
148+0.859935403, +0.999999881, +0.630496800, -0.167468503
149+0.795688748, +0.630496800, +0.999999940, -0.199580088
150-0.468120366, -0.167468503, -0.199580088, +1.00000000
152T
154T
155getMatDet(rand)
156+0.399442911E-1
157
158call setCovRand(rngf, rand, scale)
159rand
160+1.00000000, -5.46995926, +1.94590950, -0.435469180
161-5.46995926, +99.9999924, -33.1928864, -12.4557476
162+1.94590950, -33.1928864, +36.0000000, -27.3139687
163-0.435469180, -12.4557476, -27.3139687, +100.000008
165T
167T
168getMatDet(rand)
169+100682.133
170
171
172ndim = getUnifRand(3, 9)
173ndim
174+7
175call setResized(rand, [ndim, ndim])
176scale = getUnifRand(1, 10, ndim)
177scale
178+3.00000000, +1.00000000, +4.00000000, +2.00000000, +6.00000000, +5.00000000, +5.00000000
179
180call setCovRand(rngf, rand)
181rand
182+1.00000000, +0.225868165, +0.473733872, +0.146691009, -0.453077018, +0.238284275, +0.318854988
183+0.225868165, +1.00000000, +0.954240263, -0.536291957, +0.403288901, -0.504550576, -0.154705748
184+0.473733872, +0.954240263, +1.00000012, -0.355038047, +0.256449938, -0.469558716, -0.499768704E-1
185+0.146691009, -0.536291957, -0.355038047, +1.00000012, -0.136459470E-1, -0.248969764, +0.285035342
186-0.453077018, +0.403288901, +0.256449938, -0.136459470E-1, +0.999999821, -0.695822537, -0.326630354
187+0.238284275, -0.504550576, -0.469558716, -0.248969764, -0.695822537, +0.999999940, +0.113475487
188+0.318854988, -0.154705748, -0.499768704E-1, +0.285035342, -0.326630354, +0.113475487, +1.00000000
190T
192T
193getMatDet(rand)
194+0.947085994E-6
195call setCovRand(rngf, rand, scale(1))
196rand
197+9.00000000, -6.79219151, -5.17774534, +3.17225838, -2.92784262, -4.81202221, -3.20212984
198-6.79219151, +8.99999905, +4.87215614, +0.535456896, +3.96795011, +0.948569655, +2.97115159
199-5.17774534, +4.87215614, +9.00000000, -5.03748512, +6.97629642, +3.51789188, +3.23459983
200+3.17225838, +0.535456896, -5.03748512, +8.99999905, -5.91726732, -4.40062666, -4.47510433
201-2.92784262, +3.96795011, +6.97629642, -5.91726732, +9.00000000, +2.28185582, +6.06398773
202-4.81202221, +0.948569655, +3.51789188, -4.40062666, +2.28185582, +8.99999905, +3.17725801
203-3.20212984, +2.97115159, +3.23459983, -4.47510433, +6.06398773, +3.17725801, +9.00000095
205T
207T
208getMatDet(rand)
209+393.917938
210
211call setCovRand(rngf, rand, scale)
212rand
213+9.00000000, -2.43390346, +4.59102201, +0.553656518, +4.19657230, +7.46858788, -0.108692795
214-2.43390346, +1.00000000, -3.28771925, -0.103933565, -0.692990482, -2.00443339, -0.934959888
215+4.59102201, -3.28771925, +16.0000000, +1.94677031, +4.36853123, +2.08205390, +2.78367543
216+0.553656518, -0.103933565, +1.94677031, +4.00000048, +4.28786755, -4.73784924, -2.57987452
217+4.19657230, -0.692990482, +4.36853123, +4.28786755, +36.0000000, -2.45145988, -6.10959625
218+7.46858788, -2.00443339, +2.08205390, -4.73784924, -2.45145988, +25.0000000, +0.482200831
219-0.108692795, -0.934959888, +2.78367543, -2.57987452, -6.10959625, +0.482200831, +25.0000019
221T
223T
224getMatDet(rand)
225+450.903381
226
227
228ndim = getUnifRand(3, 9)
229ndim
230+4
231call setResized(rand, [ndim, ndim])
232scale = getUnifRand(1, 10, ndim)
233scale
234+4.00000000, +4.00000000, +5.00000000, +6.00000000
235
236call setCovRand(rngf, rand)
237rand
238+0.999999881, +0.136644214, -0.149096757, +0.489646077
239+0.136644214, +1.00000000, +0.599548578, -0.440750450
240-0.149096757, +0.599548578, +1.00000000, +0.448246002E-1
241+0.489646077, -0.440750450, +0.448246002E-1, +0.999999881
243T
245T
246getMatDet(rand)
247+0.975263268E-1
248call setCovRand(rngf, rand, scale(1))
249rand
250+16.0000000, +5.57715702, +4.32213497, -15.2615728
251+5.57715702, +16.0000000, -10.2397299, -2.76209116
252+4.32213497, -10.2397299, +16.0000019, -5.60533381
253-15.2615728, -2.76209116, -5.60533381, +16.0000019
255T
257T
258getMatDet(rand)
259+1006.20703
260
261call setCovRand(rngf, rand, scale)
262rand
263+16.0000000, -12.6560812, +12.2601757, -10.4007177
264-12.6560812, +16.0000038, -19.1436672, +14.0240307
265+12.2601757, -19.1436672, +25.0000000, -20.5133114
266-10.4007177, +14.0240307, -20.5133114, +36.0000076
268T
270T
271getMatDet(rand)
272+493.964905
273
274
275!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
276!Gram method for complex covariance.
277!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
278
279
280ndim = getUnifRand(3, 9)
281ndim
282+8
283call setResized(rand, [ndim, ndim])
284scale = getUnifRand(1, 10, ndim)
285scale
286+7.0000000000000000, +1.0000000000000000, +5.0000000000000000, +8.0000000000000000, +10.000000000000000, +5.0000000000000000, +2.0000000000000000, +5.0000000000000000
287
288call setCovRand(rngf, rand)
289rand
290+1.0000000000000000+0.0000000000000000i, -0.49832483866328869-0.69447989667501520i, -0.47972002420122223-0.25278382851354031E-1i, -0.53992853404372335E-1-0.51257778051012259i, -0.92321424482674219E-2+0.28110491412977828i, -0.24292094570834227+0.15185500597686385i, +0.33521229221315535E-1+0.27040602525048235i, +0.20526230740400536+0.11428534027238189i
291-0.49832483866328869+0.69447989667501520i, +1.0000000000000002+0.0000000000000000i, +0.45037142274296599-0.39477567341640390i, +0.27162568505049800+0.34678278472359092i, +0.73480786140718041E-1-0.35637139544124496i, +0.19978295535539997-0.11622916662629768E-1i, -0.30314452271131270-0.72443359474813895E-1i, -0.19465409279869628-0.54120136934354579E-1i
292-0.47972002420122223+0.25278382851354031E-1i, +0.45037142274296599+0.39477567341640390i, +1.0000000000000002+0.0000000000000000i, -0.54394560430452499+0.63726998826841597i, +0.45293024053587305-0.58937766906215372i, +0.13057382952964777-0.20506064785636452i, -0.23978566018174419-0.74026106587421592E-1i, -0.26263264488359572-0.44388202677625299i
293-0.53992853404372335E-1+0.51257778051012259i, +0.27162568505049800-0.34678278472359092i, -0.54394560430452499-0.63726998826841597i, +1.0000000000000002+0.0000000000000000i, -0.74949790603658284+0.74179737503371521E-1i, -0.96334293293891468E-1-0.48203868591324292E-1i, -0.90882879842327041E-1+0.69199766857724679E-2i, -0.14666654332803933+0.50214428042729975i
294-0.92321424482674219E-2-0.28110491412977828i, +0.73480786140718041E-1+0.35637139544124496i, +0.45293024053587305+0.58937766906215372i, -0.74949790603658284-0.74179737503371521E-1i, +0.99999999999999956+0.0000000000000000i, +0.13251080006424554+0.17457757700008039i, -0.11398602463933540+0.26245745466804554E-1i, +0.26172855995687933-0.36951618262055330i
295-0.24292094570834227-0.15185500597686385i, +0.19978295535539997+0.11622916662629768E-1i, +0.13057382952964777+0.20506064785636452i, -0.96334293293891468E-1+0.48203868591324292E-1i, +0.13251080006424554-0.17457757700008039i, +1.0000000000000000+0.0000000000000000i, -0.16269451299502213-0.24434748541414592i, -0.19823304207800527-0.33681520940505671i
296+0.33521229221315535E-1-0.27040602525048235i, -0.30314452271131270+0.72443359474813895E-1i, -0.23978566018174419+0.74026106587421592E-1i, -0.90882879842327041E-1-0.69199766857724679E-2i, -0.11398602463933540-0.26245745466804554E-1i, -0.16269451299502213+0.24434748541414592i, +1.0000000000000002+0.0000000000000000i, +0.43609467885263070E-1-0.90703647025568795E-1i
297+0.20526230740400536-0.11428534027238189i, -0.19465409279869628+0.54120136934354579E-1i, -0.26263264488359572+0.44388202677625299i, -0.14666654332803933-0.50214428042729975i, +0.26172855995687933+0.36951618262055330i, -0.19823304207800527+0.33681520940505671i, +0.43609467885263070E-1+0.90703647025568795E-1i, +1.0000000000000000+0.0000000000000000i
299T
301T
302getMatDet(rand)
303+0.82399401886198851E-5-0.72844833463869829E-19i
304call setCovRand(rngf, rand, scale(1))
305rand
306+49.000000000000000+0.0000000000000000i, -4.5174804212418893-16.560085840717001i, +8.3600036428003737+5.5151436149788973i, -25.200402790443800+19.134777030244816i, -7.0022733877475867+24.643040667601944i, -6.1701016171056358+10.930932151132515i, +0.62071493896050745+19.843955916149049i, -13.769720711096404+18.317169097622067i
307-4.5174804212418893+16.560085840717001i, +49.000000000000000+0.0000000000000000i, -20.826785484140576-27.922514394669847i, -14.387897615190655-6.2383779058941622i, +14.828008380132662-3.5668237566063525i, -0.57257417561588697+2.4855103773212326i, -1.6405232157644636-19.860471709635693i, -8.6800421247209520+10.726845697064210i
308+8.3600036428003737-5.5151436149788973i, -20.826785484140576+27.922514394669847i, +49.000000000000000+0.0000000000000000i, -8.0948367601173672+1.1167928356200152i, -3.9512782479219837+6.9814797479110826i, +8.5405082421495191+13.967914642988502i, +17.981220441529441+18.475265127863040i, -6.9526446853609940+4.4284695969111763i
309-25.200402790443800-19.134777030244816i, -14.387897615190655+6.2383779058941622i, -8.0948367601173672-1.1167928356200152i, +48.999999999999993+0.0000000000000000i, +6.4403388654424258-5.8807620307221979i, -8.0523843369438275-12.284023381745669i, -6.1903642863591983-20.205428237467494i, +19.223215471338357-20.738553148365156i
310-7.0022733877475867-24.643040667601944i, +14.828008380132662+3.5668237566063525i, -3.9512782479219837-6.9814797479110826i, +6.4403388654424258+5.8807620307221979i, +49.000000000000014+0.0000000000000000i, -10.006258478171125+17.277269646371959i, +9.0995302792721020-19.655343882269822i, -2.6535601603737549+11.103063434703525i
311-6.1701016171056358-10.930932151132515i, -0.57257417561588697-2.4855103773212326i, +8.5405082421495191-13.967914642988502i, -8.0523843369438275+12.284023381745669i, -10.006258478171125-17.277269646371959i, +48.999999999999986+0.0000000000000000i, +17.529022089657385+4.4699847317828887i, +28.431876618007319+10.326940013205563i
312+0.62071493896050745-19.843955916149049i, -1.6405232157644636+19.860471709635693i, +17.981220441529441-18.475265127863040i, -6.1903642863591983+20.205428237467494i, +9.0995302792721020+19.655343882269822i, +17.529022089657385-4.4699847317828887i, +49.000000000000000+0.0000000000000000i, +14.817728086381226+4.9252500160084738i
313-13.769720711096404-18.317169097622067i, -8.6800421247209520-10.726845697064210i, -6.9526446853609940-4.4284695969111763i, +19.223215471338357+20.738553148365156i, -2.6535601603737549-11.103063434703525i, +28.431876618007319-10.326940013205563i, +14.817728086381226-4.9252500160084738i, +49.000000000000014+0.0000000000000000i
315T
317T
318getMatDet(rand)
319+1671059209.2868035-0.11974756705805618E-4i
320
321call setCovRand(rngf, rand, scale)
322rand
323+49.000000000000000+0.0000000000000000i, -6.1612860414941073+2.0012408162165567i, +2.2240426525562991+7.5546154655638489i, +3.8586660822987309+30.114708684502581i, -28.989864511941768-42.209195888284015i, -1.1378722493831885-4.0838402059374586i, +2.8740306224584269-6.1151019226214238i, +3.7497833200578716-13.681979753815369i
324-6.1612860414941073-2.0012408162165567i, +1.0000000000000002+0.0000000000000000i, -0.74338902742851332-1.7144305494668579i, +1.6484892711513575-5.0079752464818243i, +2.4492175802020539+6.1461612372403476i, -0.83266894114766343+1.3453860400162176i, -0.58692732470136544+0.53347264646785209i, -0.66377140200104745+1.8121418261847968i
325+2.2240426525562991-7.5546154655638489i, -0.74338902742851332+1.7144305494668579i, +24.999999999999996+0.0000000000000000i, +9.0881451119497054+14.583996495037482i, -10.072020975530780+19.857294917630224i, -3.3430579786446146-3.7023698708437225i, -2.5485469814848569-1.3580112870658756i, +2.4021913759205669-4.6943359843703290i
326+3.8586660822987309-30.114708684502581i, +1.6484892711513575+5.0079752464818243i, +9.0881451119497054-14.583996495037482i, +64.000000000000000+0.0000000000000000i, -16.077786058362733+26.500416474628455i, -0.93342013892648179+6.6910242635485693i, -5.8900767835004935-2.2303202911730722i, -2.4978749791202359-6.3382448954961106i
327-28.989864511941768+42.209195888284015i, +2.4492175802020539-6.1461612372403476i, -10.072020975530780-19.857294917630224i, -16.077786058362733-26.500416474628455i, +100.00000000000000+0.0000000000000000i, +8.9119508309388333+12.899656157225131i, +2.8269335531165192+3.0969899885391161i, -1.8908095855220450-1.8916326785560891i
328-1.1378722493831885+4.0838402059374586i, -0.83266894114766343-1.3453860400162176i, -3.3430579786446146+3.7023698708437225i, -0.93342013892648179-6.6910242635485693i, +8.9119508309388333-12.899656157225131i, +24.999999999999996+0.0000000000000000i, -1.0638876925864840+2.2761744335602927i, -7.8049825794282448-8.3068224205299224i
329+2.8740306224584269+6.1151019226214238i, -0.58692732470136544-0.53347264646785209i, -2.5485469814848569+1.3580112870658756i, -5.8900767835004935+2.2303202911730722i, +2.8269335531165192-3.0969899885391161i, -1.0638876925864840-2.2761744335602927i, +4.0000000000000000+0.0000000000000000i, -0.12198309460948249-0.63703794298224503i
330+3.7497833200578716+13.681979753815369i, -0.66377140200104745-1.8121418261847968i, +2.4021913759205669+4.6943359843703290i, -2.4978749791202359+6.3382448954961106i, -1.8908095855220450+1.8916326785560891i, -7.8049825794282448+8.3068224205299224i, -0.12198309460948249+0.63703794298224503i, +25.000000000000004+0.0000000000000000i
332T
334T
335getMatDet(rand)
336+1086057.1147289227+0.17113052308559418E-7i
337
338
339ndim = getUnifRand(3, 9)
340ndim
341+9
342call setResized(rand, [ndim, ndim])
343scale = getUnifRand(1, 10, ndim)
344scale
345+6.0000000000000000, +5.0000000000000000, +9.0000000000000000, +6.0000000000000000, +10.000000000000000, +4.0000000000000000, +9.0000000000000000, +2.0000000000000000, +1.0000000000000000
346
347call setCovRand(rngf, rand)
348rand
349+1.0000000000000000+0.0000000000000000i, +0.32900590699352944-0.60215539007590091i, -0.34878222193818620+0.57382377724911748i, +0.16293905335421063-0.65573657212641386E-2i, -0.14541883389026303-0.48321938761512639i, -0.47500268105014903-0.18021057147627817i, -0.37803009397803139-0.40150551543415125i, +0.18789855007849834-0.26353962103767120i, +0.24189566659532263+0.14379234458328108i
350+0.32900590699352944+0.60215539007590091i, +1.0000000000000000+0.0000000000000000i, -0.51170406346552189-0.21284472942399368i, -0.32265236111794438E-1+0.26336784042009720i, +0.37617455313999915-0.62275982071497066E-1i, -0.25626420579401843-0.31428706749708324i, +0.14305014065680322-0.53174028392126993i, +0.49360796123049644+0.28407613520212727i, -0.18960586081750747+0.19854237796870111i
351-0.34878222193818620-0.57382377724911748i, -0.51170406346552189+0.21284472942399368i, +1.0000000000000000+0.0000000000000000i, -0.13902369578998863-0.47763544720401130i, -0.56687492425499486+0.31223469307791435i, -0.24538379829727000+0.14697527918102854i, +0.16674605509478441+0.39617620500194328i, -0.47170141579206670+0.41802140358943321E-1i, +0.92931037383273060E-1-0.16908208691741061i
352+0.16293905335421063+0.65573657212641386E-2i, -0.32265236111794438E-1-0.26336784042009720i, -0.13902369578998863+0.47763544720401130i, +1.0000000000000000+0.0000000000000000i, -0.58837852051103624E-1-0.59695899004357433i, +0.33380470567800824-0.15919370972745603i, -0.65821374869493784E-1-0.28741285216199147E-1i, -0.13073184136097038-0.24863831264263264i, -0.25491833717523438-0.15626919859592545i
353-0.14541883389026303+0.48321938761512639i, +0.37617455313999915+0.62275982071497066E-1i, -0.56687492425499486-0.31223469307791435i, -0.58837852051103624E-1+0.59695899004357433i, +1.0000000000000000+0.0000000000000000i, +0.40323704722156989E-1-0.14685274532269429i, +0.33779508856120566+0.47965506313200929E-1i, +0.40255567534359560-0.14571434722310361i, +0.18707174670290055+0.12702426507958525i
354-0.47500268105014903+0.18021057147627817i, -0.25626420579401843+0.31428706749708324i, -0.24538379829727000-0.14697527918102854i, +0.33380470567800824+0.15919370972745603i, +0.40323704722156989E-1+0.14685274532269429i, +1.0000000000000000+0.0000000000000000i, -0.40377668535180232E-1+0.30879236302780994i, +0.79268772662593204E-1+0.17606108947028576i, -0.42154880556331853-0.77186227855060047E-1i
355-0.37803009397803139+0.40150551543415125i, +0.14305014065680322+0.53174028392126993i, +0.16674605509478441-0.39617620500194328i, -0.65821374869493784E-1+0.28741285216199147E-1i, +0.33779508856120566-0.47965506313200929E-1i, -0.40377668535180232E-1-0.30879236302780994i, +1.0000000000000000+0.0000000000000000i, -0.36513437792295343E-1-0.15876156076685616i, +0.18942734694968388-0.25970121603848861i
356+0.18789855007849834+0.26353962103767120i, +0.49360796123049644-0.28407613520212727i, -0.47170141579206670-0.41802140358943321E-1i, -0.13073184136097038+0.24863831264263264i, +0.40255567534359560+0.14571434722310361i, +0.79268772662593204E-1-0.17606108947028576i, -0.36513437792295343E-1+0.15876156076685616i, +0.99999999999999989+0.0000000000000000i, +0.24380083098597743+0.49039827326152130i
357+0.24189566659532263-0.14379234458328108i, -0.18960586081750747-0.19854237796870111i, +0.92931037383273060E-1+0.16908208691741061i, -0.25491833717523438+0.15626919859592545i, +0.18707174670290055-0.12702426507958525i, -0.42154880556331853+0.77186227855060047E-1i, +0.18942734694968388+0.25970121603848861i, +0.24380083098597743-0.49039827326152130i, +0.99999999999999989+0.0000000000000000i
359T
361T
362getMatDet(rand)
363+0.38269802321705537E-5-0.17999450129153882E-19i
364call setCovRand(rngf, rand, scale(1))
365rand
366+36.000000000000000+0.0000000000000000i, +19.363829600610593+17.376844187893692i, +17.206759786725925+7.4597724196007045i, -17.821842480593190-17.494392823175971i, -0.77625016634315835-11.567041412422100i, +10.401533332524508+18.312324578761611i, -5.3978274985189421-9.5615459830395366i, -15.361305345619648+2.0940482816834938i, +0.49798982146565740+6.4597558567529214i
367+19.363829600610593-17.376844187893692i, +36.000000000000000+0.0000000000000000i, +26.523070245720216+8.5029359723746776i, -26.194663018313044-7.1927095942124764i, -15.423227910174313-15.206547553577025i, +15.295731624149660+5.7558519944592534i, -4.9222984200229014+5.3181319156497269i, -13.501024493142154+8.9278351933838351i, -1.8720409568027889+1.7302137445461248i
368+17.206759786725925-7.4597724196007045i, +26.523070245720216-8.5029359723746776i, +36.000000000000007+0.0000000000000000i, -25.505262998954109+0.48042367214056847i, -13.374683512559727-5.1771658445461810i, +15.779197784070799+2.1009045328073728i, +9.3878434491285567+1.7676819804291128i, -8.4210639671366874+5.4643958829855164i, -0.12814696131537628E-1+3.3831916342645290i
369-17.821842480593190+17.494392823175971i, -26.194663018313044+7.1927095942124764i, -25.505262998954109-0.48042367214056847i, +35.999999999999993+0.0000000000000000i, +19.980201357021095+10.150718279320520i, -17.025019249096424-6.2455219343651995i, -4.6741641502609603-11.839679468994458i, +8.4249832543982919-15.615812946268463i, -2.1463271638330133-0.30431624452945671E-1i
370-0.77625016634315835+11.567041412422100i, -15.423227910174313+15.206547553577025i, -13.374683512559727+5.1771658445461810i, +19.980201357021095-10.150718279320520i, +35.999999999999986+0.0000000000000000i, -5.4024986868039715-5.8352952993906149i, -14.725955921120875-11.456141864607607i, +6.6566248311870240-23.175051306017924i, +12.242484606818131-5.3728458443223071i
371+10.401533332524508-18.312324578761611i, +15.295731624149660-5.7558519944592534i, +15.779197784070799-2.1009045328073728i, -17.025019249096424+6.2455219343651995i, -5.4024986868039715+5.8352952993906149i, +36.000000000000000+0.0000000000000000i, -2.1514048624177935+8.7155330901141550i, +2.4344150010900130+8.8223078816080154i, +14.799474063681258+0.52050176790251801i
372-5.3978274985189421+9.5615459830395366i, -4.9222984200229014-5.3181319156497269i, +9.3878434491285567-1.7676819804291128i, -4.6741641502609603+11.839679468994458i, -14.725955921120875+11.456141864607607i, -2.1514048624177935-8.7155330901141550i, +35.999999999999986+0.0000000000000000i, +7.4413342177906801+7.6980758240236096i, -7.7430518988431638+2.2523680056663120i
373-15.361305345619648-2.0940482816834938i, -13.501024493142154-8.9278351933838351i, -8.4210639671366874-5.4643958829855164i, +8.4249832543982919+15.615812946268463i, +6.6566248311870240+23.175051306017924i, +2.4344150010900130-8.8223078816080154i, +7.4413342177906801-7.6980758240236096i, +36.000000000000000+0.0000000000000000i, +11.379112849886294-3.7025697897616570i
374+0.49798982146565740-6.4597558567529214i, -1.8720409568027889-1.7302137445461248i, -0.12814696131537628E-1-3.3831916342645290i, -2.1463271638330133+0.30431624452945671E-1i, +12.242484606818131+5.3728458443223071i, +14.799474063681258-0.52050176790251801i, -7.7430518988431638-2.2523680056663120i, +11.379112849886294+3.7025697897616570i, +36.000000000000000+0.0000000000000000i
376T
378T
379getMatDet(rand)
380+795804.00248007581+0.37316349335014820E-5i
381
382call setCovRand(rngf, rand, scale)
383rand
384+36.000000000000000+0.0000000000000000i, -6.8309581967550459-16.617582989313341i, -10.757543521582237-26.500265081753078i, -9.9392440871018675-22.085332126479525i, -24.511210007979894-2.1100018652956272i, +8.3116542489792611-3.4156797905070055i, +19.331168570763460+14.220211343248582i, -3.9924039996229146+2.1424202240379349i, +1.7922184859218482+0.61522483571171083i
385-6.8309581967550459+16.617582989313341i, +24.999999999999993+0.0000000000000000i, +16.365765831997262-3.5257407520548423i, +13.308693485204632+1.4209449352188597i, -14.873580987687731+2.3969389475399407i, -2.4831251470734532+6.8420853497323346i, +2.1222489474627744-2.9833551383408539i, -2.7790981138722297-5.2470926471078148i, +0.41688465537608077+1.9937034845182209i
386-10.757543521582237+26.500265081753078i, +16.365765831997262+3.5257407520548423i, +81.000000000000000+0.0000000000000000i, +35.610646743606765-11.858125215641920i, -30.516827569718640-15.452957481361352i, -6.9701974764560894+6.0291232144766784i, -19.608499508369306-3.7167293882003576i, -4.6569448671223537+0.61491497959616803i, +0.57088909447911473+2.5215675171168543i
387-9.9392440871018675+22.085332126479525i, +13.308693485204632-1.4209449352188597i, +35.610646743606765+11.858125215641920i, +36.000000000000014+0.0000000000000000i, -2.9834782754172888-11.243676137523659i, -8.3661026832886112+7.3510684830782287i, -19.689610856224789-3.2748627678056144i, -2.3933022343003820-2.5317627784068790i, +0.45686480238621741+1.9209745369723792i
388-24.511210007979894+2.1100018652956272i, -14.873580987687731-2.3969389475399407i, -30.516827569718640+15.452957481361352i, -2.9834782754172888+11.243676137523659i, +100.00000000000000+0.0000000000000000i, -5.4266345627782879+6.1482695700877121i, -20.261320408910684-0.58830073079258116i, +9.1763899286220187-1.1059075733631047i, -0.67737930974502847-4.7606732422380551i
389+8.3116542489792611+3.4156797905070055i, -2.4831251470734532-6.8420853497323346i, -6.9701974764560894-6.0291232144766784i, -8.3661026832886112-7.3510684830782287i, -5.4266345627782879-6.1482695700877121i, +15.999999999999998+0.0000000000000000i, +2.8746493371091622+14.787992798200060i, -3.8229747616666803+1.5509372341744034i, -0.70302257820246905-0.85123490298598625i
390+19.331168570763460-14.220211343248582i, +2.1222489474627744+2.9833551383408539i, -19.608499508369306+3.7167293882003576i, -19.689610856224789+3.2748627678056144i, -20.261320408910684+0.58830073079258116i, +2.8746493371091622-14.787992798200060i, +81.000000000000000+0.0000000000000000i, +0.43104170801621977-0.56045597906944389i, +2.8704865011386747-0.82501458903096681i
391-3.9924039996229146-2.1424202240379349i, -2.7790981138722297+5.2470926471078148i, -4.6569448671223537-0.61491497959616803i, -2.3933022343003820+2.5317627784068790i, +9.1763899286220187+1.1059075733631047i, -3.8229747616666803-1.5509372341744034i, +0.43104170801621977+0.56045597906944389i, +4.0000000000000009+0.0000000000000000i, -0.23515284397748801-0.31842297652601936i
392+1.7922184859218482-0.61522483571171083i, +0.41688465537608077-1.9937034845182209i, +0.57088909447911473-2.5215675171168543i, +0.45686480238621741-1.9209745369723792i, -0.67737930974502847+4.7606732422380551i, -0.70302257820246905+0.85123490298598625i, +2.8704865011386747+0.82501458903096681i, -0.23515284397748801+0.31842297652601936i, +1.0000000000000000+0.0000000000000000i
394T
396T
397getMatDet(rand)
398+1345255.5204133287+0.18646824173629284E-6i
399
400
401ndim = getUnifRand(3, 9)
402ndim
403+3
404call setResized(rand, [ndim, ndim])
405scale = getUnifRand(1, 10, ndim)
406scale
407+3.0000000000000000, +4.0000000000000000, +7.0000000000000000
408
409call setCovRand(rngf, rand)
410rand
411+1.0000000000000000+0.0000000000000000i, -0.24278577640165641+0.11412668922510906i, -0.68113494530342114E-1-0.66949292838016006i
412-0.24278577640165641-0.11412668922510906i, +1.0000000000000002+0.0000000000000000i, -0.43537733056428018-0.39330188025678225i
413-0.68113494530342114E-1+0.66949292838016006i, -0.43537733056428018+0.39330188025678225i, +1.0000000000000002+0.0000000000000000i
415T
417T
418getMatDet(rand)
419+0.49090313459925518E-1+0.0000000000000000i
420call setCovRand(rngf, rand, scale(1))
421rand
422+8.9999999999999982+0.0000000000000000i, -1.2268775348493075+4.7727421885624093i, -0.74955307356949619+1.9033028376614143i
423-1.2268775348493075-4.7727421885624093i, +9.0000000000000000+0.0000000000000000i, +5.3730932415273074+4.5469180126274109i
424-0.74955307356949619-1.9033028376614143i, +5.3730932415273074-4.5469180126274109i, +9.0000000000000000+0.0000000000000000i
426T
428T
429getMatDet(rand)
430+145.67814282034360-0.27753346990484932E-13i
431
432call setCovRand(rngf, rand, scale)
433rand
434+9.0000000000000000+0.0000000000000000i, +9.9383472134967903+4.5964928768958284i, +10.039673505559854+6.2027135898001040i
435+9.9383472134967903-4.5964928768958284i, +15.999999999999998+0.0000000000000000i, +21.121835772349232+5.1862711230163949i
436+10.039673505559854-6.2027135898001040i, +21.121835772349232-5.1862711230163949i, +48.999999999999993+0.0000000000000000i
438T
440T
441getMatDet(rand)
442+275.53511697642091-0.85265128291212022E-13i
443
444
445ndim = getUnifRand(3, 9)
446ndim
447+3
448call setResized(rand, [ndim, ndim])
449scale = getUnifRand(1, 10, ndim)
450scale
451+8.0000000000000000, +3.0000000000000000, +9.0000000000000000
452
453call setCovRand(rngf, rand)
454rand
455+1.0000000000000000+0.0000000000000000i, +0.51121891038361800-0.51620505590457177i, -0.47939384571359805-0.15744238857535998E-1i
456+0.51121891038361800+0.51620505590457177i, +1.0000000000000000+0.0000000000000000i, +0.17730207904467796-0.21296023539568115i
457-0.47939384571359805+0.15744238857535998E-1i, +0.17730207904467796+0.21296023539568115i, +1.0000000000000000+0.0000000000000000i
459T
461T
462getMatDet(rand)
463+0.19013917260571075-0.32764594226960008E-17i
464call setCovRand(rngf, rand, scale(1))
465rand
466+64.000000000000000+0.0000000000000000i, +22.346598233177698-30.318600807566849i, +43.570237114928226+2.3034959954035599i
467+22.346598233177698+30.318600807566849i, +64.000000000000014+0.0000000000000000i, +10.084786313213627-14.371918753463500i
468+43.570237114928226-2.3034959954035599i, +10.084786313213627+14.371918753463500i, +64.000000000000000+0.0000000000000000i
470T
472T
473getMatDet(rand)
474+8570.6010207660802+0.0000000000000000i
475
476call setCovRand(rngf, rand, scale)
477rand
478+64.000000000000000+0.0000000000000000i, +10.021961814385254-2.8173184592597456i, -35.200311334302036+3.0042895326997621i
479+10.021961814385254+2.8173184592597456i, +8.9999999999999964+0.0000000000000000i, +6.6304434088630675-11.967363751042312i
480-35.200311334302036-3.0042895326997621i, +6.6304434088630675+11.967363751042312i, +80.999999999999986+0.0000000000000000i
482T
484T
485getMatDet(rand)
486+11527.722224339237+0.0000000000000000i
487
488
489ndim = getUnifRand(3, 9)
490ndim
491+4
492call setResized(rand, [ndim, ndim])
493scale = getUnifRand(1, 10, ndim)
494scale
495+6.0000000000000000, +9.0000000000000000, +7.0000000000000000, +3.0000000000000000
496
497call setCovRand(rngf, rand)
498rand
499+1.0000000000000000+0.0000000000000000i, +0.59297559022310786+0.42443352917943344i, +0.23709311248019593E-2-0.12137061455651738i, +0.48713819195425179-0.90410769403250213E-1i
500+0.59297559022310786-0.42443352917943344i, +1.0000000000000000+0.0000000000000000i, +0.37493260256898620+0.79577368433677034E-1i, +0.14769101945066382-0.14876152195656939i
501+0.23709311248019593E-2+0.12137061455651738i, +0.37493260256898620-0.79577368433677034E-1i, +0.99999999999999978+0.0000000000000000i, +0.37169767109295143+0.14564057923451978i
502+0.48713819195425179+0.90410769403250213E-1i, +0.14769101945066382+0.14876152195656939i, +0.37169767109295143-0.14564057923451978i, +0.99999999999999978+0.0000000000000000i
504T
506T
507getMatDet(rand)
508+0.99244257229233801E-1+0.89304965079723902E-18i
509call setCovRand(rngf, rand, scale(1))
510rand
511+36.000000000000000+0.0000000000000000i, +12.437853375015122-0.27470772479900601i, -7.3601861805125992+23.786179791268921i, +3.2225116255443047+16.543099187280301i
512+12.437853375015122+0.27470772479900601i, +36.000000000000000+0.0000000000000000i, -1.9102351621555393-8.6689413483006632i, +10.614309002522168+13.894111930952757i
513-7.3601861805125992-23.786179791268921i, -1.9102351621555393+8.6689413483006632i, +36.000000000000007+0.0000000000000000i, +6.1073375154922722+9.8274126166904452i
514+3.2225116255443047-16.543099187280301i, +10.614309002522168-13.894111930952757i, +6.1073375154922722-9.8274126166904452i, +36.000000000000007+0.0000000000000000i
516T
518T
519getMatDet(rand)
520+144262.32779837449-0.13642420526593924E-11i
521
522call setCovRand(rngf, rand, scale)
523rand
524+36.000000000000000+0.0000000000000000i, +37.926476006079497+23.145741734263112i, +3.2889569612965284+18.520683921642735i, -5.3482961053538620-7.4392123483353796i
525+37.926476006079497-23.145741734263112i, +81.000000000000000+0.0000000000000000i, +36.392647399961405+35.647760423015654i, -12.045109205269821-7.4756209357306478i
526+3.2889569612965284-18.520683921642735i, +36.392647399961405-35.647760423015654i, +49.000000000000000+0.0000000000000000i, -10.901211901269006+5.2427253273340533i
527-5.3482961053538620+7.4392123483353796i, -12.045109205269821+7.4756209357306478i, -10.901211901269006-5.2427253273340533i, +9.0000000000000000+0.0000000000000000i
529T
531T
532getMatDet(rand)
533+27451.734264602517-0.14551915228366852E-10i
534
535
536!%%%%%%%%%%%%%%%%%%%%%%%%
537!Dvine and Onion methods.
538!%%%%%%%%%%%%%%%%%%%%%%%%
539
540
541eta = getUnifRand(1, 10)
542ndim = getUnifRand(2, 5)
543call setResized(rand, [ndim, ndim])
544scale = getUnifRand(1, 10)
545
546call setCovRand(rngf, rand, dvine, eta)
547onion%info
548+0
549rand
550+1.00000000, -0.768363476E-2, +0.528398871, -0.128390014
551-0.768363476E-2, +1.00000000, -0.180511743, -0.142434061
552+0.528398871, -0.180511743, +1.00000000, -0.337057561
553-0.128390014, -0.142434061, -0.337057561, +1.00000000
555T
557T
558getMatDet(rand)
559+0.577323616
560call setCovRand(rngf, rand, onion, eta)
561onion%info
562+0
563rand
564+1.00000000, -0.393439054, -0.133276591E-2, +0.211140782
565-0.393439054, +1.00000000, +0.107175291, +0.223668545
566-0.133276591E-2, +0.107175291, +1.00000000, +0.270004451
567+0.211140782, +0.223668545, +0.270004451, +1.00000000
569T
571T
572getMatDet(rand)
573+0.658506513
574call setCovRand(rngf, rand, dvine, eta, scale)
575rand
576+9.00000000, -2.02085018, +2.40207481, -1.56879616
577-2.02085018, +9.00000000, -1.13287055, -0.606026709
578+2.40207481, -1.13287055, +9.00000000, -4.06085205
579-1.56879616, -0.606026709, -4.06085205, +9.00000000
581T
583T
584getMatDet(rand)
585+4456.24512
586
587call setCovRand(rngf, rand, onion, eta, scale)
588onion%info
589+0
590rand
591+9.00000000, +0.433134913, +0.244668350E-1, +0.519786358
592+0.433134913, +9.00000000, -2.47723293, -4.54324818
593+0.244668350E-1, -2.47723293, +9.00000000, +8.13309383
594+0.519786358, -4.54324818, +8.13309383, +9.00000000
596T
598T
599getMatDet(rand)
600+651.209534
601
602
603eta = getUnifRand(1, 10)
604ndim = getUnifRand(2, 5)
605call setResized(rand, [ndim, ndim])
606scale = getUnifRand(1, 10)
607
608call setCovRand(rngf, rand, dvine, eta)
609onion%info
610+0
611rand
612+1.00000000, -0.148721457, -0.146218479, -0.997059345E-1
613-0.148721457, +1.00000000, +0.279486716, +0.295421571
614-0.146218479, +0.279486716, +1.00000000, -0.583430529E-1
615-0.997059345E-1, +0.295421571, -0.583430529E-1, +1.00000000
617T
619T
620getMatDet(rand)
621+0.788894773
622call setCovRand(rngf, rand, onion, eta)
623onion%info
624+0
625rand
626+1.00000000, +0.349811316E-1, -0.164117634, -0.409955252E-2
627+0.349811316E-1, +1.00000000, +0.354991198, -0.668283701E-1
628-0.164117634, +0.354991198, +1.00000000, +0.199157342
629-0.409955252E-2, -0.668283701E-1, +0.199157342, +1.00000000
631T
633T
634getMatDet(rand)
635+0.788508475
636call setCovRand(rngf, rand, dvine, eta, scale)
637rand
638+4.00000000, -0.242325544, +0.304498672E-1, +0.166408539
639-0.242325544, +4.00000000, -1.07088590, -1.09515822
640+0.304498672E-1, -1.07088590, +4.00000000, -0.323414981
641+0.166408539, -1.09515822, -0.323414981, +4.00000000
643T
645T
646getMatDet(rand)
647+212.819626
648
649call setCovRand(rngf, rand, onion, eta, scale)
650onion%info
651+0
652rand
653+4.00000000, -0.232072830, -1.69760776, -0.412477776E-1
654-0.232072830, +4.00000000, -0.959614992, +2.03312755
655-1.69760776, -0.959614992, +4.00000000, -3.27844453
656-0.412477776E-1, +2.03312755, -3.27844453, +4.00000000
658T
660T
661getMatDet(rand)
662+20.5815678
663
664
665eta = getUnifRand(1, 10)
666ndim = getUnifRand(2, 5)
667call setResized(rand, [ndim, ndim])
668scale = getUnifRand(1, 10)
669
670call setCovRand(rngf, rand, dvine, eta)
671onion%info
672+0
673rand
674+1.00000000, +0.321197748
675+0.321197748, +1.00000000
677T
679T
680getMatDet(rand)
681+0.896831989
682call setCovRand(rngf, rand, onion, eta)
683onion%info
684+0
685rand
686+1.00000000, +0.594866276E-2
687+0.594866276E-2, +1.00000000
689T
691T
692getMatDet(rand)
693+0.999964595
694call setCovRand(rngf, rand, dvine, eta, scale)
695rand
696+100.000000, +11.2097855
697+11.2097855, +100.000000
699T
701T
702getMatDet(rand)
703+9874.34082
704
705call setCovRand(rngf, rand, onion, eta, scale)
706onion%info
707+0
708rand
709+100.000000, -22.2780285
710-22.2780285, +100.000000
712T
714T
715getMatDet(rand)
716+9503.68945
717
718
719eta = getUnifRand(1, 10)
720ndim = getUnifRand(2, 5)
721call setResized(rand, [ndim, ndim])
722scale = getUnifRand(1, 10)
723
724call setCovRand(rngf, rand, dvine, eta)
725onion%info
726+0
727rand
728+1.00000000, +0.223436594
729+0.223436594, +1.00000000
731T
733T
734getMatDet(rand)
735+0.950076103
736call setCovRand(rngf, rand, onion, eta)
737onion%info
738+0
739rand
740+1.00000000, +0.309052467E-1
741+0.309052467E-1, +1.00000000
743T
745T
746getMatDet(rand)
747+0.999044895
748call setCovRand(rngf, rand, dvine, eta, scale)
749rand
750+9.00000000, -2.93303156
751-2.93303156, +9.00000000
753T
755T
756getMatDet(rand)
757+72.3973236
758
759call setCovRand(rngf, rand, onion, eta, scale)
760onion%info
761+0
762rand
763+9.00000000, -4.36158466
764-4.36158466, +9.00000000
766T
768T
769getMatDet(rand)
770+61.9765816
771
772
773eta = getUnifRand(1, 10)
774ndim = getUnifRand(2, 5)
775call setResized(rand, [ndim, ndim])
776scale = getUnifRand(1, 10)
777
778call setCovRand(rngf, rand, dvine, eta)
779onion%info
780+0
781rand
782+1.00000000, +0.190069318, +0.118904114
783+0.190069318, +1.00000000, +0.336788595
784+0.118904114, +0.336788595, +1.00000000
786T
788T
789getMatDet(rand)
790+0.851531744
791call setCovRand(rngf, rand, onion, eta)
792onion%info
793+0
794rand
795+1.00000000, -0.805101991E-1, -0.118841931
796-0.805101991E-1, +1.00000000, -0.342192024
797-0.118841931, -0.342192024, +1.00000000
799T
801T
802getMatDet(rand)
803+0.855751216
804call setCovRand(rngf, rand, dvine, eta, scale)
805rand
806+4.00000000, -0.874838829, +0.259295940
807-0.874838829, +4.00000000, +0.555085242
808+0.259295940, +0.555085242, +4.00000000
810T
812T
813getMatDet(rand)
814+59.1853790
815
816call setCovRand(rngf, rand, onion, eta, scale)
817onion%info
818+0
819rand
820+4.00000000, -0.600556374, +0.115844794
821-0.600556374, +4.00000000, +0.723906040
822+0.115844794, +0.723906040, +4.00000000
824T
826T
827getMatDet(rand)
828+60.3067627
829
830
831eta = getUnifRand(1, 10)
832ndim = getUnifRand(2, 5)
833call setResized(rand, [ndim, ndim])
834scale = getUnifRand(1, 10)
835
836call setCovRand(rngf, rand, dvine, eta)
837onion%info
838+0
839rand
840+1.00000000, +0.550081730E-1, +0.151851177E-1, -0.227699697, -0.105471849
841+0.550081730E-1, +1.00000000, +0.115507282, +0.203081131, -0.233924866
842+0.151851177E-1, +0.115507282, +1.00000000, -0.190082490, -0.179092839
843-0.227699697, +0.203081131, -0.190082490, +1.00000000, -0.400130272
844-0.105471849, -0.233924866, -0.179092839, -0.400130272, +1.00000000
846T
848T
849getMatDet(rand)
850+0.606201828
851call setCovRand(rngf, rand, onion, eta)
852onion%info
853+0
854rand
855+1.00000000, -0.230266452, -0.650553871E-2, +0.604735166E-1, -0.903723761E-2
856-0.230266452, +1.00000000, +0.162601739, -0.342647910, +0.219620438E-2
857-0.650553871E-2, +0.162601739, +1.00000000, +0.204351634, -0.124991737
858+0.604735166E-1, -0.342647910, +0.204351634, +1.00000000, -0.263652056E-1
859-0.903723761E-2, +0.219620438E-2, -0.124991737, -0.263652056E-1, +1.00000000
861T
863T
864getMatDet(rand)
865+0.736026824
866call setCovRand(rngf, rand, dvine, eta, scale)
867rand
868+4.00000000, -0.296523333, +1.41577721, -1.29601502, +0.621545315
869-0.296523333, +4.00000000, -1.25008404, +0.283018649, +0.800302744
870+1.41577721, -1.25008404, +4.00000000, -1.72528327, -1.07294893
871-1.29601502, +0.283018649, -1.72528327, +4.00000000, -1.20609570
872+0.621545315, +0.800302744, -1.07294893, -1.20609570, +4.00000000
874T
876T
877getMatDet(rand)
878+423.547577
879
880call setCovRand(rngf, rand, onion, eta, scale)
881onion%info
882+0
883rand
884+4.00000000, -0.557724714, -1.15070164, -0.795428693, +0.881383121E-1
885-0.557724714, +4.00000000, -0.236993670, +0.746449351, -0.653660595
886-1.15070164, -0.236993670, +4.00000000, -1.86961102, +1.98973548
887-0.795428693, +0.746449351, -1.86961102, +4.00000000, -3.70241642
888+0.881383121E-1, -0.653660595, +1.98973548, -3.70241642, +4.00000000
890T
892T
893getMatDet(rand)
894+68.5610352
895
896
897eta = getUnifRand(1, 10)
898ndim = getUnifRand(2, 5)
899call setResized(rand, [ndim, ndim])
900scale = getUnifRand(1, 10)
901
902call setCovRand(rngf, rand, dvine, eta)
903onion%info
904+0
905rand
906+1.00000000, +0.140548348, -0.374145448, +0.225747824E-1, +0.713772774E-1
907+0.140548348, +1.00000000, +0.362978317E-1, -0.241430983, +0.270006925
908-0.374145448, +0.362978317E-1, +1.00000000, -0.217398573E-1, +0.191654682
909+0.225747824E-1, -0.241430983, -0.217398573E-1, +1.00000000, -0.309428275
910+0.713772774E-1, +0.270006925, +0.191654682, -0.309428275, +1.00000000
912T
914T
915getMatDet(rand)
916+0.639279664
917call setCovRand(rngf, rand, onion, eta)
918onion%info
919+0
920rand
921+1.00000000, +0.197078705, -0.118305072, -0.421691686, +0.366986990
922+0.197078705, +1.00000000, +0.252950966, +0.343438923, +0.241906606E-1
923-0.118305072, +0.252950966, +1.00000000, +0.989452228E-1, +0.267319381E-1
924-0.421691686, +0.343438923, +0.989452228E-1, +1.00000000, -0.272463337E-1
925+0.366986990, +0.241906606E-1, +0.267319381E-1, -0.272463337E-1, +1.00000000
927T
929T
930getMatDet(rand)
931+0.446025252
932call setCovRand(rngf, rand, dvine, eta, scale)
933rand
934+81.0000000, -6.87514448, +10.4259281, -13.9465799, +17.5207844
935-6.87514448, +81.0000000, +1.04446101, +12.8645124, -7.18192768
936+10.4259281, +1.04446101, +81.0000000, +10.2399712, -6.84720087
937-13.9465799, +12.8645124, +10.2399712, +81.0000000, -19.8356571
938+17.5207844, -7.18192768, -6.84720087, -19.8356571, +81.0000000
940T
942T
943getMatDet(rand)
944+0.284292762E+10
945
946call setCovRand(rngf, rand, onion, eta, scale)
947onion%info
948+3
949rand
950+1.00000000, +0.353475094, +0.134735674, -13.9465799, +17.5207844
951-6.87514448, +1.00000000, -1.08041525, +12.8645124, -7.18192768
952+10.4259281, +1.04446101, +1.00000000, +10.2399712, -6.84720087
953-13.9465799, +12.8645124, +10.2399712, +81.0000000, -19.8356571
954+17.5207844, -7.18192768, -6.84720087, -19.8356571, +81.0000000
956F
958F
959getMatDet(rand)
960-456371.531
961
962
963eta = getUnifRand(1, 10)
964ndim = getUnifRand(2, 5)
965call setResized(rand, [ndim, ndim])
966scale = getUnifRand(1, 10)
967
968call setCovRand(rngf, rand, dvine, eta)
969onion%info
970+3
971rand
972+1.00000000, +0.789076090E-1, +0.110274076
973+0.789076090E-1, +1.00000000, -0.266308516
974+0.110274076, -0.266308516, +1.00000000
976T
978T
979getMatDet(rand)
980+0.906058431
981call setCovRand(rngf, rand, onion, eta)
982onion%info
983+0
984rand
985+1.00000000, +0.230741024, +0.121401235
986+0.230741024, +1.00000000, -0.285949763E-1
987+0.121401235, -0.285949763E-1, +1.00000000
989T
991T
992getMatDet(rand)
993+0.929600596
994call setCovRand(rngf, rand, dvine, eta, scale)
995rand
996+64.0000000, -6.30038452, +12.8439178
997-6.30038452, +64.0000000, -23.9671268
998+12.8439178, -23.9671268, +64.0000000
1000T
1001isMatClass(rand, hermitian)
1002T
1003getMatDet(rand)
1004+216161.531
1005
1006call setCovRand(rngf, rand, onion, eta, scale)
1007onion%info
1008+0
1009rand
1010+64.0000000, +15.6495132, -30.1945972
1011+15.6495132, +64.0000000, +197.622726
1012-30.1945972, +197.622726, +64.0000000
1013isMatClass(rand, posdefmat)
1014F
1015isMatClass(rand, hermitian)
1016T
1017getMatDet(rand)
1018-2498148.50
1019
1020
1021eta = getUnifRand(1, 10)
1022ndim = getUnifRand(2, 5)
1023call setResized(rand, [ndim, ndim])
1024scale = getUnifRand(1, 10)
1025
1026call setCovRand(rngf, rand, dvine, eta)
1027onion%info
1028+0
1029rand
1030+1.00000000, -0.453473151
1031-0.453473151, +1.00000000
1032isMatClass(rand, posdefmat)
1033T
1034isMatClass(rand, hermitian)
1035T
1036getMatDet(rand)
1037+0.794362068
1038call setCovRand(rngf, rand, onion, eta)
1039onion%info
1040+0
1041rand
1042+1.00000000, -0.244986296
1043-0.244986296, +1.00000000
1044isMatClass(rand, posdefmat)
1045T
1046isMatClass(rand, hermitian)
1047T
1048getMatDet(rand)
1049+0.939981699
1050call setCovRand(rngf, rand, dvine, eta, scale)
1051rand
1052+100.000000, -50.2151375
1053-50.2151375, +100.000000
1054isMatClass(rand, posdefmat)
1055T
1056isMatClass(rand, hermitian)
1057T
1058getMatDet(rand)
1059+7478.43994
1060
1061call setCovRand(rngf, rand, onion, eta, scale)
1062onion%info
1063+0
1064rand
1065+100.000000, -58.1804390
1066-58.1804390, +100.000000
1067isMatClass(rand, posdefmat)
1068T
1069isMatClass(rand, hermitian)
1070T
1071getMatDet(rand)
1072+6615.03613
1073
1074
1075eta = getUnifRand(1, 10)
1076ndim = getUnifRand(2, 5)
1077call setResized(rand, [ndim, ndim])
1078scale = getUnifRand(1, 10)
1079
1080call setCovRand(rngf, rand, dvine, eta)
1081onion%info
1082+0
1083rand
1084+1.00000000, +0.504769921, +0.194406152
1085+0.504769921, +1.00000000, +0.290999293
1086+0.194406152, +0.290999293, +1.00000000
1087isMatClass(rand, posdefmat)
1088T
1089isMatClass(rand, hermitian)
1090T
1091getMatDet(rand)
1092+0.679844737
1093call setCovRand(rngf, rand, onion, eta)
1094onion%info
1095+0
1096rand
1097+1.00000000, +0.221363544, -0.530458167E-1
1098+0.221363544, +1.00000000, -0.197407514
1099-0.530458167E-1, -0.197407514, +1.00000000
1100isMatClass(rand, posdefmat)
1101T
1102isMatClass(rand, hermitian)
1103T
1104getMatDet(rand)
1105+0.913850665
1106call setCovRand(rngf, rand, dvine, eta, scale)
1107rand
1108+64.0000000, -10.3690300, -1.84080505
1109-10.3690300, +64.0000000, +7.63591480
1110-1.84080505, +7.63591480, +64.0000000
1111isMatClass(rand, posdefmat)
1112T
1113isMatClass(rand, hermitian)
1114T
1115getMatDet(rand)
1116+251605.875
1117
1118call setCovRand(rngf, rand, onion, eta, scale)
1119onion%info
1120+0
1121rand
1122+64.0000000, -9.24347687, -14.1854877
1123-9.24347687, +64.0000000, +144.050735
1124-14.1854877, +144.050735, +64.0000000
1125isMatClass(rand, posdefmat)
1126F
1127isMatClass(rand, hermitian)
1128T
1129getMatDet(rand)
1130-1046465.44
1131
1132
1133ndim = getUnifRand(2, 10)
1134call setResized(rand, [ndim, ndim])
1135
1136call setCovRand(rngf, rand, dvine, eta = 0._TKG, scale = [(real(itry, TKG), itry = 1, ndim)])
1137rand
1138+1.00000000, +0.511458397, +1.91149914
1139+0.511458397, +4.00000000, -0.781309605
1140+1.91149914, -0.781309605, +9.00000000
1141isMatClass(rand, posdefmat)
1142T
1143
1144call setCovRand(rngf, rand, onion, eta = 0._TKG, scale = [(real(itry, TKG), itry = 1, ndim)])
1145onion%info
1146+0
1147rand
1148+1.00000000, -0.571775317, -2.37024045
1149-0.571775317, +4.00000000, -1.08830500
1150-2.37024045, -1.08830500, +9.00000000
1151isMatClass(rand, posdefmat)
1152T
1153
1154
Test:
test_pm_distCov
Todo:
High Priority: The current implementation of this generic interface uses a naive method of computing the Cholesky factorization with a default matrix packing for the Onion method.
The RFP packing format must be also implemented for this generic interface.
Todo:
High Priority: The current implementation of the Gram method can be significantly improved, both computationally and functionally.


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, Monday March 6, 2017, 3:22 pm, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin.

Definition at line 787 of file pm_distCov.F90.


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