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

Generate and return a random upper and lower Cholesky factorization.
More...

Detailed Description

Generate and return a random upper and lower Cholesky factorization.

See the documentation of pm_distChol for details.
See also setCholRand for generating only the upper or lower subset of the Cholesky factor.

Parameters
[in]mold: The input scalar 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),
whose type and kind determines the type and kind of the output rand.
The value of mold is ignored entirely within the algorithm.
[in]ndim: The input positive scalar of type integer of default kind IK, representing the rank of the matrix (the number of dimensions) of shape (ndim, ndim).
[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's responsibility to ensure it.
The generic interface getCholRand can be used to generate random Cholesky factors whose complementary subset elements are all set to zero.
Returns
rand : The output matrix of shape (1:ndim, 1:ndim) of the same type and kind as the input argument mold, containing an upper or lower Cholesky factorization as specified by the input subset.


Possible calling interfaces

rand(1:ndim, 1:ndim) = getCholRand(mold, ndim, subset)
Generate and return a random upper and lower Cholesky factorization.
This module contains classes and procedures for generating random upper or lower Cholesky factor tria...
Definition: pm_distChol.F90:43
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
8 use pm_distChol, only: getCholRand, uppDia, lowDia, subset_type
9
10 implicit none
11
12 integer(IK) :: isub, itry, ndim
13
14 type(csp_type) :: subset(2)
15
16 type(display_type) :: disp
17 disp = display_type(file = "main.out.F90")
18
19 subset = [csp_type(uppDia), csp_type(lowDia)]
20 do isub = 1, size(subset)
21
22 block
23 use pm_kind, only: TKG => RKS ! all real kinds are supported.
24 real(TKG), allocatable :: rand(:,:), cov(:,:)
25 real(TKG) :: mold
26 do itry = 1, 5
27 call disp%skip()
28 call disp%show("ndim = getUnifRand(2, 5)")
29 ndim = getUnifRand(2, 5)
30 call disp%show("[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]")
31 call disp%show( [same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)] )
32 call disp%show("rand = getCholRand(mold, ndim, subset(isub)%val)")
33 rand = getCholRand(mold, ndim, subset(isub)%val)
34 call disp%show("rand")
35 call disp%show( rand )
36 call disp%show("cov = matmul(rand, transpose(rand))")
37 cov = matmul(rand, transpose(rand))
38 call disp%show("isMatClass(cov, posdefmat)")
39 call disp%show( isMatClass(cov, posdefmat) )
40 call disp%skip()
41 end do
42 end block
43
44 block
45 use pm_kind, only: TKG => CKS ! all complex kinds are supported.
46 complex(TKG), allocatable :: rand(:,:), cov(:,:)
47 complex(TKG) :: mold
48 do itry = 1, 5
49 call disp%skip()
50 call disp%show("ndim = getUnifRand(2, 5)")
51 ndim = getUnifRand(2, 5)
52 call disp%show("[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]")
53 call disp%show( [same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)] )
54 call disp%show("rand = getCholRand(mold, ndim, subset(isub)%val)")
55 rand = getCholRand(mold, ndim, subset(isub)%val)
56 call disp%show("rand")
57 call disp%show( rand )
58 call disp%show("cov = matmul(rand, conjg(transpose(rand)))")
59 cov = matmul(rand, conjg(transpose(rand)))
60 call disp%show("isMatClass(cov, posdefmat)")
61 call disp%show( isMatClass(cov, posdefmat) )
62 call disp%skip()
63 end do
64 end block
65
66 end do
67
68end program example
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 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...
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)
3[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
4T, F
5rand = getCholRand(mold, ndim, subset(isub)%val)
6rand
7+1.00000000, -0.881023169
8+0.00000000, +0.473073095
9cov = matmul(rand, transpose(rand))
11T
12
13
14ndim = getUnifRand(2, 5)
15[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
16T, F
17rand = getCholRand(mold, ndim, subset(isub)%val)
18rand
19+0.999999940, +0.817078769, +0.714485884
20+0.00000000, +0.576526105, -0.610435367
21+0.00000000, +0.00000000, +0.341875136
22cov = matmul(rand, transpose(rand))
24T
25
26
27ndim = getUnifRand(2, 5)
28[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
29T, F
30rand = getCholRand(mold, ndim, subset(isub)%val)
31rand
32+0.999999940, -0.124576844, +0.163119212, +0.563983440
33+0.00000000, +0.992209971, -0.714186847, -0.632050633
34+0.00000000, +0.00000000, +0.680683017, +0.429951817
35+0.00000000, +0.00000000, +0.00000000, +0.312371969
36cov = matmul(rand, transpose(rand))
38T
39
40
41ndim = getUnifRand(2, 5)
42[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
43T, F
44rand = getCholRand(mold, ndim, subset(isub)%val)
45rand
46+1.00000000, +0.743592322, +0.577270150
47+0.00000000, +0.668633282, -0.591106474
48+0.00000000, +0.00000000, +0.563340366
49cov = matmul(rand, transpose(rand))
51T
52
53
54ndim = getUnifRand(2, 5)
55[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
56T, F
57rand = getCholRand(mold, ndim, subset(isub)%val)
58rand
59+1.00000000, +0.828519404
60+0.00000000, +0.559960365
61cov = matmul(rand, transpose(rand))
63T
64
65
66ndim = getUnifRand(2, 5)
67[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
68T, F
69rand = getCholRand(mold, ndim, subset(isub)%val)
70rand
71(+0.999999940, +0.00000000), (-0.520446181, +0.477501661), (-0.382391781, -0.138695091)
72(+0.00000000, +0.00000000), (+0.707903922, +0.00000000), (+0.682239234, +0.607210934)
73(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.196124483E-1, +0.00000000)
74cov = matmul(rand, conjg(transpose(rand)))
76T
77
78
79ndim = getUnifRand(2, 5)
80[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
81T, F
82rand = getCholRand(mold, ndim, subset(isub)%val)
83rand
84(+1.00000000, +0.00000000), (+0.611339629, -0.231690139), (-0.287896872, +0.634199977), (-0.919439420E-1, +0.278666586)
85(+0.00000000, +0.00000000), (+0.756692469, +0.00000000), (-0.151369691, -0.415325433), (-0.242105961, -0.617669582)
86(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.565241218, +0.00000000), (-0.142228976, +0.507636487)
87(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.442534059, +0.00000000)
88cov = matmul(rand, conjg(transpose(rand)))
90T
91
92
93ndim = getUnifRand(2, 5)
94[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
95T, F
96rand = getCholRand(mold, ndim, subset(isub)%val)
97rand
98(+1.00000000, +0.00000000), (+0.496055275, +0.749841630), (+0.668322742, +0.692601919)
99(+0.00000000, +0.00000000), (+0.437797576, +0.00000000), (-0.885681659E-1, -0.130309075)
100(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.220958278, +0.00000000)
101cov = matmul(rand, conjg(transpose(rand)))
103T
104
105
106ndim = getUnifRand(2, 5)
107[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
108T, F
109rand = getCholRand(mold, ndim, subset(isub)%val)
110rand
111(+1.00000000, +0.00000000), (-0.714058280, +0.455060363), (-0.281533003, -0.665991127), (-0.269260675, -0.557861149)
112(+0.00000000, +0.00000000), (+0.532015741, +0.00000000), (+0.390420519E-1, -0.617607415), (-0.153113753, +0.316167861)
113(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.306971997, +0.00000000), (-0.695827603E-1, +0.626029849)
114(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.310045928, +0.00000000)
115cov = matmul(rand, conjg(transpose(rand)))
117T
118
119
120ndim = getUnifRand(2, 5)
121[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
122T, F
123rand = getCholRand(mold, ndim, subset(isub)%val)
124rand
125(+1.00000000, +0.00000000), (-0.676116884, +0.316900462), (-0.443609744, -0.597368777E-1), (+0.129860222, +0.233461484)
126(+0.00000000, +0.00000000), (+0.665161729, +0.00000000), (+0.613545835, -0.554184437), (+0.577609539, +0.477946430)
127(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.340709537, +0.00000000), (+0.494210988, -0.248116046)
128(+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.246496215, +0.00000000)
129cov = matmul(rand, conjg(transpose(rand)))
131T
132
133
134ndim = getUnifRand(2, 5)
135[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
136F, T
137rand = getCholRand(mold, ndim, subset(isub)%val)
138rand
139+0.999999940, +0.00000000, +0.00000000, +0.00000000
140-0.680799782, +0.732469678, +0.00000000, +0.00000000
141+0.487655913E-2, -0.992172182, +0.124781907, +0.00000000
142-0.617822766, -0.621624708, -0.465882331, +0.121783979
143cov = matmul(rand, transpose(rand))
145T
146
147
148ndim = getUnifRand(2, 5)
149[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
150F, T
151rand = getCholRand(mold, ndim, subset(isub)%val)
152rand
153+1.00000000, +0.00000000, +0.00000000
154-0.715250671, +0.698867977, +0.00000000
155+0.590358675E-1, +0.722674131, +0.688663065
156cov = matmul(rand, transpose(rand))
158T
159
160
161ndim = getUnifRand(2, 5)
162[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
163F, T
164rand = getCholRand(mold, ndim, subset(isub)%val)
165rand
166+1.00000000, +0.00000000, +0.00000000, +0.00000000, +0.00000000
167+0.567623377, +0.823288262, +0.00000000, +0.00000000, +0.00000000
168-0.681979418, +0.705422103, +0.193089724, +0.00000000, +0.00000000
169-0.615866482, -0.369536993E-2, -0.622364283, +0.483070850, +0.00000000
170-0.364779942E-1, +0.718818843, -0.474043071, -0.497395724, +0.992440805E-1
171cov = matmul(rand, transpose(rand))
173T
174
175
176ndim = getUnifRand(2, 5)
177[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
178F, T
179rand = getCholRand(mold, ndim, subset(isub)%val)
180rand
181+1.00000000, +0.00000000
182+0.824148953, +0.566373169
183cov = matmul(rand, transpose(rand))
185T
186
187
188ndim = getUnifRand(2, 5)
189[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
190F, T
191rand = getCholRand(mold, ndim, subset(isub)%val)
192rand
193+1.00000000, +0.00000000
194-0.567899168, +0.823098123
195cov = matmul(rand, transpose(rand))
197T
198
199
200ndim = getUnifRand(2, 5)
201[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
202F, T
203rand = getCholRand(mold, ndim, subset(isub)%val)
204rand
205(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
206(-0.786486566, +0.611937344), (+0.834957361E-1, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
207(-0.186907396, +0.594342649), (+0.381694376, -0.608589947), (+0.309435040, +0.00000000), (+0.00000000, +0.00000000)
208(+0.561687946E-1, +0.587372124), (-0.135160536, -0.510476828), (-0.522009373, -0.278237522), (+0.151902333, +0.00000000)
209cov = matmul(rand, conjg(transpose(rand)))
211T
212
213
214ndim = getUnifRand(2, 5)
215[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
216F, T
217rand = getCholRand(mold, ndim, subset(isub)%val)
218rand
219(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
220(+0.896895766, -0.366310775), (+0.247779086, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
221(-0.493892729, +0.344038874), (+0.506113350, -0.548750281), (+0.283601135, +0.00000000), (+0.00000000, +0.00000000)
222(-0.812109947, +0.357648045), (+0.169093058, +0.121172853E-1), (-0.201017782, +0.311365485), (+0.215567470, +0.00000000)
223cov = matmul(rand, conjg(transpose(rand)))
225T
226
227
228ndim = getUnifRand(2, 5)
229[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
230F, T
231rand = getCholRand(mold, ndim, subset(isub)%val)
232rand
233(+0.999999940, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
234(-0.457254231, +0.240940955E-1), (+0.889009595, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
235(-0.924564749E-1, -0.186851099), (+0.610363960, +0.212777480), (+0.733975410, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
236(+0.463988870, -0.134514228), (-0.789870858, -0.660432354E-1), (+0.325407758E-1, +0.361292481), (+0.822893009E-1, +0.00000000), (+0.00000000, +0.00000000)
237(-0.414330006, -0.219350517), (-0.279329479, +0.417469144), (-0.252942741, -0.554553509), (-0.105947994, -0.278328896), (+0.260209382, +0.00000000)
238cov = matmul(rand, conjg(transpose(rand)))
240T
241
242
243ndim = getUnifRand(2, 5)
244[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
245F, T
246rand = getCholRand(mold, ndim, subset(isub)%val)
247rand
248(+1.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
249(+0.702952445, +0.234808549), (+0.671358943, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
250(-0.470153242, -0.697452128), (+0.279388547, +0.461006135), (+0.439535193E-1, +0.00000000), (+0.00000000, +0.00000000)
251(+0.508090436, -0.413282245), (+0.547055066, -0.206316561), (-0.109405583E-1, -0.475189298), (+0.572854355E-1, +0.00000000)
252cov = matmul(rand, conjg(transpose(rand)))
254T
255
256
257ndim = getUnifRand(2, 5)
258[same_type_as(subset(isub)%val, uppDia), same_type_as(subset(isub)%val, lowDia)]
259F, T
260rand = getCholRand(mold, ndim, subset(isub)%val)
261rand
262(+0.999999940, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
263(-0.580488443, -0.754666865), (+0.305795848, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
264(+0.357133225E-1, +0.501023114), (+0.422885835, -0.338838935), (+0.673836887, +0.00000000), (+0.00000000, +0.00000000), (+0.00000000, +0.00000000)
265(+0.440433145, +0.243053511), (+0.611537874, -0.316398859), (+0.416826010E-1, -0.411656708), (+0.318838805, +0.00000000), (+0.00000000, +0.00000000)
266(-0.329752773, +0.272081554), (-0.218738332, +0.407385826), (-0.300155897E-1, +0.506982386), (+0.470104843, +0.352798492), (+0.525293825E-2, +0.00000000)
267cov = matmul(rand, conjg(transpose(rand)))
269T
270
271
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 115 of file pm_distChol.F90.


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