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+0
86std = getUnifRand(1, ndim, ndim)
87std
88
89cov = getCovRand(1._TKG, std)
90cov
91
92
93cor = cov
94call setCor(cor, uppDia, cor, uppDia)
95cor
96
97cor = cov
98call setCor(cor, uppDia, cor, lowDia)
99cor
100
101cor = cov
102call setCor(cor, upp, cor, uppDia)
103cor
104
105cor = cov
106call setCor(cor, upp, cor, lowDia)
107cor
108
109cor = cov
110call setCor(cor, lowDia, cor, uppDia)
111cor
112
113cor = cov
114call setCor(cor, lowDia, cor, lowDia)
115cor
116
117cor = cov
118call setCor(cor, low, cor, uppDia)
119cor
120
121cor = cov
122call setCor(cor, low, cor, lowDia)
123cor
124
125
126cor = cov
127call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
128cor
129
130cor = cov
131call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
132cor
133
134cor = cov
135call setCor(cor, upp, cor, upp, stdinv = 1 / std)
136cor
137
138cor = cov
139call setCor(cor, upp, cor, low, stdinv = 1 / std)
140cor
141
142cor = cov
143call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
144cor
145
146cor = cov
147call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
148cor
149
150cor = cov
151call setCor(cor, low, cor, upp, stdinv = 1 / std)
152cor
153
154cor = cov
155call setCor(cor, low, cor, low, stdinv = 1 / std)
156cor
157
158
159ndim = getUnifRand(0, 4)
160ndim
161+3
162std = getUnifRand(1, ndim, ndim)
163std
164+3.000000E+00, +3.000000E+00, +3.000000E+00
165cov = getCovRand(1._TKG, std)
166cov
167+9.000000E+00, -8.644159E+00, -6.755726E+00
168-8.644159E+00, +9.000002E+00, +5.978867E+00
169-6.755726E+00, +5.978867E+00, +9.000000E+00
170
171
172cor = cov
173call setCor(cor, uppDia, cor, uppDia)
174cor
175+1.000000E+00, -9.604621E-01, -7.506363E-01
176-8.644159E+00, +1.000000E+00, +6.643185E-01
177-6.755726E+00, +5.978867E+00, +1.000000E+00
178
179cor = cov
180call setCor(cor, uppDia, cor, lowDia)
181cor
182+1.000000E+00, -9.604621E-01, -7.506363E-01
183-8.644159E+00, +1.000000E+00, +6.643185E-01
184-6.755726E+00, +5.978867E+00, +1.000000E+00
185
186cor = cov
187call setCor(cor, upp, cor, uppDia)
188cor
189+9.000000E+00, -9.604621E-01, -7.506363E-01
190-8.644159E+00, +9.000002E+00, +6.643185E-01
191-6.755726E+00, +5.978867E+00, +9.000000E+00
192
193cor = cov
194call setCor(cor, upp, cor, lowDia)
195cor
196+9.000000E+00, -9.604621E-01, -7.506363E-01
197-8.644159E+00, +9.000002E+00, +6.643185E-01
198-6.755726E+00, +5.978867E+00, +9.000000E+00
199
200cor = cov
201call setCor(cor, lowDia, cor, uppDia)
202cor
203+1.000000E+00, -8.644159E+00, -6.755726E+00
204-8.644159E+00, +1.000000E+00, +5.978867E+00
205-6.755726E+00, +5.978867E+00, +1.000000E+00
206
207cor = cov
208call setCor(cor, lowDia, cor, lowDia)
209cor
210+1.000000E+00, -8.644159E+00, -6.755726E+00
211-9.604621E-01, +1.000000E+00, +5.978867E+00
212-7.506363E-01, +6.643185E-01, +1.000000E+00
213
214cor = cov
215call setCor(cor, low, cor, uppDia)
216cor
217+9.000000E+00, -8.644159E+00, -6.755726E+00
218-9.604621E-01, +9.000002E+00, +5.978867E+00
219-7.506363E-01, +6.643185E-01, +9.000000E+00
220
221cor = cov
222call setCor(cor, low, cor, lowDia)
223cor
224+9.000000E+00, -8.644159E+00, -6.755726E+00
225-9.604621E-01, +9.000002E+00, +5.978867E+00
226-7.506363E-01, +6.643185E-01, +9.000000E+00
227
228
229cor = cov
230call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
231cor
232+1.000000E+00, -9.604622E-01, -7.506363E-01
233-8.644159E+00, +1.000000E+00, +6.643186E-01
234-6.755726E+00, +5.978867E+00, +1.000000E+00
235
236cor = cov
237call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
238cor
239+1.000000E+00, -9.604622E-01, -7.506363E-01
240-8.644159E+00, +1.000000E+00, +6.643186E-01
241-6.755726E+00, +5.978867E+00, +1.000000E+00
242
243cor = cov
244call setCor(cor, upp, cor, upp, stdinv = 1 / std)
245cor
246+9.00000000, -0.960462213, -0.750636280
247-8.64415932, +9.00000191, +0.664318621
248-6.75572586, +5.97886705, +9.00000000
249
250cor = cov
251call setCor(cor, upp, cor, low, stdinv = 1 / std)
252cor
253+9.000000E+00, -9.604622E-01, -7.506363E-01
254-8.644159E+00, +9.000002E+00, +6.643186E-01
255-6.755726E+00, +5.978867E+00, +9.000000E+00
256
257cor = cov
258call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
259cor
260+1.000000E+00, -8.644159E+00, -6.755726E+00
261-9.604622E-01, +1.000000E+00, +5.978867E+00
262-7.506363E-01, +6.643186E-01, +1.000000E+00
263
264cor = cov
265call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
266cor
267+1.000000E+00, -8.644159E+00, -6.755726E+00
268-9.604622E-01, +1.000000E+00, +5.978867E+00
269-7.506363E-01, +6.643186E-01, +1.000000E+00
270
271cor = cov
272call setCor(cor, low, cor, upp, stdinv = 1 / std)
273cor
274+9.000000E+00, -8.644159E+00, -6.755726E+00
275-9.604622E-01, +9.000002E+00, +5.978867E+00
276-7.506363E-01, +6.643186E-01, +9.000000E+00
277
278cor = cov
279call setCor(cor, low, cor, low, stdinv = 1 / std)
280cor
281+9.000000E+00, -8.644159E+00, -6.755726E+00
282-9.604622E-01, +9.000002E+00, +5.978867E+00
283-7.506363E-01, +6.643186E-01, +9.000000E+00
284
285
286ndim = getUnifRand(0, 4)
287ndim
288+4
289std = getUnifRand(1, ndim, ndim)
290std
291+3.000000E+00, +1.000000E+00, +2.000000E+00, +3.000000E+00
292cov = getCovRand(1._TKG, std)
293cov
294+9.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
295+1.945879E+00, +1.000000E+00, +1.062772E+00, +8.429271E-01
296+2.330843E+00, +1.062772E+00, +4.000000E+00, -3.917000E+00
297+2.180150E+00, +8.429271E-01, -3.917000E+00, +9.000000E+00
298
299
300cor = cov
301call setCor(cor, uppDia, cor, uppDia)
302cor
303+1.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
304+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
305+2.330843E+00, +1.062772E+00, +1.000000E+00, -6.528333E-01
306+2.180150E+00, +8.429271E-01, -3.917000E+00, +1.000000E+00
307
308cor = cov
309call setCor(cor, uppDia, cor, lowDia)
310cor
311+1.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
312+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
313+2.330843E+00, +1.062772E+00, +1.000000E+00, -6.528333E-01
314+2.180150E+00, +8.429271E-01, -3.917000E+00, +1.000000E+00
315
316cor = cov
317call setCor(cor, upp, cor, uppDia)
318cor
319+9.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
320+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
321+2.330843E+00, +1.062772E+00, +4.000000E+00, -6.528333E-01
322+2.180150E+00, +8.429271E-01, -3.917000E+00, +9.000000E+00
323
324cor = cov
325call setCor(cor, upp, cor, lowDia)
326cor
327+9.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
328+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
329+2.330843E+00, +1.062772E+00, +4.000000E+00, -6.528333E-01
330+2.180150E+00, +8.429271E-01, -3.917000E+00, +9.000000E+00
331
332cor = cov
333call setCor(cor, lowDia, cor, uppDia)
334cor
335+1.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
336+1.945879E+00, +1.000000E+00, +1.062772E+00, +8.429271E-01
337+2.330843E+00, +1.062772E+00, +1.000000E+00, -3.917000E+00
338+2.180150E+00, +8.429271E-01, -3.917000E+00, +1.000000E+00
339
340cor = cov
341call setCor(cor, lowDia, cor, lowDia)
342cor
343+1.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
344+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
345+3.884738E-01, +5.313858E-01, +1.000000E+00, -3.917000E+00
346+2.422389E-01, +2.809757E-01, -6.528333E-01, +1.000000E+00
347
348cor = cov
349call setCor(cor, low, cor, uppDia)
350cor
351+9.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
352+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
353+3.884738E-01, +5.313858E-01, +4.000000E+00, -3.917000E+00
354+2.422389E-01, +2.809757E-01, -6.528333E-01, +9.000000E+00
355
356cor = cov
357call setCor(cor, low, cor, lowDia)
358cor
359+9.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
360+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
361+3.884738E-01, +5.313858E-01, +4.000000E+00, -3.917000E+00
362+2.422389E-01, +2.809757E-01, -6.528333E-01, +9.000000E+00
363
364
365cor = cov
366call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
367cor
368+1.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
369+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
370+2.330843E+00, +1.062772E+00, +1.000000E+00, -6.528333E-01
371+2.180150E+00, +8.429271E-01, -3.917000E+00, +1.000000E+00
372
373cor = cov
374call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
375cor
376+1.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
377+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
378+2.330843E+00, +1.062772E+00, +1.000000E+00, -6.528333E-01
379+2.180150E+00, +8.429271E-01, -3.917000E+00, +1.000000E+00
380
381cor = cov
382call setCor(cor, upp, cor, upp, stdinv = 1 / std)
383cor
384+9.00000000, +0.648626328, +0.388473839, +0.242238864
385+1.94587898, +1.00000000, +0.531385839, +0.280975699
386+2.33084297, +1.06277168, +4.00000048, -0.652833343
387+2.18014956, +0.842927098, -3.91699982, +9.00000000
388
389cor = cov
390call setCor(cor, upp, cor, low, stdinv = 1 / std)
391cor
392+9.000000E+00, +6.486263E-01, +3.884738E-01, +2.422389E-01
393+1.945879E+00, +1.000000E+00, +5.313858E-01, +2.809757E-01
394+2.330843E+00, +1.062772E+00, +4.000000E+00, -6.528333E-01
395+2.180150E+00, +8.429271E-01, -3.917000E+00, +9.000000E+00
396
397cor = cov
398call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
399cor
400+1.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
401+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
402+3.884738E-01, +5.313858E-01, +1.000000E+00, -3.917000E+00
403+2.422389E-01, +2.809757E-01, -6.528333E-01, +1.000000E+00
404
405cor = cov
406call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
407cor
408+1.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
409+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
410+3.884738E-01, +5.313858E-01, +1.000000E+00, -3.917000E+00
411+2.422389E-01, +2.809757E-01, -6.528333E-01, +1.000000E+00
412
413cor = cov
414call setCor(cor, low, cor, upp, stdinv = 1 / std)
415cor
416+9.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
417+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
418+3.884738E-01, +5.313858E-01, +4.000000E+00, -3.917000E+00
419+2.422389E-01, +2.809757E-01, -6.528333E-01, +9.000000E+00
420
421cor = cov
422call setCor(cor, low, cor, low, stdinv = 1 / std)
423cor
424+9.000000E+00, +1.945879E+00, +2.330843E+00, +2.180150E+00
425+6.486263E-01, +1.000000E+00, +1.062772E+00, +8.429271E-01
426+3.884738E-01, +5.313858E-01, +4.000000E+00, -3.917000E+00
427+2.422389E-01, +2.809757E-01, -6.528333E-01, +9.000000E+00
428
429
430ndim = getUnifRand(0, 4)
431ndim
432+0
433std = getUnifRand(1, ndim, ndim)
434std
435
436cov = getCovRand(1._TKG, std)
437cov
438
439
440cor = cov
441call setCor(cor, uppDia, cor, uppDia)
442cor
443
444cor = cov
445call setCor(cor, uppDia, cor, lowDia)
446cor
447
448cor = cov
449call setCor(cor, upp, cor, uppDia)
450cor
451
452cor = cov
453call setCor(cor, upp, cor, lowDia)
454cor
455
456cor = cov
457call setCor(cor, lowDia, cor, uppDia)
458cor
459
460cor = cov
461call setCor(cor, lowDia, cor, lowDia)
462cor
463
464cor = cov
465call setCor(cor, low, cor, uppDia)
466cor
467
468cor = cov
469call setCor(cor, low, cor, lowDia)
470cor
471
472
473cor = cov
474call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
475cor
476
477cor = cov
478call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
479cor
480
481cor = cov
482call setCor(cor, upp, cor, upp, stdinv = 1 / std)
483cor
484
485cor = cov
486call setCor(cor, upp, cor, low, stdinv = 1 / std)
487cor
488
489cor = cov
490call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
491cor
492
493cor = cov
494call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
495cor
496
497cor = cov
498call setCor(cor, low, cor, upp, stdinv = 1 / std)
499cor
500
501cor = cov
502call setCor(cor, low, cor, low, stdinv = 1 / std)
503cor
504
505
506ndim = getUnifRand(0, 4)
507ndim
508+0
509std = getUnifRand(1, ndim, ndim)
510std
511
512cov = getCovRand(1._TKG, std)
513cov
514
515
516cor = cov
517call setCor(cor, uppDia, cor, uppDia)
518cor
519
520cor = cov
521call setCor(cor, uppDia, cor, lowDia)
522cor
523
524cor = cov
525call setCor(cor, upp, cor, uppDia)
526cor
527
528cor = cov
529call setCor(cor, upp, cor, lowDia)
530cor
531
532cor = cov
533call setCor(cor, lowDia, cor, uppDia)
534cor
535
536cor = cov
537call setCor(cor, lowDia, cor, lowDia)
538cor
539
540cor = cov
541call setCor(cor, low, cor, uppDia)
542cor
543
544cor = cov
545call setCor(cor, low, cor, lowDia)
546cor
547
548
549cor = cov
550call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
551cor
552
553cor = cov
554call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
555cor
556
557cor = cov
558call setCor(cor, upp, cor, upp, stdinv = 1 / std)
559cor
560
561cor = cov
562call setCor(cor, upp, cor, low, stdinv = 1 / std)
563cor
564
565cor = cov
566call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
567cor
568
569cor = cov
570call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
571cor
572
573cor = cov
574call setCor(cor, low, cor, upp, stdinv = 1 / std)
575cor
576
577cor = cov
578call setCor(cor, low, cor, low, stdinv = 1 / std)
579cor
580
581
582ndim = getUnifRand(0, 4)
583ndim
584+4
585std = getUnifRand(1, ndim, ndim)
586std
587+4.000000E+00, +2.000000E+00, +3.000000E+00, +4.000000E+00
588cov = getCovRand(1._TKG, std)
589cov
590+1.600000E+01, +7.412047E+00, -1.942240E-01, +2.439008E+00
591+7.412047E+00, +4.000000E+00, -2.190776E+00, +1.851153E+00
592-1.942240E-01, -2.190776E+00, +9.000000E+00, +1.033868E+00
593+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.600000E+01
594
595
596cor = cov
597call setCor(cor, uppDia, cor, uppDia)
598cor
599+1.000000E+00, +9.265060E-01, -1.618534E-02, +1.524380E-01
600+7.412047E+00, +1.000000E+00, -3.651294E-01, +2.313942E-01
601-1.942240E-01, -2.190776E+00, +1.000000E+00, +8.615566E-02
602+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.000000E+00
603
604cor = cov
605call setCor(cor, uppDia, cor, lowDia)
606cor
607+1.000000E+00, +9.265060E-01, -1.618534E-02, +1.524380E-01
608+7.412047E+00, +1.000000E+00, -3.651294E-01, +2.313942E-01
609-1.942240E-01, -2.190776E+00, +1.000000E+00, +8.615566E-02
610+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.000000E+00
611
612cor = cov
613call setCor(cor, upp, cor, uppDia)
614cor
615+1.600000E+01, +9.265060E-01, -1.618534E-02, +1.524380E-01
616+7.412047E+00, +4.000000E+00, -3.651294E-01, +2.313942E-01
617-1.942240E-01, -2.190776E+00, +9.000000E+00, +8.615566E-02
618+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.600000E+01
619
620cor = cov
621call setCor(cor, upp, cor, lowDia)
622cor
623+1.600000E+01, +9.265060E-01, -1.618534E-02, +1.524380E-01
624+7.412047E+00, +4.000000E+00, -3.651294E-01, +2.313942E-01
625-1.942240E-01, -2.190776E+00, +9.000000E+00, +8.615566E-02
626+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.600000E+01
627
628cor = cov
629call setCor(cor, lowDia, cor, uppDia)
630cor
631+1.000000E+00, +7.412047E+00, -1.942240E-01, +2.439008E+00
632+7.412047E+00, +1.000000E+00, -2.190776E+00, +1.851153E+00
633-1.942240E-01, -2.190776E+00, +1.000000E+00, +1.033868E+00
634+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.000000E+00
635
636cor = cov
637call setCor(cor, lowDia, cor, lowDia)
638cor
639+1.000000E+00, +7.412047E+00, -1.942240E-01, +2.439008E+00
640+9.265060E-01, +1.000000E+00, -2.190776E+00, +1.851153E+00
641-1.618534E-02, -3.651294E-01, +1.000000E+00, +1.033868E+00
642+1.524380E-01, +2.313942E-01, +8.615566E-02, +1.000000E+00
643
644cor = cov
645call setCor(cor, low, cor, uppDia)
646cor
647+1.600000E+01, +7.412047E+00, -1.942240E-01, +2.439008E+00
648+9.265060E-01, +4.000000E+00, -2.190776E+00, +1.851153E+00
649-1.618534E-02, -3.651294E-01, +9.000000E+00, +1.033868E+00
650+1.524380E-01, +2.313942E-01, +8.615566E-02, +1.600000E+01
651
652cor = cov
653call setCor(cor, low, cor, lowDia)
654cor
655+1.600000E+01, +7.412047E+00, -1.942240E-01, +2.439008E+00
656+9.265060E-01, +4.000000E+00, -2.190776E+00, +1.851153E+00
657-1.618534E-02, -3.651294E-01, +9.000000E+00, +1.033868E+00
658+1.524380E-01, +2.313942E-01, +8.615566E-02, +1.600000E+01
659
660
661cor = cov
662call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
663cor
664+1.000000E+00, +9.265059E-01, -1.618534E-02, +1.524380E-01
665+7.412047E+00, +1.000000E+00, -3.651293E-01, +2.313941E-01
666-1.942240E-01, -2.190776E+00, +1.000000E+00, +8.615565E-02
667+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.000000E+00
668
669cor = cov
670call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
671cor
672+1.000000E+00, +9.265059E-01, -1.618534E-02, +1.524380E-01
673+7.412047E+00, +1.000000E+00, -3.651293E-01, +2.313941E-01
674-1.942240E-01, -2.190776E+00, +1.000000E+00, +8.615565E-02
675+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.000000E+00
676
677cor = cov
678call setCor(cor, upp, cor, upp, stdinv = 1 / std)
679cor
680+16.0000000, +0.926505864, -0.161853358E-1, +0.152438030
681+7.41204691, +3.99999976, -0.365129322, +0.231394112
682-0.194224030, -2.19077587, +9.00000000, +0.861556530E-1
683+2.43900847, +1.85115290, +1.03386784, +15.9999981
684
685cor = cov
686call setCor(cor, upp, cor, low, stdinv = 1 / std)
687cor
688+1.600000E+01, +9.265059E-01, -1.618534E-02, +1.524380E-01
689+7.412047E+00, +4.000000E+00, -3.651293E-01, +2.313941E-01
690-1.942240E-01, -2.190776E+00, +9.000000E+00, +8.615565E-02
691+2.439008E+00, +1.851153E+00, +1.033868E+00, +1.600000E+01
692
693cor = cov
694call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
695cor
696+1.000000E+00, +7.412047E+00, -1.942240E-01, +2.439008E+00
697+9.265059E-01, +1.000000E+00, -2.190776E+00, +1.851153E+00
698-1.618534E-02, -3.651293E-01, +1.000000E+00, +1.033868E+00
699+1.524380E-01, +2.313941E-01, +8.615565E-02, +1.000000E+00
700
701cor = cov
702call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
703cor
704+1.000000E+00, +7.412047E+00, -1.942240E-01, +2.439008E+00
705+9.265059E-01, +1.000000E+00, -2.190776E+00, +1.851153E+00
706-1.618534E-02, -3.651293E-01, +1.000000E+00, +1.033868E+00
707+1.524380E-01, +2.313941E-01, +8.615565E-02, +1.000000E+00
708
709cor = cov
710call setCor(cor, low, cor, upp, stdinv = 1 / std)
711cor
712+1.600000E+01, +7.412047E+00, -1.942240E-01, +2.439008E+00
713+9.265059E-01, +4.000000E+00, -2.190776E+00, +1.851153E+00
714-1.618534E-02, -3.651293E-01, +9.000000E+00, +1.033868E+00
715+1.524380E-01, +2.313941E-01, +8.615565E-02, +1.600000E+01
716
717cor = cov
718call setCor(cor, low, cor, low, stdinv = 1 / std)
719cor
720+1.600000E+01, +7.412047E+00, -1.942240E-01, +2.439008E+00
721+9.265059E-01, +4.000000E+00, -2.190776E+00, +1.851153E+00
722-1.618534E-02, -3.651293E-01, +9.000000E+00, +1.033868E+00
723+1.524380E-01, +2.313941E-01, +8.615565E-02, +1.600000E+01
724
725
726ndim = getUnifRand(0, 4)
727ndim
728+2
729std = getUnifRand(1, ndim, ndim)
730std
731+1.000000E+00, +2.000000E+00
732cov = getCovRand(1._TKG, std)
733cov
734+1.000000E+00, +1.768927E+00
735+1.768927E+00, +4.000000E+00
736
737
738cor = cov
739call setCor(cor, uppDia, cor, uppDia)
740cor
741+1.000000E+00, +8.844635E-01
742+1.768927E+00, +1.000000E+00
743
744cor = cov
745call setCor(cor, uppDia, cor, lowDia)
746cor
747+1.000000E+00, +8.844635E-01
748+1.768927E+00, +1.000000E+00
749
750cor = cov
751call setCor(cor, upp, cor, uppDia)
752cor
753+1.000000E+00, +8.844635E-01
754+1.768927E+00, +4.000000E+00
755
756cor = cov
757call setCor(cor, upp, cor, lowDia)
758cor
759+1.000000E+00, +8.844635E-01
760+1.768927E+00, +4.000000E+00
761
762cor = cov
763call setCor(cor, lowDia, cor, uppDia)
764cor
765+1.000000E+00, +1.768927E+00
766+1.768927E+00, +1.000000E+00
767
768cor = cov
769call setCor(cor, lowDia, cor, lowDia)
770cor
771+1.000000E+00, +1.768927E+00
772+8.844635E-01, +1.000000E+00
773
774cor = cov
775call setCor(cor, low, cor, uppDia)
776cor
777+1.000000E+00, +1.768927E+00
778+8.844635E-01, +4.000000E+00
779
780cor = cov
781call setCor(cor, low, cor, lowDia)
782cor
783+1.000000E+00, +1.768927E+00
784+8.844635E-01, +4.000000E+00
785
786
787cor = cov
788call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
789cor
790+1.000000E+00, +8.844634E-01
791+1.768927E+00, +1.000000E+00
792
793cor = cov
794call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
795cor
796+1.000000E+00, +8.844634E-01
797+1.768927E+00, +1.000000E+00
798
799cor = cov
800call setCor(cor, upp, cor, upp, stdinv = 1 / std)
801cor
802+1.00000000, +0.884463370
803+1.76892674, +3.99999976
804
805cor = cov
806call setCor(cor, upp, cor, low, stdinv = 1 / std)
807cor
808+1.000000E+00, +8.844634E-01
809+1.768927E+00, +4.000000E+00
810
811cor = cov
812call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
813cor
814+1.000000E+00, +1.768927E+00
815+8.844634E-01, +1.000000E+00
816
817cor = cov
818call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
819cor
820+1.000000E+00, +1.768927E+00
821+8.844634E-01, +1.000000E+00
822
823cor = cov
824call setCor(cor, low, cor, upp, stdinv = 1 / std)
825cor
826+1.000000E+00, +1.768927E+00
827+8.844634E-01, +4.000000E+00
828
829cor = cov
830call setCor(cor, low, cor, low, stdinv = 1 / std)
831cor
832+1.000000E+00, +1.768927E+00
833+8.844634E-01, +4.000000E+00
834
835
836ndim = getUnifRand(0, 4)
837ndim
838+4
839std = getUnifRand(1, ndim, ndim)
840std
841+3.000000E+00, +2.000000E+00, +4.000000E+00, +2.000000E+00
842cov = getCovRand(1._TKG, std)
843cov
844+9.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
845+5.530821E+00, +4.000000E+00, -6.607532E+00, -3.342036E+00
846-1.050892E+01, -6.607532E+00, +1.600000E+01, +4.399990E+00
847-3.728874E+00, -3.342036E+00, +4.399990E+00, +4.000000E+00
848
849
850cor = cov
851call setCor(cor, uppDia, cor, uppDia)
852cor
853+1.000000E+00, +9.218036E-01, -8.757436E-01, -6.214790E-01
854+5.530821E+00, +1.000000E+00, -8.259416E-01, -8.355091E-01
855-1.050892E+01, -6.607532E+00, +1.000000E+00, +5.499988E-01
856-3.728874E+00, -3.342036E+00, +4.399990E+00, +1.000000E+00
857
858cor = cov
859call setCor(cor, uppDia, cor, lowDia)
860cor
861+1.000000E+00, +9.218036E-01, -8.757436E-01, -6.214790E-01
862+5.530821E+00, +1.000000E+00, -8.259416E-01, -8.355091E-01
863-1.050892E+01, -6.607532E+00, +1.000000E+00, +5.499988E-01
864-3.728874E+00, -3.342036E+00, +4.399990E+00, +1.000000E+00
865
866cor = cov
867call setCor(cor, upp, cor, uppDia)
868cor
869+9.000000E+00, +9.218036E-01, -8.757436E-01, -6.214790E-01
870+5.530821E+00, +4.000000E+00, -8.259416E-01, -8.355091E-01
871-1.050892E+01, -6.607532E+00, +1.600000E+01, +5.499988E-01
872-3.728874E+00, -3.342036E+00, +4.399990E+00, +4.000000E+00
873
874cor = cov
875call setCor(cor, upp, cor, lowDia)
876cor
877+9.000000E+00, +9.218036E-01, -8.757436E-01, -6.214790E-01
878+5.530821E+00, +4.000000E+00, -8.259416E-01, -8.355091E-01
879-1.050892E+01, -6.607532E+00, +1.600000E+01, +5.499988E-01
880-3.728874E+00, -3.342036E+00, +4.399990E+00, +4.000000E+00
881
882cor = cov
883call setCor(cor, lowDia, cor, uppDia)
884cor
885+1.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
886+5.530821E+00, +1.000000E+00, -6.607532E+00, -3.342036E+00
887-1.050892E+01, -6.607532E+00, +1.000000E+00, +4.399990E+00
888-3.728874E+00, -3.342036E+00, +4.399990E+00, +1.000000E+00
889
890cor = cov
891call setCor(cor, lowDia, cor, lowDia)
892cor
893+1.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
894+9.218036E-01, +1.000000E+00, -6.607532E+00, -3.342036E+00
895-8.757436E-01, -8.259416E-01, +1.000000E+00, +4.399990E+00
896-6.214790E-01, -8.355091E-01, +5.499988E-01, +1.000000E+00
897
898cor = cov
899call setCor(cor, low, cor, uppDia)
900cor
901+9.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
902+9.218036E-01, +4.000000E+00, -6.607532E+00, -3.342036E+00
903-8.757436E-01, -8.259416E-01, +1.600000E+01, +4.399990E+00
904-6.214790E-01, -8.355091E-01, +5.499988E-01, +4.000000E+00
905
906cor = cov
907call setCor(cor, low, cor, lowDia)
908cor
909+9.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
910+9.218036E-01, +4.000000E+00, -6.607532E+00, -3.342036E+00
911-8.757436E-01, -8.259416E-01, +1.600000E+01, +4.399990E+00
912-6.214790E-01, -8.355091E-01, +5.499988E-01, +4.000000E+00
913
914
915cor = cov
916call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
917cor
918+1.000000E+00, +9.218036E-01, -8.757435E-01, -6.214790E-01
919+5.530821E+00, +1.000000E+00, -8.259415E-01, -8.355091E-01
920-1.050892E+01, -6.607532E+00, +1.000000E+00, +5.499988E-01
921-3.728874E+00, -3.342036E+00, +4.399990E+00, +1.000000E+00
922
923cor = cov
924call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
925cor
926+1.000000E+00, +9.218036E-01, -8.757435E-01, -6.214790E-01
927+5.530821E+00, +1.000000E+00, -8.259415E-01, -8.355091E-01
928-1.050892E+01, -6.607532E+00, +1.000000E+00, +5.499988E-01
929-3.728874E+00, -3.342036E+00, +4.399990E+00, +1.000000E+00
930
931cor = cov
932call setCor(cor, upp, cor, upp, stdinv = 1 / std)
933cor
934+9.00000000, +0.921803594, -0.875743508, -0.621478975
935+5.53082132, +4.00000048, -0.825941503, -0.835509062
936-10.5089216, -6.60753202, +15.9999971, +0.549998760
937-3.72887373, -3.34203625, +4.39999008, +4.00000048
938
939cor = cov
940call setCor(cor, upp, cor, low, stdinv = 1 / std)
941cor
942+9.000000E+00, +9.218036E-01, -8.757435E-01, -6.214790E-01
943+5.530821E+00, +4.000000E+00, -8.259415E-01, -8.355091E-01
944-1.050892E+01, -6.607532E+00, +1.600000E+01, +5.499988E-01
945-3.728874E+00, -3.342036E+00, +4.399990E+00, +4.000000E+00
946
947cor = cov
948call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
949cor
950+1.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
951+9.218036E-01, +1.000000E+00, -6.607532E+00, -3.342036E+00
952-8.757435E-01, -8.259415E-01, +1.000000E+00, +4.399990E+00
953-6.214790E-01, -8.355091E-01, +5.499988E-01, +1.000000E+00
954
955cor = cov
956call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
957cor
958+1.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
959+9.218036E-01, +1.000000E+00, -6.607532E+00, -3.342036E+00
960-8.757435E-01, -8.259415E-01, +1.000000E+00, +4.399990E+00
961-6.214790E-01, -8.355091E-01, +5.499988E-01, +1.000000E+00
962
963cor = cov
964call setCor(cor, low, cor, upp, stdinv = 1 / std)
965cor
966+9.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
967+9.218036E-01, +4.000000E+00, -6.607532E+00, -3.342036E+00
968-8.757435E-01, -8.259415E-01, +1.600000E+01, +4.399990E+00
969-6.214790E-01, -8.355091E-01, +5.499988E-01, +4.000000E+00
970
971cor = cov
972call setCor(cor, low, cor, low, stdinv = 1 / std)
973cor
974+9.000000E+00, +5.530821E+00, -1.050892E+01, -3.728874E+00
975+9.218036E-01, +4.000000E+00, -6.607532E+00, -3.342036E+00
976-8.757435E-01, -8.259415E-01, +1.600000E+01, +4.399990E+00
977-6.214790E-01, -8.355091E-01, +5.499988E-01, +4.000000E+00
978
979
980ndim = getUnifRand(0, 4)
981ndim
982+0
983std = getUnifRand(1, ndim, ndim)
984std
985
986cov = getCovRand(1._TKG, std)
987cov
988
989
990cor = cov
991call setCor(cor, uppDia, cor, uppDia)
992cor
993
994cor = cov
995call setCor(cor, uppDia, cor, lowDia)
996cor
997
998cor = cov
999call setCor(cor, upp, cor, uppDia)
1000cor
1001
1002cor = cov
1003call setCor(cor, upp, cor, lowDia)
1004cor
1005
1006cor = cov
1007call setCor(cor, lowDia, cor, uppDia)
1008cor
1009
1010cor = cov
1011call setCor(cor, lowDia, cor, lowDia)
1012cor
1013
1014cor = cov
1015call setCor(cor, low, cor, uppDia)
1016cor
1017
1018cor = cov
1019call setCor(cor, low, cor, lowDia)
1020cor
1021
1022
1023cor = cov
1024call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
1025cor
1026
1027cor = cov
1028call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
1029cor
1030
1031cor = cov
1032call setCor(cor, upp, cor, upp, stdinv = 1 / std)
1033cor
1034
1035cor = cov
1036call setCor(cor, upp, cor, low, stdinv = 1 / std)
1037cor
1038
1039cor = cov
1040call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
1041cor
1042
1043cor = cov
1044call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
1045cor
1046
1047cor = cov
1048call setCor(cor, low, cor, upp, stdinv = 1 / std)
1049cor
1050
1051cor = cov
1052call setCor(cor, low, cor, low, stdinv = 1 / std)
1053cor
1054
1055
1056!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1057!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1058!Compute the Pearson correlation matrix for a pair of time series.
1059!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1060!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1061
1062ndim = 2; nsam = 10; dim = 2
1063sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
1064sample
1065+1.300000E+01, +1.000000E+00, +1.600000E+01, +1.000000E+01, +1.100000E+01, +1.000000E+01, +1.000000E+00, +5.000000E+00, +1.500000E+01, +1.200000E+01
1066+1.600000E+01, +5.000000E+00, +3.000000E+00, +2.000000E+01, +1.800000E+01, +1.700000E+01, +6.000000E+00, +1.000000E+01, +7.000000E+00, +1.800000E+01
1067mean = getMean(sample, dim)
1068mean
1069+9.400001E+00, +1.200000E+01
1070samShifted = getShifted(sample, dim, -mean)
1071samShifted
1072+3.599999E+00, -8.400001E+00, +6.599999E+00, +5.999994E-01, +1.599999E+00, +5.999994E-01, -8.400001E+00, -4.400001E+00, +5.599999E+00, +2.599999E+00
1073+4.000000E+00, -7.000000E+00, -9.000000E+00, +8.000000E+00, +6.000000E+00, +5.000000E+00, -6.000000E+00, -2.000000E+00, -5.000000E+00, +6.000000E+00
1074cor = getFilled(0., ndim, ndim)
1075call setCor(cor, uppDia, mean, sample, dim)
1076cor
1077+1.000000E+00, +2.515804E-01
1078+0.000000E+00, +1.000000E+00
1079cor = getFilled(0., ndim, ndim)
1080call setCor(cor, uppDia, samShifted, dim) ! same result as above.
1081cor
1082+1.000000E+00, +2.515804E-01
1083+0.000000E+00, +1.000000E+00
1084
1085cor = getFilled(0., ndim, ndim)
1086call setCor(cor, lowDia, mean, sample, dim)
1087cor
1088+1.000000E+00, +0.000000E+00
1089+2.515804E-01, +1.000000E+00
1090cor = getFilled(0., ndim, ndim)
1091call setCor(cor, lowDia, samShifted, dim) ! same result as above.
1092cor
1093+1.000000E+00, +0.000000E+00
1094+2.515804E-01, +1.000000E+00
1095
1096'Compute the sample correlation along the first dimension.'
1097
1098dim = 1
1099cor = getFilled(0., ndim, ndim)
1100call setCor(cor, uppDia, mean, transpose(sample), dim)
1101cor
1102+1.000000E+00, +2.515804E-01
1103+0.000000E+00, +1.000000E+00
1104cor = getFilled(0., ndim, ndim)
1105call setCor(cor, uppDia, transpose(samShifted), dim) ! same result as above.
1106cor
1107+1.000000E+00, +2.515804E-01
1108+0.000000E+00, +1.000000E+00
1109
1110cor = getFilled(0., ndim, ndim)
1111call setCor(cor, lowDia, mean, transpose(sample), dim)
1112cor
1113+1.000000E+00, +0.000000E+00
1114+2.515804E-01, +1.000000E+00
1115cor = getFilled(0., ndim, ndim)
1116call setCor(cor, lowDia, transpose(samShifted), dim) ! same result as above.
1117cor
1118+1.000000E+00, +0.000000E+00
1119+2.515804E-01, +1.000000E+00
1120
1121'Compute the full sample correlation for a pair of time series.'
1122
1123call setCor(cor(1,1), mean, sample(1,:), sample(2,:))
1124cor(1,1)
1125+2.515804E-01
1126call setCor(cor(1,1), samShifted(1,:), samShifted(2,:)) ! same result as above.
1127cor(1,1)
1128+2.515804E-01
1129
1130
1131!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1132!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1133!Compute the Pearson correlation matrix for a weighted pair of time series.
1134!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1135!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1136
1137ndim = 2; nsam = 10; dim = 2
1138sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
1139sample
1140+1.500000E+01, +4.000000E+00, +1.000000E+01, +1.000000E+00, +1.300000E+01, +1.600000E+01, +1.400000E+01, +1.300000E+01, +9.000000E+00, +5.000000E+00
1141+1.300000E+01, +1.400000E+01, +1.700000E+01, +3.000000E+00, +1.600000E+01, +1.000000E+00, +8.000000E+00, +1.500000E+01, +1.900000E+01, +5.000000E+00
1142call setResized(mean, ndim)
1143iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
1144iweight
1145+2, +3, +10, +1, +9, +10, +4, +9, +5, +10
1146call setMean(mean, sample, dim, iweight, iweisum)
1147mean
1148+1.092064E+01, +1.122222E+01
1149iweisum
1150+63
1151rweight = iweight ! or real-valued weights.
1152iweight
1153+2, +3, +10, +1, +9, +10, +4, +9, +5, +10
1154call setMean(mean, sample, dim, rweight, rweisum)
1155mean
1156+1.092064E+01, +1.122222E+01
1157rweisum
1158+6.300000E+01
1159
1160!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1161!Compute the correlation matrix with integer weights.
1162!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1163
1164samShifted = getShifted(sample, dim, -mean)
1165samShifted
1166+4.079365E+00, -6.920635E+00, -9.206352E-01, -9.920635E+00, +2.079365E+00, +5.079365E+00, +3.079365E+00, +2.079365E+00, -1.920635E+00, -5.920635E+00
1167+1.777777E+00, +2.777777E+00, +5.777777E+00, -8.222223E+00, +4.777777E+00, -1.022222E+01, -3.222223E+00, +3.777777E+00, +7.777777E+00, -6.222223E+00
1168cor = getFilled(0., ndim, ndim)
1169call setCor(cor, uppDia, mean, sample, dim, iweight, iweisum)
1170cor
1171+1.000000E+00, -7.423090E-02
1172+0.000000E+00, +1.000000E+00
1173cor = getFilled(0., ndim, ndim)
1174call setCor(cor, uppDia, samShifted, dim, iweight, iweisum) ! same result as above.
1175cor
1176+1.000000E+00, -7.423090E-02
1177+0.000000E+00, +1.000000E+00
1178
1179cor = getFilled(0., ndim, ndim)
1180call setCor(cor, lowDia, mean, sample, dim, iweight, iweisum)
1181cor
1182+1.000000E+00, +0.000000E+00
1183-7.423091E-02, +1.000000E+00
1184cor = getFilled(0., ndim, ndim)
1185call setCor(cor, lowDia, samShifted, dim, iweight, iweisum) ! same result as above.
1186cor
1187+1.000000E+00, +0.000000E+00
1188-7.423091E-02, +1.000000E+00
1189cor = getFilled(0., ndim, ndim)
1190call setCor(cor, lowDia, getVerbose(samShifted, iweight, sum(iweight), dim), dim) ! same result as above.
1191cor
1192+1.000000E+00, +0.000000E+00
1193-7.423089E-02, +1.000000E+00
1194
1195'Compute the sample correlation along the first dimension.'
1196
1197dim = 1
1198cor = getFilled(0., ndim, ndim)
1199call setCor(cor, uppDia, mean, transpose(sample), dim, iweight, iweisum)
1200cor
1201+1.000000E+00, -7.423090E-02
1202+0.000000E+00, +1.000000E+00
1203cor = getFilled(0., ndim, ndim)
1204call setCor(cor, uppDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.
1205cor
1206+1.000000E+00, -7.423090E-02
1207+0.000000E+00, +1.000000E+00
1208
1209cor = getFilled(0., ndim, ndim)
1210call setCor(cor, lowDia, mean, transpose(sample), dim, iweight, iweisum)
1211cor
1212+1.000000E+00, +0.000000E+00
1213-7.423091E-02, +1.000000E+00
1214cor = getFilled(0., ndim, ndim)
1215call setCor(cor, lowDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.
1216cor
1217+1.000000E+00, +0.000000E+00
1218-7.423091E-02, +1.000000E+00
1219
1220'Compute the full sample correlation for a pair of time series.'
1221
1222call setCor(cor(1,1), mean, sample(1,:), sample(2,:), iweight, iweisum)
1223cor(1,1)
1224-7.423089E-02
1225call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), iweight, iweisum) ! same result as above.
1226cor(1,1)
1227-7.423089E-02
1228
1229
1230!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1231!Compute the correlation matrix with real weights.
1232!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1233
1234dim = 2
1235samShifted = getShifted(sample, dim, -mean)
1236samShifted
1237+4.07936478, -6.92063522, -0.920635223, -9.92063522, +2.07936478, +5.07936478, +3.07936478, +2.07936478, -1.92063522, -5.92063522
1238+1.77777672, +2.77777672, +5.77777672, -8.22222328, +4.77777672, -10.2222233, -3.22222328, +3.77777672, +7.77777672, -6.22222328
1239cor = getFilled(0., ndim, ndim)
1240call setCor(cor, uppDia, mean, sample, dim, rweight, rweisum)
1241cor
1242+1.000000E+00, -7.423090E-02
1243+0.000000E+00, +1.000000E+00
1244cor = getFilled(0., ndim, ndim)
1245call setCor(cor, uppDia, samShifted, dim, rweight, rweisum) ! same result as above.
1246cor
1247+1.000000E+00, -7.423090E-02
1248+0.000000E+00, +1.000000E+00
1249
1250cor = getFilled(0., ndim, ndim)
1251call setCor(cor, lowDia, mean, sample, dim, rweight, rweisum)
1252cor
1253+1.000000E+00, +0.000000E+00
1254-7.423091E-02, +1.000000E+00
1255cor = getFilled(0., ndim, ndim)
1256call setCor(cor, lowDia, samShifted, dim, rweight, rweisum) ! same result as above.
1257cor
1258+1.000000E+00, +0.000000E+00
1259-7.423091E-02, +1.000000E+00
1260
1261'Compute the sample correlation along the first dimension.'
1262
1263dim = 1
1264cor = getFilled(0., ndim, ndim)
1265call setCor(cor, uppDia, mean, transpose(sample), dim, rweight, rweisum)
1266cor
1267+1.000000E+00, -7.423090E-02
1268+0.000000E+00, +1.000000E+00
1269cor = getFilled(0., ndim, ndim)
1270call setCor(cor, uppDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.
1271cor
1272+1.000000E+00, -7.423090E-02
1273+0.000000E+00, +1.000000E+00
1274
1275cor = getFilled(0., ndim, ndim)
1276call setCor(cor, lowDia, mean, transpose(sample), dim, rweight, rweisum)
1277cor
1278+1.000000E+00, +0.000000E+00
1279-7.423091E-02, +1.000000E+00
1280cor = getFilled(0., ndim, ndim)
1281call setCor(cor, lowDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.
1282cor
1283+1.000000E+00, +0.000000E+00
1284-7.423091E-02, +1.000000E+00
1285
1286'Compute the full sample correlation for a pair of time series.'
1287
1288call setCor(cor(1,1), mean, sample(1,:), sample(2,:), rweight, rweisum)
1289cor(1,1)
1290-7.423089E-02
1291call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), rweight, rweisum) ! same result as above.
1292cor(1,1)
1293-7.423089E-02
1294
1295
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: