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

Generate and return the (optionally unbiased) covariance matrix of a pair of (potentially weighted) time series x(1:nsam) and y(1:nsam) or of an input (potentially weighted) array of shape (ndim, nsam) or (nsam, ndim) where ndim is the number of data dimensions (the number of data attributes) and nsam is the number of data points.
More...

Detailed Description

Generate and return the (optionally unbiased) covariance matrix of a pair of (potentially weighted) time series x(1:nsam) and y(1:nsam) or of an input (potentially weighted) array of shape (ndim, nsam) or (nsam, ndim) 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 covariance matrix corresponding to an input correlation matrix and vector of standard deviations in arbitrary ndim dimensions.
  2. compute the sample covariance matrix of a random sample of nsam observations in arbitrary ndim dimensional space.
  3. compute the sample covariance matrix of a pair of nsam observations in two separated data vectors x and y.

See the documentation of the parent module pm_sampleCov for algorithmic details and sample covariance matrix definition.

Parameters
[in]cor: The input positive semi-definite square matrix of shape (1:ndim, 1:ndim) of the same type and kind as the input cov, representing the correlation matrix based upon which the output covariance matrix cov is to be computed.
(optional. It must be present if and only if the input arguments std and subsetr are present and the rest are missing.)
[in]subsetr: 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 correlation matrix must be used.
  2. The constant uppDia, implying that only the upper-diagonal subset of the input correlation matrix must be used.
This input argument is merely serves to resolve the different procedures of this generic interface from each other at compile-time.
Although the allowed subset constants imply the use of the diagonal elements, they are by definition assumed to be 1 and therefore are not referenced in the algorithm. (optional. If missing, only the upper triangle will be used. It must be present if and only if the input arguments cor is present.)
[in]std: The input positive vector of shape (1:ndim) of type real of the same kind as the input cov, containing the standard deviations of the ndim data attributes based upon which the output covariance matrix cov is to be computed.
(optional. It must be present if and only if the input argument cor is present.)
[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 sample is 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 sample is missing.)
[in]sample: The input contiguous array of shape (ndim, nsam) or (nsam, ndim) of the same type and kind as the output cov, containing the sample comprised of nsam observations each with ndim attributes.
If sample is a matrix, then the direction along which the covariance matrix is computed is dictated by the input argument dim.
(optional. It must be present if and only if the input argument dim is present and x and y are missing.)
[in]dim: The input scalar integer of default kind IK indicating the dimension of sample along which the covariance 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.)
[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 output cov,
containing the corresponding weights of individual nsam observations in sample.
(optional. default = getFilled(1, nsam). It can be present only if the input arguments cor, subsetr, and std are missing.)
[in]correction: The input scalar object that can be the following:
  1. The constant fweight or an object of type fweight implying a bias correction based on the assumption of frequency weights for the sample observations, even if the weight argument is missing.
    This is the most popular type of correction, also known as the Bessel correction.
  2. The constant rweight or an object of type rweight_type implying a bias correction based on the assumption of reliability weights for the sample observations.
(optional. If missing, no bias-correction will be applied to the output cov.
It can be present only if the input arguments cor, subsetr, and std are missing.)
Returns
cov : The output positive semi-definite 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 computed covariance matrix.
When the input arguments x and y are present, then ndim == 2 holds by definition and the output covariance is of the form,

\begin{equation} \ms{cov} = \begin{bmatrix} \sigma_{xx} && \sigma_{yx} \\ \sigma_{xy} && \sigma_{yy} \end{bmatrix} \end{equation}


Possible calling interfaces

use pm_sampleCov, only: getCov
! correlation matrix to covariance matrix.
cov(1:ndim, 1:ndim) = getCov(cor(1:ndim, 1:ndim), subsetr, std(1:ndim))
! XY time series covariance matrix.
cov(1:2, 1:2) = getCov(x(1:nsam), y(1:nsam) , correction = correction) ! full covariance matrix.
cov(1:2, 1:2) = getCov(x(1:nsam), y(1:nsam), weight(1:nsam), correction = correction) ! full covariance matrix.
! sample covariance matrix.
cov(1:ndim, 1:ndim) = getCov(sample(:,:), dim , correction = correction) ! full covariance matrix.
cov(1:ndim, 1:ndim) = getCov(sample(:,:), dim, weight(1:nsam), correction = correction) ! full covariance matrix.
!
Generate and return the (optionally unbiased) covariance matrix of a pair of (potentially weighted) t...
This module contains classes and procedures for computing the properties related to the covariance ma...
Warning
This generic interface is merely a flexible wrapper around the generic subroutine interface setCov.
As such, all conditions and warnings associated with setCov equally hold for this generic interface.
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.
Note
Note the effects of bias-correction in computing the covariance become noticeable only for very sample sample sizes (i.e., when nsam is small).
See also
getCor
setCor
getCov
setCov
getVar
setVar
getMean
setMean
getCovMerged
setCovMerged
setCovMeanMerged
getVarCorrection
getMeanMerged
setMeanMerged
getVarMerged
setVarMerged


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK
4 use pm_sampleMean, only: getMean
5 use pm_sampleMean, only: setMean
9 use pm_distUnif, only: getUnifRand
10 use pm_arrayFill, only: getFilled
11 use pm_sampleCov, only: uppDia
12 use pm_sampleCov, only: lowDia
13 use pm_sampleCov, only: getCov
14 use pm_sampleCor, only: getCor
15 use pm_io, only: display_type
16 use pm_io, only: getFormat
17
18 implicit none
19
20 integer(IK) :: itry, ntry = 10
21 type(display_type) :: disp
22 character(:), allocatable :: format
23 disp = display_type(file = "main.out.F90")
24
25 call disp%skip()
26 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
27 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
28 call disp%show("!Convert correlation matrix and standard deviation to covariance matrix.")
29 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%skip()
32
33 block
34 use pm_kind, only: TKG => RKS ! All other real types are also supported.
35 use pm_matrixCopy, only: setMatCopy, rdpack
36 use pm_distCov, only: getCovRand
37 integer(IK) :: ndim
38 real(TKG), allocatable :: cov(:,:), cor(:,:), std(:)
39 format = getFormat(mold = [0._TKG], ed = SK_"es", signed = .true._LK)
40 do itry = 1, ntry
41 call disp%skip()
42 call disp%show("ndim = getUnifRand(1, 7)")
43 ndim = getUnifRand(1, 7)
44 call disp%show("ndim")
45 call disp%show( ndim )
46 call disp%show("std = getUnifRand(1, 10, ndim)")
47 std = getUnifRand(1, 10, ndim)
48 call disp%show("std")
49 call disp%show( std , format = format )
50 call disp%show("call setResized(cov, [ndim, ndim])")
51 call setResized(cov, [ndim, ndim])
52 call disp%show("cor = getCovRand(1., ndim)")
53 cor = getCovRand(1., ndim)
54 call disp%show("cor")
55 call disp%show( cor , format = format )
56 call disp%show("cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.")
57 cov = getCov(cor, uppDia, std)
58 call disp%show("cov")
59 call disp%show( cov , format = format )
60 call disp%show("cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.")
61 cov = getCov(cor, lowDia, std)
62 call disp%show("cov")
63 call disp%show( cov , format = format )
64 call disp%show("getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.")
65 call disp%show( getCor(getCov(cor, lowDia, std), lowDia) , format = format )
66 call disp%skip()
67 end do
68 end block
69
70 call disp%skip()
71 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
72 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
73 call disp%show("!Compute the covariance matrix of a 2-D sample.")
74 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
75 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
76 call disp%skip()
77
78 block
79 use pm_kind, only: TKG => RKS ! All other real types are also supported.
80 real(TKG), allocatable :: sample(:,:), cov(:,:), mean(:)
81 integer(IK) :: ndim, nsam
82 call disp%show("ndim = 2; nsam = 10")
83 ndim = 2; nsam = 10
84 call disp%show("sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
85 sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
86 call disp%show("sample")
87 call disp%show( sample )
88 call disp%skip()
89 call disp%show("Compute the sample covariance along the second dimension.", deliml = SK_'''')
90 call disp%skip()
91 call disp%show("cov = getCov(sample, dim = 2_IK)")
92 cov = getCov(sample, dim = 2_IK)
93 call disp%show("cov")
94 call disp%show( cov )
95 call disp%skip()
96 call disp%show("Compute the sample covariance along the first dimension.", deliml = SK_'''')
97 call disp%skip()
98 call disp%show("cov = getCov(transpose(sample), dim = 1_IK)")
99 cov = getCov(transpose(sample), dim = 1_IK)
100 call disp%show("cov")
101 call disp%show( cov )
102 call disp%skip()
103 call disp%show("Compute the full sample covariance for a pair of time series.", deliml = SK_'''')
104 call disp%skip()
105 call disp%show("cov = getCov(sample(1,:), sample(2,:))")
106 cov = getCov(sample(1,:), sample(2,:))
107 call disp%show("cov")
108 call disp%show( cov )
109 call disp%skip()
110 end block
111
112 block
113 use pm_kind, only: TKG => RKS ! All other real types are also supported.
114 complex(TKG), allocatable :: sample(:,:), cov(:,:), mean(:)
115 integer(IK) :: ndim, nsam
116 call disp%show("ndim = 2; nsam = 10")
117 ndim = 2; nsam = 10
118 call disp%show("sample = reshape(cmplx(getUnifRand(1, 20, ndim * nsam), -getUnifRand(1, 20, ndim * nsam), TKG), shape = [ndim, nsam], order = [2, 1])")
119 sample = reshape(cmplx(getUnifRand(1, 20, ndim * nsam), -getUnifRand(1, 20, ndim * nsam), TKG), shape = [ndim, nsam], order = [2, 1])
120 call disp%show("sample")
121 call disp%show( sample )
122 call disp%skip()
123 call disp%show("Compute the sample covariance along the second dimension.", deliml = SK_'''')
124 call disp%skip()
125 call disp%show("cov = getCov(sample, dim = 2_IK)")
126 cov = getCov(sample, dim = 2_IK)
127 call disp%show("cov")
128 call disp%show( cov )
129 call disp%skip()
130 call disp%show("Compute the sample covariance along the first dimension.", deliml = SK_'''')
131 call disp%skip()
132 call disp%show("cov = getCov(transpose(sample), dim = 1_IK)")
133 cov = getCov(transpose(sample), dim = 1_IK)
134 call disp%show("cov")
135 call disp%show( cov )
136 call disp%skip()
137 call disp%show("Compute the full sample covariance for a pair of time series.", deliml = SK_'''')
138 call disp%skip()
139 call disp%show("cov = getCov(sample(1,:), sample(2,:))")
140 cov = getCov(sample(1,:), sample(2,:))
141 call disp%show("cov")
142 call disp%show( cov )
143 call disp%skip()
144 end block
145
146 call disp%skip()
147 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
148 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
149 call disp%show("!Compute the biased covariance matrix of a weighted 2-D sample.")
150 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
151 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
152 call disp%skip()
153
154 block
155 use pm_kind, only: TKG => RKS ! All other real types are also supported.
156 real(TKG) :: rweisum
157 integer(IK) :: iweisum
158 real(TKG), allocatable :: rweight(:)
159 integer(IK), allocatable :: iweight(:)
160 real(TKG), allocatable :: sample(:,:), cov(:,:), mean(:)
161 integer(IK) :: ndim, nsam
162 call disp%show("ndim = 2; nsam = 10")
163 ndim = 2; nsam = 10
164 call disp%show("sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
165 sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
166 call disp%show("sample")
167 call disp%show( sample )
168 call disp%show("call setResized(mean, ndim)")
169 call setResized(mean, ndim)
170 call disp%show("iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.")
171 iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
172 call disp%show("iweight")
173 call disp%show( iweight )
174 call disp%show("call setMean(mean, sample, 2_IK, iweight, iweisum)")
175 call setMean(mean, sample, 2_IK, iweight, iweisum)
176 call disp%show("mean")
177 call disp%show( mean )
178 call disp%show("iweisum")
179 call disp%show( iweisum )
180 call disp%show("rweight = iweight ! or real-valued weights.")
181 rweight = iweight ! or real-valued weights.
182 call disp%show("iweight")
183 call disp%show( iweight )
184 call disp%show("call setMean(mean, sample, 2_IK, rweight, rweisum)")
185 call setMean(mean, sample, 2_IK, rweight, rweisum)
186 call disp%show("mean")
187 call disp%show( mean )
188 call disp%show("rweisum")
189 call disp%show( rweisum )
190
191 call disp%skip()
192 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
193 call disp%show("!Compute the covariance matrix with integer weights.")
194 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
195 call disp%skip()
196
197 call disp%show("cov = getCov(sample, 2_IK, iweight)")
198 cov = getCov(sample, 2_IK, iweight)
199 call disp%show("cov")
200 call disp%show( cov )
201 call disp%skip()
202 call disp%show("Compute the sample covariance along the first dimension.", deliml = SK_'''')
203 call disp%skip()
204 call disp%show("cov = getCov(transpose(sample), 1_IK, iweight)")
205 cov = getCov(transpose(sample), 1_IK, iweight)
206 call disp%show("cov")
207 call disp%show( cov )
208 call disp%skip()
209 call disp%show("Compute the full sample covariance for a pair of time series.", deliml = SK_'''')
210 call disp%skip()
211 call disp%show("cov = getCov(sample(1,:), sample(2,:), weight = iweight)")
212 cov = getCov(sample(1,:), sample(2,:), weight = iweight)
213 call disp%show("cov")
214 call disp%show( cov )
215 call disp%skip()
216
217 call disp%skip()
218 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
219 call disp%show("!Compute the covariance matrix with real weights.")
220 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
221 call disp%skip()
222
223 call disp%show("cov = getCov(sample, 2_IK, rweight)")
224 cov = getCov(sample, 2_IK, rweight)
225 call disp%show("cov")
226 call disp%show( cov )
227 call disp%skip()
228 call disp%show("Compute the sample covariance along the first dimension.", deliml = SK_'''')
229 call disp%skip()
230 call disp%show("cov = getCov(transpose(sample), 1_IK, rweight)")
231 cov = getCov(transpose(sample), 1_IK, rweight)
232 call disp%show("cov")
233 call disp%show( cov )
234 call disp%skip()
235 call disp%show("Compute the full sample covariance for a pair of time series.", deliml = SK_'''')
236 call disp%skip()
237 call disp%show("cov = getCov(sample(1,:), sample(2,:), weight = rweight)")
238 cov = getCov(sample(1,:), sample(2,:), weight = rweight)
239 call disp%show("cov")
240 call disp%show( cov )
241 call disp%skip()
242 end block
243
244 block
245 use pm_kind, only: TKG => RKS ! All other real types are also supported.
246 real(TKG) :: rweisum
247 integer(IK) :: iweisum
248 real(TKG), allocatable :: rweight(:)
249 integer(IK), allocatable :: iweight(:)
250 complex(TKG), allocatable :: sample(:,:), cov(:,:), mean(:)
251 integer(IK) :: ndim, nsam
252 call disp%show("ndim = 2; nsam = 10")
253 ndim = 2; nsam = 10
254 call disp%show("sample = reshape(cmplx(getUnifRand(1, 20, ndim * nsam), -getUnifRand(1, 20, ndim * nsam), TKG), shape = [ndim, nsam], order = [2, 1])")
255 sample = reshape(cmplx(getUnifRand(1, 20, ndim * nsam), -getUnifRand(1, 20, ndim * nsam), TKG), shape = [ndim, nsam], order = [2, 1])
256 call disp%show("sample")
257 call disp%show( sample )
258 call disp%show("call setResized(mean, ndim)")
259 call setResized(mean, ndim)
260 call disp%show("iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.")
261 iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
262 call disp%show("iweight")
263 call disp%show( iweight )
264 call disp%show("call setMean(mean, sample, 2_IK, iweight, iweisum)")
265 call setMean(mean, sample, 2_IK, iweight, iweisum)
266 call disp%show("mean")
267 call disp%show( mean )
268 call disp%show("iweisum")
269 call disp%show( iweisum )
270 call disp%show("rweight = iweight ! or real-valued weights.")
271 rweight = iweight ! or real-valued weights.
272 call disp%show("iweight")
273 call disp%show( iweight )
274 call disp%show("call setMean(mean, sample, 2_IK, rweight, rweisum)")
275 call setMean(mean, sample, 2_IK, rweight, rweisum)
276 call disp%show("mean")
277 call disp%show( mean )
278 call disp%show("rweisum")
279 call disp%show( rweisum )
280
281 call disp%skip()
282 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
283 call disp%show("!Compute the covariance matrix with integer weights.")
284 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
285 call disp%skip()
286
287 call disp%show("cov = getCov(sample, 2_IK, iweight)")
288 cov = getCov(sample, 2_IK, iweight)
289 call disp%show("cov")
290 call disp%show( cov )
291 call disp%skip()
292 call disp%show("Compute the sample covariance along the first dimension.", deliml = SK_'''')
293 call disp%skip()
294 call disp%show("cov = getCov(transpose(sample), 1_IK, iweight)")
295 cov = getCov(transpose(sample), 1_IK, iweight)
296 call disp%show("cov")
297 call disp%show( cov )
298 call disp%skip()
299 call disp%show("Compute the full sample covariance for a pair of time series.", deliml = SK_'''')
300 call disp%skip()
301 call disp%show("cov = getCov(sample(1,:), sample(2,:), weight = iweight)")
302 cov = getCov(sample(1,:), sample(2,:), weight = iweight)
303 call disp%show("cov")
304 call disp%show( cov )
305 call disp%skip()
306
307 call disp%skip()
308 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
309 call disp%show("!Compute the covariance matrix with real weights.")
310 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
311 call disp%skip()
312
313 call disp%show("cov = getCov(sample, 2_IK, rweight)")
314 cov = getCov(sample, 2_IK, rweight)
315 call disp%show("cov")
316 call disp%show( cov )
317 call disp%skip()
318 call disp%show("Compute the sample covariance along the first dimension.", deliml = SK_'''')
319 call disp%skip()
320 call disp%show("cov = getCov(transpose(sample), 1_IK, rweight)")
321 cov = getCov(transpose(sample), 1_IK, rweight)
322 call disp%show("cov")
323 call disp%show( cov )
324 call disp%skip()
325 call disp%show("Compute the full sample covariance for a pair of time series.", deliml = SK_'''')
326 call disp%skip()
327 call disp%show("cov = getCov(sample(1,:), sample(2,:), weight = rweight)")
328 cov = getCov(sample(1,:), sample(2,:), weight = rweight)
329 call disp%show("cov")
330 call disp%show( cov )
331 call disp%skip()
332 end block
333
334end program example
Generate and return an array of the specified rank and shape of arbitrary intrinsic type and kind wit...
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 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
Copy a desired subset of the input source matrix of arbitrary shape (:) or (:,:) to the target subset...
Generate and return the (Pearson) correlation coefficient or matrix of a pair of (weighted) time seri...
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 resizing allocatable arrays of various typ...
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
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 procedures and generic interfaces relevant to copying (diagonal or upper/lower t...
This module contains classes and procedures for computing properties related to the correlation matri...
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!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
4!Convert correlation matrix and standard deviation to covariance matrix.
5!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
6!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
7
8
9ndim = getUnifRand(1, 7)
10ndim
11+1
12std = getUnifRand(1, 10, ndim)
13std
14+8.000000E+00
15call setResized(cov, [ndim, ndim])
16cor = getCovRand(1., ndim)
17cor
18+1.000000E+00
19cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
20cov
21+6.400000E+01
22cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
23cov
24+6.400000E+01
25getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
26+1.000000E+00
27
28
29ndim = getUnifRand(1, 7)
30ndim
31+5
32std = getUnifRand(1, 10, ndim)
33std
34+8.000000E+00, +5.000000E+00, +8.000000E+00, +8.000000E+00, +7.000000E+00
35call setResized(cov, [ndim, ndim])
36cor = getCovRand(1., ndim)
37cor
38+1.000000E+00, -9.370748E-01, +8.963406E-01, +2.980443E-01, +5.642448E-01
39-9.370748E-01, +1.000000E+00, -9.789037E-01, -5.025175E-01, -3.917778E-01
40+8.963406E-01, -9.789037E-01, +1.000000E+00, +6.118757E-01, +3.180515E-01
41+2.980443E-01, -5.025175E-01, +6.118757E-01, +1.000000E+00, +4.018869E-02
42+5.642448E-01, -3.917778E-01, +3.180515E-01, +4.018869E-02, +1.000000E+00
43cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
44cov
45+6.400000E+01, -3.748299E+01, +5.736580E+01, +1.907484E+01, +3.159771E+01
46-3.748299E+01, +2.500000E+01, -3.915615E+01, -2.010070E+01, -1.371222E+01
47+5.736580E+01, -3.915615E+01, +6.400000E+01, +3.916005E+01, +1.781088E+01
48+1.907484E+01, -2.010070E+01, +3.916005E+01, +6.400000E+01, +2.250566E+00
49+3.159771E+01, -1.371222E+01, +1.781088E+01, +2.250566E+00, +4.900000E+01
50cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
51cov
52+6.400000E+01, -3.748299E+01, +5.736580E+01, +1.907484E+01, +3.159771E+01
53-3.748299E+01, +2.500000E+01, -3.915615E+01, -2.010070E+01, -1.371222E+01
54+5.736580E+01, -3.915615E+01, +6.400000E+01, +3.916005E+01, +1.781088E+01
55+1.907484E+01, -2.010070E+01, +3.916005E+01, +6.400000E+01, +2.250566E+00
56+3.159771E+01, -1.371222E+01, +1.781088E+01, +2.250566E+00, +4.900000E+01
57getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
58+1.000000E+00, -9.370748E-01, +8.963406E-01, +2.980443E-01, +5.642449E-01
59-9.370748E-01, +1.000000E+00, -9.789037E-01, -5.025175E-01, -3.917778E-01
60+8.963406E-01, -9.789037E-01, +1.000000E+00, +6.118757E-01, +3.180515E-01
61+2.980443E-01, -5.025175E-01, +6.118757E-01, +1.000000E+00, +4.018869E-02
62+5.642449E-01, -3.917778E-01, +3.180515E-01, +4.018869E-02, +1.000000E+00
63
64
65ndim = getUnifRand(1, 7)
66ndim
67+6
68std = getUnifRand(1, 10, ndim)
69std
70+6.000000E+00, +5.000000E+00, +6.000000E+00, +5.000000E+00, +2.000000E+00, +5.000000E+00
71call setResized(cov, [ndim, ndim])
72cor = getCovRand(1., ndim)
73cor
74+1.000000E+00, +7.692201E-01, +5.479082E-01, -3.191917E-01, +1.859854E-01, +6.699873E-01
75+7.692201E-01, +1.000000E+00, +8.324206E-02, +9.449770E-02, -2.210099E-01, +7.964961E-01
76+5.479082E-01, +8.324206E-02, +1.000000E+00, -2.154084E-01, +6.956096E-01, -1.365258E-01
77-3.191917E-01, +9.449770E-02, -2.154084E-01, +1.000000E+00, -4.653472E-01, -7.786135E-02
78+1.859854E-01, -2.210099E-01, +6.956096E-01, -4.653472E-01, +1.000000E+00, -1.932095E-01
79+6.699873E-01, +7.964961E-01, -1.365258E-01, -7.786135E-02, -1.932095E-01, +1.000000E+00
80cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
81cov
82+3.600000E+01, +2.307660E+01, +1.972470E+01, -9.575750E+00, +2.231825E+00, +2.009962E+01
83+2.307660E+01, +2.500000E+01, +2.497262E+00, +2.362442E+00, -2.210099E+00, +1.991240E+01
84+1.972470E+01, +2.497262E+00, +3.600000E+01, -6.462253E+00, +8.347315E+00, -4.095774E+00
85-9.575750E+00, +2.362442E+00, -6.462253E+00, +2.500000E+01, -4.653472E+00, -1.946534E+00
86+2.231825E+00, -2.210099E+00, +8.347315E+00, -4.653472E+00, +4.000000E+00, -1.932095E+00
87+2.009962E+01, +1.991240E+01, -4.095774E+00, -1.946534E+00, -1.932095E+00, +2.500000E+01
88cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
89cov
90+3.600000E+01, +2.307660E+01, +1.972470E+01, -9.575750E+00, +2.231825E+00, +2.009962E+01
91+2.307660E+01, +2.500000E+01, +2.497262E+00, +2.362442E+00, -2.210099E+00, +1.991240E+01
92+1.972470E+01, +2.497262E+00, +3.600000E+01, -6.462253E+00, +8.347315E+00, -4.095774E+00
93-9.575750E+00, +2.362442E+00, -6.462253E+00, +2.500000E+01, -4.653472E+00, -1.946534E+00
94+2.231825E+00, -2.210099E+00, +8.347315E+00, -4.653472E+00, +4.000000E+00, -1.932095E+00
95+2.009962E+01, +1.991240E+01, -4.095774E+00, -1.946534E+00, -1.932095E+00, +2.500000E+01
96getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
97+1.000000E+00, +7.692201E-01, +5.479083E-01, -3.191917E-01, +1.859854E-01, +6.699873E-01
98+7.692201E-01, +1.000000E+00, +8.324207E-02, +9.449770E-02, -2.210099E-01, +7.964962E-01
99+5.479083E-01, +8.324207E-02, +1.000000E+00, -2.154084E-01, +6.956096E-01, -1.365258E-01
100-3.191917E-01, +9.449770E-02, -2.154084E-01, +1.000000E+00, -4.653473E-01, -7.786135E-02
101+1.859854E-01, -2.210099E-01, +6.956096E-01, -4.653473E-01, +1.000000E+00, -1.932095E-01
102+6.699873E-01, +7.964962E-01, -1.365258E-01, -7.786135E-02, -1.932095E-01, +1.000000E+00
103
104
105ndim = getUnifRand(1, 7)
106ndim
107+5
108std = getUnifRand(1, 10, ndim)
109std
110+1.000000E+01, +2.000000E+00, +2.000000E+00, +6.000000E+00, +2.000000E+00
111call setResized(cov, [ndim, ndim])
112cor = getCovRand(1., ndim)
113cor
114+1.000000E+00, +2.443270E-03, -8.852635E-01, -1.088063E-01, +2.729474E-01
115+2.443270E-03, +1.000000E+00, +3.190553E-01, -9.789006E-01, -5.788258E-01
116-8.852635E-01, +3.190553E-01, +1.000000E+00, -1.799237E-01, -3.007606E-01
117-1.088063E-01, -9.789006E-01, -1.799237E-01, +1.000000E+00, +6.657165E-01
118+2.729474E-01, -5.788258E-01, -3.007606E-01, +6.657165E-01, +1.000000E+00
119cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
120cov
121+1.000000E+02, +4.886541E-02, -1.770527E+01, -6.528380E+00, +5.458947E+00
122+4.886541E-02, +4.000000E+00, +1.276221E+00, -1.174681E+01, -2.315303E+00
123-1.770527E+01, +1.276221E+00, +4.000000E+00, -2.159085E+00, -1.203043E+00
124-6.528380E+00, -1.174681E+01, -2.159085E+00, +3.600000E+01, +7.988598E+00
125+5.458947E+00, -2.315303E+00, -1.203043E+00, +7.988598E+00, +4.000000E+00
126cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
127cov
128+1.000000E+02, +4.886541E-02, -1.770527E+01, -6.528380E+00, +5.458947E+00
129+4.886541E-02, +4.000000E+00, +1.276221E+00, -1.174681E+01, -2.315303E+00
130-1.770527E+01, +1.276221E+00, +4.000000E+00, -2.159085E+00, -1.203043E+00
131-6.528380E+00, -1.174681E+01, -2.159085E+00, +3.600000E+01, +7.988598E+00
132+5.458947E+00, -2.315303E+00, -1.203043E+00, +7.988598E+00, +4.000000E+00
133getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
134+1.000000E+00, +2.443271E-03, -8.852636E-01, -1.088063E-01, +2.729474E-01
135+2.443271E-03, +1.000000E+00, +3.190553E-01, -9.789006E-01, -5.788258E-01
136-8.852636E-01, +3.190553E-01, +1.000000E+00, -1.799237E-01, -3.007606E-01
137-1.088063E-01, -9.789006E-01, -1.799237E-01, +1.000000E+00, +6.657165E-01
138+2.729474E-01, -5.788258E-01, -3.007606E-01, +6.657165E-01, +1.000000E+00
139
140
141ndim = getUnifRand(1, 7)
142ndim
143+2
144std = getUnifRand(1, 10, ndim)
145std
146+9.000000E+00, +7.000000E+00
147call setResized(cov, [ndim, ndim])
148cor = getCovRand(1., ndim)
149cor
150+1.000000E+00, +5.869347E-01
151+5.869347E-01, +1.000000E+00
152cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
153cov
154+8.100000E+01, +3.697689E+01
155+3.697689E+01, +4.900000E+01
156cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
157cov
158+8.100000E+01, +3.697689E+01
159+3.697689E+01, +4.900000E+01
160getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
161+1.000000E+00, +5.869348E-01
162+5.869348E-01, +1.000000E+00
163
164
165ndim = getUnifRand(1, 7)
166ndim
167+6
168std = getUnifRand(1, 10, ndim)
169std
170+2.000000E+00, +9.000000E+00, +1.000000E+00, +2.000000E+00, +2.000000E+00, +5.000000E+00
171call setResized(cov, [ndim, ndim])
172cor = getCovRand(1., ndim)
173cor
174+1.000000E+00, -5.708815E-01, +7.001553E-01, -5.786425E-01, -5.864874E-01, +7.140691E-01
175-5.708815E-01, +1.000000E+00, -7.327796E-01, -1.113626E-01, +5.871239E-01, -4.990965E-01
176+7.001553E-01, -7.327796E-01, +1.000000E+00, -2.631226E-01, -3.553498E-01, +6.805356E-01
177-5.786425E-01, -1.113626E-01, -2.631226E-01, +1.000000E+00, -1.560526E-01, -7.037156E-01
178-5.864874E-01, +5.871239E-01, -3.553498E-01, -1.560526E-01, +1.000000E+00, +5.629885E-02
179+7.140691E-01, -4.990965E-01, +6.805356E-01, -7.037156E-01, +5.629885E-02, +1.000000E+00
180cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
181cov
182+4.000000E+00, -1.027587E+01, +1.400311E+00, -2.314570E+00, -2.345950E+00, +7.140691E+00
183-1.027587E+01, +8.100000E+01, -6.595016E+00, -2.004526E+00, +1.056823E+01, -2.245934E+01
184+1.400311E+00, -6.595016E+00, +1.000000E+00, -5.262451E-01, -7.106996E-01, +3.402678E+00
185-2.314570E+00, -2.004526E+00, -5.262451E-01, +4.000000E+00, -6.242104E-01, -7.037156E+00
186-2.345950E+00, +1.056823E+01, -7.106996E-01, -6.242104E-01, +4.000000E+00, +5.629885E-01
187+7.140691E+00, -2.245934E+01, +3.402678E+00, -7.037156E+00, +5.629885E-01, +2.500000E+01
188cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
189cov
190+4.000000E+00, -1.027587E+01, +1.400311E+00, -2.314570E+00, -2.345950E+00, +7.140691E+00
191-1.027587E+01, +8.100000E+01, -6.595016E+00, -2.004526E+00, +1.056823E+01, -2.245934E+01
192+1.400311E+00, -6.595016E+00, +1.000000E+00, -5.262451E-01, -7.106996E-01, +3.402678E+00
193-2.314570E+00, -2.004526E+00, -5.262451E-01, +4.000000E+00, -6.242104E-01, -7.037156E+00
194-2.345950E+00, +1.056823E+01, -7.106996E-01, -6.242104E-01, +4.000000E+00, +5.629885E-01
195+7.140691E+00, -2.245934E+01, +3.402678E+00, -7.037156E+00, +5.629885E-01, +2.500000E+01
196getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
197+1.000000E+00, -5.708815E-01, +7.001553E-01, -5.786425E-01, -5.864874E-01, +7.140691E-01
198-5.708815E-01, +1.000000E+00, -7.327796E-01, -1.113626E-01, +5.871239E-01, -4.990965E-01
199+7.001553E-01, -7.327796E-01, +1.000000E+00, -2.631226E-01, -3.553498E-01, +6.805356E-01
200-5.786425E-01, -1.113626E-01, -2.631226E-01, +1.000000E+00, -1.560526E-01, -7.037156E-01
201-5.864874E-01, +5.871239E-01, -3.553498E-01, -1.560526E-01, +1.000000E+00, +5.629885E-02
202+7.140691E-01, -4.990965E-01, +6.805356E-01, -7.037156E-01, +5.629885E-02, +1.000000E+00
203
204
205ndim = getUnifRand(1, 7)
206ndim
207+2
208std = getUnifRand(1, 10, ndim)
209std
210+7.000000E+00, +7.000000E+00
211call setResized(cov, [ndim, ndim])
212cor = getCovRand(1., ndim)
213cor
214+1.000000E+00, -6.439612E-01
215-6.439612E-01, +1.000000E+00
216cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
217cov
218+4.900000E+01, -3.155410E+01
219-3.155410E+01, +4.900000E+01
220cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
221cov
222+4.900000E+01, -3.155410E+01
223-3.155410E+01, +4.900000E+01
224getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
225+1.000000E+00, -6.439613E-01
226-6.439613E-01, +1.000000E+00
227
228
229ndim = getUnifRand(1, 7)
230ndim
231+7
232std = getUnifRand(1, 10, ndim)
233std
234+5.000000E+00, +1.000000E+00, +1.000000E+00, +4.000000E+00, +1.000000E+00, +3.000000E+00, +1.000000E+00
235call setResized(cov, [ndim, ndim])
236cor = getCovRand(1., ndim)
237cor
238+1.000000E+00, -6.139759E-01, +9.018455E-01, +2.172033E-01, -7.290598E-02, +3.643507E-01, +5.489259E-01
239-6.139759E-01, +1.000000E+00, -6.272851E-01, -7.969132E-02, +2.896488E-01, -7.127306E-01, -7.307767E-01
240+9.018455E-01, -6.272851E-01, +1.000000E+00, +2.123346E-01, -2.867150E-01, +6.079513E-01, +5.876704E-01
241+2.172033E-01, -7.969132E-02, +2.123346E-01, +1.000000E+00, +7.303581E-01, +4.086391E-01, +5.376651E-01
242-7.290598E-02, +2.896488E-01, -2.867150E-01, +7.303581E-01, +1.000000E+00, -2.314590E-01, +1.470055E-01
243+3.643507E-01, -7.127306E-01, +6.079513E-01, +4.086391E-01, -2.314590E-01, +1.000000E+00, +5.986668E-01
244+5.489259E-01, -7.307767E-01, +5.876704E-01, +5.376651E-01, +1.470055E-01, +5.986668E-01, +1.000000E+00
245cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
246cov
247+2.500000E+01, -3.069880E+00, +4.509227E+00, +4.344065E+00, -3.645299E-01, +5.465260E+00, +2.744630E+00
248-3.069880E+00, +1.000000E+00, -6.272851E-01, -3.187653E-01, +2.896488E-01, -2.138192E+00, -7.307767E-01
249+4.509227E+00, -6.272851E-01, +1.000000E+00, +8.493386E-01, -2.867150E-01, +1.823854E+00, +5.876704E-01
250+4.344065E+00, -3.187653E-01, +8.493386E-01, +1.600000E+01, +2.921432E+00, +4.903669E+00, +2.150661E+00
251-3.645299E-01, +2.896488E-01, -2.867150E-01, +2.921432E+00, +1.000000E+00, -6.943769E-01, +1.470055E-01
252+5.465260E+00, -2.138192E+00, +1.823854E+00, +4.903669E+00, -6.943769E-01, +9.000000E+00, +1.796000E+00
253+2.744630E+00, -7.307767E-01, +5.876704E-01, +2.150661E+00, +1.470055E-01, +1.796000E+00, +1.000000E+00
254cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
255cov
256+2.500000E+01, -3.069880E+00, +4.509227E+00, +4.344065E+00, -3.645299E-01, +5.465260E+00, +2.744630E+00
257-3.069880E+00, +1.000000E+00, -6.272851E-01, -3.187653E-01, +2.896488E-01, -2.138192E+00, -7.307767E-01
258+4.509227E+00, -6.272851E-01, +1.000000E+00, +8.493386E-01, -2.867150E-01, +1.823854E+00, +5.876704E-01
259+4.344065E+00, -3.187653E-01, +8.493386E-01, +1.600000E+01, +2.921432E+00, +4.903669E+00, +2.150661E+00
260-3.645299E-01, +2.896488E-01, -2.867150E-01, +2.921432E+00, +1.000000E+00, -6.943769E-01, +1.470055E-01
261+5.465260E+00, -2.138192E+00, +1.823854E+00, +4.903669E+00, -6.943769E-01, +9.000000E+00, +1.796000E+00
262+2.744630E+00, -7.307767E-01, +5.876704E-01, +2.150661E+00, +1.470055E-01, +1.796000E+00, +1.000000E+00
263getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
264+1.000000E+00, -6.139759E-01, +9.018455E-01, +2.172033E-01, -7.290598E-02, +3.643507E-01, +5.489259E-01
265-6.139759E-01, +1.000000E+00, -6.272851E-01, -7.969132E-02, +2.896488E-01, -7.127306E-01, -7.307767E-01
266+9.018455E-01, -6.272851E-01, +1.000000E+00, +2.123346E-01, -2.867150E-01, +6.079513E-01, +5.876704E-01
267+2.172033E-01, -7.969132E-02, +2.123346E-01, +1.000000E+00, +7.303581E-01, +4.086391E-01, +5.376651E-01
268-7.290598E-02, +2.896488E-01, -2.867150E-01, +7.303581E-01, +1.000000E+00, -2.314590E-01, +1.470055E-01
269+3.643507E-01, -7.127306E-01, +6.079513E-01, +4.086391E-01, -2.314590E-01, +1.000000E+00, +5.986668E-01
270+5.489259E-01, -7.307767E-01, +5.876704E-01, +5.376651E-01, +1.470055E-01, +5.986668E-01, +1.000000E+00
271
272
273ndim = getUnifRand(1, 7)
274ndim
275+6
276std = getUnifRand(1, 10, ndim)
277std
278+1.000000E+01, +5.000000E+00, +5.000000E+00, +7.000000E+00, +5.000000E+00, +6.000000E+00
279call setResized(cov, [ndim, ndim])
280cor = getCovRand(1., ndim)
281cor
282+1.000000E+00, +5.646743E-01, +7.850605E-01, -2.208656E-02, +9.809051E-02, +1.253671E-01
283+5.646743E-01, +1.000000E+00, -5.891359E-02, -4.176201E-01, -5.934682E-01, +4.898605E-01
284+7.850605E-01, -5.891359E-02, +1.000000E+00, +3.819137E-01, +5.283641E-01, -1.520058E-01
285-2.208656E-02, -4.176201E-01, +3.819137E-01, +1.000000E+00, +1.547616E-01, +1.914386E-01
286+9.809051E-02, -5.934682E-01, +5.283641E-01, +1.547616E-01, +1.000000E+00, -3.745260E-01
287+1.253671E-01, +4.898605E-01, -1.520058E-01, +1.914386E-01, -3.745260E-01, +1.000000E+00
288cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
289cov
290+1.000000E+02, +2.823372E+01, +3.925303E+01, -1.546059E+00, +4.904525E+00, +7.522026E+00
291+2.823372E+01, +2.500000E+01, -1.472840E+00, -1.461670E+01, -1.483671E+01, +1.469582E+01
292+3.925303E+01, -1.472840E+00, +2.500000E+01, +1.336698E+01, +1.320910E+01, -4.560174E+00
293-1.546059E+00, -1.461670E+01, +1.336698E+01, +4.900000E+01, +5.416655E+00, +8.040422E+00
294+4.904525E+00, -1.483671E+01, +1.320910E+01, +5.416655E+00, +2.500000E+01, -1.123578E+01
295+7.522026E+00, +1.469582E+01, -4.560174E+00, +8.040422E+00, -1.123578E+01, +3.600000E+01
296cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
297cov
298+1.000000E+02, +2.823372E+01, +3.925303E+01, -1.546059E+00, +4.904525E+00, +7.522026E+00
299+2.823372E+01, +2.500000E+01, -1.472840E+00, -1.461670E+01, -1.483671E+01, +1.469582E+01
300+3.925303E+01, -1.472840E+00, +2.500000E+01, +1.336698E+01, +1.320910E+01, -4.560174E+00
301-1.546059E+00, -1.461670E+01, +1.336698E+01, +4.900000E+01, +5.416655E+00, +8.040422E+00
302+4.904525E+00, -1.483671E+01, +1.320910E+01, +5.416655E+00, +2.500000E+01, -1.123578E+01
303+7.522026E+00, +1.469582E+01, -4.560174E+00, +8.040422E+00, -1.123578E+01, +3.600000E+01
304getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
305+1.000000E+00, +5.646743E-01, +7.850606E-01, -2.208656E-02, +9.809051E-02, +1.253671E-01
306+5.646743E-01, +1.000000E+00, -5.891359E-02, -4.176201E-01, -5.934683E-01, +4.898606E-01
307+7.850606E-01, -5.891359E-02, +1.000000E+00, +3.819138E-01, +5.283641E-01, -1.520058E-01
308-2.208656E-02, -4.176201E-01, +3.819138E-01, +1.000000E+00, +1.547616E-01, +1.914386E-01
309+9.809051E-02, -5.934683E-01, +5.283641E-01, +1.547616E-01, +1.000000E+00, -3.745261E-01
310+1.253671E-01, +4.898606E-01, -1.520058E-01, +1.914386E-01, -3.745261E-01, +1.000000E+00
311
312
313ndim = getUnifRand(1, 7)
314ndim
315+2
316std = getUnifRand(1, 10, ndim)
317std
318+1.000000E+01, +4.000000E+00
319call setResized(cov, [ndim, ndim])
320cor = getCovRand(1., ndim)
321cor
322+1.000000E+00, -5.784964E-01
323-5.784964E-01, +1.000000E+00
324cov = getCov(cor, uppDia, std) ! convert upper correlation matrix to full covariance matrix.
325cov
326+1.000000E+02, -2.313986E+01
327-2.313986E+01, +1.600000E+01
328cov = getCov(cor, lowDia, std) ! convert upper correlation matrix to full covariance matrix.
329cov
330+1.000000E+02, -2.313986E+01
331-2.313986E+01, +1.600000E+01
332getCor(getCov(cor, lowDia, std), lowDia) ! reconstruct the original correlation matrix.
333+1.000000E+00, -5.784964E-01
334-5.784964E-01, +1.000000E+00
335
336
337!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
338!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
339!Compute the covariance matrix of a 2-D sample.
340!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
341!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
342
343ndim = 2; nsam = 10
344sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
345sample
346+3.00000000, +9.00000000, +8.00000000, +9.00000000, +14.0000000, +4.00000000, +6.00000000, +2.00000000, +8.00000000, +5.00000000
347+15.0000000, +18.0000000, +6.00000000, +6.00000000, +15.0000000, +6.00000000, +13.0000000, +9.00000000, +19.0000000, +10.0000000
348
349'Compute the sample covariance along the second dimension.'
350
351cov = getCov(sample, dim = 2_IK)
352cov
353+11.3600006, +4.54000044
354+4.54000044, +22.4100018
355
356'Compute the sample covariance along the first dimension.'
357
358cov = getCov(transpose(sample), dim = 1_IK)
359cov
360+11.3599997, +4.54000044
361+4.54000044, +22.4100018
362
363'Compute the full sample covariance for a pair of time series.'
364
365cov = getCov(sample(1,:), sample(2,:))
366cov
367+11.3600006, +4.54000044
368+4.54000044, +22.4100018
369
370ndim = 2; nsam = 10
371sample = reshape(cmplx(getUnifRand(1, 20, ndim * nsam), -getUnifRand(1, 20, ndim * nsam), TKG), shape = [ndim, nsam], order = [2, 1])
372sample
373(+12.0000000, -9.00000000), (+19.0000000, -20.0000000), (+17.0000000, -14.0000000), (+3.00000000, -19.0000000), (+19.0000000, -16.0000000), (+3.00000000, -14.0000000), (+5.00000000, -19.0000000), (+5.00000000, -4.00000000), (+13.0000000, -17.0000000), (+12.0000000, -1.00000000)
374(+8.00000000, -18.0000000), (+20.0000000, -10.0000000), (+20.0000000, -19.0000000), (+1.00000000, -14.0000000), (+4.00000000, -15.0000000), (+1.00000000, -10.0000000), (+2.00000000, -6.00000000), (+5.00000000, -2.00000000), (+11.0000000, -7.00000000), (+19.0000000, -18.0000000)
375
376'Compute the sample covariance along the second dimension.'
377
378cov = getCov(sample, dim = 2_IK)
379cov
380(+75.7699966, +0.00000000), (+27.1499996, +24.5100002)
381(+27.1499996, -24.5100002), (+86.7800064, +0.00000000)
382
383'Compute the sample covariance along the first dimension.'
384
385cov = getCov(transpose(sample), dim = 1_IK)
386cov
387(+75.7699966, +0.00000000), (+27.1499996, +24.5100002)
388(+27.1499996, -24.5100002), (+86.7800064, +0.00000000)
389
390'Compute the full sample covariance for a pair of time series.'
391
392cov = getCov(sample(1,:), sample(2,:))
393cov
394(+75.7699966, +0.00000000), (+27.1499996, +24.5100002)
395(+27.1499996, -24.5100002), (+86.7800064, +0.00000000)
396
397
398!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
399!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
400!Compute the biased covariance matrix of a weighted 2-D sample.
401!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
402!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
403
404ndim = 2; nsam = 10
405sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
406sample
407+19.0000000, +18.0000000, +10.0000000, +3.00000000, +1.00000000, +6.00000000, +3.00000000, +1.00000000, +20.0000000, +10.0000000
408+18.0000000, +3.00000000, +6.00000000, +13.0000000, +4.00000000, +2.00000000, +12.0000000, +14.0000000, +10.0000000, +15.0000000
409call setResized(mean, ndim)
410iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
411iweight
412+9, +3, +5, +8, +7, +7, +10, +6, +7, +7
413call setMean(mean, sample, 2_IK, iweight, iweisum)
414mean
415+8.60869598, +10.5217390
416iweisum
417+69
418rweight = iweight ! or real-valued weights.
419iweight
420+9, +3, +5, +8, +7, +7, +10, +6, +7, +7
421call setMean(mean, sample, 2_IK, rweight, rweisum)
422mean
423+8.60869598, +10.5217390
424rweisum
425+69.0000000
426
427!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
428!Compute the covariance matrix with integer weights.
429!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
430
431cov = getCov(sample, 2_IK, iweight)
432cov
433+51.2236824, +8.81284714
434+8.81284714, +27.0611210
435
436'Compute the sample covariance along the first dimension.'
437
438cov = getCov(transpose(sample), 1_IK, iweight)
439cov
440+51.2236938, +8.81284714
441+8.81284714, +27.0611210
442
443'Compute the full sample covariance for a pair of time series.'
444
445cov = getCov(sample(1,:), sample(2,:), weight = iweight)
446cov
447+51.2236824, +8.81284714
448+8.81284714, +27.0611210
449
450
451!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
452!Compute the covariance matrix with real weights.
453!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
454
455cov = getCov(sample, 2_IK, rweight)
456cov
457+51.2236824, +8.81284714
458+8.81284714, +27.0611210
459
460'Compute the sample covariance along the first dimension.'
461
462cov = getCov(transpose(sample), 1_IK, rweight)
463cov
464+51.2236938, +8.81284714
465+8.81284714, +27.0611210
466
467'Compute the full sample covariance for a pair of time series.'
468
469cov = getCov(sample(1,:), sample(2,:), weight = rweight)
470cov
471+51.2236824, +8.81284714
472+8.81284714, +27.0611210
473
474ndim = 2; nsam = 10
475sample = reshape(cmplx(getUnifRand(1, 20, ndim * nsam), -getUnifRand(1, 20, ndim * nsam), TKG), shape = [ndim, nsam], order = [2, 1])
476sample
477(+15.0000000, -1.00000000), (+18.0000000, -16.0000000), (+1.00000000, -16.0000000), (+12.0000000, -10.0000000), (+18.0000000, -9.00000000), (+12.0000000, -19.0000000), (+3.00000000, -2.00000000), (+20.0000000, -17.0000000), (+15.0000000, -1.00000000), (+17.0000000, -14.0000000)
478(+4.00000000, -8.00000000), (+6.00000000, -10.0000000), (+2.00000000, -5.00000000), (+2.00000000, -18.0000000), (+16.0000000, -14.0000000), (+17.0000000, -2.00000000), (+12.0000000, -2.00000000), (+18.0000000, -2.00000000), (+5.00000000, -14.0000000), (+12.0000000, -1.00000000)
479call setResized(mean, ndim)
480iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
481iweight
482+6, +3, +4, +1, +7, +5, +2, +4, +5, +8
483call setMean(mean, sample, 2_IK, iweight, iweisum)
484mean
485(+14.2888889, -10.5555553), (+10.3555555, -6.97777796)
486iweisum
487+45
488rweight = iweight ! or real-valued weights.
489iweight
490+6, +3, +4, +1, +7, +5, +2, +4, +5, +8
491call setMean(mean, sample, 2_IK, rweight, rweisum)
492mean
493(+14.2888889, -10.5555553), (+10.3555555, -6.97777796)
494rweisum
495+45.0000000
496
497!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
498!Compute the covariance matrix with integer weights.
499!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
500
501cov = getCov(sample, 2_IK, iweight)
502cov
503(+73.7412491, +0.00000000), (-6.13481474, -12.3071613)
504(-6.13481474, +12.3071613), (+62.6064262, +0.00000000)
505
506'Compute the sample covariance along the first dimension.'
507
508cov = getCov(transpose(sample), 1_IK, iweight)
509cov
510(+73.7412338, +0.00000000), (-6.13481474, -12.3071613)
511(-6.13481474, +12.3071613), (+62.6064262, +0.00000000)
512
513'Compute the full sample covariance for a pair of time series.'
514
515cov = getCov(sample(1,:), sample(2,:), weight = iweight)
516cov
517(+73.7412491, +0.00000000), (-6.13481474, -12.3071613)
518(-6.13481474, +12.3071613), (+62.6064262, +0.00000000)
519
520
521!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
522!Compute the covariance matrix with real weights.
523!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
524
525cov = getCov(sample, 2_IK, rweight)
526cov
527(+73.7412491, +0.00000000), (-6.13481474, -12.3071613)
528(-6.13481474, +12.3071613), (+62.6064262, +0.00000000)
529
530'Compute the sample covariance along the first dimension.'
531
532cov = getCov(transpose(sample), 1_IK, rweight)
533cov
534(+73.7412338, +0.00000000), (-6.13481474, -12.3071613)
535(-6.13481474, +12.3071613), (+62.6064262, +0.00000000)
536
537'Compute the full sample covariance for a pair of time series.'
538
539cov = getCov(sample(1,:), sample(2,:), weight = rweight)
540cov
541(+73.7412491, +0.00000000), (-6.13481474, -12.3071613)
542(-6.13481474, +12.3071613), (+62.6064262, +0.00000000)
543
544
Test:
test_pm_sampleCov
Bug:

Status: Unresolved
Source: Intel Classic Fortran Compiler ifort version 2021.8
Description: Intel ifort (possibly only with heap memory allocation) yields a segfault error in OpenMP parallel code at pm_sampleCov@routines.inc.F90:238.
The root cause of the segmentation fault was determined to be due to the use of do concurrent construct in the implementation.

Remedy (as of ParaMonte Library version 2.0.0): For now, all do concurrent statements are converted back to normal do loops.


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

Definition at line 344 of file pm_sampleCov.F90.


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