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

Return the (weighted) correlation matrix corresponding to the input (weighted) covariance matrix or return the (weighted) sample Pearson correlation matrix of the input array of shape (ndim, nsam) or (nsam, ndim) or the Pearson correlation coefficient a pair of (weighted) time series x(1:nsam) and y(1:nsam) where ndim is the number of data dimensions (the number of data attributes) and nsam is the number of data points.
More...

Detailed Description

Return the (weighted) correlation matrix corresponding to the input (weighted) covariance matrix or return the (weighted) sample Pearson correlation matrix of the input array of shape (ndim, nsam) or (nsam, ndim) or the Pearson correlation coefficient a pair of (weighted) time series x(1:nsam) and y(1:nsam) where ndim is the number of data dimensions (the number of data attributes) and nsam is the number of data points.

This generic interface performs one of the following computational tasks:

  1. Compute the correlation matrix corresponding to an input covariance matrix in arbitrary ndim dimensions.
  2. Compute the Pearson correlation coefficient corresponding to an input pair of time series x and y of nsam observations.
  3. Compute the Pearson correlation matrix corresponding to an input multivariate sample of nsam observations each with ndim attributes.

See the documentation of the parent module pm_sampleCor for algorithmic details and sample correlation matrix definition.

Note
The sample (Pearson) correlation matrix can also be readily computed from the sample covariance matrix as in the following example,
! unweighted sample.
call setCov(cor, subset, .false., sample, dim)
call setCor(cor, subset, cor, subset)
! weighted sample.
call setCov(cor, subset, .false., sample, dim, weight)
call setCor(cor, subset, cor, subset)
!
Parameters
[in,out]cor: The output or input/output positive semi-definite scalar or square 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 the (Pearson) correlation coefficient or matrix corresponding to the input sample or time series x and y or correlation matrix cov, whichever is present.
  1. If x(:) and y(:) are present, then cor shall be a scalar, corresponding to the element \(r_xy\) (cor(1,2)) of the computed correlation matrix.
    Since the correlation matrix is Hermitian, cor(2,1) = conjg(cor(1,2) while the diagonal elements are unit values.
  2. If cov is present, then cor shall be a square matrix of shape shape(cov).
    On output, only the specified input subset will be overwritten with the correlation matrix.
    Any elements not in the specified input subset remains intact.
  3. If sample is present, then cor shall be a square matrix of shape [size(sample, 3 - dim), size(sample, 3 - dim)].
    On output, only the specified input subset will be overwritten with the correlation matrix.
    Any elements not in the specified input subset remains intact.
[in]subset: The input scalar constant argument that can be any of the following:
  1. The constant low, implying that only the lower subset of the output correlation matrix must be computed.
    Specifying this value is useful when the contents of the diagonal elements of the input cor must be preserved because they contain valuable information, e.g., sample variance.
    This option is currently available only if the input argument cov is present.
    The diagonal elements of cor will remain untouched upon return.
  2. The constant upp, implying that only the upper subset of the output correlation matrix must be computed.
    Specifying this value is useful when the contents of the diagonal elements of the input cor must be preserved because they contain valuable information, e.g., sample variance.
    This option is currently available only if the input argument cov is present.
    The diagonal elements of cor will remain untouched upon return.
  3. The constant lowDia, implying that only the lower-diagonal subset of the output correlation matrix must be computed.
    This option is available only if either of the input arguments cov or sample are present.
    By definition, all diagonal elements of cor will be set to 1.
  4. The constant uppDia, implying that only the upper-diagonal subset of the output correlation matrix must be computed.
    This option is available only if either of the input arguments cov or sample are present.
    By definition, all diagonal elements of cor will be set to 1.
This input argument is merely serves to resolve the different procedures of this generic interface from each other at compile-time.
[in]cov: The input positive semi-definite square matrix of shape (1 : ndim, 1 : ndim) of the same type and kind as the input cor, representing the covariance matrix based upon which the output correlation matrix cor is to be computed.
(optional. It must be present if and only if the input arguments sample, x, and y are missing.)
[in]subsetv: The input scalar constant argument that can be any of the following:
  1. The constant lowDia, implying that only the lower-diagonal subset of the input covariance matrix must be used.
  2. The constant uppDia, implying that only the upper-diagonal subset of the input covariance matrix must be used.
  3. The constant low, implying that only the lower subset of the input covariance matrix must be used.
    This option is available if and only if the input argument stdinv is present (otherwise, how do we know the variances?).
  4. The constant upp, implying that only the upper subset of the input covariance matrix must be used.
    This option is available if and only if the input argument stdinv is present (otherwise, how do we know the variances?).
This input argument is merely serves to resolve the different procedures of this generic interface from each other at compile-time.
(optional. It must be present if and only if the input argument cov is present.)
[in]stdinv: The input positive vector of shape (1 : ndim) of type real of the same kind as the input cor, containing the inverse of the standard deviations of the ndim data attributes based upon which the output correlation matrix cor is to be computed.

\begin{equation} \ms{stdinv}(i) = \frac{1} {\sqrt{\ms{cov}(i,i)}} ~, \end{equation}

(optional. It must be present if and only if the input argument subsetv is set to either low or upp.)
[in]mean: The input scalar or contiguous vector of shape (1:ndim) of the same type and kind as the output cor containing the sample mean.
  1. If the input sample is a 1D array, then mean must be a scalar.
  2. If the input sample is array, then mean must be a vector of size ndim = size(sample, 3 - dim) (i.e., computed along the specified input dimension dim).
The mean vector can be readily computed via getMean or setMean.
(optional. default = getFilled(0., ndim) or [0., 0.] depending on the rank of the input sample and dim or the presence of x and y.)
[in]x: The input contiguous vector of shape (nsam) of the same type and kind as the output cov, containing the first attribute x of the observational sample, where nsam is the number of observations in the sample.
(optional. It must be present if and only if the input argument y is present and cov and sample are missing.)
[in]y: The input contiguous vector of shape (nsam) of the same type and kind as the output cov, containing the second attribute x of the observational sample, where nsam is the number of observations in the sample.
(optional. It must be present if and only if the input argument x is present and cov and sample are missing.)
[in]sample: The input contiguous array of shape (ndim, nsam) or (nsam, ndim) of the same type and kind as the output cor, containing the sample comprised of nsam observations each with ndim attributes.
If sample is a matrix, then the input argument dim dictates the direction along which the correlation matrix cor must be computed (i.e., the direction of individual observations).
(optional. It must be present if and only if the input arguments cov, x, and y are missing.)
[in]dim: The input scalar integer of default kind IK indicating the dimension of sample along which the correlation matrix must be computed.
  1. If dim = 1, the input sample is assumed to have the shape (nsam, ndim).
  2. If dim = 2, the input sample is assumed to have the shape (ndim, nsam).
(optional. It must be present if and only if the input argument sample is present and is of rank 2.)
[in]weight: The contiguous vector of length nsam of,
  1. type integer of default kind IK, or
  2. type real of the same kind as the kind of the output cor,
containing the corresponding weights of individual nsam observations in sample or the pair of vectors x and y.
(optional. default = getFilled(1, nsam). It can be present if and only if the input arguments sample or x and y are present.)
[in]weisum: The input scalar of the same type and kind as the input weight containing sum(weight).
This quantity is a byproduct of computing the mean of a sample and is automatically returned by setMean.
(optional. It must be present if and only if both mean and weight arguments are present and the input argument cov is missing.)


Possible calling interfaces

use pm_sampleCor, only: setCor
! correlation matrix to covariance matrix.
call setCor(cor(1:ndim, 1:ndim), subset, cov(1:ndim, 1:ndim), subsetv) ! subsetv = uppDia, lowDia
call setCor(cor(1:ndim, 1:ndim), subset, cov(1:ndim, 1:ndim), subsetv, stdinv(1:ndim)) ! subsetv = upp, low
! XY time series Pearson correlation coefficient.
call setCor(cor, x(1:nsam), y(1:nsam)) ! Pearson correlation coefficient.
call setCor(cor, x(1:nsam), y(1:nsam), weight(1:nsam)) ! Pearson correlation coefficient.
call setCor(cor, mean(1:2), x(1:nsam), y(1:nsam)) ! Pearson correlation coefficient.
call setCor(cor, mean(1:2), x(1:nsam), y(1:nsam), weight(1:nsam), weisum) ! Pearson correlation coefficient.
! sample correlation matrix.
call setCor(cor(1:ndim, 1:ndim), subset, sample(:,:), dim)
call setCor(cor(1:ndim, 1:ndim), subset, sample(:,:), dim, weight(1:nsam))
call setCor(cor(1:ndim, 1:ndim), subset, mean(1:ndim), sample(:,:), dim)
call setCor(cor(1:ndim, 1:ndim), subset, mean(1:ndim), sample(:,:), dim, weight(1:nsam), weisum)
Return the (weighted) correlation matrix corresponding to the input (weighted) covariance matrix or r...
This module contains classes and procedures for computing properties related to the correlation matri...
Warning
The condition any(x(1) /= x) must hold for the corresponding input arguments.
The condition any(y(1) /= x) must hold for the corresponding input arguments.
The input sample must contain at least two unique values per sample attribute.
The condition 0 < sum(weight) must hold for the corresponding input arguments.
The condition all(0. <= weight) must hold for the corresponding input arguments.
The condition 1 < size(sample, dim) must hold for the corresponding input arguments.
The condition 0 < size(sample, 3 - dim) must hold for the corresponding input arguments.
The condition all(shape(cov) == shape(cor)) must hold for the corresponding input arguments.
The condition 1 <= dim .and. dim <= rank(sample) must hold for the corresponding input arguments.
The condition size(sample, 3 - dim) == size(variance) must hold for the corresponding input arguments.
The condition size(sample, 3 - dim) == size(mean) must hold for the corresponding input arguments.
The condition size(sample, dim) == size(weight) must hold for the corresponding input arguments.
The condition size(x) == size(weight) must hold for the corresponding input arguments.
The condition size(x) == size(y) must hold for the corresponding input arguments.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1.
The pure procedure(s) documented herein become impure when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1.
By default, these procedures are pure in release build and impure in debug and testing builds.
See also
getCor
setCor
getRho
setRho
getCov
setCov
setECDF
getMean
setMean
getShifted
setShifted
getVar
setVar


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_distCov, only: getCovRand
7 use pm_arrayRange, only: getRange
8 use pm_sampleCor, only: setCor
9 use pm_sampleCor, only: uppDia
10 use pm_sampleCor, only: lowDia
11 use pm_sampleCor, only: upp
12 use pm_sampleCor, only: low
13 use pm_sampleMean, only: getMean
14 use pm_sampleMean, only: setMean
15 use pm_sampleShift, only: getShifted
16 use pm_arraySpace, only: getLinSpace
17 use pm_arrayResize, only: setResized
18 use pm_arrayFill, only: getFilled
19 use pm_io, only: getFormat
20
21 implicit none
22
23 type(display_type) :: disp
24 integer(IK) :: itry, ntry = 10
25 character(:), allocatable :: format
26 disp = display_type(file = "main.out.F90")
27
28 call disp%skip()
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%show("!Compute the correlation matrix from correlation matrix.")
31 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
32 call disp%skip()
33
34 block
35 use pm_kind, only: TKG => RKS ! All processor kinds are supported.
36 real(TKG), allocatable :: cov(:,:), cor(:,:), std(:)
37 integer(IK) :: ndim
38 format = getFormat(mold = [0._TKG], ed = SK_"es", signed = .true._LK)
39 do itry = 1, ntry
40
41 call disp%skip()
42 call disp%show("ndim = getUnifRand(0, 4)")
43 ndim = getUnifRand(0, 4)
44 call disp%show("ndim")
45 call disp%show( ndim )
46 call disp%show("std = getUnifRand(1, ndim, ndim)")
47 std = getUnifRand(1, ndim, ndim)
48 call disp%show("std")
49 call disp%show( std , format = format )
50 call disp%show("cov = getCovRand(1._TKG, std)")
51 cov = getCovRand(1._TKG, std)
52 call disp%show("cov")
53 call disp%show( cov , format = format )
54 call disp%skip()
55
56
57 call disp%skip()
58 call disp%show("cor = cov")
59 cor = cov
60 call disp%show("call setCor(cor, uppDia, cor, uppDia)")
61 call setCor(cor, uppDia, cor, uppDia)
62 call disp%show("cor")
63 call disp%show( cor , format = format )
64 call disp%skip()
65 call disp%show("cor = cov")
66 cor = cov
67 call disp%show("call setCor(cor, uppDia, cor, lowDia)")
68 call setCor(cor, uppDia, cor, lowDia)
69 call disp%show("cor")
70 call disp%show( cor , format = format )
71 call disp%skip()
72 call disp%show("cor = cov")
73 cor = cov
74 call disp%show("call setCor(cor, upp, cor, uppDia)")
75 call setCor(cor, upp, cor, uppDia)
76 call disp%show("cor")
77 call disp%show( cor , format = format )
78 call disp%skip()
79 call disp%show("cor = cov")
80 cor = cov
81 call disp%show("call setCor(cor, upp, cor, lowDia)")
82 call setCor(cor, upp, cor, lowDia)
83 call disp%show("cor")
84 call disp%show( cor , format = format )
85 call disp%skip()
86 call disp%show("cor = cov")
87 cor = cov
88 call disp%show("call setCor(cor, lowDia, cor, uppDia)")
89 call setCor(cor, lowDia, cor, uppDia)
90 call disp%show("cor")
91 call disp%show( cor , format = format )
92 call disp%skip()
93 call disp%show("cor = cov")
94 cor = cov
95 call disp%show("call setCor(cor, lowDia, cor, lowDia)")
96 call setCor(cor, lowDia, cor, lowDia)
97 call disp%show("cor")
98 call disp%show( cor , format = format )
99 call disp%skip()
100 call disp%show("cor = cov")
101 cor = cov
102 call disp%show("call setCor(cor, low, cor, uppDia)")
103 call setCor(cor, low, cor, uppDia)
104 call disp%show("cor")
105 call disp%show( cor , format = format )
106 call disp%skip()
107 call disp%show("cor = cov")
108 cor = cov
109 call disp%show("call setCor(cor, low, cor, lowDia)")
110 call setCor(cor, low, cor, lowDia)
111 call disp%show("cor")
112 call disp%show( cor , format = format )
113 call disp%skip()
114
115
116 call disp%skip()
117 call disp%show("cor = cov")
118 cor = cov
119 call disp%show("call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)")
120 call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
121 call disp%show("cor")
122 call disp%show( cor , format = format )
123 call disp%skip()
124 call disp%show("cor = cov")
125 cor = cov
126 call disp%show("call setCor(cor, uppDia, cor, low, stdinv = 1 / std)")
127 call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
128 call disp%show("cor")
129 call disp%show( cor , format = format )
130 call disp%skip()
131 call disp%show("cor = cov")
132 cor = cov
133 call disp%show("call setCor(cor, upp, cor, upp, stdinv = 1 / std)")
134 call setCor(cor, upp, cor, upp, stdinv = 1 / std)
135 call disp%show("cor")
136 call disp%show( cor )
137 call disp%skip()
138 call disp%show("cor = cov")
139 cor = cov
140 call disp%show("call setCor(cor, upp, cor, low, stdinv = 1 / std)")
141 call setCor(cor, upp, cor, low, stdinv = 1 / std)
142 call disp%show("cor")
143 call disp%show( cor , format = format )
144 call disp%skip()
145 call disp%show("cor = cov")
146 cor = cov
147 call disp%show("call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)")
148 call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
149 call disp%show("cor")
150 call disp%show( cor , format = format )
151 call disp%skip()
152 call disp%show("cor = cov")
153 cor = cov
154 call disp%show("call setCor(cor, lowDia, cor, low, stdinv = 1 / std)")
155 call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
156 call disp%show("cor")
157 call disp%show( cor , format = format )
158 call disp%skip()
159 call disp%show("cor = cov")
160 cor = cov
161 call disp%show("call setCor(cor, low, cor, upp, stdinv = 1 / std)")
162 call setCor(cor, low, cor, upp, stdinv = 1 / std)
163 call disp%show("cor")
164 call disp%show( cor , format = format )
165 call disp%skip()
166 call disp%show("cor = cov")
167 cor = cov
168 call disp%show("call setCor(cor, low, cor, low, stdinv = 1 / std)")
169 call setCor(cor, low, cor, low, stdinv = 1 / std)
170 call disp%show("cor")
171 call disp%show( cor , format = format )
172 call disp%skip()
173
174 end do
175 end block
176
177 call disp%skip()
178 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
179 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
180 call disp%show("!Compute the Pearson correlation matrix for a pair of time series.")
181 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
182 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
183 call disp%skip()
184
185 block
186 use pm_kind, only: TKG => RKS ! All other real types are also supported.
187 integer(IK) :: ndim, nsam, dim
188 real(TKG), allocatable :: sample(:,:), samShifted(:,:), cor(:,:), mean(:)
189 format = getFormat(mold = [0._TKG], ed = SK_"es", signed = .true._LK)
190 call disp%show("ndim = 2; nsam = 10; dim = 2")
191 ndim = 2; nsam = 10; dim = 2
192 call disp%show("sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
193 sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
194 call disp%show("sample")
195 call disp%show( sample , format = format )
196 call disp%show("mean = getMean(sample, dim)")
197 mean = getMean(sample, dim)
198 call disp%show("mean")
199 call disp%show( mean , format = format )
200 call disp%show("samShifted = getShifted(sample, dim, -mean)")
201 samShifted = getShifted(sample, dim, -mean)
202 call disp%show("samShifted")
203 call disp%show( samShifted , format = format )
204 call disp%show("cor = getFilled(0., ndim, ndim)")
205 cor = getFilled(0., ndim, ndim)
206 call disp%show("call setCor(cor, uppDia, mean, sample, dim)")
207 call setCor(cor, uppDia, mean, sample, dim)
208 call disp%show("cor")
209 call disp%show( cor , format = format )
210 call disp%show("cor = getFilled(0., ndim, ndim)")
211 cor = getFilled(0., ndim, ndim)
212 call disp%show("call setCor(cor, uppDia, samShifted, dim) ! same result as above.")
213 call setCor(cor, uppDia, samShifted, dim)
214 call disp%show("cor")
215 call disp%show( cor , format = format )
216 call disp%skip()
217 call disp%show("cor = getFilled(0., ndim, ndim)")
218 cor = getFilled(0., ndim, ndim)
219 call disp%show("call setCor(cor, lowDia, mean, sample, dim)")
220 call setCor(cor, lowDia, mean, sample, dim)
221 call disp%show("cor")
222 call disp%show( cor , format = format )
223 call disp%show("cor = getFilled(0., ndim, ndim)")
224 cor = getFilled(0., ndim, ndim)
225 call disp%show("call setCor(cor, lowDia, samShifted, dim) ! same result as above.")
226 call setCor(cor, lowDia, samShifted, dim)
227 call disp%show("cor")
228 call disp%show( cor , format = format )
229 call disp%skip()
230 call disp%show("Compute the sample correlation along the first dimension.", deliml = SK_'''')
231 call disp%skip()
232 call disp%show("dim = 1")
233 dim = 1
234 call disp%show("cor = getFilled(0., ndim, ndim)")
235 cor = getFilled(0., ndim, ndim)
236 call disp%show("call setCor(cor, uppDia, mean, transpose(sample), dim)")
237 call setCor(cor, uppDia, mean, transpose(sample), dim)
238 call disp%show("cor")
239 call disp%show( cor , format = format )
240 call disp%show("cor = getFilled(0., ndim, ndim)")
241 cor = getFilled(0., ndim, ndim)
242 call disp%show("call setCor(cor, uppDia, transpose(samShifted), dim) ! same result as above.")
243 call setCor(cor, uppDia, transpose(samShifted), dim)
244 call disp%show("cor")
245 call disp%show( cor , format = format )
246 call disp%skip()
247 call disp%show("cor = getFilled(0., ndim, ndim)")
248 cor = getFilled(0., ndim, ndim)
249 call disp%show("call setCor(cor, lowDia, mean, transpose(sample), dim)")
250 call setCor(cor, lowDia, mean, transpose(sample), dim)
251 call disp%show("cor")
252 call disp%show( cor , format = format )
253 call disp%show("cor = getFilled(0., ndim, ndim)")
254 cor = getFilled(0., ndim, ndim)
255 call disp%show("call setCor(cor, lowDia, transpose(samShifted), dim) ! same result as above.")
256 call setCor(cor, lowDia, transpose(samShifted), dim)
257 call disp%show("cor")
258 call disp%show( cor , format = format )
259 call disp%skip()
260 call disp%show("Compute the full sample correlation for a pair of time series.", deliml = SK_'''')
261 call disp%skip()
262 call disp%show("call setCor(cor(1,1), mean, sample(1,:), sample(2,:))")
263 call setCor(cor(1,1), mean, sample(1,:), sample(2,:))
264 call disp%show("cor(1,1)")
265 call disp%show( cor(1,1) , format = format )
266 call disp%show("call setCor(cor(1,1), samShifted(1,:), samShifted(2,:)) ! same result as above.")
267 call setCor(cor(1,1), samShifted(1,:), samShifted(2,:))
268 call disp%show("cor(1,1)")
269 call disp%show( cor(1,1) , format = format )
270 call disp%skip()
271 end block
272
273 call disp%skip()
274 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
275 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
276 call disp%show("!Compute the Pearson correlation matrix for a weighted pair of time series.")
277 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
278 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
279 call disp%skip()
280
281 block
282 use pm_arrayVerbose, only: getVerbose
283 use pm_kind, only: TKG => RKS ! All other real types are also supported.
284 integer(IK), allocatable :: iweight(:)
285 real(TKG), allocatable :: rweight(:)
286 real(TKG) :: rweisum
287 integer(IK) :: dim
288 integer(IK) :: iweisum
289 integer(IK) :: ndim, nsam
290 real(TKG), allocatable :: sample(:,:), samShifted(:,:), cor(:,:), mean(:)
291 format = getFormat(mold = [0._TKG], ed = SK_"es", signed = .true._LK)
292 call disp%show("ndim = 2; nsam = 10; dim = 2")
293 ndim = 2; nsam = 10; dim = 2
294 call disp%show("sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
295 sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
296 call disp%show("sample")
297 call disp%show( sample , format = format )
298 call disp%show("call setResized(mean, ndim)")
299 call setResized(mean, ndim)
300 call disp%show("iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.")
301 iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
302 call disp%show("iweight")
303 call disp%show( iweight )
304 call disp%show("call setMean(mean, sample, dim, iweight, iweisum)")
305 call setMean(mean, sample, dim, iweight, iweisum)
306 call disp%show("mean")
307 call disp%show( mean , format = format )
308 call disp%show("iweisum")
309 call disp%show( iweisum )
310 call disp%show("rweight = iweight ! or real-valued weights.")
311 rweight = iweight ! or real-valued weights.
312 call disp%show("iweight")
313 call disp%show( iweight )
314 call disp%show("call setMean(mean, sample, dim, rweight, rweisum)")
315 call setMean(mean, sample, dim, rweight, rweisum)
316 call disp%show("mean")
317 call disp%show( mean , format = format )
318 call disp%show("rweisum")
319 call disp%show( rweisum , format = format )
320
321 call disp%skip()
322 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
323 call disp%show("!Compute the correlation matrix with integer weights.")
324 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
325 call disp%skip()
326
327 call disp%show("samShifted = getShifted(sample, dim, -mean)")
328 samShifted = getShifted(sample, dim, -mean)
329 call disp%show("samShifted")
330 call disp%show( samShifted , format = format )
331 call disp%show("cor = getFilled(0., ndim, ndim)")
332 cor = getFilled(0., ndim, ndim)
333 call disp%show("call setCor(cor, uppDia, mean, sample, dim, iweight, iweisum)")
334 call setCor(cor, uppDia, mean, sample, dim, iweight, iweisum)
335 call disp%show("cor")
336 call disp%show( cor , format = format )
337 call disp%show("cor = getFilled(0., ndim, ndim)")
338 cor = getFilled(0., ndim, ndim)
339 call disp%show("call setCor(cor, uppDia, samShifted, dim, iweight, iweisum) ! same result as above.")
340 call setCor(cor, uppDia, samShifted, dim, iweight, iweisum)
341 call disp%show("cor")
342 call disp%show( cor , format = format )
343 call disp%skip()
344 call disp%show("cor = getFilled(0., ndim, ndim)")
345 cor = getFilled(0., ndim, ndim)
346 call disp%show("call setCor(cor, lowDia, mean, sample, dim, iweight, iweisum)")
347 call setCor(cor, lowDia, mean, sample, dim, iweight, iweisum)
348 call disp%show("cor")
349 call disp%show( cor , format = format )
350 call disp%show("cor = getFilled(0., ndim, ndim)")
351 cor = getFilled(0., ndim, ndim)
352 call disp%show("call setCor(cor, lowDia, samShifted, dim, iweight, iweisum) ! same result as above.")
353 call setCor(cor, lowDia, samShifted, dim, iweight, iweisum)
354 call disp%show("cor")
355 call disp%show( cor , format = format )
356 call disp%show("cor = getFilled(0., ndim, ndim)")
357 cor = getFilled(0., ndim, ndim)
358 call disp%show("call setCor(cor, lowDia, getVerbose(samShifted, iweight, sum(iweight), dim), dim) ! same result as above.")
359 call setCor(cor, lowDia, getVerbose(samShifted, iweight, sum(iweight), dim), dim)
360 call disp%show("cor")
361 call disp%show( cor , format = format )
362 call disp%skip()
363 call disp%show("Compute the sample correlation along the first dimension.", deliml = SK_'''')
364 call disp%skip()
365 call disp%show("dim = 1")
366 dim = 1
367 call disp%show("cor = getFilled(0., ndim, ndim)")
368 cor = getFilled(0., ndim, ndim)
369 call disp%show("call setCor(cor, uppDia, mean, transpose(sample), dim, iweight, iweisum)")
370 call setCor(cor, uppDia, mean, transpose(sample), dim, iweight, iweisum)
371 call disp%show("cor")
372 call disp%show( cor , format = format )
373 call disp%show("cor = getFilled(0., ndim, ndim)")
374 cor = getFilled(0., ndim, ndim)
375 call disp%show("call setCor(cor, uppDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.")
376 call setCor(cor, uppDia, transpose(samShifted), dim, iweight, iweisum)
377 call disp%show("cor")
378 call disp%show( cor , format = format )
379 call disp%skip()
380 call disp%show("cor = getFilled(0., ndim, ndim)")
381 cor = getFilled(0., ndim, ndim)
382 call disp%show("call setCor(cor, lowDia, mean, transpose(sample), dim, iweight, iweisum)")
383 call setCor(cor, lowDia, mean, transpose(sample), dim, iweight, iweisum)
384 call disp%show("cor")
385 call disp%show( cor , format = format )
386 call disp%show("cor = getFilled(0., ndim, ndim)")
387 cor = getFilled(0., ndim, ndim)
388 call disp%show("call setCor(cor, lowDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.")
389 call setCor(cor, lowDia, transpose(samShifted), dim, iweight, iweisum)
390 call disp%show("cor")
391 call disp%show( cor , format = format )
392 call disp%skip()
393 call disp%show("Compute the full sample correlation for a pair of time series.", deliml = SK_'''')
394 call disp%skip()
395 call disp%show("call setCor(cor(1,1), mean, sample(1,:), sample(2,:), iweight, iweisum)")
396 call setCor(cor(1,1), mean, sample(1,:), sample(2,:), iweight, iweisum)
397 call disp%show("cor(1,1)")
398 call disp%show( cor(1,1) , format = format )
399 call disp%show("call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), iweight, iweisum) ! same result as above.")
400 call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), iweight, iweisum)
401 call disp%show("cor(1,1)")
402 call disp%show( cor(1,1) , format = format )
403 call disp%skip()
404
405 call disp%skip()
406 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
407 call disp%show("!Compute the correlation matrix with real weights.")
408 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
409 call disp%skip()
410
411 call disp%show("dim = 2")
412 dim = 2
413 call disp%show("samShifted = getShifted(sample, dim, -mean)")
414 samShifted = getShifted(sample, dim, -mean)
415 call disp%show("samShifted")
416 call disp%show( samShifted )
417 call disp%show("cor = getFilled(0., ndim, ndim)")
418 cor = getFilled(0., ndim, ndim)
419 call disp%show("call setCor(cor, uppDia, mean, sample, dim, rweight, rweisum)")
420 call setCor(cor, uppDia, mean, sample, dim, rweight, rweisum)
421 call disp%show("cor")
422 call disp%show( cor , format = format )
423 call disp%show("cor = getFilled(0., ndim, ndim)")
424 cor = getFilled(0., ndim, ndim)
425 call disp%show("call setCor(cor, uppDia, samShifted, dim, rweight, rweisum) ! same result as above.")
426 call setCor(cor, uppDia, samShifted, dim, rweight, rweisum)
427 call disp%show("cor")
428 call disp%show( cor , format = format )
429 call disp%skip()
430 call disp%show("cor = getFilled(0., ndim, ndim)")
431 cor = getFilled(0., ndim, ndim)
432 call disp%show("call setCor(cor, lowDia, mean, sample, dim, rweight, rweisum)")
433 call setCor(cor, lowDia, mean, sample, dim, rweight, rweisum)
434 call disp%show("cor")
435 call disp%show( cor , format = format )
436 call disp%show("cor = getFilled(0., ndim, ndim)")
437 cor = getFilled(0., ndim, ndim)
438 call disp%show("call setCor(cor, lowDia, samShifted, dim, rweight, rweisum) ! same result as above.")
439 call setCor(cor, lowDia, samShifted, dim, rweight, rweisum)
440 call disp%show("cor")
441 call disp%show( cor , format = format )
442 call disp%skip()
443 call disp%show("Compute the sample correlation along the first dimension.", deliml = SK_'''')
444 call disp%skip()
445 call disp%show("dim = 1")
446 dim = 1
447 call disp%show("cor = getFilled(0., ndim, ndim)")
448 cor = getFilled(0., ndim, ndim)
449 call disp%show("call setCor(cor, uppDia, mean, transpose(sample), dim, rweight, rweisum)")
450 call setCor(cor, uppDia, mean, transpose(sample), dim, rweight, rweisum)
451 call disp%show("cor")
452 call disp%show( cor , format = format )
453 call disp%show("cor = getFilled(0., ndim, ndim)")
454 cor = getFilled(0., ndim, ndim)
455 call disp%show("call setCor(cor, uppDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.")
456 call setCor(cor, uppDia, transpose(samShifted), dim, rweight, rweisum)
457 call disp%show("cor")
458 call disp%show( cor , format = format )
459 call disp%skip()
460 call disp%show("cor = getFilled(0., ndim, ndim)")
461 cor = getFilled(0., ndim, ndim)
462 call disp%show("call setCor(cor, lowDia, mean, transpose(sample), dim, rweight, rweisum)")
463 call setCor(cor, lowDia, mean, transpose(sample), dim, rweight, rweisum)
464 call disp%show("cor")
465 call disp%show( cor , format = format )
466 call disp%show("cor = getFilled(0., ndim, ndim)")
467 cor = getFilled(0., ndim, ndim)
468 call disp%show("call setCor(cor, lowDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.")
469 call setCor(cor, lowDia, transpose(samShifted), dim, rweight, rweisum)
470 call disp%show("cor")
471 call disp%show( cor , format = format )
472 call disp%skip()
473 call disp%show("Compute the full sample correlation for a pair of time series.", deliml = SK_'''')
474 call disp%skip()
475 call disp%show("call setCor(cor(1,1), mean, sample(1,:), sample(2,:), rweight, rweisum)")
476 call setCor(cor(1,1), mean, sample(1,:), sample(2,:), rweight, rweisum)
477 call disp%show("cor(1,1)")
478 call disp%show( cor(1,1) , format = format )
479 call disp%show("call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), rweight, rweisum) ! same result as above.")
480 call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), rweight, rweisum)
481 call disp%show("cor(1,1)")
482 call disp%show( cor(1,1) , format = format )
483 call disp%skip()
484 end block
485
486end program example
Generate and return an array of the specified rank and shape of arbitrary intrinsic type and kind wit...
Generate minimally-spaced character, integer, real sequences or sequences at fixed intervals of size ...
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate an equally-weighted (verbose or flattened) array of the input weighted array of rank 1 or 2.
Generate and return a random positive-definite (correlation or covariance) matrix using the Gram meth...
Definition: pm_distCov.F90:394
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Generate and return a generic or type/kind-specific IO format with the requested specifications that ...
Definition: pm_io.F90:18485
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 the (weighted) mean of an input sample of nsam observations with ndim = 1 or 2 at...
Return the (weighted) mean of a pair of time series or of an input sample of nsam observations with n...
Generate a sample of shape (nsam), or (ndim, nsam) or (nsam, ndim) that is shifted by the specified i...
This module contains procedures and generic interfaces for convenient allocation and filling of array...
This module contains procedures and generic interfaces for generating ranges of discrete character,...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
This module contains procedures and generic interfaces for flattening (duplicating the elements of) a...
This module contains classes and procedures for generating random matrices distributed on the space o...
Definition: pm_distCov.F90:72
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 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 classes and procedures for computing the first moment (i.e., the statistical mea...
This module contains classes and procedures for shifting univariate or multivariate samples by arbitr...
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
2!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
3!Compute the correlation matrix from correlation matrix.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = getUnifRand(0, 4)
8ndim
9+0
10std = getUnifRand(1, ndim, ndim)
11std
12
13cov = getCovRand(1._TKG, std)
14cov
15
16
17cor = cov
18call setCor(cor, uppDia, cor, uppDia)
19cor
20
21cor = cov
22call setCor(cor, uppDia, cor, lowDia)
23cor
24
25cor = cov
26call setCor(cor, upp, cor, uppDia)
27cor
28
29cor = cov
30call setCor(cor, upp, cor, lowDia)
31cor
32
33cor = cov
34call setCor(cor, lowDia, cor, uppDia)
35cor
36
37cor = cov
38call setCor(cor, lowDia, cor, lowDia)
39cor
40
41cor = cov
42call setCor(cor, low, cor, uppDia)
43cor
44
45cor = cov
46call setCor(cor, low, cor, lowDia)
47cor
48
49
50cor = cov
51call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
52cor
53
54cor = cov
55call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
56cor
57
58cor = cov
59call setCor(cor, upp, cor, upp, stdinv = 1 / std)
60cor
61
62cor = cov
63call setCor(cor, upp, cor, low, stdinv = 1 / std)
64cor
65
66cor = cov
67call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
68cor
69
70cor = cov
71call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
72cor
73
74cor = cov
75call setCor(cor, low, cor, upp, stdinv = 1 / std)
76cor
77
78cor = cov
79call setCor(cor, low, cor, low, stdinv = 1 / std)
80cor
81
82
83ndim = getUnifRand(0, 4)
84ndim
85+4
86std = getUnifRand(1, ndim, ndim)
87std
88+2.000000E+00, +1.000000E+00, +3.000000E+00, +1.000000E+00
89cov = getCovRand(1._TKG, std)
90cov
91+4.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
92-1.033379E+00, +9.999999E-01, +1.892676E+00, -3.661638E-01
93-5.841490E+00, +1.892676E+00, +9.000002E+00, -1.619798E+00
94+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
95
96
97cor = cov
98call setCor(cor, uppDia, cor, uppDia)
99cor
100+1.000000E+00, -5.166894E-01, -9.735816E-01, +6.258060E-01
101-1.033379E+00, +1.000000E+00, +6.308920E-01, -3.661639E-01
102-5.841490E+00, +1.892676E+00, +1.000000E+00, -5.399325E-01
103+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
104
105cor = cov
106call setCor(cor, uppDia, cor, lowDia)
107cor
108+1.000000E+00, -5.166894E-01, -9.735816E-01, +6.258060E-01
109-1.033379E+00, +1.000000E+00, +6.308920E-01, -3.661639E-01
110-5.841490E+00, +1.892676E+00, +1.000000E+00, -5.399325E-01
111+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
112
113cor = cov
114call setCor(cor, upp, cor, uppDia)
115cor
116+4.000000E+00, -5.166894E-01, -9.735816E-01, +6.258060E-01
117-1.033379E+00, +9.999999E-01, +6.308920E-01, -3.661639E-01
118-5.841490E+00, +1.892676E+00, +9.000002E+00, -5.399325E-01
119+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
120
121cor = cov
122call setCor(cor, upp, cor, lowDia)
123cor
124+4.000000E+00, -5.166894E-01, -9.735816E-01, +6.258060E-01
125-1.033379E+00, +9.999999E-01, +6.308920E-01, -3.661639E-01
126-5.841490E+00, +1.892676E+00, +9.000002E+00, -5.399325E-01
127+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
128
129cor = cov
130call setCor(cor, lowDia, cor, uppDia)
131cor
132+1.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
133-1.033379E+00, +1.000000E+00, +1.892676E+00, -3.661638E-01
134-5.841490E+00, +1.892676E+00, +1.000000E+00, -1.619798E+00
135+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
136
137cor = cov
138call setCor(cor, lowDia, cor, lowDia)
139cor
140+1.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
141-5.166894E-01, +1.000000E+00, +1.892676E+00, -3.661638E-01
142-9.735816E-01, +6.308920E-01, +1.000000E+00, -1.619798E+00
143+6.258060E-01, -3.661639E-01, -5.399325E-01, +1.000000E+00
144
145cor = cov
146call setCor(cor, low, cor, uppDia)
147cor
148+4.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
149-5.166894E-01, +9.999999E-01, +1.892676E+00, -3.661638E-01
150-9.735816E-01, +6.308920E-01, +9.000002E+00, -1.619798E+00
151+6.258060E-01, -3.661639E-01, -5.399325E-01, +1.000000E+00
152
153cor = cov
154call setCor(cor, low, cor, lowDia)
155cor
156+4.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
157-5.166894E-01, +9.999999E-01, +1.892676E+00, -3.661638E-01
158-9.735816E-01, +6.308920E-01, +9.000002E+00, -1.619798E+00
159+6.258060E-01, -3.661639E-01, -5.399325E-01, +1.000000E+00
160
161
162cor = cov
163call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
164cor
165+1.000000E+00, -5.166894E-01, -9.735817E-01, +6.258060E-01
166-1.033379E+00, +1.000000E+00, +6.308920E-01, -3.661638E-01
167-5.841490E+00, +1.892676E+00, +1.000000E+00, -5.399326E-01
168+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
169
170cor = cov
171call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
172cor
173+1.000000E+00, -5.166894E-01, -9.735817E-01, +6.258060E-01
174-1.033379E+00, +1.000000E+00, +6.308920E-01, -3.661638E-01
175-5.841490E+00, +1.892676E+00, +1.000000E+00, -5.399326E-01
176+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
177
178cor = cov
179call setCor(cor, upp, cor, upp, stdinv = 1 / std)
180cor
181+4.00000000, -0.516689360, -0.973581672, +0.625806034
182-1.03337872, +0.999999940, +0.630891979, -0.366163850
183-5.84148979, +1.89267588, +9.00000191, -0.539932609
184+1.25161207, -0.366163850, -1.61979783, +1.00000000
185
186cor = cov
187call setCor(cor, upp, cor, low, stdinv = 1 / std)
188cor
189+4.000000E+00, -5.166894E-01, -9.735817E-01, +6.258060E-01
190-1.033379E+00, +9.999999E-01, +6.308920E-01, -3.661638E-01
191-5.841490E+00, +1.892676E+00, +9.000002E+00, -5.399326E-01
192+1.251612E+00, -3.661638E-01, -1.619798E+00, +1.000000E+00
193
194cor = cov
195call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
196cor
197+1.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
198-5.166894E-01, +1.000000E+00, +1.892676E+00, -3.661638E-01
199-9.735817E-01, +6.308920E-01, +1.000000E+00, -1.619798E+00
200+6.258060E-01, -3.661638E-01, -5.399326E-01, +1.000000E+00
201
202cor = cov
203call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
204cor
205+1.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
206-5.166894E-01, +1.000000E+00, +1.892676E+00, -3.661638E-01
207-9.735817E-01, +6.308920E-01, +1.000000E+00, -1.619798E+00
208+6.258060E-01, -3.661638E-01, -5.399326E-01, +1.000000E+00
209
210cor = cov
211call setCor(cor, low, cor, upp, stdinv = 1 / std)
212cor
213+4.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
214-5.166894E-01, +9.999999E-01, +1.892676E+00, -3.661638E-01
215-9.735817E-01, +6.308920E-01, +9.000002E+00, -1.619798E+00
216+6.258060E-01, -3.661638E-01, -5.399326E-01, +1.000000E+00
217
218cor = cov
219call setCor(cor, low, cor, low, stdinv = 1 / std)
220cor
221+4.000000E+00, -1.033379E+00, -5.841490E+00, +1.251612E+00
222-5.166894E-01, +9.999999E-01, +1.892676E+00, -3.661638E-01
223-9.735817E-01, +6.308920E-01, +9.000002E+00, -1.619798E+00
224+6.258060E-01, -3.661638E-01, -5.399326E-01, +1.000000E+00
225
226
227ndim = getUnifRand(0, 4)
228ndim
229+0
230std = getUnifRand(1, ndim, ndim)
231std
232
233cov = getCovRand(1._TKG, std)
234cov
235
236
237cor = cov
238call setCor(cor, uppDia, cor, uppDia)
239cor
240
241cor = cov
242call setCor(cor, uppDia, cor, lowDia)
243cor
244
245cor = cov
246call setCor(cor, upp, cor, uppDia)
247cor
248
249cor = cov
250call setCor(cor, upp, cor, lowDia)
251cor
252
253cor = cov
254call setCor(cor, lowDia, cor, uppDia)
255cor
256
257cor = cov
258call setCor(cor, lowDia, cor, lowDia)
259cor
260
261cor = cov
262call setCor(cor, low, cor, uppDia)
263cor
264
265cor = cov
266call setCor(cor, low, cor, lowDia)
267cor
268
269
270cor = cov
271call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
272cor
273
274cor = cov
275call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
276cor
277
278cor = cov
279call setCor(cor, upp, cor, upp, stdinv = 1 / std)
280cor
281
282cor = cov
283call setCor(cor, upp, cor, low, stdinv = 1 / std)
284cor
285
286cor = cov
287call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
288cor
289
290cor = cov
291call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
292cor
293
294cor = cov
295call setCor(cor, low, cor, upp, stdinv = 1 / std)
296cor
297
298cor = cov
299call setCor(cor, low, cor, low, stdinv = 1 / std)
300cor
301
302
303ndim = getUnifRand(0, 4)
304ndim
305+4
306std = getUnifRand(1, ndim, ndim)
307std
308+2.000000E+00, +3.000000E+00, +4.000000E+00, +4.000000E+00
309cov = getCovRand(1._TKG, std)
310cov
311+4.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
312-2.730345E+00, +9.000000E+00, +8.524410E+00, +2.967734E+00
313-2.237197E+00, +8.524410E+00, +1.600000E+01, -3.488663E+00
314-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.600000E+01
315
316
317cor = cov
318call setCor(cor, uppDia, cor, uppDia)
319cor
320+1.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
321-2.730345E+00, +1.000000E+00, +7.103676E-01, +2.473111E-01
322-2.237197E+00, +8.524410E+00, +1.000000E+00, -2.180414E-01
323-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.000000E+00
324
325cor = cov
326call setCor(cor, uppDia, cor, lowDia)
327cor
328+1.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
329-2.730345E+00, +1.000000E+00, +7.103676E-01, +2.473111E-01
330-2.237197E+00, +8.524410E+00, +1.000000E+00, -2.180414E-01
331-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.000000E+00
332
333cor = cov
334call setCor(cor, upp, cor, uppDia)
335cor
336+4.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
337-2.730345E+00, +9.000000E+00, +7.103676E-01, +2.473111E-01
338-2.237197E+00, +8.524410E+00, +1.600000E+01, -2.180414E-01
339-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.600000E+01
340
341cor = cov
342call setCor(cor, upp, cor, lowDia)
343cor
344+4.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
345-2.730345E+00, +9.000000E+00, +7.103676E-01, +2.473111E-01
346-2.237197E+00, +8.524410E+00, +1.600000E+01, -2.180414E-01
347-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.600000E+01
348
349cor = cov
350call setCor(cor, lowDia, cor, uppDia)
351cor
352+1.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
353-2.730345E+00, +1.000000E+00, +8.524410E+00, +2.967734E+00
354-2.237197E+00, +8.524410E+00, +1.000000E+00, -3.488663E+00
355-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.000000E+00
356
357cor = cov
358call setCor(cor, lowDia, cor, lowDia)
359cor
360+1.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
361-4.550575E-01, +1.000000E+00, +8.524410E+00, +2.967734E+00
362-2.796496E-01, +7.103676E-01, +1.000000E+00, -3.488663E+00
363-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.000000E+00
364
365cor = cov
366call setCor(cor, low, cor, uppDia)
367cor
368+4.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
369-4.550575E-01, +9.000000E+00, +8.524410E+00, +2.967734E+00
370-2.796496E-01, +7.103676E-01, +1.600000E+01, -3.488663E+00
371-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.600000E+01
372
373cor = cov
374call setCor(cor, low, cor, lowDia)
375cor
376+4.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
377-4.550575E-01, +9.000000E+00, +8.524410E+00, +2.967734E+00
378-2.796496E-01, +7.103676E-01, +1.600000E+01, -3.488663E+00
379-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.600000E+01
380
381
382cor = cov
383call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
384cor
385+1.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
386-2.730345E+00, +1.000000E+00, +7.103676E-01, +2.473111E-01
387-2.237197E+00, +8.524410E+00, +1.000000E+00, -2.180414E-01
388-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.000000E+00
389
390cor = cov
391call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
392cor
393+1.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
394-2.730345E+00, +1.000000E+00, +7.103676E-01, +2.473111E-01
395-2.237197E+00, +8.524410E+00, +1.000000E+00, -2.180414E-01
396-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.000000E+00
397
398cor = cov
399call setCor(cor, upp, cor, upp, stdinv = 1 / std)
400cor
401+4.00000000, -0.455057472, -0.279649615, -0.423663676
402-2.73034477, +9.00000000, +0.710367560, +0.247311145
403-2.23719692, +8.52441025, +15.9999990, -0.218041420
404-3.38930941, +2.96773362, -3.48866272, +16.0000000
405
406cor = cov
407call setCor(cor, upp, cor, low, stdinv = 1 / std)
408cor
409+4.000000E+00, -4.550575E-01, -2.796496E-01, -4.236637E-01
410-2.730345E+00, +9.000000E+00, +7.103676E-01, +2.473111E-01
411-2.237197E+00, +8.524410E+00, +1.600000E+01, -2.180414E-01
412-3.389309E+00, +2.967734E+00, -3.488663E+00, +1.600000E+01
413
414cor = cov
415call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
416cor
417+1.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
418-4.550575E-01, +1.000000E+00, +8.524410E+00, +2.967734E+00
419-2.796496E-01, +7.103676E-01, +1.000000E+00, -3.488663E+00
420-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.000000E+00
421
422cor = cov
423call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
424cor
425+1.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
426-4.550575E-01, +1.000000E+00, +8.524410E+00, +2.967734E+00
427-2.796496E-01, +7.103676E-01, +1.000000E+00, -3.488663E+00
428-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.000000E+00
429
430cor = cov
431call setCor(cor, low, cor, upp, stdinv = 1 / std)
432cor
433+4.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
434-4.550575E-01, +9.000000E+00, +8.524410E+00, +2.967734E+00
435-2.796496E-01, +7.103676E-01, +1.600000E+01, -3.488663E+00
436-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.600000E+01
437
438cor = cov
439call setCor(cor, low, cor, low, stdinv = 1 / std)
440cor
441+4.000000E+00, -2.730345E+00, -2.237197E+00, -3.389309E+00
442-4.550575E-01, +9.000000E+00, +8.524410E+00, +2.967734E+00
443-2.796496E-01, +7.103676E-01, +1.600000E+01, -3.488663E+00
444-4.236637E-01, +2.473111E-01, -2.180414E-01, +1.600000E+01
445
446
447ndim = getUnifRand(0, 4)
448ndim
449+1
450std = getUnifRand(1, ndim, ndim)
451std
452+1.000000E+00
453cov = getCovRand(1._TKG, std)
454cov
455+1.000000E+00
456
457
458cor = cov
459call setCor(cor, uppDia, cor, uppDia)
460cor
461+1.000000E+00
462
463cor = cov
464call setCor(cor, uppDia, cor, lowDia)
465cor
466+1.000000E+00
467
468cor = cov
469call setCor(cor, upp, cor, uppDia)
470cor
471+1.000000E+00
472
473cor = cov
474call setCor(cor, upp, cor, lowDia)
475cor
476+1.000000E+00
477
478cor = cov
479call setCor(cor, lowDia, cor, uppDia)
480cor
481+1.000000E+00
482
483cor = cov
484call setCor(cor, lowDia, cor, lowDia)
485cor
486+1.000000E+00
487
488cor = cov
489call setCor(cor, low, cor, uppDia)
490cor
491+1.000000E+00
492
493cor = cov
494call setCor(cor, low, cor, lowDia)
495cor
496+1.000000E+00
497
498
499cor = cov
500call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
501cor
502+1.000000E+00
503
504cor = cov
505call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
506cor
507+1.000000E+00
508
509cor = cov
510call setCor(cor, upp, cor, upp, stdinv = 1 / std)
511cor
512+1.00000000
513
514cor = cov
515call setCor(cor, upp, cor, low, stdinv = 1 / std)
516cor
517+1.000000E+00
518
519cor = cov
520call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
521cor
522+1.000000E+00
523
524cor = cov
525call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
526cor
527+1.000000E+00
528
529cor = cov
530call setCor(cor, low, cor, upp, stdinv = 1 / std)
531cor
532+1.000000E+00
533
534cor = cov
535call setCor(cor, low, cor, low, stdinv = 1 / std)
536cor
537+1.000000E+00
538
539
540ndim = getUnifRand(0, 4)
541ndim
542+1
543std = getUnifRand(1, ndim, ndim)
544std
545+1.000000E+00
546cov = getCovRand(1._TKG, std)
547cov
548+1.000000E+00
549
550
551cor = cov
552call setCor(cor, uppDia, cor, uppDia)
553cor
554+1.000000E+00
555
556cor = cov
557call setCor(cor, uppDia, cor, lowDia)
558cor
559+1.000000E+00
560
561cor = cov
562call setCor(cor, upp, cor, uppDia)
563cor
564+1.000000E+00
565
566cor = cov
567call setCor(cor, upp, cor, lowDia)
568cor
569+1.000000E+00
570
571cor = cov
572call setCor(cor, lowDia, cor, uppDia)
573cor
574+1.000000E+00
575
576cor = cov
577call setCor(cor, lowDia, cor, lowDia)
578cor
579+1.000000E+00
580
581cor = cov
582call setCor(cor, low, cor, uppDia)
583cor
584+1.000000E+00
585
586cor = cov
587call setCor(cor, low, cor, lowDia)
588cor
589+1.000000E+00
590
591
592cor = cov
593call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
594cor
595+1.000000E+00
596
597cor = cov
598call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
599cor
600+1.000000E+00
601
602cor = cov
603call setCor(cor, upp, cor, upp, stdinv = 1 / std)
604cor
605+1.00000000
606
607cor = cov
608call setCor(cor, upp, cor, low, stdinv = 1 / std)
609cor
610+1.000000E+00
611
612cor = cov
613call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
614cor
615+1.000000E+00
616
617cor = cov
618call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
619cor
620+1.000000E+00
621
622cor = cov
623call setCor(cor, low, cor, upp, stdinv = 1 / std)
624cor
625+1.000000E+00
626
627cor = cov
628call setCor(cor, low, cor, low, stdinv = 1 / std)
629cor
630+1.000000E+00
631
632
633ndim = getUnifRand(0, 4)
634ndim
635+2
636std = getUnifRand(1, ndim, ndim)
637std
638+1.000000E+00, +2.000000E+00
639cov = getCovRand(1._TKG, std)
640cov
641+1.000000E+00, -1.999785E+00
642-1.999785E+00, +4.000000E+00
643
644
645cor = cov
646call setCor(cor, uppDia, cor, uppDia)
647cor
648+1.000000E+00, -9.998927E-01
649-1.999785E+00, +1.000000E+00
650
651cor = cov
652call setCor(cor, uppDia, cor, lowDia)
653cor
654+1.000000E+00, -9.998927E-01
655-1.999785E+00, +1.000000E+00
656
657cor = cov
658call setCor(cor, upp, cor, uppDia)
659cor
660+1.000000E+00, -9.998927E-01
661-1.999785E+00, +4.000000E+00
662
663cor = cov
664call setCor(cor, upp, cor, lowDia)
665cor
666+1.000000E+00, -9.998927E-01
667-1.999785E+00, +4.000000E+00
668
669cor = cov
670call setCor(cor, lowDia, cor, uppDia)
671cor
672+1.000000E+00, -1.999785E+00
673-1.999785E+00, +1.000000E+00
674
675cor = cov
676call setCor(cor, lowDia, cor, lowDia)
677cor
678+1.000000E+00, -1.999785E+00
679-9.998927E-01, +1.000000E+00
680
681cor = cov
682call setCor(cor, low, cor, uppDia)
683cor
684+1.000000E+00, -1.999785E+00
685-9.998927E-01, +4.000000E+00
686
687cor = cov
688call setCor(cor, low, cor, lowDia)
689cor
690+1.000000E+00, -1.999785E+00
691-9.998927E-01, +4.000000E+00
692
693
694cor = cov
695call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
696cor
697+1.000000E+00, -9.998927E-01
698-1.999785E+00, +1.000000E+00
699
700cor = cov
701call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
702cor
703+1.000000E+00, -9.998927E-01
704-1.999785E+00, +1.000000E+00
705
706cor = cov
707call setCor(cor, upp, cor, upp, stdinv = 1 / std)
708cor
709+1.00000000, -0.999892652
710-1.99978530, +4.00000000
711
712cor = cov
713call setCor(cor, upp, cor, low, stdinv = 1 / std)
714cor
715+1.000000E+00, -9.998927E-01
716-1.999785E+00, +4.000000E+00
717
718cor = cov
719call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
720cor
721+1.000000E+00, -1.999785E+00
722-9.998927E-01, +1.000000E+00
723
724cor = cov
725call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
726cor
727+1.000000E+00, -1.999785E+00
728-9.998927E-01, +1.000000E+00
729
730cor = cov
731call setCor(cor, low, cor, upp, stdinv = 1 / std)
732cor
733+1.000000E+00, -1.999785E+00
734-9.998927E-01, +4.000000E+00
735
736cor = cov
737call setCor(cor, low, cor, low, stdinv = 1 / std)
738cor
739+1.000000E+00, -1.999785E+00
740-9.998927E-01, +4.000000E+00
741
742
743ndim = getUnifRand(0, 4)
744ndim
745+2
746std = getUnifRand(1, ndim, ndim)
747std
748+2.000000E+00, +2.000000E+00
749cov = getCovRand(1._TKG, std)
750cov
751+4.000000E+00, +2.328817E+00
752+2.328817E+00, +4.000000E+00
753
754
755cor = cov
756call setCor(cor, uppDia, cor, uppDia)
757cor
758+1.000000E+00, +5.822042E-01
759+2.328817E+00, +1.000000E+00
760
761cor = cov
762call setCor(cor, uppDia, cor, lowDia)
763cor
764+1.000000E+00, +5.822042E-01
765+2.328817E+00, +1.000000E+00
766
767cor = cov
768call setCor(cor, upp, cor, uppDia)
769cor
770+4.000000E+00, +5.822042E-01
771+2.328817E+00, +4.000000E+00
772
773cor = cov
774call setCor(cor, upp, cor, lowDia)
775cor
776+4.000000E+00, +5.822042E-01
777+2.328817E+00, +4.000000E+00
778
779cor = cov
780call setCor(cor, lowDia, cor, uppDia)
781cor
782+1.000000E+00, +2.328817E+00
783+2.328817E+00, +1.000000E+00
784
785cor = cov
786call setCor(cor, lowDia, cor, lowDia)
787cor
788+1.000000E+00, +2.328817E+00
789+5.822042E-01, +1.000000E+00
790
791cor = cov
792call setCor(cor, low, cor, uppDia)
793cor
794+4.000000E+00, +2.328817E+00
795+5.822042E-01, +4.000000E+00
796
797cor = cov
798call setCor(cor, low, cor, lowDia)
799cor
800+4.000000E+00, +2.328817E+00
801+5.822042E-01, +4.000000E+00
802
803
804cor = cov
805call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
806cor
807+1.000000E+00, +5.822042E-01
808+2.328817E+00, +1.000000E+00
809
810cor = cov
811call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
812cor
813+1.000000E+00, +5.822042E-01
814+2.328817E+00, +1.000000E+00
815
816cor = cov
817call setCor(cor, upp, cor, upp, stdinv = 1 / std)
818cor
819+4.00000000, +0.582204163
820+2.32881665, +4.00000048
821
822cor = cov
823call setCor(cor, upp, cor, low, stdinv = 1 / std)
824cor
825+4.000000E+00, +5.822042E-01
826+2.328817E+00, +4.000000E+00
827
828cor = cov
829call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
830cor
831+1.000000E+00, +2.328817E+00
832+5.822042E-01, +1.000000E+00
833
834cor = cov
835call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
836cor
837+1.000000E+00, +2.328817E+00
838+5.822042E-01, +1.000000E+00
839
840cor = cov
841call setCor(cor, low, cor, upp, stdinv = 1 / std)
842cor
843+4.000000E+00, +2.328817E+00
844+5.822042E-01, +4.000000E+00
845
846cor = cov
847call setCor(cor, low, cor, low, stdinv = 1 / std)
848cor
849+4.000000E+00, +2.328817E+00
850+5.822042E-01, +4.000000E+00
851
852
853ndim = getUnifRand(0, 4)
854ndim
855+2
856std = getUnifRand(1, ndim, ndim)
857std
858+1.000000E+00, +1.000000E+00
859cov = getCovRand(1._TKG, std)
860cov
861+9.999999E-01, +7.193822E-01
862+7.193822E-01, +1.000000E+00
863
864
865cor = cov
866call setCor(cor, uppDia, cor, uppDia)
867cor
868+1.000000E+00, +7.193822E-01
869+7.193822E-01, +1.000000E+00
870
871cor = cov
872call setCor(cor, uppDia, cor, lowDia)
873cor
874+1.000000E+00, +7.193822E-01
875+7.193822E-01, +1.000000E+00
876
877cor = cov
878call setCor(cor, upp, cor, uppDia)
879cor
880+9.999999E-01, +7.193822E-01
881+7.193822E-01, +1.000000E+00
882
883cor = cov
884call setCor(cor, upp, cor, lowDia)
885cor
886+9.999999E-01, +7.193822E-01
887+7.193822E-01, +1.000000E+00
888
889cor = cov
890call setCor(cor, lowDia, cor, uppDia)
891cor
892+1.000000E+00, +7.193822E-01
893+7.193822E-01, +1.000000E+00
894
895cor = cov
896call setCor(cor, lowDia, cor, lowDia)
897cor
898+1.000000E+00, +7.193822E-01
899+7.193822E-01, +1.000000E+00
900
901cor = cov
902call setCor(cor, low, cor, uppDia)
903cor
904+9.999999E-01, +7.193822E-01
905+7.193822E-01, +1.000000E+00
906
907cor = cov
908call setCor(cor, low, cor, lowDia)
909cor
910+9.999999E-01, +7.193822E-01
911+7.193822E-01, +1.000000E+00
912
913
914cor = cov
915call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
916cor
917+1.000000E+00, +7.193822E-01
918+7.193822E-01, +1.000000E+00
919
920cor = cov
921call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
922cor
923+1.000000E+00, +7.193822E-01
924+7.193822E-01, +1.000000E+00
925
926cor = cov
927call setCor(cor, upp, cor, upp, stdinv = 1 / std)
928cor
929+0.999999881, +0.719382167
930+0.719382167, +1.00000000
931
932cor = cov
933call setCor(cor, upp, cor, low, stdinv = 1 / std)
934cor
935+9.999999E-01, +7.193822E-01
936+7.193822E-01, +1.000000E+00
937
938cor = cov
939call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
940cor
941+1.000000E+00, +7.193822E-01
942+7.193822E-01, +1.000000E+00
943
944cor = cov
945call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
946cor
947+1.000000E+00, +7.193822E-01
948+7.193822E-01, +1.000000E+00
949
950cor = cov
951call setCor(cor, low, cor, upp, stdinv = 1 / std)
952cor
953+9.999999E-01, +7.193822E-01
954+7.193822E-01, +1.000000E+00
955
956cor = cov
957call setCor(cor, low, cor, low, stdinv = 1 / std)
958cor
959+9.999999E-01, +7.193822E-01
960+7.193822E-01, +1.000000E+00
961
962
963ndim = getUnifRand(0, 4)
964ndim
965+4
966std = getUnifRand(1, ndim, ndim)
967std
968+1.000000E+00, +4.000000E+00, +3.000000E+00, +1.000000E+00
969cov = getCovRand(1._TKG, std)
970cov
971+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
972-1.485865E+00, +1.600000E+01, -5.481150E+00, +3.242330E+00
973+2.808865E+00, -5.481150E+00, +9.000002E+00, +1.088034E-01
974-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
975
976
977cor = cov
978call setCor(cor, uppDia, cor, uppDia)
979cor
980+1.000000E+00, -3.714662E-01, +9.362884E-01, -3.269074E-02
981-1.485865E+00, +1.000000E+00, -4.567625E-01, +8.105825E-01
982+2.808865E+00, -5.481150E+00, +1.000000E+00, +3.626780E-02
983-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
984
985cor = cov
986call setCor(cor, uppDia, cor, lowDia)
987cor
988+1.000000E+00, -3.714662E-01, +9.362884E-01, -3.269074E-02
989-1.485865E+00, +1.000000E+00, -4.567625E-01, +8.105825E-01
990+2.808865E+00, -5.481150E+00, +1.000000E+00, +3.626780E-02
991-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
992
993cor = cov
994call setCor(cor, upp, cor, uppDia)
995cor
996+1.000000E+00, -3.714662E-01, +9.362884E-01, -3.269074E-02
997-1.485865E+00, +1.600000E+01, -4.567625E-01, +8.105825E-01
998+2.808865E+00, -5.481150E+00, +9.000002E+00, +3.626780E-02
999-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
1000
1001cor = cov
1002call setCor(cor, upp, cor, lowDia)
1003cor
1004+1.000000E+00, -3.714662E-01, +9.362884E-01, -3.269074E-02
1005-1.485865E+00, +1.600000E+01, -4.567625E-01, +8.105825E-01
1006+2.808865E+00, -5.481150E+00, +9.000002E+00, +3.626780E-02
1007-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
1008
1009cor = cov
1010call setCor(cor, lowDia, cor, uppDia)
1011cor
1012+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1013-1.485865E+00, +1.000000E+00, -5.481150E+00, +3.242330E+00
1014+2.808865E+00, -5.481150E+00, +1.000000E+00, +1.088034E-01
1015-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
1016
1017cor = cov
1018call setCor(cor, lowDia, cor, lowDia)
1019cor
1020+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1021-3.714662E-01, +1.000000E+00, -5.481150E+00, +3.242330E+00
1022+9.362884E-01, -4.567625E-01, +1.000000E+00, +1.088034E-01
1023-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1024
1025cor = cov
1026call setCor(cor, low, cor, uppDia)
1027cor
1028+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1029-3.714662E-01, +1.600000E+01, -5.481150E+00, +3.242330E+00
1030+9.362884E-01, -4.567625E-01, +9.000002E+00, +1.088034E-01
1031-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1032
1033cor = cov
1034call setCor(cor, low, cor, lowDia)
1035cor
1036+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1037-3.714662E-01, +1.600000E+01, -5.481150E+00, +3.242330E+00
1038+9.362884E-01, -4.567625E-01, +9.000002E+00, +1.088034E-01
1039-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1040
1041
1042cor = cov
1043call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
1044cor
1045+1.000000E+00, -3.714662E-01, +9.362885E-01, -3.269074E-02
1046-1.485865E+00, +1.000000E+00, -4.567625E-01, +8.105825E-01
1047+2.808865E+00, -5.481150E+00, +1.000000E+00, +3.626780E-02
1048-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
1049
1050cor = cov
1051call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
1052cor
1053+1.000000E+00, -3.714662E-01, +9.362885E-01, -3.269074E-02
1054-1.485865E+00, +1.000000E+00, -4.567625E-01, +8.105825E-01
1055+2.808865E+00, -5.481150E+00, +1.000000E+00, +3.626780E-02
1056-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
1057
1058cor = cov
1059call setCor(cor, upp, cor, upp, stdinv = 1 / std)
1060cor
1061+1.00000000, -0.371466160, +0.936288476, -0.326907411E-1
1062-1.48586464, +16.0000000, -0.456762522, +0.810582519
1063+2.80886531, -5.48115015, +9.00000191, +0.362678021E-1
1064-0.326907411E-1, +3.24233007, +0.108803406, +1.00000000
1065
1066cor = cov
1067call setCor(cor, upp, cor, low, stdinv = 1 / std)
1068cor
1069+1.000000E+00, -3.714662E-01, +9.362885E-01, -3.269074E-02
1070-1.485865E+00, +1.600000E+01, -4.567625E-01, +8.105825E-01
1071+2.808865E+00, -5.481150E+00, +9.000002E+00, +3.626780E-02
1072-3.269074E-02, +3.242330E+00, +1.088034E-01, +1.000000E+00
1073
1074cor = cov
1075call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
1076cor
1077+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1078-3.714662E-01, +1.000000E+00, -5.481150E+00, +3.242330E+00
1079+9.362885E-01, -4.567625E-01, +1.000000E+00, +1.088034E-01
1080-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1081
1082cor = cov
1083call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
1084cor
1085+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1086-3.714662E-01, +1.000000E+00, -5.481150E+00, +3.242330E+00
1087+9.362885E-01, -4.567625E-01, +1.000000E+00, +1.088034E-01
1088-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1089
1090cor = cov
1091call setCor(cor, low, cor, upp, stdinv = 1 / std)
1092cor
1093+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1094-3.714662E-01, +1.600000E+01, -5.481150E+00, +3.242330E+00
1095+9.362885E-01, -4.567625E-01, +9.000002E+00, +1.088034E-01
1096-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1097
1098cor = cov
1099call setCor(cor, low, cor, low, stdinv = 1 / std)
1100cor
1101+1.000000E+00, -1.485865E+00, +2.808865E+00, -3.269074E-02
1102-3.714662E-01, +1.600000E+01, -5.481150E+00, +3.242330E+00
1103+9.362885E-01, -4.567625E-01, +9.000002E+00, +1.088034E-01
1104-3.269074E-02, +8.105825E-01, +3.626780E-02, +1.000000E+00
1105
1106
1107!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1108!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1109!Compute the Pearson correlation matrix for a pair of time series.
1110!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1111!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1112
1113ndim = 2; nsam = 10; dim = 2
1114sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
1115sample
1116+7.000000E+00, +9.000000E+00, +1.800000E+01, +1.000000E+00, +4.000000E+00, +1.600000E+01, +5.000000E+00, +1.300000E+01, +2.000000E+01, +1.400000E+01
1117+4.000000E+00, +7.000000E+00, +9.000000E+00, +1.000000E+00, +1.900000E+01, +3.000000E+00, +1.600000E+01, +1.500000E+01, +1.200000E+01, +1.400000E+01
1118mean = getMean(sample, dim)
1119mean
1120+1.070000E+01, +1.000000E+01
1121samShifted = getShifted(sample, dim, -mean)
1122samShifted
1123-3.700000E+00, -1.700000E+00, +7.300000E+00, -9.700000E+00, -6.700000E+00, +5.300000E+00, -5.700000E+00, +2.300000E+00, +9.300000E+00, +3.300000E+00
1124-6.000000E+00, -3.000000E+00, -1.000000E+00, -9.000000E+00, +9.000000E+00, -7.000000E+00, +6.000000E+00, +5.000000E+00, +2.000000E+00, +4.000000E+00
1125cor = getFilled(0., ndim, ndim)
1126call setCor(cor, uppDia, mean, sample, dim)
1127cor
1128+1.000000E+00, +5.357539E-02
1129+0.000000E+00, +1.000000E+00
1130cor = getFilled(0., ndim, ndim)
1131call setCor(cor, uppDia, samShifted, dim) ! same result as above.
1132cor
1133+1.000000E+00, +5.357539E-02
1134+0.000000E+00, +1.000000E+00
1135
1136cor = getFilled(0., ndim, ndim)
1137call setCor(cor, lowDia, mean, sample, dim)
1138cor
1139+1.000000E+00, +0.000000E+00
1140+5.357540E-02, +1.000000E+00
1141cor = getFilled(0., ndim, ndim)
1142call setCor(cor, lowDia, samShifted, dim) ! same result as above.
1143cor
1144+1.000000E+00, +0.000000E+00
1145+5.357540E-02, +1.000000E+00
1146
1147'Compute the sample correlation along the first dimension.'
1148
1149dim = 1
1150cor = getFilled(0., ndim, ndim)
1151call setCor(cor, uppDia, mean, transpose(sample), dim)
1152cor
1153+1.000000E+00, +5.357539E-02
1154+0.000000E+00, +1.000000E+00
1155cor = getFilled(0., ndim, ndim)
1156call setCor(cor, uppDia, transpose(samShifted), dim) ! same result as above.
1157cor
1158+1.000000E+00, +5.357539E-02
1159+0.000000E+00, +1.000000E+00
1160
1161cor = getFilled(0., ndim, ndim)
1162call setCor(cor, lowDia, mean, transpose(sample), dim)
1163cor
1164+1.000000E+00, +0.000000E+00
1165+5.357540E-02, +1.000000E+00
1166cor = getFilled(0., ndim, ndim)
1167call setCor(cor, lowDia, transpose(samShifted), dim) ! same result as above.
1168cor
1169+1.000000E+00, +0.000000E+00
1170+5.357540E-02, +1.000000E+00
1171
1172'Compute the full sample correlation for a pair of time series.'
1173
1174call setCor(cor(1,1), mean, sample(1,:), sample(2,:))
1175cor(1,1)
1176+5.357540E-02
1177call setCor(cor(1,1), samShifted(1,:), samShifted(2,:)) ! same result as above.
1178cor(1,1)
1179+5.357540E-02
1180
1181
1182!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1183!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1184!Compute the Pearson correlation matrix for a weighted pair of time series.
1185!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1186!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1187
1188ndim = 2; nsam = 10; dim = 2
1189sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
1190sample
1191+1.900000E+01, +1.700000E+01, +3.000000E+00, +1.800000E+01, +1.200000E+01, +1.400000E+01, +8.000000E+00, +1.000000E+00, +1.300000E+01, +1.000000E+00
1192+1.000000E+00, +1.600000E+01, +9.000000E+00, +2.000000E+00, +1.000000E+01, +2.000000E+00, +1.200000E+01, +1.000000E+00, +9.000000E+00, +1.000000E+00
1193call setResized(mean, ndim)
1194iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
1195iweight
1196+4, +1, +7, +7, +3, +7, +1, +6, +1, +9
1197call setMean(mean, sample, dim, iweight, iweisum)
1198mean
1199+8.913044E+00, +3.847826E+00
1200iweisum
1201+46
1202rweight = iweight ! or real-valued weights.
1203iweight
1204+4, +1, +7, +7, +3, +7, +1, +6, +1, +9
1205call setMean(mean, sample, dim, rweight, rweisum)
1206mean
1207+8.913044E+00, +3.847826E+00
1208rweisum
1209+4.600000E+01
1210
1211!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1212!Compute the correlation matrix with integer weights.
1213!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1214
1215samShifted = getShifted(sample, dim, -mean)
1216samShifted
1217+1.008696E+01, +8.086956E+00, -5.913044E+00, +9.086956E+00, +3.086956E+00, +5.086956E+00, -9.130440E-01, -7.913044E+00, +4.086956E+00, -7.913044E+00
1218-2.847826E+00, +1.215217E+01, +5.152174E+00, -1.847826E+00, +6.152174E+00, -1.847826E+00, +8.152174E+00, -2.847826E+00, +5.152174E+00, -2.847826E+00
1219cor = getFilled(0., ndim, ndim)
1220call setCor(cor, uppDia, mean, sample, dim, iweight, iweisum)
1221cor
1222+1.000000E+00, -3.410210E-03
1223+0.000000E+00, +1.000000E+00
1224cor = getFilled(0., ndim, ndim)
1225call setCor(cor, uppDia, samShifted, dim, iweight, iweisum) ! same result as above.
1226cor
1227+1.000000E+00, -3.410210E-03
1228+0.000000E+00, +1.000000E+00
1229
1230cor = getFilled(0., ndim, ndim)
1231call setCor(cor, lowDia, mean, sample, dim, iweight, iweisum)
1232cor
1233+1.000000E+00, +0.000000E+00
1234-3.410244E-03, +1.000000E+00
1235cor = getFilled(0., ndim, ndim)
1236call setCor(cor, lowDia, samShifted, dim, iweight, iweisum) ! same result as above.
1237cor
1238+1.000000E+00, +0.000000E+00
1239-3.410244E-03, +1.000000E+00
1240cor = getFilled(0., ndim, ndim)
1241call setCor(cor, lowDia, getVerbose(samShifted, iweight, sum(iweight), dim), dim) ! same result as above.
1242cor
1243+1.000000E+00, +0.000000E+00
1244-3.410219E-03, +1.000000E+00
1245
1246'Compute the sample correlation along the first dimension.'
1247
1248dim = 1
1249cor = getFilled(0., ndim, ndim)
1250call setCor(cor, uppDia, mean, transpose(sample), dim, iweight, iweisum)
1251cor
1252+1.000000E+00, -3.410210E-03
1253+0.000000E+00, +1.000000E+00
1254cor = getFilled(0., ndim, ndim)
1255call setCor(cor, uppDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.
1256cor
1257+1.000000E+00, -3.410210E-03
1258+0.000000E+00, +1.000000E+00
1259
1260cor = getFilled(0., ndim, ndim)
1261call setCor(cor, lowDia, mean, transpose(sample), dim, iweight, iweisum)
1262cor
1263+1.000000E+00, +0.000000E+00
1264-3.410210E-03, +1.000000E+00
1265cor = getFilled(0., ndim, ndim)
1266call setCor(cor, lowDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.
1267cor
1268+1.000000E+00, +0.000000E+00
1269-3.410210E-03, +1.000000E+00
1270
1271'Compute the full sample correlation for a pair of time series.'
1272
1273call setCor(cor(1,1), mean, sample(1,:), sample(2,:), iweight, iweisum)
1274cor(1,1)
1275-3.410210E-03
1276call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), iweight, iweisum) ! same result as above.
1277cor(1,1)
1278-3.410210E-03
1279
1280
1281!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1282!Compute the correlation matrix with real weights.
1283!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1284
1285dim = 2
1286samShifted = getShifted(sample, dim, -mean)
1287samShifted
1288+10.0869560, +8.08695602, -5.91304398, +9.08695602, +3.08695602, +5.08695602, -0.913043976, -7.91304398, +4.08695602, -7.91304398
1289-2.84782624, +12.1521740, +5.15217400, -1.84782624, +6.15217400, -1.84782624, +8.15217400, -2.84782624, +5.15217400, -2.84782624
1290cor = getFilled(0., ndim, ndim)
1291call setCor(cor, uppDia, mean, sample, dim, rweight, rweisum)
1292cor
1293+1.000000E+00, -3.410210E-03
1294+0.000000E+00, +1.000000E+00
1295cor = getFilled(0., ndim, ndim)
1296call setCor(cor, uppDia, samShifted, dim, rweight, rweisum) ! same result as above.
1297cor
1298+1.000000E+00, -3.410210E-03
1299+0.000000E+00, +1.000000E+00
1300
1301cor = getFilled(0., ndim, ndim)
1302call setCor(cor, lowDia, mean, sample, dim, rweight, rweisum)
1303cor
1304+1.000000E+00, +0.000000E+00
1305-3.410244E-03, +1.000000E+00
1306cor = getFilled(0., ndim, ndim)
1307call setCor(cor, lowDia, samShifted, dim, rweight, rweisum) ! same result as above.
1308cor
1309+1.000000E+00, +0.000000E+00
1310-3.410244E-03, +1.000000E+00
1311
1312'Compute the sample correlation along the first dimension.'
1313
1314dim = 1
1315cor = getFilled(0., ndim, ndim)
1316call setCor(cor, uppDia, mean, transpose(sample), dim, rweight, rweisum)
1317cor
1318+1.000000E+00, -3.410210E-03
1319+0.000000E+00, +1.000000E+00
1320cor = getFilled(0., ndim, ndim)
1321call setCor(cor, uppDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.
1322cor
1323+1.000000E+00, -3.410210E-03
1324+0.000000E+00, +1.000000E+00
1325
1326cor = getFilled(0., ndim, ndim)
1327call setCor(cor, lowDia, mean, transpose(sample), dim, rweight, rweisum)
1328cor
1329+1.000000E+00, +0.000000E+00
1330-3.410210E-03, +1.000000E+00
1331cor = getFilled(0., ndim, ndim)
1332call setCor(cor, lowDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.
1333cor
1334+1.000000E+00, +0.000000E+00
1335-3.410210E-03, +1.000000E+00
1336
1337'Compute the full sample correlation for a pair of time series.'
1338
1339call setCor(cor(1,1), mean, sample(1,:), sample(2,:), rweight, rweisum)
1340cor(1,1)
1341-3.410210E-03
1342call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), rweight, rweisum) ! same result as above.
1343cor(1,1)
1344-3.410210E-03
1345
1346
Test:
test_pm_sampleCor


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:
Fatemeh Bagheri, Monday 02:15 AM, September 27, 2021, Dallas, TX
Amir Shahmoradi, Wednesday 4:13 AM, August 13, 2016, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 2323 of file pm_sampleCor.F90.


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