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+1.00000000, -0.739765838E-1, -0.663845301, +0.694602549, +0.347815424
8+0.00000000, +0.997259974, +0.185639009, +0.656342387, +0.643875837
9+0.00000000, +0.00000000, +0.724463582, -0.181059301, +0.330593616
10+0.00000000, +0.00000000, +0.00000000, +0.232291520, -0.159986481
11+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.574073553
12cov = matmul(rand, transpose(rand))
14T
15
16
17ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
18[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
19T, F
20call setCholRand(rngf_type(), rand, subset(isub)%val)
21rand
22+1.00000000, +0.999812543
23+0.00000000, +0.193599667E-1
24cov = matmul(rand, transpose(rand))
26T
27
28
29ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
30[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
31T, F
32call setCholRand(rngf_type(), rand, subset(isub)%val)
33rand
34+1.00000000, -0.580451727
35+0.00000000, +0.814294577
36cov = matmul(rand, transpose(rand))
38T
39
40
41ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
42[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
43T, F
44call setCholRand(rngf_type(), rand, subset(isub)%val)
45rand
46+0.999999940, -0.941551685, -0.533253670, -0.493680149
47+0.00000000, +0.336868227, -0.142660692, -0.239754394
48+0.00000000, +0.00000000, +0.833839536, +0.739271641
49+0.00000000, +0.00000000, +0.00000000, +0.390224397
50cov = matmul(rand, transpose(rand))
52T
53
54
55ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
56[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
57T, F
58call setCholRand(rngf_type(), rand, subset(isub)%val)
59rand
60+1.00000000, -0.175480679, -0.720953822, -0.631816030
61+0.00000000, +0.984482825, -0.669096053, +0.480709612
62+0.00000000, +0.00000000, +0.180377632, -0.163244173
63+0.00000000, +0.00000000, +0.00000000, +0.585728705
64cov = matmul(rand, transpose(rand))
66T
67
68
69ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
70[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
71T, F
72call setCholRand(rngf_type(), rand, subset(isub)%val)
73rand
74(+1.00000000, +0.00000000), (+0.683468640, -0.640034854), (+0.163929924, -0.801959217), (+0.948456079E-1, +0.319611907)
75(+0.00000000, +0.00000000), (+0.351035774, +0.00000000), (+0.233231410, +0.235806003), (-0.527997971, -0.289958388)
76(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.469027847, +0.00000000), (+0.515764534, +0.190313280)
77(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.473035574, +0.00000000)
78cov = matmul(rand, conjg(transpose(rand)))
80T
81
82
83ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
84[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
85T, F
86call setCholRand(rngf_type(), rand, subset(isub)%val)
87rand
88(+1.00000000, +0.00000000), (-0.707968116, -0.612926343E-2)
89(+0.00000000, +0.00000000), (+0.706217706, +0.00000000)
90cov = matmul(rand, conjg(transpose(rand)))
92T
93
94
95ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
96[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
97T, F
98call setCholRand(rngf_type(), rand, subset(isub)%val)
99rand
100(+1.00000000, +0.00000000), (+0.836065412, -0.282887012), (+0.217805460, +0.542472482), (-0.441470653, +0.110357068)
101(+0.00000000, +0.00000000), (+0.470074028, +0.00000000), (-0.623225629, -0.977301523E-1), (+0.428905427, +0.623144031)
102(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.510218501, +0.00000000), (-0.392975032, -0.251899987)
103(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.526662990E-1, +0.00000000)
104cov = matmul(rand, conjg(transpose(rand)))
106T
107
108
109ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
110[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
111T, F
112call setCholRand(rngf_type(), rand, subset(isub)%val)
113rand
114(+1.00000000, +0.00000000), (-0.847072303, +0.488581508), (-0.612751722, +0.108179614)
115(+0.00000000, +0.00000000), (+0.209180802, +0.00000000), (-0.260590881, -0.616063237)
116(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.406682670, +0.00000000)
117cov = matmul(rand, conjg(transpose(rand)))
119T
120
121
122ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
123[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
124T, F
125call setCholRand(rngf_type(), rand, subset(isub)%val)
126rand
127(+1.00000000, +0.00000000), (-0.617743313, +0.685746253), (+0.765011385E-1, +0.198219433)
128(+0.00000000, +0.00000000), (+0.384896427, +0.00000000), (+0.699901998, +0.568517409E-1)
129(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.679530501, +0.00000000)
130cov = matmul(rand, conjg(transpose(rand)))
132T
133
134
135ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
136[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
137F, T
138call setCholRand(rngf_type(), rand, subset(isub)%val)
139rand
140+1.00000000, +0.00000000, +0.00000000, +0.00000000
141+0.957042992, +0.289946169, +0.00000000, +0.00000000
142-0.777617753, -0.225902706, +0.586752594, +0.00000000
143+0.629367307E-1, +0.405321062, +0.622657895, +0.666371465
144cov = matmul(rand, transpose(rand))
146T
147
148
149ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
150[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
151F, T
152call setCholRand(rngf_type(), rand, subset(isub)%val)
153rand
154+0.999999940, +0.00000000, +0.00000000, +0.00000000
155+0.861693799, +0.507428646, +0.00000000, +0.00000000
156+0.133946091E-1, -0.992824137, +0.118830971, +0.00000000
157-0.165356889, -0.370114148, -0.782581091, +0.472482204
158cov = matmul(rand, transpose(rand))
160T
161
162
163ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
164[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
165F, T
166call setCholRand(rngf_type(), rand, subset(isub)%val)
167rand
168+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
169+0.510989726, +0.859586895, +0.00000000, +0.00000000, +0.00000000
170+0.143555447, +0.780935705, +0.607890904, +0.00000000, +0.00000000
171+0.446217269, +0.544218242, -0.497979462, +0.506688297, +0.00000000
172-0.209729761, -0.174255729, -0.487680316, -0.521542549, +0.644832969
173cov = matmul(rand, transpose(rand))
175T
176
177
178ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
179[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
180F, T
181call setCholRand(rngf_type(), rand, subset(isub)%val)
182rand
183+1.00000000, +0.00000000
184+0.504680872, +0.863305926
185cov = matmul(rand, transpose(rand))
187T
188
189
190ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
191[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
192F, T
193call setCholRand(rngf_type(), rand, subset(isub)%val)
194rand
195+0.999999940, +0.00000000
196+0.217102304, +0.976148903
197cov = matmul(rand, transpose(rand))
199T
200
201
202ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
203[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
204F, T
205call setCholRand(rngf_type(), rand, subset(isub)%val)
206rand
207(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
208(+0.592895091, -0.576422691), (+0.562327504, +0.00000000), (+0.00000000, +0.00000000)
209(-0.717085958, -0.520690799), (+0.415858775, -0.202744886), (+0.249974467E-1, +0.00000000)
210cov = matmul(rand, conjg(transpose(rand)))
212T
213
214
215ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
216[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
217F, T
218call setCholRand(rngf_type(), rand, subset(isub)%val)
219rand
220(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
221(-0.713647187, +0.524680376), (+0.464131683, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
222(-0.409436882, -0.500087082), (+0.389481157, -0.637860477), (+0.153989613, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
223(-0.267093211, +0.660974145), (+0.228817374, +0.459787011), (-0.157180950, +0.285677612), (+0.348848701, +0.00000000), (+0.00000000, +0.00000000)
224(+0.740294531E-1, +0.402090669), (+0.503392518, -0.259241164), (+0.275559932, -0.950609520E-1), (-0.330883950, +0.284531057), (+0.486642212, +0.00000000)
225cov = matmul(rand, conjg(transpose(rand)))
227T
228
229
230ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
231[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
232F, T
233call setCholRand(rngf_type(), rand, subset(isub)%val)
234rand
235(+1.00000000, +0.00000000), (+0.00000000, +0.00000000)
236(+0.946827710, -0.318956882), (+0.422331989E-1, +0.00000000)
237cov = matmul(rand, conjg(transpose(rand)))
239T
240
241
242ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
243[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
244F, T
245call setCholRand(rngf_type(), rand, subset(isub)%val)
246rand
247(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
248(-0.524434507, -0.511814952), (+0.680451214, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
249(-0.512067556, +0.550266445), (-0.110606804, -0.379096180), (+0.528248072, +0.00000000), (+0.00000000, +0.00000000)
250(-0.562052071, +0.933619067E-1), (+0.279817462, -0.582947373), (+0.224070121E-1, +0.505284905), (+0.379551612E-1, +0.00000000)
251cov = matmul(rand, conjg(transpose(rand)))
253T
254
255
256ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
257[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
258F, T
259call setCholRand(rngf_type(), rand, subset(isub)%val)
260rand
261(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
262(+0.545466959, +0.489833951), (+0.680094481, +0.00000000), (+0.00000000, +0.00000000)
263(-0.715112507, -0.558512583E-1), (+0.300503499E-2, +0.691222250), (+0.877348855E-1, +0.00000000)
264cov = matmul(rand, conjg(transpose(rand)))
266T
267
268
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: