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.425585151
8+0.00000000, +0.904918373
9cov = matmul(rand, transpose(rand))
11T
12
13
14ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
15[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
16T, F
17call setCholRand(rngf_type(), rand, subset(isub)%val)
18rand
19+1.00000000, -0.662149906, +0.666454554
20+0.00000000, +0.749371350, -0.686088085
21+0.00000000, +0.00000000, +0.291755944
22cov = matmul(rand, transpose(rand))
24T
25
26
27ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
28[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
29T, F
30call setCholRand(rngf_type(), rand, subset(isub)%val)
31rand
32+1.00000000, +0.914508760
33+0.00000000, +0.404565990
34cov = matmul(rand, transpose(rand))
36T
37
38
39ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
40[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
41T, F
42call setCholRand(rngf_type(), rand, subset(isub)%val)
43rand
44+1.00000000, +0.652417004, -0.465997756
45+0.00000000, +0.757860184, -0.564494371
46+0.00000000, +0.00000000, +0.681316555
47cov = matmul(rand, transpose(rand))
49T
50
51
52ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
53[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
54T, F
55call setCholRand(rngf_type(), rand, subset(isub)%val)
56rand
57+1.00000000, -0.919773340, -0.524613142, -0.591849744, +0.615483634E-1
58+0.00000000, +0.392450124, -0.758221209, -0.591785192, +0.590273857
59+0.00000000, +0.00000000, +0.387145370, -0.237296417, -0.559922040
60+0.00000000, +0.00000000, +0.00000000, +0.493147492, +0.518885195
61+0.00000000, +0.00000000, +0.00000000, +0.00000000, +0.255017906
62cov = matmul(rand, transpose(rand))
64T
65
66
67ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
68[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
69T, F
70call setCholRand(rngf_type(), rand, subset(isub)%val)
71rand
72(+1.00000000, +0.00000000), (+0.958259225, +0.115046669E-1), (-0.142230794, -0.833119005E-1), (-0.159731761, -0.427446187), (-0.475065522E-1, -0.358440697)
73(+0.00000000, +0.00000000), (+0.285669208, +0.00000000), (-0.600978553, -0.777166784), (-0.395708233, -0.306872666), (-0.268540740, -0.266245663)
74(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.875554606E-1, +0.00000000), (-0.415262729, -0.544336081), (-0.508443356, +0.474356681)
75(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.268840104, +0.00000000), (-0.996307805E-1, -0.475130171)
76(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.840154290E-1, +0.00000000)
77cov = matmul(rand, conjg(transpose(rand)))
79T
80
81
82ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
83[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
84T, F
85call setCholRand(rngf_type(), rand, subset(isub)%val)
86rand
87(+1.00000000, +0.00000000), (-0.349207632E-1, +0.805915356), (+0.429173142, +0.210942030), (+0.370629251, +0.231197536)
88(+0.00000000, +0.00000000), (+0.591000021, +0.00000000), (+0.510435402, -0.319310158), (-0.595185578, +0.325467557)
89(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.639383018, +0.00000000), (+0.517262638, -0.169772238)
90(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.229398057, +0.00000000)
91cov = matmul(rand, conjg(transpose(rand)))
93T
94
95
96ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
97[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
98T, F
99call setCholRand(rngf_type(), rand, subset(isub)%val)
100rand
101(+1.00000000, +0.00000000), (+0.844284773, +0.240977585), (+0.269895315, -0.491238922)
102(+0.00000000, +0.00000000), (+0.478657454, +0.00000000), (-0.791895807, +0.123631977)
103(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.208463624, +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.546078980, -0.511225283), (+0.360939860, -0.832776576E-1), (+0.595321119, -0.523544252), (-0.365521461, +0.351287991)
115(+0.00000000, +0.00000000), (+0.663661420, +0.00000000), (+0.791177332, +0.465306759), (-0.474569947, -0.198456898), (+0.472455084, +0.319023192)
116(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.142532006, +0.00000000), (-0.704541877E-1, -0.791039541E-1), (-0.423498988, -0.113260699E-1)
117(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.309307754, +0.00000000), (+0.832130238E-1, +0.413113952)
118(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.246848002, +0.00000000)
119cov = matmul(rand, conjg(transpose(rand)))
121T
122
123
124ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
125[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
126T, F
127call setCholRand(rngf_type(), rand, subset(isub)%val)
128rand
129(+0.999999940, +0.00000000), (+0.665597796, +0.744406939)
130(+0.00000000, +0.00000000), (+0.532718636E-1, +0.00000000)
131cov = matmul(rand, conjg(transpose(rand)))
133T
134
135
136ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
137[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
138F, T
139call setCholRand(rngf_type(), rand, subset(isub)%val)
140rand
141+1.00000000, +0.00000000, +0.00000000
142+0.233753592, +0.972295880, +0.00000000
143+0.831458628, +0.191351444, +0.521594822
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+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
155-0.965386093, +0.260825127, +0.00000000, +0.00000000, +0.00000000
156+0.967681646, +0.178959623, +0.177667439, +0.00000000, +0.00000000
157-0.590900302, -0.219072565, +0.714601517, +0.303625673, +0.00000000
158+0.528525233, +0.111198761, +0.708474457, +0.284662515, +0.354015648
159cov = matmul(rand, transpose(rand))
161T
162
163
164ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
165[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
166F, T
167call setCholRand(rngf_type(), rand, subset(isub)%val)
168rand
169+1.00000000, +0.00000000
170+0.617926300, +0.786235988
171cov = matmul(rand, transpose(rand))
173T
174
175
176ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
177[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
178F, T
179call setCholRand(rngf_type(), rand, subset(isub)%val)
180rand
181+1.00000000, +0.00000000
182-0.965245068, +0.261346161
183cov = matmul(rand, transpose(rand))
185T
186
187
188ndim = getUnifRand(2, 5); call setRefilled(rand, 0._TKG, [ndim, ndim])
189[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
190F, T
191call setCholRand(rngf_type(), rand, subset(isub)%val)
192rand
193+1.00000000, +0.00000000, +0.00000000
194-0.855897307, +0.517145753, +0.00000000
195+0.900649250, +0.676939264E-1, +0.429241717
196cov = matmul(rand, transpose(rand))
198T
199
200
201ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
202[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
203F, T
204call setCholRand(rngf_type(), rand, subset(isub)%val)
205rand
206(+1.00000000, +0.00000000), (+0.00000000, +0.00000000)
207(-0.179492041, +0.434517443), (+0.882596731, +0.00000000)
208cov = matmul(rand, conjg(transpose(rand)))
210T
211
212
213ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
214[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
215F, T
216call setCholRand(rngf_type(), rand, subset(isub)%val)
217rand
218(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
219(-0.588443458, -0.596750498), (+0.545548618, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
220(-0.195323572, +0.579052925), (+0.849264339E-1, -0.503517270), (+0.604817569, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
221(-0.285100907, -0.206219733), (+0.382310718, -0.461558312), (+0.606806874, +0.304751098), (+0.236443251, +0.00000000), (+0.00000000, +0.00000000)
222(-0.567227341E-1, +0.448033959), (+0.417219698, -0.282653451), (-0.278970748, -0.521920860), (-0.251283675, -0.328902185), (+0.143305883, +0.00000000)
223cov = matmul(rand, conjg(transpose(rand)))
225T
226
227
228ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
229[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
230F, T
231call setCholRand(rngf_type(), rand, subset(isub)%val)
232rand
233(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
234(+0.897334039, -0.182382748), (+0.401905507, +0.00000000), (+0.00000000, +0.00000000)
235(-0.135075852, +0.587263823), (+0.377067894, +0.471984863), (+0.521465182, +0.00000000)
236cov = matmul(rand, conjg(transpose(rand)))
238T
239
240
241ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
242[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
243F, T
244call setCholRand(rngf_type(), rand, subset(isub)%val)
245rand
246(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
247(-0.501975477, +0.860542715), (+0.865257531E-1, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
248(+0.498893529, -0.249854982), (+0.447148412, +0.486257166), (+0.502284825, +0.00000000), (+0.00000000, +0.00000000)
249(+0.141714454, -0.579694092), (-0.368996203, +0.364295363), (-0.582066536, -0.175259858), (+0.740604252E-1, +0.00000000)
250cov = matmul(rand, conjg(transpose(rand)))
252T
253
254
255ndim = getUnifRand(2, 5); call setRefilled(rand, cmplx(0, 0, TKG), [ndim, ndim])
256[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
257F, T
258call setCholRand(rngf_type(), rand, subset(isub)%val)
259rand
260(+1.00000000, +0.00000000), (+0.00000000, +0.00000000)
261(-0.742168576E-1, +0.686697185), (+0.723145127, +0.00000000)
262cov = matmul(rand, conjg(transpose(rand)))
264T
265
266
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: