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+3
10std = getUnifRand(1, ndim, ndim)
11std
12+2.000000E+00, +2.000000E+00, +3.000000E+00
13cov = getCovRand(1._TKG, std)
14cov
15+4.000000E+00, +2.731849E+00, -3.342749E+00
16+2.731849E+00, +4.000000E+00, -2.704294E+00
17-3.342749E+00, -2.704294E+00, +8.999998E+00
18
19
20cor = cov
21call setCor(cor, uppDia, cor, uppDia)
22cor
23+1.000000E+00, +6.829623E-01, -5.571249E-01
24+2.731849E+00, +1.000000E+00, -4.507158E-01
25-3.342749E+00, -2.704294E+00, +1.000000E+00
26
27cor = cov
28call setCor(cor, uppDia, cor, lowDia)
29cor
30+1.000000E+00, +6.829623E-01, -5.571249E-01
31+2.731849E+00, +1.000000E+00, -4.507158E-01
32-3.342749E+00, -2.704294E+00, +1.000000E+00
33
34cor = cov
35call setCor(cor, upp, cor, uppDia)
36cor
37+4.000000E+00, +6.829623E-01, -5.571249E-01
38+2.731849E+00, +4.000000E+00, -4.507158E-01
39-3.342749E+00, -2.704294E+00, +8.999998E+00
40
41cor = cov
42call setCor(cor, upp, cor, lowDia)
43cor
44+4.000000E+00, +6.829623E-01, -5.571249E-01
45+2.731849E+00, +4.000000E+00, -4.507158E-01
46-3.342749E+00, -2.704294E+00, +8.999998E+00
47
48cor = cov
49call setCor(cor, lowDia, cor, uppDia)
50cor
51+1.000000E+00, +2.731849E+00, -3.342749E+00
52+2.731849E+00, +1.000000E+00, -2.704294E+00
53-3.342749E+00, -2.704294E+00, +1.000000E+00
54
55cor = cov
56call setCor(cor, lowDia, cor, lowDia)
57cor
58+1.000000E+00, +2.731849E+00, -3.342749E+00
59+6.829623E-01, +1.000000E+00, -2.704294E+00
60-5.571249E-01, -4.507158E-01, +1.000000E+00
61
62cor = cov
63call setCor(cor, low, cor, uppDia)
64cor
65+4.000000E+00, +2.731849E+00, -3.342749E+00
66+6.829623E-01, +4.000000E+00, -2.704294E+00
67-5.571249E-01, -4.507158E-01, +8.999998E+00
68
69cor = cov
70call setCor(cor, low, cor, lowDia)
71cor
72+4.000000E+00, +2.731849E+00, -3.342749E+00
73+6.829623E-01, +4.000000E+00, -2.704294E+00
74-5.571249E-01, -4.507158E-01, +8.999998E+00
75
76
77cor = cov
78call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
79cor
80+1.000000E+00, +6.829623E-01, -5.571248E-01
81+2.731849E+00, +1.000000E+00, -4.507157E-01
82-3.342749E+00, -2.704294E+00, +1.000000E+00
83
84cor = cov
85call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
86cor
87+1.000000E+00, +6.829623E-01, -5.571248E-01
88+2.731849E+00, +1.000000E+00, -4.507157E-01
89-3.342749E+00, -2.704294E+00, +1.000000E+00
90
91cor = cov
92call setCor(cor, upp, cor, upp, stdinv = 1 / std)
93cor
94+4.00000000, +0.682962298, -0.557124794
95+2.73184919, +4.00000000, -0.450715721
96-3.34274864, -2.70429420, +8.99999809
97
98cor = cov
99call setCor(cor, upp, cor, low, stdinv = 1 / std)
100cor
101+4.000000E+00, +6.829623E-01, -5.571248E-01
102+2.731849E+00, +4.000000E+00, -4.507157E-01
103-3.342749E+00, -2.704294E+00, +8.999998E+00
104
105cor = cov
106call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
107cor
108+1.000000E+00, +2.731849E+00, -3.342749E+00
109+6.829623E-01, +1.000000E+00, -2.704294E+00
110-5.571248E-01, -4.507157E-01, +1.000000E+00
111
112cor = cov
113call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
114cor
115+1.000000E+00, +2.731849E+00, -3.342749E+00
116+6.829623E-01, +1.000000E+00, -2.704294E+00
117-5.571248E-01, -4.507157E-01, +1.000000E+00
118
119cor = cov
120call setCor(cor, low, cor, upp, stdinv = 1 / std)
121cor
122+4.000000E+00, +2.731849E+00, -3.342749E+00
123+6.829623E-01, +4.000000E+00, -2.704294E+00
124-5.571248E-01, -4.507157E-01, +8.999998E+00
125
126cor = cov
127call setCor(cor, low, cor, low, stdinv = 1 / std)
128cor
129+4.000000E+00, +2.731849E+00, -3.342749E+00
130+6.829623E-01, +4.000000E+00, -2.704294E+00
131-5.571248E-01, -4.507157E-01, +8.999998E+00
132
133
134ndim = getUnifRand(0, 4)
135ndim
136+0
137std = getUnifRand(1, ndim, ndim)
138std
139
140cov = getCovRand(1._TKG, std)
141cov
142
143
144cor = cov
145call setCor(cor, uppDia, cor, uppDia)
146cor
147
148cor = cov
149call setCor(cor, uppDia, cor, lowDia)
150cor
151
152cor = cov
153call setCor(cor, upp, cor, uppDia)
154cor
155
156cor = cov
157call setCor(cor, upp, cor, lowDia)
158cor
159
160cor = cov
161call setCor(cor, lowDia, cor, uppDia)
162cor
163
164cor = cov
165call setCor(cor, lowDia, cor, lowDia)
166cor
167
168cor = cov
169call setCor(cor, low, cor, uppDia)
170cor
171
172cor = cov
173call setCor(cor, low, cor, lowDia)
174cor
175
176
177cor = cov
178call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
179cor
180
181cor = cov
182call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
183cor
184
185cor = cov
186call setCor(cor, upp, cor, upp, stdinv = 1 / std)
187cor
188
189cor = cov
190call setCor(cor, upp, cor, low, stdinv = 1 / std)
191cor
192
193cor = cov
194call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
195cor
196
197cor = cov
198call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
199cor
200
201cor = cov
202call setCor(cor, low, cor, upp, stdinv = 1 / std)
203cor
204
205cor = cov
206call setCor(cor, low, cor, low, stdinv = 1 / std)
207cor
208
209
210ndim = getUnifRand(0, 4)
211ndim
212+4
213std = getUnifRand(1, ndim, ndim)
214std
215+2.000000E+00, +4.000000E+00, +4.000000E+00, +4.000000E+00
216cov = getCovRand(1._TKG, std)
217cov
218+4.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
219+6.281831E+00, +1.600000E+01, +1.212170E+01, -1.554213E+01
220+6.033614E+00, +1.212170E+01, +1.600000E+01, -1.003023E+01
221-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.600000E+01
222
223
224cor = cov
225call setCor(cor, uppDia, cor, uppDia)
226cor
227+1.000000E+00, +7.852287E-01, +7.542017E-01, -6.940523E-01
228+6.281831E+00, +1.000000E+00, +7.576063E-01, -9.713830E-01
229+6.033614E+00, +1.212170E+01, +1.000000E+00, -6.268895E-01
230-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.000000E+00
231
232cor = cov
233call setCor(cor, uppDia, cor, lowDia)
234cor
235+1.000000E+00, +7.852287E-01, +7.542017E-01, -6.940523E-01
236+6.281831E+00, +1.000000E+00, +7.576063E-01, -9.713830E-01
237+6.033614E+00, +1.212170E+01, +1.000000E+00, -6.268895E-01
238-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.000000E+00
239
240cor = cov
241call setCor(cor, upp, cor, uppDia)
242cor
243+4.000000E+00, +7.852287E-01, +7.542017E-01, -6.940523E-01
244+6.281831E+00, +1.600000E+01, +7.576063E-01, -9.713830E-01
245+6.033614E+00, +1.212170E+01, +1.600000E+01, -6.268895E-01
246-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.600000E+01
247
248cor = cov
249call setCor(cor, upp, cor, lowDia)
250cor
251+4.000000E+00, +7.852287E-01, +7.542017E-01, -6.940523E-01
252+6.281831E+00, +1.600000E+01, +7.576063E-01, -9.713830E-01
253+6.033614E+00, +1.212170E+01, +1.600000E+01, -6.268895E-01
254-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.600000E+01
255
256cor = cov
257call setCor(cor, lowDia, cor, uppDia)
258cor
259+1.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
260+6.281831E+00, +1.000000E+00, +1.212170E+01, -1.554213E+01
261+6.033614E+00, +1.212170E+01, +1.000000E+00, -1.003023E+01
262-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.000000E+00
263
264cor = cov
265call setCor(cor, lowDia, cor, lowDia)
266cor
267+1.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
268+7.852287E-01, +1.000000E+00, +1.212170E+01, -1.554213E+01
269+7.542017E-01, +7.576063E-01, +1.000000E+00, -1.003023E+01
270-6.940523E-01, -9.713830E-01, -6.268895E-01, +1.000000E+00
271
272cor = cov
273call setCor(cor, low, cor, uppDia)
274cor
275+4.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
276+7.852287E-01, +1.600000E+01, +1.212170E+01, -1.554213E+01
277+7.542017E-01, +7.576063E-01, +1.600000E+01, -1.003023E+01
278-6.940523E-01, -9.713830E-01, -6.268895E-01, +1.600000E+01
279
280cor = cov
281call setCor(cor, low, cor, lowDia)
282cor
283+4.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
284+7.852287E-01, +1.600000E+01, +1.212170E+01, -1.554213E+01
285+7.542017E-01, +7.576063E-01, +1.600000E+01, -1.003023E+01
286-6.940523E-01, -9.713830E-01, -6.268895E-01, +1.600000E+01
287
288
289cor = cov
290call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
291cor
292+1.000000E+00, +7.852288E-01, +7.542017E-01, -6.940523E-01
293+6.281831E+00, +1.000000E+00, +7.576064E-01, -9.713831E-01
294+6.033614E+00, +1.212170E+01, +1.000000E+00, -6.268895E-01
295-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.000000E+00
296
297cor = cov
298call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
299cor
300+1.000000E+00, +7.852288E-01, +7.542017E-01, -6.940523E-01
301+6.281831E+00, +1.000000E+00, +7.576064E-01, -9.713831E-01
302+6.033614E+00, +1.212170E+01, +1.000000E+00, -6.268895E-01
303-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.000000E+00
304
305cor = cov
306call setCor(cor, upp, cor, upp, stdinv = 1 / std)
307cor
308+4.00000000, +0.785228848, +0.754201710, -0.694052339
309+6.28183079, +16.0000038, +0.757606387, -0.971383095
310+6.03361368, +12.1217022, +16.0000019, -0.626889467
311-5.55241871, -15.5421295, -10.0302315, +16.0000019
312
313cor = cov
314call setCor(cor, upp, cor, low, stdinv = 1 / std)
315cor
316+4.000000E+00, +7.852288E-01, +7.542017E-01, -6.940523E-01
317+6.281831E+00, +1.600000E+01, +7.576064E-01, -9.713831E-01
318+6.033614E+00, +1.212170E+01, +1.600000E+01, -6.268895E-01
319-5.552419E+00, -1.554213E+01, -1.003023E+01, +1.600000E+01
320
321cor = cov
322call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
323cor
324+1.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
325+7.852288E-01, +1.000000E+00, +1.212170E+01, -1.554213E+01
326+7.542017E-01, +7.576064E-01, +1.000000E+00, -1.003023E+01
327-6.940523E-01, -9.713831E-01, -6.268895E-01, +1.000000E+00
328
329cor = cov
330call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
331cor
332+1.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
333+7.852288E-01, +1.000000E+00, +1.212170E+01, -1.554213E+01
334+7.542017E-01, +7.576064E-01, +1.000000E+00, -1.003023E+01
335-6.940523E-01, -9.713831E-01, -6.268895E-01, +1.000000E+00
336
337cor = cov
338call setCor(cor, low, cor, upp, stdinv = 1 / std)
339cor
340+4.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
341+7.852288E-01, +1.600000E+01, +1.212170E+01, -1.554213E+01
342+7.542017E-01, +7.576064E-01, +1.600000E+01, -1.003023E+01
343-6.940523E-01, -9.713831E-01, -6.268895E-01, +1.600000E+01
344
345cor = cov
346call setCor(cor, low, cor, low, stdinv = 1 / std)
347cor
348+4.000000E+00, +6.281831E+00, +6.033614E+00, -5.552419E+00
349+7.852288E-01, +1.600000E+01, +1.212170E+01, -1.554213E+01
350+7.542017E-01, +7.576064E-01, +1.600000E+01, -1.003023E+01
351-6.940523E-01, -9.713831E-01, -6.268895E-01, +1.600000E+01
352
353
354ndim = getUnifRand(0, 4)
355ndim
356+4
357std = getUnifRand(1, ndim, ndim)
358std
359+3.000000E+00, +4.000000E+00, +2.000000E+00, +4.000000E+00
360cov = getCovRand(1._TKG, std)
361cov
362+9.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
363+8.984529E+00, +1.600000E+01, -3.727107E+00, +1.758899E+00
364+6.489483E-01, -3.727107E+00, +4.000000E+00, +4.784607E+00
365+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.600000E+01
366
367
368cor = cov
369call setCor(cor, uppDia, cor, uppDia)
370cor
371+1.000000E+00, +7.487109E-01, +1.081581E-01, +3.671595E-01
372+8.984529E+00, +1.000000E+00, -4.658885E-01, +1.099312E-01
373+6.489483E-01, -3.727107E+00, +1.000000E+00, +5.980760E-01
374+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.000000E+00
375
376cor = cov
377call setCor(cor, uppDia, cor, lowDia)
378cor
379+1.000000E+00, +7.487109E-01, +1.081581E-01, +3.671595E-01
380+8.984529E+00, +1.000000E+00, -4.658885E-01, +1.099312E-01
381+6.489483E-01, -3.727107E+00, +1.000000E+00, +5.980760E-01
382+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.000000E+00
383
384cor = cov
385call setCor(cor, upp, cor, uppDia)
386cor
387+9.000000E+00, +7.487109E-01, +1.081581E-01, +3.671595E-01
388+8.984529E+00, +1.600000E+01, -4.658885E-01, +1.099312E-01
389+6.489483E-01, -3.727107E+00, +4.000000E+00, +5.980760E-01
390+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.600000E+01
391
392cor = cov
393call setCor(cor, upp, cor, lowDia)
394cor
395+9.000000E+00, +7.487109E-01, +1.081581E-01, +3.671595E-01
396+8.984529E+00, +1.600000E+01, -4.658885E-01, +1.099312E-01
397+6.489483E-01, -3.727107E+00, +4.000000E+00, +5.980760E-01
398+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.600000E+01
399
400cor = cov
401call setCor(cor, lowDia, cor, uppDia)
402cor
403+1.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
404+8.984529E+00, +1.000000E+00, -3.727107E+00, +1.758899E+00
405+6.489483E-01, -3.727107E+00, +1.000000E+00, +4.784607E+00
406+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.000000E+00
407
408cor = cov
409call setCor(cor, lowDia, cor, lowDia)
410cor
411+1.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
412+7.487109E-01, +1.000000E+00, -3.727107E+00, +1.758899E+00
413+1.081581E-01, -4.658885E-01, +1.000000E+00, +4.784607E+00
414+3.671595E-01, +1.099312E-01, +5.980760E-01, +1.000000E+00
415
416cor = cov
417call setCor(cor, low, cor, uppDia)
418cor
419+9.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
420+7.487109E-01, +1.600000E+01, -3.727107E+00, +1.758899E+00
421+1.081581E-01, -4.658885E-01, +4.000000E+00, +4.784607E+00
422+3.671595E-01, +1.099312E-01, +5.980760E-01, +1.600000E+01
423
424cor = cov
425call setCor(cor, low, cor, lowDia)
426cor
427+9.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
428+7.487109E-01, +1.600000E+01, -3.727107E+00, +1.758899E+00
429+1.081581E-01, -4.658885E-01, +4.000000E+00, +4.784607E+00
430+3.671595E-01, +1.099312E-01, +5.980760E-01, +1.600000E+01
431
432
433cor = cov
434call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
435cor
436+1.000000E+00, +7.487108E-01, +1.081581E-01, +3.671595E-01
437+8.984529E+00, +1.000000E+00, -4.658884E-01, +1.099312E-01
438+6.489483E-01, -3.727107E+00, +1.000000E+00, +5.980759E-01
439+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.000000E+00
440
441cor = cov
442call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
443cor
444+1.000000E+00, +7.487108E-01, +1.081581E-01, +3.671595E-01
445+8.984529E+00, +1.000000E+00, -4.658884E-01, +1.099312E-01
446+6.489483E-01, -3.727107E+00, +1.000000E+00, +5.980759E-01
447+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.000000E+00
448
449cor = cov
450call setCor(cor, upp, cor, upp, stdinv = 1 / std)
451cor
452+9.00000000, +0.748710811, +0.108158052, +0.367159456
453+8.98452950, +15.9999981, -0.465888381, +0.109931186
454+0.648948312, -3.72710705, +3.99999952, +0.598075867
455+4.40591335, +1.75889897, +4.78460693, +15.9999990
456
457cor = cov
458call setCor(cor, upp, cor, low, stdinv = 1 / std)
459cor
460+9.000000E+00, +7.487108E-01, +1.081581E-01, +3.671595E-01
461+8.984529E+00, +1.600000E+01, -4.658884E-01, +1.099312E-01
462+6.489483E-01, -3.727107E+00, +4.000000E+00, +5.980759E-01
463+4.405913E+00, +1.758899E+00, +4.784607E+00, +1.600000E+01
464
465cor = cov
466call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
467cor
468+1.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
469+7.487108E-01, +1.000000E+00, -3.727107E+00, +1.758899E+00
470+1.081581E-01, -4.658884E-01, +1.000000E+00, +4.784607E+00
471+3.671595E-01, +1.099312E-01, +5.980759E-01, +1.000000E+00
472
473cor = cov
474call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
475cor
476+1.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
477+7.487108E-01, +1.000000E+00, -3.727107E+00, +1.758899E+00
478+1.081581E-01, -4.658884E-01, +1.000000E+00, +4.784607E+00
479+3.671595E-01, +1.099312E-01, +5.980759E-01, +1.000000E+00
480
481cor = cov
482call setCor(cor, low, cor, upp, stdinv = 1 / std)
483cor
484+9.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
485+7.487108E-01, +1.600000E+01, -3.727107E+00, +1.758899E+00
486+1.081581E-01, -4.658884E-01, +4.000000E+00, +4.784607E+00
487+3.671595E-01, +1.099312E-01, +5.980759E-01, +1.600000E+01
488
489cor = cov
490call setCor(cor, low, cor, low, stdinv = 1 / std)
491cor
492+9.000000E+00, +8.984529E+00, +6.489483E-01, +4.405913E+00
493+7.487108E-01, +1.600000E+01, -3.727107E+00, +1.758899E+00
494+1.081581E-01, -4.658884E-01, +4.000000E+00, +4.784607E+00
495+3.671595E-01, +1.099312E-01, +5.980759E-01, +1.600000E+01
496
497
498ndim = getUnifRand(0, 4)
499ndim
500+3
501std = getUnifRand(1, ndim, ndim)
502std
503+2.000000E+00, +2.000000E+00, +1.000000E+00
504cov = getCovRand(1._TKG, std)
505cov
506+4.000000E+00, -2.356627E+00, +6.051202E-01
507-2.356627E+00, +4.000000E+00, +1.020138E+00
508+6.051202E-01, +1.020138E+00, +1.000000E+00
509
510
511cor = cov
512call setCor(cor, uppDia, cor, uppDia)
513cor
514+1.000000E+00, -5.891567E-01, +3.025601E-01
515-2.356627E+00, +1.000000E+00, +5.100691E-01
516+6.051202E-01, +1.020138E+00, +1.000000E+00
517
518cor = cov
519call setCor(cor, uppDia, cor, lowDia)
520cor
521+1.000000E+00, -5.891567E-01, +3.025601E-01
522-2.356627E+00, +1.000000E+00, +5.100691E-01
523+6.051202E-01, +1.020138E+00, +1.000000E+00
524
525cor = cov
526call setCor(cor, upp, cor, uppDia)
527cor
528+4.000000E+00, -5.891567E-01, +3.025601E-01
529-2.356627E+00, +4.000000E+00, +5.100691E-01
530+6.051202E-01, +1.020138E+00, +1.000000E+00
531
532cor = cov
533call setCor(cor, upp, cor, lowDia)
534cor
535+4.000000E+00, -5.891567E-01, +3.025601E-01
536-2.356627E+00, +4.000000E+00, +5.100691E-01
537+6.051202E-01, +1.020138E+00, +1.000000E+00
538
539cor = cov
540call setCor(cor, lowDia, cor, uppDia)
541cor
542+1.000000E+00, -2.356627E+00, +6.051202E-01
543-2.356627E+00, +1.000000E+00, +1.020138E+00
544+6.051202E-01, +1.020138E+00, +1.000000E+00
545
546cor = cov
547call setCor(cor, lowDia, cor, lowDia)
548cor
549+1.000000E+00, -2.356627E+00, +6.051202E-01
550-5.891567E-01, +1.000000E+00, +1.020138E+00
551+3.025601E-01, +5.100691E-01, +1.000000E+00
552
553cor = cov
554call setCor(cor, low, cor, uppDia)
555cor
556+4.000000E+00, -2.356627E+00, +6.051202E-01
557-5.891567E-01, +4.000000E+00, +1.020138E+00
558+3.025601E-01, +5.100691E-01, +1.000000E+00
559
560cor = cov
561call setCor(cor, low, cor, lowDia)
562cor
563+4.000000E+00, -2.356627E+00, +6.051202E-01
564-5.891567E-01, +4.000000E+00, +1.020138E+00
565+3.025601E-01, +5.100691E-01, +1.000000E+00
566
567
568cor = cov
569call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
570cor
571+1.000000E+00, -5.891566E-01, +3.025601E-01
572-2.356627E+00, +1.000000E+00, +5.100690E-01
573+6.051202E-01, +1.020138E+00, +1.000000E+00
574
575cor = cov
576call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
577cor
578+1.000000E+00, -5.891566E-01, +3.025601E-01
579-2.356627E+00, +1.000000E+00, +5.100690E-01
580+6.051202E-01, +1.020138E+00, +1.000000E+00
581
582cor = cov
583call setCor(cor, upp, cor, upp, stdinv = 1 / std)
584cor
585+4.00000000, -0.589156628, +0.302560091
586-2.35662651, +3.99999952, +0.510069013
587+0.605120182, +1.02013803, +1.00000000
588
589cor = cov
590call setCor(cor, upp, cor, low, stdinv = 1 / std)
591cor
592+4.000000E+00, -5.891566E-01, +3.025601E-01
593-2.356627E+00, +4.000000E+00, +5.100690E-01
594+6.051202E-01, +1.020138E+00, +1.000000E+00
595
596cor = cov
597call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
598cor
599+1.000000E+00, -2.356627E+00, +6.051202E-01
600-5.891566E-01, +1.000000E+00, +1.020138E+00
601+3.025601E-01, +5.100690E-01, +1.000000E+00
602
603cor = cov
604call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
605cor
606+1.000000E+00, -2.356627E+00, +6.051202E-01
607-5.891566E-01, +1.000000E+00, +1.020138E+00
608+3.025601E-01, +5.100690E-01, +1.000000E+00
609
610cor = cov
611call setCor(cor, low, cor, upp, stdinv = 1 / std)
612cor
613+4.000000E+00, -2.356627E+00, +6.051202E-01
614-5.891566E-01, +4.000000E+00, +1.020138E+00
615+3.025601E-01, +5.100690E-01, +1.000000E+00
616
617cor = cov
618call setCor(cor, low, cor, low, stdinv = 1 / std)
619cor
620+4.000000E+00, -2.356627E+00, +6.051202E-01
621-5.891566E-01, +4.000000E+00, +1.020138E+00
622+3.025601E-01, +5.100690E-01, +1.000000E+00
623
624
625ndim = getUnifRand(0, 4)
626ndim
627+3
628std = getUnifRand(1, ndim, ndim)
629std
630+3.000000E+00, +2.000000E+00, +2.000000E+00
631cov = getCovRand(1._TKG, std)
632cov
633+9.000000E+00, +5.503739E+00, +3.622772E+00
634+5.503739E+00, +3.999999E+00, +1.997328E+00
635+3.622772E+00, +1.997328E+00, +3.999999E+00
636
637
638cor = cov
639call setCor(cor, uppDia, cor, uppDia)
640cor
641+1.000000E+00, +9.172900E-01, +6.037955E-01
642+5.503739E+00, +1.000000E+00, +4.993322E-01
643+3.622772E+00, +1.997328E+00, +1.000000E+00
644
645cor = cov
646call setCor(cor, uppDia, cor, lowDia)
647cor
648+1.000000E+00, +9.172900E-01, +6.037955E-01
649+5.503739E+00, +1.000000E+00, +4.993322E-01
650+3.622772E+00, +1.997328E+00, +1.000000E+00
651
652cor = cov
653call setCor(cor, upp, cor, uppDia)
654cor
655+9.000000E+00, +9.172900E-01, +6.037955E-01
656+5.503739E+00, +3.999999E+00, +4.993322E-01
657+3.622772E+00, +1.997328E+00, +3.999999E+00
658
659cor = cov
660call setCor(cor, upp, cor, lowDia)
661cor
662+9.000000E+00, +9.172900E-01, +6.037955E-01
663+5.503739E+00, +3.999999E+00, +4.993322E-01
664+3.622772E+00, +1.997328E+00, +3.999999E+00
665
666cor = cov
667call setCor(cor, lowDia, cor, uppDia)
668cor
669+1.000000E+00, +5.503739E+00, +3.622772E+00
670+5.503739E+00, +1.000000E+00, +1.997328E+00
671+3.622772E+00, +1.997328E+00, +1.000000E+00
672
673cor = cov
674call setCor(cor, lowDia, cor, lowDia)
675cor
676+1.000000E+00, +5.503739E+00, +3.622772E+00
677+9.172900E-01, +1.000000E+00, +1.997328E+00
678+6.037955E-01, +4.993322E-01, +1.000000E+00
679
680cor = cov
681call setCor(cor, low, cor, uppDia)
682cor
683+9.000000E+00, +5.503739E+00, +3.622772E+00
684+9.172900E-01, +3.999999E+00, +1.997328E+00
685+6.037955E-01, +4.993322E-01, +3.999999E+00
686
687cor = cov
688call setCor(cor, low, cor, lowDia)
689cor
690+9.000000E+00, +5.503739E+00, +3.622772E+00
691+9.172900E-01, +3.999999E+00, +1.997328E+00
692+6.037955E-01, +4.993322E-01, +3.999999E+00
693
694
695cor = cov
696call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
697cor
698+1.000000E+00, +9.172899E-01, +6.037954E-01
699+5.503739E+00, +1.000000E+00, +4.993321E-01
700+3.622772E+00, +1.997328E+00, +1.000000E+00
701
702cor = cov
703call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
704cor
705+1.000000E+00, +9.172899E-01, +6.037954E-01
706+5.503739E+00, +1.000000E+00, +4.993321E-01
707+3.622772E+00, +1.997328E+00, +1.000000E+00
708
709cor = cov
710call setCor(cor, upp, cor, upp, stdinv = 1 / std)
711cor
712+9.00000000, +0.917289913, +0.603795409
713+5.50373936, +3.99999928, +0.499332070
714+3.62277222, +1.99732828, +3.99999905
715
716cor = cov
717call setCor(cor, upp, cor, low, stdinv = 1 / std)
718cor
719+9.000000E+00, +9.172899E-01, +6.037954E-01
720+5.503739E+00, +3.999999E+00, +4.993321E-01
721+3.622772E+00, +1.997328E+00, +3.999999E+00
722
723cor = cov
724call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
725cor
726+1.000000E+00, +5.503739E+00, +3.622772E+00
727+9.172899E-01, +1.000000E+00, +1.997328E+00
728+6.037954E-01, +4.993321E-01, +1.000000E+00
729
730cor = cov
731call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
732cor
733+1.000000E+00, +5.503739E+00, +3.622772E+00
734+9.172899E-01, +1.000000E+00, +1.997328E+00
735+6.037954E-01, +4.993321E-01, +1.000000E+00
736
737cor = cov
738call setCor(cor, low, cor, upp, stdinv = 1 / std)
739cor
740+9.000000E+00, +5.503739E+00, +3.622772E+00
741+9.172899E-01, +3.999999E+00, +1.997328E+00
742+6.037954E-01, +4.993321E-01, +3.999999E+00
743
744cor = cov
745call setCor(cor, low, cor, low, stdinv = 1 / std)
746cor
747+9.000000E+00, +5.503739E+00, +3.622772E+00
748+9.172899E-01, +3.999999E+00, +1.997328E+00
749+6.037954E-01, +4.993321E-01, +3.999999E+00
750
751
752ndim = getUnifRand(0, 4)
753ndim
754+3
755std = getUnifRand(1, ndim, ndim)
756std
757+2.000000E+00, +1.000000E+00, +2.000000E+00
758cov = getCovRand(1._TKG, std)
759cov
760+4.000000E+00, +7.658734E-01, +2.311712E+00
761+7.658734E-01, +1.000000E+00, -9.103572E-02
762+2.311712E+00, -9.103572E-02, +4.000000E+00
763
764
765cor = cov
766call setCor(cor, uppDia, cor, uppDia)
767cor
768+1.000000E+00, +3.829367E-01, +5.779281E-01
769+7.658734E-01, +1.000000E+00, -4.551786E-02
770+2.311712E+00, -9.103572E-02, +1.000000E+00
771
772cor = cov
773call setCor(cor, uppDia, cor, lowDia)
774cor
775+1.000000E+00, +3.829367E-01, +5.779281E-01
776+7.658734E-01, +1.000000E+00, -4.551786E-02
777+2.311712E+00, -9.103572E-02, +1.000000E+00
778
779cor = cov
780call setCor(cor, upp, cor, uppDia)
781cor
782+4.000000E+00, +3.829367E-01, +5.779281E-01
783+7.658734E-01, +1.000000E+00, -4.551786E-02
784+2.311712E+00, -9.103572E-02, +4.000000E+00
785
786cor = cov
787call setCor(cor, upp, cor, lowDia)
788cor
789+4.000000E+00, +3.829367E-01, +5.779281E-01
790+7.658734E-01, +1.000000E+00, -4.551786E-02
791+2.311712E+00, -9.103572E-02, +4.000000E+00
792
793cor = cov
794call setCor(cor, lowDia, cor, uppDia)
795cor
796+1.000000E+00, +7.658734E-01, +2.311712E+00
797+7.658734E-01, +1.000000E+00, -9.103572E-02
798+2.311712E+00, -9.103572E-02, +1.000000E+00
799
800cor = cov
801call setCor(cor, lowDia, cor, lowDia)
802cor
803+1.000000E+00, +7.658734E-01, +2.311712E+00
804+3.829367E-01, +1.000000E+00, -9.103572E-02
805+5.779281E-01, -4.551786E-02, +1.000000E+00
806
807cor = cov
808call setCor(cor, low, cor, uppDia)
809cor
810+4.000000E+00, +7.658734E-01, +2.311712E+00
811+3.829367E-01, +1.000000E+00, -9.103572E-02
812+5.779281E-01, -4.551786E-02, +4.000000E+00
813
814cor = cov
815call setCor(cor, low, cor, lowDia)
816cor
817+4.000000E+00, +7.658734E-01, +2.311712E+00
818+3.829367E-01, +1.000000E+00, -9.103572E-02
819+5.779281E-01, -4.551786E-02, +4.000000E+00
820
821
822cor = cov
823call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
824cor
825+1.000000E+00, +3.829367E-01, +5.779281E-01
826+7.658734E-01, +1.000000E+00, -4.551786E-02
827+2.311712E+00, -9.103572E-02, +1.000000E+00
828
829cor = cov
830call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
831cor
832+1.000000E+00, +3.829367E-01, +5.779281E-01
833+7.658734E-01, +1.000000E+00, -4.551786E-02
834+2.311712E+00, -9.103572E-02, +1.000000E+00
835
836cor = cov
837call setCor(cor, upp, cor, upp, stdinv = 1 / std)
838cor
839+4.00000000, +0.382936716, +0.577928066
840+0.765873432, +1.00000000, -0.455178618E-1
841+2.31171227, -0.910357237E-1, +4.00000000
842
843cor = cov
844call setCor(cor, upp, cor, low, stdinv = 1 / std)
845cor
846+4.000000E+00, +3.829367E-01, +5.779281E-01
847+7.658734E-01, +1.000000E+00, -4.551786E-02
848+2.311712E+00, -9.103572E-02, +4.000000E+00
849
850cor = cov
851call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
852cor
853+1.000000E+00, +7.658734E-01, +2.311712E+00
854+3.829367E-01, +1.000000E+00, -9.103572E-02
855+5.779281E-01, -4.551786E-02, +1.000000E+00
856
857cor = cov
858call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
859cor
860+1.000000E+00, +7.658734E-01, +2.311712E+00
861+3.829367E-01, +1.000000E+00, -9.103572E-02
862+5.779281E-01, -4.551786E-02, +1.000000E+00
863
864cor = cov
865call setCor(cor, low, cor, upp, stdinv = 1 / std)
866cor
867+4.000000E+00, +7.658734E-01, +2.311712E+00
868+3.829367E-01, +1.000000E+00, -9.103572E-02
869+5.779281E-01, -4.551786E-02, +4.000000E+00
870
871cor = cov
872call setCor(cor, low, cor, low, stdinv = 1 / std)
873cor
874+4.000000E+00, +7.658734E-01, +2.311712E+00
875+3.829367E-01, +1.000000E+00, -9.103572E-02
876+5.779281E-01, -4.551786E-02, +4.000000E+00
877
878
879ndim = getUnifRand(0, 4)
880ndim
881+1
882std = getUnifRand(1, ndim, ndim)
883std
884+1.000000E+00
885cov = getCovRand(1._TKG, std)
886cov
887+1.000000E+00
888
889
890cor = cov
891call setCor(cor, uppDia, cor, uppDia)
892cor
893+1.000000E+00
894
895cor = cov
896call setCor(cor, uppDia, cor, lowDia)
897cor
898+1.000000E+00
899
900cor = cov
901call setCor(cor, upp, cor, uppDia)
902cor
903+1.000000E+00
904
905cor = cov
906call setCor(cor, upp, cor, lowDia)
907cor
908+1.000000E+00
909
910cor = cov
911call setCor(cor, lowDia, cor, uppDia)
912cor
913+1.000000E+00
914
915cor = cov
916call setCor(cor, lowDia, cor, lowDia)
917cor
918+1.000000E+00
919
920cor = cov
921call setCor(cor, low, cor, uppDia)
922cor
923+1.000000E+00
924
925cor = cov
926call setCor(cor, low, cor, lowDia)
927cor
928+1.000000E+00
929
930
931cor = cov
932call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
933cor
934+1.000000E+00
935
936cor = cov
937call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
938cor
939+1.000000E+00
940
941cor = cov
942call setCor(cor, upp, cor, upp, stdinv = 1 / std)
943cor
944+1.00000000
945
946cor = cov
947call setCor(cor, upp, cor, low, stdinv = 1 / std)
948cor
949+1.000000E+00
950
951cor = cov
952call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
953cor
954+1.000000E+00
955
956cor = cov
957call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
958cor
959+1.000000E+00
960
961cor = cov
962call setCor(cor, low, cor, upp, stdinv = 1 / std)
963cor
964+1.000000E+00
965
966cor = cov
967call setCor(cor, low, cor, low, stdinv = 1 / std)
968cor
969+1.000000E+00
970
971
972ndim = getUnifRand(0, 4)
973ndim
974+2
975std = getUnifRand(1, ndim, ndim)
976std
977+2.000000E+00, +2.000000E+00
978cov = getCovRand(1._TKG, std)
979cov
980+4.000000E+00, -1.403646E-01
981-1.403646E-01, +3.999999E+00
982
983
984cor = cov
985call setCor(cor, uppDia, cor, uppDia)
986cor
987+1.000000E+00, -3.509115E-02
988-1.403646E-01, +1.000000E+00
989
990cor = cov
991call setCor(cor, uppDia, cor, lowDia)
992cor
993+1.000000E+00, -3.509115E-02
994-1.403646E-01, +1.000000E+00
995
996cor = cov
997call setCor(cor, upp, cor, uppDia)
998cor
999+4.000000E+00, -3.509115E-02
1000-1.403646E-01, +3.999999E+00
1001
1002cor = cov
1003call setCor(cor, upp, cor, lowDia)
1004cor
1005+4.000000E+00, -3.509115E-02
1006-1.403646E-01, +3.999999E+00
1007
1008cor = cov
1009call setCor(cor, lowDia, cor, uppDia)
1010cor
1011+1.000000E+00, -1.403646E-01
1012-1.403646E-01, +1.000000E+00
1013
1014cor = cov
1015call setCor(cor, lowDia, cor, lowDia)
1016cor
1017+1.000000E+00, -1.403646E-01
1018-3.509115E-02, +1.000000E+00
1019
1020cor = cov
1021call setCor(cor, low, cor, uppDia)
1022cor
1023+4.000000E+00, -1.403646E-01
1024-3.509115E-02, +3.999999E+00
1025
1026cor = cov
1027call setCor(cor, low, cor, lowDia)
1028cor
1029+4.000000E+00, -1.403646E-01
1030-3.509115E-02, +3.999999E+00
1031
1032
1033cor = cov
1034call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
1035cor
1036+1.000000E+00, -3.509115E-02
1037-1.403646E-01, +1.000000E+00
1038
1039cor = cov
1040call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
1041cor
1042+1.000000E+00, -3.509115E-02
1043-1.403646E-01, +1.000000E+00
1044
1045cor = cov
1046call setCor(cor, upp, cor, upp, stdinv = 1 / std)
1047cor
1048+4.00000000, -0.350911506E-1
1049-0.140364602, +3.99999928
1050
1051cor = cov
1052call setCor(cor, upp, cor, low, stdinv = 1 / std)
1053cor
1054+4.000000E+00, -3.509115E-02
1055-1.403646E-01, +3.999999E+00
1056
1057cor = cov
1058call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
1059cor
1060+1.000000E+00, -1.403646E-01
1061-3.509115E-02, +1.000000E+00
1062
1063cor = cov
1064call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
1065cor
1066+1.000000E+00, -1.403646E-01
1067-3.509115E-02, +1.000000E+00
1068
1069cor = cov
1070call setCor(cor, low, cor, upp, stdinv = 1 / std)
1071cor
1072+4.000000E+00, -1.403646E-01
1073-3.509115E-02, +3.999999E+00
1074
1075cor = cov
1076call setCor(cor, low, cor, low, stdinv = 1 / std)
1077cor
1078+4.000000E+00, -1.403646E-01
1079-3.509115E-02, +3.999999E+00
1080
1081
1082ndim = getUnifRand(0, 4)
1083ndim
1084+0
1085std = getUnifRand(1, ndim, ndim)
1086std
1087
1088cov = getCovRand(1._TKG, std)
1089cov
1090
1091
1092cor = cov
1093call setCor(cor, uppDia, cor, uppDia)
1094cor
1095
1096cor = cov
1097call setCor(cor, uppDia, cor, lowDia)
1098cor
1099
1100cor = cov
1101call setCor(cor, upp, cor, uppDia)
1102cor
1103
1104cor = cov
1105call setCor(cor, upp, cor, lowDia)
1106cor
1107
1108cor = cov
1109call setCor(cor, lowDia, cor, uppDia)
1110cor
1111
1112cor = cov
1113call setCor(cor, lowDia, cor, lowDia)
1114cor
1115
1116cor = cov
1117call setCor(cor, low, cor, uppDia)
1118cor
1119
1120cor = cov
1121call setCor(cor, low, cor, lowDia)
1122cor
1123
1124
1125cor = cov
1126call setCor(cor, uppDia, cor, upp, stdinv = 1 / std)
1127cor
1128
1129cor = cov
1130call setCor(cor, uppDia, cor, low, stdinv = 1 / std)
1131cor
1132
1133cor = cov
1134call setCor(cor, upp, cor, upp, stdinv = 1 / std)
1135cor
1136
1137cor = cov
1138call setCor(cor, upp, cor, low, stdinv = 1 / std)
1139cor
1140
1141cor = cov
1142call setCor(cor, lowDia, cor, upp, stdinv = 1 / std)
1143cor
1144
1145cor = cov
1146call setCor(cor, lowDia, cor, low, stdinv = 1 / std)
1147cor
1148
1149cor = cov
1150call setCor(cor, low, cor, upp, stdinv = 1 / std)
1151cor
1152
1153cor = cov
1154call setCor(cor, low, cor, low, stdinv = 1 / std)
1155cor
1156
1157
1158!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1159!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1160!Compute the Pearson correlation matrix for a pair of time series.
1161!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1162!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1163
1164ndim = 2; nsam = 10; dim = 2
1165sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
1166sample
1167+1.300000E+01, +3.000000E+00, +1.700000E+01, +1.400000E+01, +8.000000E+00, +2.000000E+00, +2.000000E+00, +8.000000E+00, +1.800000E+01, +2.000000E+01
1168+1.000000E+01, +9.000000E+00, +1.500000E+01, +1.700000E+01, +1.700000E+01, +9.000000E+00, +6.000000E+00, +5.000000E+00, +1.000000E+01, +1.000000E+01
1169mean = getMean(sample, dim)
1170mean
1171+1.050000E+01, +1.080000E+01
1172samShifted = getShifted(sample, dim, -mean)
1173samShifted
1174+2.500000E+00, -7.500000E+00, +6.500000E+00, +3.500000E+00, -2.500000E+00, -8.500000E+00, -8.500000E+00, -2.500000E+00, +7.500000E+00, +9.500000E+00
1175-8.000002E-01, -1.800000E+00, +4.200000E+00, +6.200000E+00, +6.200000E+00, -1.800000E+00, -4.800000E+00, -5.800000E+00, -8.000002E-01, -8.000002E-01
1176cor = getFilled(0., ndim, ndim)
1177call setCor(cor, uppDia, mean, sample, dim)
1178cor
1179+1.000000E+00, +3.937320E-01
1180+0.000000E+00, +1.000000E+00
1181cor = getFilled(0., ndim, ndim)
1182call setCor(cor, uppDia, samShifted, dim) ! same result as above.
1183cor
1184+1.000000E+00, +3.937320E-01
1185+0.000000E+00, +1.000000E+00
1186
1187cor = getFilled(0., ndim, ndim)
1188call setCor(cor, lowDia, mean, sample, dim)
1189cor
1190+1.000000E+00, +0.000000E+00
1191+3.937320E-01, +1.000000E+00
1192cor = getFilled(0., ndim, ndim)
1193call setCor(cor, lowDia, samShifted, dim) ! same result as above.
1194cor
1195+1.000000E+00, +0.000000E+00
1196+3.937320E-01, +1.000000E+00
1197
1198'Compute the sample correlation along the first dimension.'
1199
1200dim = 1
1201cor = getFilled(0., ndim, ndim)
1202call setCor(cor, uppDia, mean, transpose(sample), dim)
1203cor
1204+1.000000E+00, +3.937320E-01
1205+0.000000E+00, +1.000000E+00
1206cor = getFilled(0., ndim, ndim)
1207call setCor(cor, uppDia, transpose(samShifted), dim) ! same result as above.
1208cor
1209+1.000000E+00, +3.937320E-01
1210+0.000000E+00, +1.000000E+00
1211
1212cor = getFilled(0., ndim, ndim)
1213call setCor(cor, lowDia, mean, transpose(sample), dim)
1214cor
1215+1.000000E+00, +0.000000E+00
1216+3.937320E-01, +1.000000E+00
1217cor = getFilled(0., ndim, ndim)
1218call setCor(cor, lowDia, transpose(samShifted), dim) ! same result as above.
1219cor
1220+1.000000E+00, +0.000000E+00
1221+3.937320E-01, +1.000000E+00
1222
1223'Compute the full sample correlation for a pair of time series.'
1224
1225call setCor(cor(1,1), mean, sample(1,:), sample(2,:))
1226cor(1,1)
1227+3.937321E-01
1228call setCor(cor(1,1), samShifted(1,:), samShifted(2,:)) ! same result as above.
1229cor(1,1)
1230+3.937321E-01
1231
1232
1233!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1234!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1235!Compute the Pearson correlation matrix for a weighted pair of time series.
1236!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1237!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1238
1239ndim = 2; nsam = 10; dim = 2
1240sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
1241sample
1242+1.300000E+01, +5.000000E+00, +2.000000E+01, +1.400000E+01, +4.000000E+00, +9.000000E+00, +1.700000E+01, +1.000000E+00, +6.000000E+00, +8.000000E+00
1243+2.000000E+00, +6.000000E+00, +1.000000E+01, +1.700000E+01, +2.000000E+00, +1.200000E+01, +1.300000E+01, +4.000000E+00, +1.700000E+01, +1.200000E+01
1244call setResized(mean, ndim)
1245iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
1246iweight
1247+7, +3, +7, +6, +1, +1, +8, +1, +4, +10
1248call setMean(mean, sample, dim, iweight, iweisum)
1249mean
1250+1.216667E+01, +1.070833E+01
1251iweisum
1252+48
1253rweight = iweight ! or real-valued weights.
1254iweight
1255+7, +3, +7, +6, +1, +1, +8, +1, +4, +10
1256call setMean(mean, sample, dim, rweight, rweisum)
1257mean
1258+1.216667E+01, +1.070833E+01
1259rweisum
1260+4.800000E+01
1261
1262!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1263!Compute the correlation matrix with integer weights.
1264!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1265
1266samShifted = getShifted(sample, dim, -mean)
1267samShifted
1268+8.333330E-01, -7.166667E+00, +7.833333E+00, +1.833333E+00, -8.166667E+00, -3.166667E+00, +4.833333E+00, -1.116667E+01, -6.166667E+00, -4.166667E+00
1269-8.708334E+00, -4.708334E+00, -7.083340E-01, +6.291666E+00, -8.708334E+00, +1.291666E+00, +2.291666E+00, -6.708334E+00, +6.291666E+00, +1.291666E+00
1270cor = getFilled(0., ndim, ndim)
1271call setCor(cor, uppDia, mean, sample, dim, iweight, iweisum)
1272cor
1273+1.000000E+00, +8.198504E-02
1274+0.000000E+00, +1.000000E+00
1275cor = getFilled(0., ndim, ndim)
1276call setCor(cor, uppDia, samShifted, dim, iweight, iweisum) ! same result as above.
1277cor
1278+1.000000E+00, +8.198504E-02
1279+0.000000E+00, +1.000000E+00
1280
1281cor = getFilled(0., ndim, ndim)
1282call setCor(cor, lowDia, mean, sample, dim, iweight, iweisum)
1283cor
1284+1.000000E+00, +0.000000E+00
1285+8.198504E-02, +1.000000E+00
1286cor = getFilled(0., ndim, ndim)
1287call setCor(cor, lowDia, samShifted, dim, iweight, iweisum) ! same result as above.
1288cor
1289+1.000000E+00, +0.000000E+00
1290+8.198504E-02, +1.000000E+00
1291cor = getFilled(0., ndim, ndim)
1292call setCor(cor, lowDia, getVerbose(samShifted, iweight, sum(iweight), dim), dim) ! same result as above.
1293cor
1294+1.000000E+00, +0.000000E+00
1295+8.198504E-02, +1.000000E+00
1296
1297'Compute the sample correlation along the first dimension.'
1298
1299dim = 1
1300cor = getFilled(0., ndim, ndim)
1301call setCor(cor, uppDia, mean, transpose(sample), dim, iweight, iweisum)
1302cor
1303+1.000000E+00, +8.198504E-02
1304+0.000000E+00, +1.000000E+00
1305cor = getFilled(0., ndim, ndim)
1306call setCor(cor, uppDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.
1307cor
1308+1.000000E+00, +8.198504E-02
1309+0.000000E+00, +1.000000E+00
1310
1311cor = getFilled(0., ndim, ndim)
1312call setCor(cor, lowDia, mean, transpose(sample), dim, iweight, iweisum)
1313cor
1314+1.000000E+00, +0.000000E+00
1315+8.198504E-02, +1.000000E+00
1316cor = getFilled(0., ndim, ndim)
1317call setCor(cor, lowDia, transpose(samShifted), dim, iweight, iweisum) ! same result as above.
1318cor
1319+1.000000E+00, +0.000000E+00
1320+8.198504E-02, +1.000000E+00
1321
1322'Compute the full sample correlation for a pair of time series.'
1323
1324call setCor(cor(1,1), mean, sample(1,:), sample(2,:), iweight, iweisum)
1325cor(1,1)
1326+8.198504E-02
1327call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), iweight, iweisum) ! same result as above.
1328cor(1,1)
1329+8.198504E-02
1330
1331
1332!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1333!Compute the correlation matrix with real weights.
1334!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
1335
1336dim = 2
1337samShifted = getShifted(sample, dim, -mean)
1338samShifted
1339+0.833333015, -7.16666698, +7.83333302, +1.83333302, -8.16666698, -3.16666698, +4.83333302, -11.1666670, -6.16666698, -4.16666698
1340-8.70833397, -4.70833397, -0.708333969, +6.29166603, -8.70833397, +1.29166603, +2.29166603, -6.70833397, +6.29166603, +1.29166603
1341cor = getFilled(0., ndim, ndim)
1342call setCor(cor, uppDia, mean, sample, dim, rweight, rweisum)
1343cor
1344+1.000000E+00, +8.198504E-02
1345+0.000000E+00, +1.000000E+00
1346cor = getFilled(0., ndim, ndim)
1347call setCor(cor, uppDia, samShifted, dim, rweight, rweisum) ! same result as above.
1348cor
1349+1.000000E+00, +8.198504E-02
1350+0.000000E+00, +1.000000E+00
1351
1352cor = getFilled(0., ndim, ndim)
1353call setCor(cor, lowDia, mean, sample, dim, rweight, rweisum)
1354cor
1355+1.000000E+00, +0.000000E+00
1356+8.198504E-02, +1.000000E+00
1357cor = getFilled(0., ndim, ndim)
1358call setCor(cor, lowDia, samShifted, dim, rweight, rweisum) ! same result as above.
1359cor
1360+1.000000E+00, +0.000000E+00
1361+8.198504E-02, +1.000000E+00
1362
1363'Compute the sample correlation along the first dimension.'
1364
1365dim = 1
1366cor = getFilled(0., ndim, ndim)
1367call setCor(cor, uppDia, mean, transpose(sample), dim, rweight, rweisum)
1368cor
1369+1.000000E+00, +8.198504E-02
1370+0.000000E+00, +1.000000E+00
1371cor = getFilled(0., ndim, ndim)
1372call setCor(cor, uppDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.
1373cor
1374+1.000000E+00, +8.198504E-02
1375+0.000000E+00, +1.000000E+00
1376
1377cor = getFilled(0., ndim, ndim)
1378call setCor(cor, lowDia, mean, transpose(sample), dim, rweight, rweisum)
1379cor
1380+1.000000E+00, +0.000000E+00
1381+8.198504E-02, +1.000000E+00
1382cor = getFilled(0., ndim, ndim)
1383call setCor(cor, lowDia, transpose(samShifted), dim, rweight, rweisum) ! same result as above.
1384cor
1385+1.000000E+00, +0.000000E+00
1386+8.198504E-02, +1.000000E+00
1387
1388'Compute the full sample correlation for a pair of time series.'
1389
1390call setCor(cor(1,1), mean, sample(1,:), sample(2,:), rweight, rweisum)
1391cor(1,1)
1392+8.198504E-02
1393call setCor(cor(1,1), samShifted(1,:), samShifted(2,:), rweight, rweisum) ! same result as above.
1394cor(1,1)
1395+8.198504E-02
1396
1397
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 at Austin

Definition at line 2323 of file pm_sampleCor.F90.


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