Generate and return a random positive-definite (correlation or covariance) matrix using the Gram method.
More...
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,
-
type
complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128),
-
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...
This module contains classes and procedures for generating random matrices distributed on the space o...
- 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 ⛓
13 integer(IK) :: itry, ndim
15 type(display_type) :: disp
20 real(TKG),
allocatable :: rand(:,:)
25 call disp%show(
"ndim = getUnifRand(2, 3)")
27 call disp%show(
"scale = getUnifRand(1, 10)")
29 call disp%show(
"rand = getCovRand(mold, ndim)")
33 call disp%show(
"isMatClass(rand, posdefmat)")
35 call disp%show(
"rand = getCovRand(mold, ndim, scale)")
39 call disp%show(
"isMatClass(rand, posdefmat)")
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)])
47 call disp%show(
"isMatClass(rand, posdefmat)")
55 complex(TKG),
allocatable :: rand(:,:)
60 call disp%show(
"ndim = getUnifRand(2, 3)")
62 call disp%show(
"scale = getUnifRand(1, 10)")
64 call disp%show(
"rand = getCovRand(mold, ndim)")
68 call disp%show(
"isMatClass(rand, posdefmat)")
70 call disp%show(
"rand = getCovRand(mold, ndim, scale)")
74 call disp%show(
"isMatClass(rand, posdefmat)")
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)])
82 call disp%show(
"isMatClass(rand, posdefmat)")
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.
This is a generic method of the derived type display_type with pass attribute.
[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...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter CKS
The single-precision complex kind in Fortran mode. On most platforms, this is a 32-bit real kind.
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
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.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
6+1.00000000,
-0.678253531
7-0.678253531,
+1.00000000
12+64.0000000,
-3.33755255
13-3.33755255,
+63.9999962
18rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
20+1.00000000,
+1.87044573,
-1.68540144,
+1.63426292,
-0.722899914
21+1.87044573,
+3.99999976,
-1.86545885,
+0.920164824,
+0.177879333
22-1.68540144,
-1.86545885,
+8.99999905,
-8.71189404,
+3.64769506
23+1.63426292,
+0.920164824,
-8.71189404,
+16.0000000,
-2.45268774
24-0.722899914,
+0.177879333,
+3.64769506,
-2.45268774,
+25.0000000
33+1.00000000,
+0.850945950
34+0.850945950,
+1.00000000
39+4.00000000,
-1.38679838
40-1.38679838,
+4.00000000
45rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
47+1.00000000,
+1.55580175,
+1.82904339,
-0.210393891,
+3.56267190
48+1.55580175,
+4.00000000,
+4.05974817,
+3.07578659,
+7.74551964
49+1.82904339,
+4.05974817,
+9.00000000,
+7.66835022,
+10.6735783
50-0.210393891,
+3.07578659,
+7.66835022,
+15.9999971,
+10.5352306
51+3.56267190,
+7.74551964,
+10.6735783,
+10.5352306,
+25.0000000
60+1.00000000,
-0.576553762,
+0.712023318
61-0.576553762,
+1.00000000,
-0.718852520
62+0.712023318,
-0.718852520,
+1.00000000
67+81.0000000,
+77.1657639,
-65.4065247
68+77.1657639,
+80.9999924,
-61.6757469
69-65.4065247,
-61.6757469,
+80.9999924
74rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
76+1.00000000,
-1.13147545,
-2.19191217,
-2.20564246,
+2.74008274
77-1.13147545,
+4.00000000,
+2.64308381,
+2.28928328,
-1.62833118
78-2.19191217,
+2.64308381,
+8.99999809,
+10.4355431,
-0.602020919
79-2.20564246,
+2.28928328,
+10.4355431,
+16.0000000,
+3.75626707
80+2.74008274,
-1.62833118,
-0.602020919,
+3.75626707,
+24.9999981
89+1.00000000,
-0.797513843,
-0.622656584
90-0.797513843,
+1.00000000,
+0.887649179
91-0.622656584,
+0.887649179,
+1.00000000
96+1.00000000,
-0.754223943,
-0.931917250
97-0.754223943,
+1.00000000,
+0.467446566
98-0.931917250,
+0.467446566,
+0.999999940
103rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
105+0.999999881,
-0.705764830,
+1.10615718,
-1.05442274,
-1.60689688
106-0.705764830,
+4.00000000,
-3.93263412,
+6.84730196,
+2.06657386
107+1.10615718,
-3.93263412,
+8.99999809,
-2.90267181,
+7.29608011
108-1.05442274,
+6.84730196,
-2.90267181,
+15.9999990,
+9.07913399
109-1.60689688,
+2.06657386,
+7.29608011,
+9.07913399,
+25.0000019
118+1.00000000,
-0.187521651
119-0.187521651,
+1.00000000
124+64.0000000,
+49.2762146
125+49.2762146,
+63.9999962
130rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
132+1.00000000,
-1.47697520,
-0.110004004E-1,
+1.99357927,
-1.98655176
133-1.47697520,
+4.00000048,
-3.89132237,
-5.92055607,
+7.07246113
134-0.110004004E-1,
-3.89132237,
+8.99999905,
+7.61685944,
-11.0060549
135+1.99357927,
-5.92055607,
+7.61685944,
+15.9999971,
-15.9632530
136-1.98655176,
+7.07246113,
-11.0060549,
-15.9632530,
+25.0000000
145(
+1.00000000,
+0.00000000), (
-0.686078370,
+0.434828699)
146(
-0.686078370,
-0.434828699), (
+1.00000000,
+0.00000000)
151(
+36.0000000,
+0.00000000), (
-20.9644279,
+9.64115238)
152(
-20.9644279,
-9.64115238), (
+36.0000038,
+0.00000000)
157rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
159(
+1.00000000,
+0.00000000), (
-0.457986832,
+0.00000000), (
-1.64626336,
+0.00000000), (
+2.20828080,
+0.00000000), (
+1.63638210,
+0.00000000)
160(
-0.457986832,
+0.00000000), (
+4.00000048,
+0.00000000), (
-1.64370024,
+0.00000000), (
+4.61983347,
+0.00000000), (
-4.17207336,
+0.00000000)
161(
-1.64626336,
+0.00000000), (
-1.64370024,
+0.00000000), (
+8.99999809,
+0.00000000), (
-10.5807819,
+0.00000000), (
+7.99656963,
+0.00000000)
162(
+2.20828080,
+0.00000000), (
+4.61983347,
+0.00000000), (
-10.5807819,
+0.00000000), (
+16.0000000,
+0.00000000), (
-6.74015951,
+0.00000000)
163(
+1.63638210,
+0.00000000), (
-4.17207336,
+0.00000000), (
+7.99656963,
+0.00000000), (
-6.74015951,
+0.00000000), (
+25.0000000,
+0.00000000)
172(
+1.00000000,
+0.00000000), (
-0.587565482,
-0.661676705), (
-0.222535193,
+0.158793837)
173(
-0.587565482,
+0.661676705), (
+1.00000000,
+0.00000000), (
+0.700510964E-1,
-0.586710453)
174(
-0.222535193,
-0.158793837), (
+0.700510964E-1,
+0.586710453), (
+1.00000000,
+0.00000000)
179(
+48.9999924,
+0.00000000), (
-31.0752697,
-32.9097061), (
+23.7876492,
-23.1091728)
180(
-31.0752697,
+32.9097061), (
+48.9999924,
+0.00000000), (
+0.157827437,
+26.0913677)
181(
+23.7876492,
+23.1091728), (
+0.157827437,
-26.0913677), (
+48.9999924,
+0.00000000)
186rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
188(
+1.00000000,
+0.00000000), (
-1.54050994,
+0.00000000), (
+0.333717972,
+0.00000000), (
-1.19415343,
+0.00000000), (
+0.736712933,
+0.00000000)
189(
-1.54050994,
+0.00000000), (
+3.99999976,
+0.00000000), (
-3.44799280,
+0.00000000), (
+4.57718182,
+0.00000000), (
+2.07279730,
+0.00000000)
190(
+0.333717972,
+0.00000000), (
-3.44799280,
+0.00000000), (
+8.99999809,
+0.00000000), (
-1.40387714,
+0.00000000), (
+1.37851501,
+0.00000000)
191(
-1.19415343,
+0.00000000), (
+4.57718182,
+0.00000000), (
-1.40387714,
+0.00000000), (
+16.0000019,
+0.00000000), (
+13.5879316,
+0.00000000)
192(
+0.736712933,
+0.00000000), (
+2.07279730,
+0.00000000), (
+1.37851501,
+0.00000000), (
+13.5879316,
+0.00000000), (
+25.0000000,
+0.00000000)
201(
+1.00000000,
+0.00000000), (
-0.517126858,
+0.663563609)
202(
-0.517126858,
-0.663563609), (
+1.00000000,
+0.00000000)
207(
+81.0000000,
+0.00000000), (
+61.4656792,
-26.1785831)
208(
+61.4656792,
+26.1785831), (
+81.0000000,
+0.00000000)
213rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
215(
+1.00000000,
+0.00000000), (
+1.54776490,
+0.00000000), (
+1.45657814,
+0.00000000), (
+3.43243718,
+0.00000000), (
+2.93747044,
+0.00000000)
216(
+1.54776490,
+0.00000000), (
+4.00000048,
+0.00000000), (
+5.57432985,
+0.00000000), (
+6.66049480,
+0.00000000), (
+1.95156074,
+0.00000000)
217(
+1.45657814,
+0.00000000), (
+5.57432985,
+0.00000000), (
+9.00000191,
+0.00000000), (
+7.62362909,
+0.00000000), (
-1.02295768,
+0.00000000)
218(
+3.43243718,
+0.00000000), (
+6.66049480,
+0.00000000), (
+7.62362909,
+0.00000000), (
+15.9999981,
+0.00000000), (
+6.46011639,
+0.00000000)
219(
+2.93747044,
+0.00000000), (
+1.95156074,
+0.00000000), (
-1.02295768,
+0.00000000), (
+6.46011639,
+0.00000000), (
+25.0000019,
+0.00000000)
228(
+1.00000000,
+0.00000000), (
-0.762964785E-1,
+0.261936903), (
+0.400371492,
+0.284627825)
229(
-0.762964785E-1,
-0.261936903), (
+1.00000000,
+0.00000000), (
-0.363773525,
+0.224467188)
230(
+0.400371492,
-0.284627825), (
-0.363773525,
-0.224467188), (
+1.00000000,
+0.00000000)
235(
+16.0000000,
+0.00000000), (
-0.574790597,
-10.3225441), (
-7.97072029,
+8.25149918)
236(
-0.574790597,
+10.3225441), (
+16.0000038,
+0.00000000), (
-4.97169638,
-9.81996250)
237(
-7.97072029,
-8.25149918), (
-4.97169638,
+9.81996250), (
+16.0000019,
+0.00000000)
242rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
244(
+1.00000000,
+0.00000000), (
-1.95634520,
+0.00000000), (
-0.360667944,
+0.00000000), (
+0.541187823,
+0.00000000), (
+0.291445613,
+0.00000000)
245(
-1.95634520,
+0.00000000), (
+3.99999976,
+0.00000000), (
-0.454427004E-1,
+0.00000000), (
+0.185179353,
+0.00000000), (
-1.27388930,
+0.00000000)
246(
-0.360667944,
+0.00000000), (
-0.454427004E-1,
+0.00000000), (
+8.99999905,
+0.00000000), (
-6.97992039,
+0.00000000), (
-1.30858243,
+0.00000000)
247(
+0.541187823,
+0.00000000), (
+0.185179353,
+0.00000000), (
-6.97992039,
+0.00000000), (
+15.9999990,
+0.00000000), (
-12.0468855,
+0.00000000)
248(
+0.291445613,
+0.00000000), (
-1.27388930,
+0.00000000), (
-1.30858243,
+0.00000000), (
-12.0468855,
+0.00000000), (
+25.0000000,
+0.00000000)
257(
+1.00000000,
+0.00000000), (
-0.375887960,
+0.998196900E-1), (
+0.960442871E-1,
-0.499969721)
258(
-0.375887960,
-0.998196900E-1), (
+1.00000000,
+0.00000000), (
+0.401803136,
-0.917571783E-2)
259(
+0.960442871E-1,
+0.499969721), (
+0.401803136,
+0.917571783E-2), (
+1.00000000,
+0.00000000)
264(
+36.0000000,
+0.00000000), (
+7.46848583,
+0.224649966), (
+19.4333019,
+14.4216156)
265(
+7.46848583,
-0.224649966), (
+35.9999962,
+0.00000000), (
+4.68500662,
+11.6969585)
266(
+19.4333019,
-14.4216156), (
+4.68500662,
-11.6969585), (
+36.0000000,
+0.00000000)
271rand
= getCovRand(
mold = 0._TKG, scale
= [(
real(itry, TKG), itry
= 1,
5)])
273(
+1.00000000,
+0.00000000), (
-1.31374502,
+0.00000000), (
-0.987170577,
+0.00000000), (
-0.516361177,
+0.00000000), (
+1.40397084,
+0.00000000)
274(
-1.31374502,
+0.00000000), (
+4.00000000,
+0.00000000), (
+2.34914804,
+0.00000000), (
+3.43394709,
+0.00000000), (
-5.70366859,
+0.00000000)
275(
-0.987170577,
+0.00000000), (
+2.34914804,
+0.00000000), (
+9.00000000,
+0.00000000), (
+10.5729828,
+0.00000000), (
+3.94087410,
+0.00000000)
276(
-0.516361177,
+0.00000000), (
+3.43394709,
+0.00000000), (
+10.5729828,
+0.00000000), (
+16.0000000,
+0.00000000), (
+4.61187267,
+0.00000000)
277(
+1.40397084,
+0.00000000), (
-5.70366859,
+0.00000000), (
+3.94087410,
+0.00000000), (
+4.61187267,
+0.00000000), (
+25.0000000,
+0.00000000)
- 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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- 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.