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

Return a random upper or lower Cholesky factorization.
More...

Detailed Description

Return a random upper or lower Cholesky factorization.

See the documentation of pm_distChol 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]subset: The input scalar constant that can be any of the following:
  1. the constant uppDia or an object of type uppDia_type implying that the upper-diagonal triangular block of the argument rand must be used while the lower subset is not referenced.
  2. the constant lowDia or an object of type lowDia_type implying that the lower-diagonal triangular block of the argument rand must be used while the upper subset is not referenced.
This argument is merely a convenience to differentiate the different procedure functionalities within this generic interface.
Beware that the oppsite subset will not be set to zero explicitly and it is the user responsibility to ensure it.
If needed, the generic interface setMatInit can be used to set the complementary subset to zero.
The generic interface getCholRand can be used to generate random Cholesky factors whose complementary subset elements are all set to zero.


Possible calling interfaces

call setCholRand(rng, rand(1:ndim, 1:ndim), subset)
Return a random upper or lower Cholesky factorization.
This module contains classes and procedures for generating random upper or lower Cholesky factor tria...
Definition: pm_distChol.F90:43
Warning
The condition size(rand, 1) == size(rand, 2) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
See also
getCovRand
setCovRand
getCholRand
setCholRand


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_io, only: display_type
5 use pm_container, only: csp_type
6 use pm_distUnif, only: getUnifRand
9 use pm_distChol, only: setCholRand, uppDia, lowDia, subset_type
11
12 implicit none
13
14 integer(IK) :: isub, itry, ndim
15
16 type(csp_type) :: subset(2)
17
18 type(display_type) :: disp
19 disp = display_type(file = "main.out.F90")
20
21 subset = [csp_type(uppDia), csp_type(lowDia)]
22 do isub = 1, size(subset)
23
24 block
25 use pm_kind, only: TKG => RKS ! all real kinds are supported.
26 real(TKG), allocatable :: rand(:,:), cov(:,:)
27 real(TKG) :: mold
28 do itry = 1, 5
29 call disp%skip()
30 call disp%show("ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])")
31 ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
32 call disp%show("[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]")
33 call disp%show( [same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)] )
34 call disp%show("call setCholRand(rngf_type(), rand, subset(isub)%val)")
35 call setCholRand(rngf_type(), rand, subset(isub)%val)
36 call disp%show("rand")
37 call disp%show( rand )
38 call disp%show("cov = matmul(rand, transpose(rand))")
39 cov = matmul(rand, transpose(rand))
40 call disp%show("isMatClass(cov, posdefmat)")
41 call disp%show( isMatClass(cov, posdefmat) )
42 call disp%skip()
43 end do
44 end block
45
46 block
47 use pm_kind, only: TKG => CKS ! all complex kinds are supported.
48 complex(TKG), allocatable :: rand(:,:), cov(:,:)
49 complex(TKG) :: mold
50 do itry = 1, 5
51 call disp%skip()
52 call disp%show("ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])")
53 ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
54 call disp%show("[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]")
55 call disp%show( [same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)] )
56 call disp%show("call setCholRand(rngf_type(), rand, subset(isub)%val)")
57 call setCholRand(rngf_type(), rand, subset(isub)%val)
58 call disp%show("rand")
59 call disp%show( rand )
60 call disp%show("cov = matmul(rand, conjg(transpose(rand)))")
61 cov = matmul(rand, conjg(transpose(rand)))
62 call disp%show("isMatClass(cov, posdefmat)")
63 call disp%show( isMatClass(cov, posdefmat) )
64 call disp%skip()
65 end do
66 end block
67
68 end do
69
70end program example
Allocate or resize (shrink or expand) and refill an input allocatable scalar string or array of rank ...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11726
This is a generic method of the derived type display_type with pass attribute.
Definition: pm_io.F90:11508
Generate and return .true. if and only if the input matrix is of the specified input class.
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains the derived types for generating allocatable containers of scalar,...
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 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...
This is the csp_type type for generating instances of container of scalar unlimited polymorphic objec...
This is a concrete derived type whose instances can be used to define/request the default uniform ran...
This is the derived type for declaring and generating objects of type xoshiro256ssw_type containing a...
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
2ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
3[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
4T, F
5call setCholRand(rngf_type(), rand, subset(isub)%val)
6rand
7+0.999999940, +0.741289437, -0.874833941, +0.348693103
8+0.00000000, +0.671185493, -0.224786177, +0.774243116
9+0.00000000, +0.00000000, +0.429111779, -0.527429879
10+0.00000000, +0.00000000, +0.00000000, +0.279006436E-1
11cov = matmul(rand, transpose(rand))
13T
14
15
16ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
17[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
18T, F
19call setCholRand(rngf_type(), rand, subset(isub)%val)
20rand
21+0.999999940, -0.980810642, -0.720449030, -0.565279305
22+0.00000000, +0.194962636, +0.675444663, -0.454110682
23+0.00000000, +0.00000000, +0.157250702, -0.176420242
24+0.00000000, +0.00000000, +0.00000000, +0.665671468
25cov = matmul(rand, transpose(rand))
27T
28
29
30ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
31[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
32T, F
33call setCholRand(rngf_type(), rand, subset(isub)%val)
34rand
35+1.00000000, -0.367873430, -0.665421247
36+0.00000000, +0.929875910, +0.844389126E-1
37+0.00000000, +0.00000000, +0.741676807
38cov = matmul(rand, transpose(rand))
40T
41
42
43ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
44[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
45T, F
46call setCholRand(rngf_type(), rand, subset(isub)%val)
47rand
48+0.999999940, +0.521290064, +0.153810203, +0.622930527, -0.174366191
49+0.00000000, +0.853379548, -0.976830840, -0.214354336, -0.766758204
50+0.00000000, +0.00000000, +0.148808330, +0.251696557, +0.104969434
51+0.00000000, +0.00000000, +0.00000000, +0.708984256, +0.243308783
52+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.558086693
53cov = matmul(rand, transpose(rand))
55T
56
57
58ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
59[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
60T, F
61call setCholRand(rngf_type(), rand, subset(isub)%val)
62rand
63+1.00000000, -0.541476905E-1, +0.506222725
64+0.00000000, +0.998532891, +0.532744527
65+0.00000000, +0.00000000, +0.678175390
66cov = matmul(rand, transpose(rand))
68T
69
70
71ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
72[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
73T, F
74call setCholRand(rngf_type(), rand, subset(isub)%val)
75rand
76(+1.00000000, +0.00000000), (+0.415614516, -0.748362839), (-0.393448710, -0.265420467), (+0.120104358, +0.403686315), (-0.181095943, +0.119459532)
77(+0.00000000, +0.00000000), (+0.516930997, +0.00000000), (+0.238732874, -0.844059348), (+0.316377282, -0.288823009), (+0.286252260, -0.434528172)
78(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.729428753E-1, +0.00000000), (-0.260494381, -0.483011931), (-0.400578022, -0.494320005)
79(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.581327081, +0.00000000), (-0.204525322, +0.372404963)
80(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.311203271, +0.00000000)
81cov = matmul(rand, conjg(transpose(rand)))
83T
84
85
86ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
87[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
88T, F
89call setCholRand(rngf_type(), rand, subset(isub)%val)
90rand
91(+0.999999940, +0.00000000), (+0.334222138, -0.448117018)
92(+0.00000000, +0.00000000), (+0.829148233, +0.00000000)
93cov = matmul(rand, conjg(transpose(rand)))
95T
96
97
98ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
99[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
100T, F
101call setCholRand(rngf_type(), rand, subset(isub)%val)
102rand
103(+0.999999940, +0.00000000), (+0.455087647E-1, +0.660011113), (-0.696164489, +0.365032077), (+0.214418247, +0.378500670), (-0.736810938E-1, +0.954043567E-1)
104(+0.00000000, +0.00000000), (+0.749876142, +0.00000000), (+0.195506305, -0.432707906), (-0.530280590, +0.227295384E-1), (-0.342573732, -0.452991724E-1)
105(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.395787328, +0.00000000), (-0.361173064, -0.442054302), (-0.400152415, +0.424528867)
106(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.450765908, +0.00000000), (-0.268666625, -0.451288968)
107(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.499870151, +0.00000000)
108cov = matmul(rand, conjg(transpose(rand)))
110T
111
112
113ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
114[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
115T, F
116call setCholRand(rngf_type(), rand, subset(isub)%val)
117rand
118(+0.999999940, +0.00000000), (-0.238311827, +0.255759537), (+0.482432693, +0.435055673), (-0.203566730, +0.435889632)
119(+0.00000000, +0.00000000), (+0.936906934, +0.00000000), (+0.577356458, +0.387609452E-1), (-0.584613323, -0.281243771)
120(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.493094623, +0.00000000), (-0.396584541, +0.118936852)
121(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.419838935, +0.00000000)
122cov = matmul(rand, conjg(transpose(rand)))
124T
125
126
127ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
128[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
129T, F
130call setCholRand(rngf_type(), rand, subset(isub)%val)
131rand
132(+1.00000000, +0.00000000), (+0.538379490, +0.646714270)
133(+0.00000000, +0.00000000), (+0.540285170, +0.00000000)
134cov = matmul(rand, conjg(transpose(rand)))
136T
137
138
139ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
140[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
141F, T
142call setCholRand(rngf_type(), rand, subset(isub)%val)
143rand
144+1.00000000, +0.00000000
145-0.711827457, +0.702354372
146cov = matmul(rand, transpose(rand))
148T
149
150
151ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
152[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
153F, T
154call setCholRand(rngf_type(), rand, subset(isub)%val)
155rand
156+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
157+0.953756988, +0.300578773, +0.00000000, +0.00000000, +0.00000000
158+0.148063794, -0.667945087, +0.729332864, +0.00000000, +0.00000000
159+0.694149494, -0.186891332, +0.399936326E-1, +0.693994701, +0.00000000
160+0.587197185, +0.158222482, +0.107513621, -0.343557686, +0.707512617
161cov = matmul(rand, transpose(rand))
163T
164
165
166ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
167[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
168F, T
169call setCholRand(rngf_type(), rand, subset(isub)%val)
170rand
171+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
172-0.925177991, +0.379533559, +0.00000000, +0.00000000, +0.00000000
173-0.243577793, -0.561038733, +0.791141808, +0.00000000, +0.00000000
174+0.102712572, +0.204309016, -0.755814612, +0.613557041, +0.00000000
175-0.300667614, +0.700518861E-1, +0.850727618, +0.863466114E-1, +0.416531503
176cov = matmul(rand, transpose(rand))
178T
179
180
181ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
182[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
183F, T
184call setCholRand(rngf_type(), rand, subset(isub)%val)
185rand
186+1.00000000, +0.00000000
187+0.895340621, +0.445381880
188cov = matmul(rand, transpose(rand))
190T
191
192
193ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
194[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
195F, T
196call setCholRand(rngf_type(), rand, subset(isub)%val)
197rand
198+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
199+0.611669183, +0.791113615, +0.00000000, +0.00000000, +0.00000000
200-0.666761160, +0.609569550, +0.428782642, +0.00000000, +0.00000000
201+0.281817079, +0.708146036, -0.594973385, +0.255176514, +0.00000000
202-0.248658031, -0.548592031, +0.242300376, -0.347132474, +0.676761031
203cov = matmul(rand, transpose(rand))
205T
206
207
208ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
209[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
210F, T
211call setCholRand(rngf_type(), rand, subset(isub)%val)
212rand
213(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
214(+0.939479172, -0.342443109), (+0.105548743E-1, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
215(-0.445869356, -0.152105018), (-0.530295670, -0.270102650), (+0.651072681, +0.00000000), (+0.00000000, +0.00000000)
216(-0.579761147, -0.359207660), (-0.262124777, +0.384900808), (-0.351156861, +0.180376336), (+0.402668685, +0.00000000)
217cov = matmul(rand, conjg(transpose(rand)))
219T
220
221
222ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
223[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
224F, T
225call setCholRand(rngf_type(), rand, subset(isub)%val)
226rand
227(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
228(+0.432393342, -0.768933713), (+0.470931828, +0.00000000), (+0.00000000, +0.00000000)
229(+0.477651507, +0.264958531), (-0.570802927, -0.421097636), (+0.445541084, +0.00000000)
230cov = matmul(rand, conjg(transpose(rand)))
232T
233
234
235ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
236[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
237F, T
238call setCholRand(rngf_type(), rand, subset(isub)%val)
239rand
240(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
241(+0.880752131E-1, +0.105935186), (+0.990464807, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
242(-0.529297233, -0.333430290), (-0.288646966, -0.573459327), (+0.443278670, +0.00000000), (+0.00000000, +0.00000000)
243(-0.286162764, +0.258986533), (-0.233083770, -0.316485256), (-0.652252138, +0.472353339), (+0.219078735, +0.00000000)
244cov = matmul(rand, conjg(transpose(rand)))
246T
247
248
249ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
250[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
251F, T
252call setCholRand(rngf_type(), rand, subset(isub)%val)
253rand
254(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
255(-0.596251227E-1, -0.588439286), (+0.806339979, +0.00000000), (+0.00000000, +0.00000000)
256(+0.685182393, +0.257756650), (-0.758722574E-1, -0.616509914), (+0.279723972, +0.00000000)
257cov = matmul(rand, conjg(transpose(rand)))
259T
260
261
262ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
263[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
264F, T
265call setCholRand(rngf_type(), rand, subset(isub)%val)
266rand
267(+1.00000000, +0.00000000), (+0.00000000, +0.00000000)
268(+0.829807639, +0.499466062), (+0.248903215, +0.00000000)
269cov = matmul(rand, conjg(transpose(rand)))
271T
272
273
Test:
test_pm_distChol


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 333 of file pm_distChol.F90.


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