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,
-
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) .
|
[in] | subset | : The input scalar constant that can be any of the following:
-
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.
-
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 ⛓
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...
- See also
- getCovRand
setCovRand
getCholRand
setCholRand
Example usage ⛓
12 integer(IK) :: isub, itry, ndim
14 type(csp_type) :: subset(
2)
16 type(display_type) :: disp
20 do isub
= 1,
size(subset)
24 real(TKG),
allocatable :: rand(:,:), cov(:,:)
28 call disp%show(
"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)")
36 call disp%show(
"cov = matmul(rand, transpose(rand))")
37 cov
= matmul(rand,
transpose(rand))
38 call disp%show(
"isMatClass(cov, posdefmat)")
46 complex(TKG),
allocatable :: rand(:,:), cov(:,:)
50 call disp%show(
"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)")
58 call disp%show(
"cov = matmul(rand, conjg(transpose(rand)))")
59 cov
= matmul(rand,
conjg(
transpose(rand)))
60 call disp%show(
"isMatClass(cov, posdefmat)")
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.
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...
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 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.
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 ⛓
3[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
7+1.00000000,
+0.431011409,
+0.190460011,
+0.593082964
8+0.00000000,
+0.902346432,
-0.702898920,
+0.585011184
9+0.00000000,
+0.00000000,
+0.685316026,
+0.224273711
10+0.00000000,
+0.00000000,
+0.00000000,
+0.505683601
11cov
= matmul(rand,
transpose(rand))
17[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
21+1.00000000,
-0.999967039
22+0.00000000,
+0.812121667E-2
23cov
= matmul(rand,
transpose(rand))
29[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
33+1.00000000,
+0.344397962,
-0.405866913E-1,
-0.678654790,
+0.440554470
34+0.00000000,
+0.938823700,
+0.310386062,
+0.661256373,
-0.128440931
35+0.00000000,
+0.00000000,
+0.949743807,
-0.319285095,
+0.135735199
36+0.00000000,
+0.00000000,
+0.00000000,
+0.149888098E-1,
-0.221862957
37+0.00000000,
+0.00000000,
+0.00000000,
+0.00000000,
+0.849568963
38cov
= matmul(rand,
transpose(rand))
44[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
48+1.00000000,
-0.655208349
49+0.00000000,
+0.755448222
50cov
= matmul(rand,
transpose(rand))
56[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
60+1.00000000,
+0.599283457,
-0.369215429,
-0.103403524,
+0.528383136
61+0.00000000,
+0.800536871,
-0.698699415,
-0.337057449E-1,
+0.277623609E-1
62+0.00000000,
+0.00000000,
+0.612779796,
-0.702681720,
+0.707158506
63+0.00000000,
+0.00000000,
+0.00000000,
+0.703142941,
+0.451549232
64+0.00000000,
+0.00000000,
+0.00000000,
+0.00000000,
+0.126769796
65cov
= matmul(rand,
transpose(rand))
71[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
75(
+1.00000000,
+0.00000000), (
+0.501353681,
+0.517161131)
76(
+0.00000000,
+0.00000000), (
+0.693677783,
+0.00000000)
77cov
= matmul(rand,
conjg(
transpose(rand)))
83[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
87(
+1.00000000,
+0.00000000), (
-0.958445966,
-0.268557906)
88(
+0.00000000,
+0.00000000), (
+0.962181538E-1,
+0.00000000)
89cov
= matmul(rand,
conjg(
transpose(rand)))
95[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
99(
+0.999999940,
+0.00000000), (
-0.651302457,
-0.737995207), (
-0.593800768E-1,
-0.963791072)
100(
+0.00000000,
+0.00000000), (
+0.176545456,
+0.00000000), (
-0.373499021E-1,
-0.256699622)
101(
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.170607101E-1,
+0.00000000)
102cov
= matmul(rand,
conjg(
transpose(rand)))
108[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
112(
+0.999999940,
+0.00000000), (
+0.546174049,
-0.463748164E-1), (
+0.488137633,
+0.404933095)
113(
+0.00000000,
+0.00000000), (
+0.836386979,
+0.00000000), (
-0.455201685,
+0.624451339)
114(
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.245520603E-1,
+0.00000000)
115cov
= matmul(rand,
conjg(
transpose(rand)))
121[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
125(
+1.00000000,
+0.00000000), (
+0.585761905,
-0.328417152), (
-0.393315881,
+0.379337490), (
-0.372362107,
+0.350343198)
126(
+0.00000000,
+0.00000000), (
+0.740962386,
+0.00000000), (
+0.491754651,
-0.555692315), (
-0.670220554,
+0.108992733)
127(
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.388315916,
+0.00000000), (
+0.285452873,
-0.435228407)
128(
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.813885033E-1,
+0.00000000)
129cov
= matmul(rand,
conjg(
transpose(rand)))
135[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
139+1.00000000,
+0.00000000
140-0.985380352,
+0.170368791
141cov
= matmul(rand,
transpose(rand))
147[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
151+1.00000000,
+0.00000000,
+0.00000000,
+0.00000000
152+0.143819168,
+0.989603937,
+0.00000000,
+0.00000000
153-0.599744201,
+0.669297203E-1,
+0.797387779,
+0.00000000
154-0.566567898,
-0.732030153,
-0.627516657E-1,
+0.373088300
155cov
= matmul(rand,
transpose(rand))
161[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
165+0.999999940,
+0.00000000,
+0.00000000,
+0.00000000
166-0.297733128,
+0.954649091,
+0.00000000,
+0.00000000
167-0.875454582E-2,
-0.837437689,
+0.546462715,
+0.00000000
168-0.204685748,
-0.729362965,
-0.446407348,
+0.476291865
169cov
= matmul(rand,
transpose(rand))
175[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
179+1.00000000,
+0.00000000,
+0.00000000,
+0.00000000
180-0.498659104,
+0.866798222,
+0.00000000,
+0.00000000
181+0.140598938,
-0.570691943,
+0.809038162,
+0.00000000
182-0.345104933,
+0.750771224,
-0.426404960E-1,
+0.561628878
183cov
= matmul(rand,
transpose(rand))
189[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
193+1.00000000,
+0.00000000,
+0.00000000,
+0.00000000,
+0.00000000
194+0.931225538,
+0.364443421,
+0.00000000,
+0.00000000,
+0.00000000
195+0.841458321,
+0.232422918,
+0.487778008,
+0.00000000,
+0.00000000
196-0.427990258,
+0.534431398,
-0.359917581,
+0.633771837,
+0.00000000
197+0.638648152,
+0.494456530,
-0.337863237,
+0.297248363,
+0.380963385
198cov
= matmul(rand,
transpose(rand))
204[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
208(
+1.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
209(
-0.750156403,
-0.486277342), (
+0.448106885,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
210(
-0.424236029,
-0.770430744), (
-0.134289891,
-0.338273078), (
+0.306590587,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
211(
+0.422100842,
+0.298661411), (
+0.346862018,
-0.432208240), (
+0.219184563,
+0.560490787), (
+0.251641065,
+0.00000000), (
+0.00000000,
+0.00000000)
212(
-0.299449027,
+0.505881667), (
+0.339169443,
+0.460669965), (
-0.623089001E-1,
-0.572213195E-1), (
+0.278725594,
-0.212214768), (
+0.444164008,
+0.00000000)
213cov
= matmul(rand,
conjg(
transpose(rand)))
219[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
223(
+1.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
224(
-0.164456919,
+0.203668982), (
+0.965128422,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
225(
-0.613006711,
+0.344371289), (
-0.391044289,
+0.552321635E-1), (
+0.591324747,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
226(
-0.489179611,
-0.304220766), (
+0.165346071,
+0.385288179), (
+0.463985950,
-0.133975372), (
+0.509052515,
+0.00000000), (
+0.00000000,
+0.00000000)
227(
-0.308671862,
-0.138337892E-1), (
-0.178560063,
-0.311820507), (
-0.140370905,
-0.719681561), (
-0.210986391,
+0.878540277E-1), (
+0.430738330,
+0.00000000)
228cov
= matmul(rand,
conjg(
transpose(rand)))
234[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
238(
+1.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
239(
-0.205078840,
-0.974743664), (
+0.884168744E-1,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
240(
+0.437355787,
-0.800350130), (
-0.251017213,
-0.201735243), (
+0.253875852,
+0.00000000), (
+0.00000000,
+0.00000000)
241(
-0.156143472,
-0.120431624), (
-0.188329324,
+0.491471350), (
-0.579425395,
+0.499466389), (
+0.314488441,
+0.00000000)
242cov
= matmul(rand,
conjg(
transpose(rand)))
248[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
252(
+1.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
253(
-0.778123617,
+0.285814285), (
+0.559315383,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
254(
-0.907061845E-1,
+0.508100353E-1), (
+0.729050219,
+0.665900171), (
+0.119387902,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
255(
-0.551060140,
-0.142469049), (
-0.121511482,
+0.352817088), (
-0.356894732,
+0.439488471), (
+0.465044498,
+0.00000000), (
+0.00000000,
+0.00000000)
256(
-0.363008901E-1,
-0.453984082), (
-0.228594363,
+0.127662629), (
-0.549087077E-1,
-0.493640333), (
-0.452644944,
-0.488726407), (
+0.183277994,
+0.00000000)
257cov
= matmul(rand,
conjg(
transpose(rand)))
263[
same_type_as(subset(isub)
%val, uppDia),
same_type_as(subset(isub)
%val, lowDia)]
267(
+0.999999940,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
268(
+0.691514552,
-0.165196881), (
+0.703219414,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
269(
-0.603671789,
+0.432656795), (
-0.648069561,
+0.478674583E-1), (
+0.161564484,
+0.00000000), (
+0.00000000,
+0.00000000), (
+0.00000000,
+0.00000000)
270(
-0.321031511,
-0.445555627), (
+0.266406685,
-0.481722951), (
-0.144210504E-1,
-0.524895787), (
+0.345927477,
+0.00000000), (
+0.00000000,
+0.00000000)
271(
+0.826192647E-2,
-0.186399132), (
-0.459752262,
-0.769128203E-1), (
+0.502102494,
+0.241792127), (
-0.516032159,
+0.407368273), (
+0.713516623E-1,
+0.00000000)
272cov
= matmul(rand,
conjg(
transpose(rand)))
- 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.
-
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 115 of file pm_distChol.F90.