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

Generate and return a random positive-definite (correlation or covariance) matrix using the Gram method.
More...

Detailed Description

Generate and return a random positive-definite (correlation or covariance) matrix using the Gram method.

See the documentation of pm_distCov for details.
See also setCovRand for generating random covariance matrices using the method of Dvine or Onion.

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).
(optional. It must be present if and only if the input scale argument is missing or is a scalar.)
[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.. It can be present if and only if it is a scalar or, it is a vector and the input argument ndim is missing.)
Returns
rand : The output matrix of shape (1:ndim, 1:ndim) of the same type and kind as the input argument mold, containing a random positive-definite matrix.


Possible calling interfaces

rand(1:ndim, 1:ndim) = getCovRand(mold, scale(1:ndim))
rand(1:ndim, 1:ndim) = getCovRand(mold, ndim, scale = scale)
Generate and return a random positive-definite (correlation or covariance) matrix using the Gram meth...
Definition: pm_distCov.F90:394
This module contains classes and procedures for generating random matrices distributed on the space o...
Definition: pm_distCov.F90:72
Warning
The condition all([0 < scale]) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
Note
Unlike the case for setCovRand, when the input argument scale is missing, the diagonal elements of the output correlation matrix are strictly enforced to 1.


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_io, only: display_type
5 use pm_distUnif, only: getUnifRand
6 use pm_matrixChol, only: setChoLow
9 use pm_distCov, only: getCovRand
10
11 implicit none
12
13 integer(IK) :: itry, ndim
14
15 type(display_type) :: disp
16 disp = display_type(file = "main.out.F90")
17
18 block
19 use pm_kind, only: TKG => RKS ! all real kinds are supported.
20 real(TKG), allocatable :: rand(:,:)
21 real(TKG) :: mold
22 real(TKG) :: scale
23 do itry = 1, 5
24 call disp%skip()
25 call disp%show("ndim = getUnifRand(2, 3)")
26 ndim = getUnifRand(2, 3)
27 call disp%show("scale = getUnifRand(1, 10)")
28 scale = getUnifRand(1, 10)
29 call disp%show("rand = getCovRand(mold, ndim)")
30 rand = getCovRand(mold, ndim)
31 call disp%show("rand")
32 call disp%show( rand )
33 call disp%show("isMatClass(rand, posdefmat)")
34 call disp%show( isMatClass(rand, posdefmat) )
35 call disp%show("rand = getCovRand(mold, ndim, scale)")
36 rand = getCovRand(mold, ndim, scale)
37 call disp%show("rand")
38 call disp%show( rand )
39 call disp%show("isMatClass(rand, posdefmat)")
40 call disp%show( isMatClass(rand, posdefmat) )
41 call disp%skip()
42 call disp%skip()
43 call disp%show("rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])")
44 rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
45 call disp%show("rand")
46 call disp%show( rand )
47 call disp%show("isMatClass(rand, posdefmat)")
48 call disp%show( isMatClass(rand, posdefmat) )
49 call disp%skip()
50 end do
51 end block
52
53 block
54 use pm_kind, only: TKG => CKS ! all real kinds are supported.
55 complex(TKG), allocatable :: rand(:,:)
56 complex(TKG) :: mold
57 real(TKG) :: scale
58 do itry = 1, 5
59 call disp%skip()
60 call disp%show("ndim = getUnifRand(2, 3)")
61 ndim = getUnifRand(2, 3)
62 call disp%show("scale = getUnifRand(1, 10)")
63 scale = getUnifRand(1, 10)
64 call disp%show("rand = getCovRand(mold, ndim)")
65 rand = getCovRand(mold, ndim)
66 call disp%show("rand")
67 call disp%show( rand )
68 call disp%show("isMatClass(rand, posdefmat)")
69 call disp%show( isMatClass(rand, posdefmat) )
70 call disp%show("rand = getCovRand(mold, ndim, scale)")
71 rand = getCovRand(mold, ndim, scale)
72 call disp%show("rand")
73 call disp%show( rand )
74 call disp%show("isMatClass(rand, posdefmat)")
75 call disp%show( isMatClass(rand, posdefmat) )
76 call disp%skip()
77 call disp%skip()
78 call disp%show("rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])")
79 rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
80 call disp%show("rand")
81 call disp%show( rand )
82 call disp%show("isMatClass(rand, posdefmat)")
83 call disp%show( isMatClass(rand, posdefmat) )
84 call disp%skip()
85 end do
86 end block
87
88end program example
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.
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...
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 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...
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, 3)
3scale = getUnifRand(1, 10)
4rand = getCovRand(mold, ndim)
5rand
6+1.00000000, +0.747877121, +0.160131469
7+0.747877121, +1.00000000, +0.547551095
8+0.160131469, +0.547551095, +1.00000000
10T
11rand = getCovRand(mold, ndim, scale)
12rand
13+1.00000000, -0.886900425, -0.677867174
14-0.886900425, +0.999999940, +0.616202950
15-0.677867174, +0.616202950, +1.00000000
17T
18
19
20rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
21rand
22+1.00000000, -0.426910102, -2.56361914, -0.130611569, +1.83145535
23-0.426910102, +3.99999976, +2.21671295, +4.03206253, -5.38523674
24-2.56361914, +2.21671295, +8.99999905, -2.30833435, -10.3025818
25-0.130611569, +4.03206253, -2.30833435, +15.9999981, -1.73274398
26+1.83145535, -5.38523674, -10.3025818, -1.73274398, +24.9999981
28T
29
30
31ndim = getUnifRand(2, 3)
32scale = getUnifRand(1, 10)
33rand = getCovRand(mold, ndim)
34rand
35+1.00000000, +0.668615550E-1, +0.386919558
36+0.668615550E-1, +1.00000000, +0.849065125
37+0.386919558, +0.849065125, +1.00000000
39T
40rand = getCovRand(mold, ndim, scale)
41rand
42+100.000000, +88.8880768, -53.9552765
43+88.8880768, +100.000000, -83.4008560
44-53.9552765, -83.4008560, +99.9999924
46T
47
48
49rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
50rand
51+1.00000000, -1.37277853, +2.45538735, +2.68960500, -2.92528677
52-1.37277853, +4.00000048, -3.07490301, -1.12012267, +6.66605759
53+2.45538735, -3.07490301, +9.00000191, +10.9935493, -6.09769535
54+2.68960500, -1.12012267, +10.9935493, +16.0000019, -4.38999176
55-2.92528677, +6.66605759, -6.09769535, -4.38999176, +25.0000019
57T
58
59
60ndim = getUnifRand(2, 3)
61scale = getUnifRand(1, 10)
62rand = getCovRand(mold, ndim)
63rand
64+1.00000000, -0.740277469
65-0.740277469, +1.00000000
67T
68rand = getCovRand(mold, ndim, scale)
69rand
70+100.000000, -99.9758987
71-99.9758987, +100.000015
73T
74
75
76rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
77rand
78+1.00000000, +1.48557889, +0.309450291E-1, +0.289537907, +1.39877415
79+1.48557889, +4.00000048, +3.45683527, -1.21328533, -2.17302418
80+0.309450291E-1, +3.45683527, +8.99999809, -0.416234672, -12.7757549
81+0.289537907, -1.21328533, -0.416234672, +16.0000000, +5.77055836
82+1.39877415, -2.17302418, -12.7757549, +5.77055836, +25.0000019
84T
85
86
87ndim = getUnifRand(2, 3)
88scale = getUnifRand(1, 10)
89rand = getCovRand(mold, ndim)
90rand
91+1.00000000, +0.478561729, -0.924712777
92+0.478561729, +1.00000000, -0.108680874
93-0.924712777, -0.108680874, +1.00000000
95T
96rand = getCovRand(mold, ndim, scale)
97rand
98+1.00000000, -0.136964336, -0.133559898
99-0.136964336, +1.00000000, +0.507089138
100-0.133559898, +0.507089138, +1.00000012
102T
103
104
105rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
106rand
107+1.00000000, -0.352264196, -2.90221095, +2.06953335, -1.34628832
108-0.352264196, +3.99999976, +1.03480017, +5.92855883, +1.75670934
109-2.90221095, +1.03480017, +8.99999809, -5.64722633, +1.02246833
110+2.06953335, +5.92855883, -5.64722633, +15.9999990, -2.71537566
111-1.34628832, +1.75670934, +1.02246833, -2.71537566, +25.0000000
113T
114
115
116ndim = getUnifRand(2, 3)
117scale = getUnifRand(1, 10)
118rand = getCovRand(mold, ndim)
119rand
120+1.00000000, -0.535483778
121-0.535483778, +1.00000000
123T
124rand = getCovRand(mold, ndim, scale)
125rand
126+81.0000000, +13.6812983
127+13.6812983, +81.0000000
129T
130
131
132rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
133rand
134+0.999999881, +1.69643664, +1.22547245, -1.11216950, +0.643765867
135+1.69643664, +4.00000000, +2.37604666, +1.32404470, +4.66061640
136+1.22547245, +2.37604666, +9.00000191, -2.26288939, -7.84399509
137-1.11216950, +1.32404470, -2.26288939, +15.9999981, +9.65211678
138+0.643765867, +4.66061640, -7.84399509, +9.65211678, +25.0000000
140T
141
142
143ndim = getUnifRand(2, 3)
144scale = getUnifRand(1, 10)
145rand = getCovRand(mold, ndim)
146rand
147(+1.00000000, +0.00000000), (+0.558829546, +0.449122995), (+0.519566953, -0.178082741E-2)
148(+0.558829546, -0.449122995), (+1.00000000, +0.00000000), (-0.783119500E-1, -0.613214076)
149(+0.519566953, +0.178082741E-2), (-0.783119500E-1, +0.613214076), (+1.00000000, +0.00000000)
151T
152rand = getCovRand(mold, ndim, scale)
153rand
154(+9.00000000, +0.00000000), (+0.735334516, +1.42338598), (-2.23954988, +0.147947028)
155(+0.735334516, -1.42338598), (+8.99999905, +0.00000000), (+3.41638088, -5.19160604)
156(-2.23954988, -0.147947028), (+3.41638088, +5.19160604), (+9.00000191, +0.00000000)
158T
159
160
161rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
162rand
163(+1.00000000, +0.00000000), (-1.91282868, +0.00000000), (+1.98555136, +0.00000000), (+0.909577727, +0.00000000), (+2.48975229, +0.00000000)
164(-1.91282868, +0.00000000), (+4.00000000, +0.00000000), (-4.77456284, +0.00000000), (+0.155988216, +0.00000000), (-4.09898710, +0.00000000)
165(+1.98555136, +0.00000000), (-4.77456284, +0.00000000), (+9.00000191, +0.00000000), (-0.474937677, +0.00000000), (+7.30296326, +0.00000000)
166(+0.909577727, +0.00000000), (+0.155988216, +0.00000000), (-0.474937677, +0.00000000), (+16.0000038, +0.00000000), (+13.2090816, +0.00000000)
167(+2.48975229, +0.00000000), (-4.09898710, +0.00000000), (+7.30296326, +0.00000000), (+13.2090816, +0.00000000), (+25.0000019, +0.00000000)
169T
170
171
172ndim = getUnifRand(2, 3)
173scale = getUnifRand(1, 10)
174rand = getCovRand(mold, ndim)
175rand
176(+1.00000000, +0.00000000), (-0.753797948, -0.498124510)
177(-0.753797948, +0.498124510), (+1.00000000, +0.00000000)
179T
180rand = getCovRand(mold, ndim, scale)
181rand
182(+25.0000000, +0.00000000), (-2.21847892, +15.4024448)
183(-2.21847892, -15.4024448), (+25.0000038, +0.00000000)
185T
186
187
188rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
189rand
190(+1.00000000, +0.00000000), (+1.32659864, +0.00000000), (-1.65774441, +0.00000000), (+0.500658490E-1, +0.00000000), (-1.88713694, +0.00000000)
191(+1.32659864, +0.00000000), (+4.00000000, +0.00000000), (-5.59919357, +0.00000000), (-4.25570869, +0.00000000), (+0.291871727, +0.00000000)
192(-1.65774441, +0.00000000), (-5.59919357, +0.00000000), (+9.00000191, +0.00000000), (+4.36084747, +0.00000000), (+1.26752746, +0.00000000)
193(+0.500658490E-1, +0.00000000), (-4.25570869, +0.00000000), (+4.36084747, +0.00000000), (+15.9999971, +0.00000000), (-11.7452898, +0.00000000)
194(-1.88713694, +0.00000000), (+0.291871727, +0.00000000), (+1.26752746, +0.00000000), (-11.7452898, +0.00000000), (+24.9999981, +0.00000000)
196T
197
198
199ndim = getUnifRand(2, 3)
200scale = getUnifRand(1, 10)
201rand = getCovRand(mold, ndim)
202rand
203(+1.00000000, +0.00000000), (+0.624205112, -0.156064793), (+0.517979324, -0.753848493)
204(+0.624205112, +0.156064793), (+1.00000000, +0.00000000), (+0.457151234, -0.473889977)
205(+0.517979324, +0.753848493), (+0.457151234, +0.473889977), (+1.00000000, +0.00000000)
207T
208rand = getCovRand(mold, ndim, scale)
209rand
210(+100.000000, +0.00000000), (-71.9668274, -4.17559862), (-32.0589142, +16.0958099)
211(-71.9668274, +4.17559862), (+99.9999847, +0.00000000), (+43.8245125, -58.0102310)
212(-32.0589142, -16.0958099), (+43.8245125, +58.0102310), (+100.000000, +0.00000000)
214T
215
216
217rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
218rand
219(+1.00000000, +0.00000000), (+1.90523672, +0.00000000), (-2.11284924, +0.00000000), (-2.22213721, +0.00000000), (-2.89092064, +0.00000000)
220(+1.90523672, +0.00000000), (+4.00000048, +0.00000000), (-4.55563402, +0.00000000), (-5.55646610, +0.00000000), (-4.26904774, +0.00000000)
221(-2.11284924, +0.00000000), (-4.55563402, +0.00000000), (+9.00000000, +0.00000000), (+2.00845003, +0.00000000), (+8.07615852, +0.00000000)
222(-2.22213721, +0.00000000), (-5.55646610, +0.00000000), (+2.00845003, +0.00000000), (+16.0000019, +0.00000000), (-2.02138257, +0.00000000)
223(-2.89092064, +0.00000000), (-4.26904774, +0.00000000), (+8.07615852, +0.00000000), (-2.02138257, +0.00000000), (+25.0000019, +0.00000000)
225T
226
227
228ndim = getUnifRand(2, 3)
229scale = getUnifRand(1, 10)
230rand = getCovRand(mold, ndim)
231rand
232(+1.00000000, +0.00000000), (+0.255815506, +0.483223706)
233(+0.255815506, -0.483223706), (+1.00000000, +0.00000000)
235T
236rand = getCovRand(mold, ndim, scale)
237rand
238(+36.0000000, +0.00000000), (-17.0286980, +5.47461271)
239(-17.0286980, -5.47461271), (+36.0000000, +0.00000000)
241T
242
243
244rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
245rand
246(+1.00000000, +0.00000000), (-1.85860157, +0.00000000), (+1.56948590, +0.00000000), (+1.22511518, +0.00000000), (-1.96078491, +0.00000000)
247(-1.85860157, +0.00000000), (+4.00000000, +0.00000000), (-4.30787849, +0.00000000), (-0.464299560, +0.00000000), (+3.92996907, +0.00000000)
248(+1.56948590, +0.00000000), (-4.30787849, +0.00000000), (+9.00000000, +0.00000000), (-5.72658110, +0.00000000), (+0.295363516, +0.00000000)
249(+1.22511518, +0.00000000), (-0.464299560, +0.00000000), (-5.72658110, +0.00000000), (+16.0000000, +0.00000000), (-10.9778929, +0.00000000)
250(-1.96078491, +0.00000000), (+3.92996907, +0.00000000), (+0.295363516, +0.00000000), (-10.9778929, +0.00000000), (+24.9999981, +0.00000000)
252T
253
254
255ndim = getUnifRand(2, 3)
256scale = getUnifRand(1, 10)
257rand = getCovRand(mold, ndim)
258rand
259(+1.00000000, +0.00000000), (-0.144456744, +0.663907528), (-0.145413294, -0.449386686)
260(-0.144456744, -0.663907528), (+1.00000000, +0.00000000), (+0.346612930E-2, +0.347764552)
261(-0.145413294, +0.449386686), (+0.346612930E-2, -0.347764552), (+1.00000000, +0.00000000)
263T
264rand = getCovRand(mold, ndim, scale)
265rand
266(+25.0000000, +0.00000000), (+14.0006218, -11.6441727), (-1.92677546, -17.5996399)
267(+14.0006218, +11.6441727), (+24.9999962, +0.00000000), (+8.04820251, -1.62719035)
268(-1.92677546, +17.5996399), (+8.04820251, +1.62719035), (+25.0000000, +0.00000000)
270T
271
272
273rand = getCovRand(mold = 0._TKG, scale = [(real(itry, TKG), itry = 1, 5)])
274rand
275(+1.00000000, +0.00000000), (-1.74850202, +0.00000000), (+1.07058811, +0.00000000), (+2.08959079, +0.00000000), (-3.30059981, +0.00000000)
276(-1.74850202, +0.00000000), (+3.99999952, +0.00000000), (-2.31068301, +0.00000000), (-1.19517899, +0.00000000), (+7.81518221, +0.00000000)
277(+1.07058811, +0.00000000), (-2.31068301, +0.00000000), (+9.00000000, +0.00000000), (-5.03986788, +0.00000000), (+0.759493411, +0.00000000)
278(+2.08959079, +0.00000000), (-1.19517899, +0.00000000), (-5.03986788, +0.00000000), (+16.0000000, +0.00000000), (-7.05084562, +0.00000000)
279(-3.30059981, +0.00000000), (+7.81518221, +0.00000000), (+0.759493411, +0.00000000), (-7.05084562, +0.00000000), (+24.9999981, +0.00000000)
281T
282
283
Test:
test_pm_distCov


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 394 of file pm_distCov.F90.


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