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

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

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

Note
The sample correlation matrix can be also readily computed from the sample covariance matrix as in the following example,
cor(:,:) = getCor(getCov(sample, dim))
cor(:,:) = getCor(getCov(sample, dim, weight)) ! weighted sample.
!
Parameters
[in]cov: The input positive semi-definite square matrix of shape (1:ndim, 1:ndim) of the same type and kind as the output cor, representing the correlation 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 and 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]x: The input contiguous vector of shape (nsam) of the same type and kind as the output cor, 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 cor, 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 cor, containing the sample comprised of nsam observations each with ndim attributes.
If sample is a matrix, then the direction along which the correlation 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 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.)
[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.
(optional. default = getFilled(1, nsam). It can be present only if the input arguments cov, subsetv, and stdinv are missing.)
Returns
cor : 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),
representing its Pearson correlation matrix.
When the input arguments x and y are present, then ndim == 2 holds by definition and the output cor is a scalar of value \(r_{xy}\),

\begin{equation} \ms{cor} = \begin{bmatrix} 1 && r_{xy} \\ r_{yx} && 1 ~. \end{bmatrix} \end{equation}


Possible calling interfaces

use pm_sampleCor, only: getCor
! correlation matrix to correlation matrix.
cor(1:ndim, 1:ndim) = getCor(cov(1:ndim, 1:ndim), subsetv) ! subsetv = uppDia, lowDia
cor(1:ndim, 1:ndim) = getCor(cov(1:ndim, 1:ndim), subsetv, stdinv(1:ndim)) ! subsetv = upp, low
! XY time series correlation matrix.
cor = getCor(x(1:nsam), y(1:nsam)) ! correlation coefficient.
cor = getCor(x(1:nsam), y(1:nsam), weight(1:nsam)) ! correlation coefficient.
! sample correlation matrix.
cor(:,:) = getCor(sample(:,:), dim) ! full correlation matrix.
cor(:,:) = getCor(sample(:,:), dim, weight(1:nsam)) ! upper/lower correlation matrix.
!
Generate and return the (Pearson) correlation coefficient or matrix of a pair of (weighted) time seri...
This module contains classes and procedures for computing properties related to the correlation matri...
Warning
This generic interface is merely a flexible wrapper around the generic subroutine interface setCor.
As such, all conditions and warnings associated with setCor 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.
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: getCor
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_sampleCov, only: getCov
14 use pm_sampleMean, only: getMean
15 use pm_sampleMean, only: setMean
16 use pm_sampleShift, only: getShifted
17 use pm_arraySpace, only: getLinSpace
18 use pm_arrayResize, only: setResized
19 use pm_arrayFill, only: getFilled
20 use pm_io, only: getFormat
21
22 implicit none
23
24 type(display_type) :: disp
25 integer(IK) :: itry, ntry = 10
26 character(:), allocatable :: format
27 disp = display_type(file = "main.out.F90")
28
29 call disp%skip()
30 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%show("!Compute the correlation matrix from covariance matrix.")
32 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
33 call disp%skip()
34
35 block
36 use pm_kind, only: TKG => RKS ! All other real types are also supported.
37 real(TKG), allocatable :: cov(:,:), cor(:,:), std(:)
38 integer(IK) :: ndim
39 do itry = 1, 10
40 call disp%skip()
41 call disp%show("ndim = getUnifRand(0, 7)")
42 ndim = getUnifRand(0, 7)
43 call disp%show("ndim")
44 call disp%show( ndim )
45 call disp%show("std = getUnifRand(1, ndim, ndim)")
46 std = getUnifRand(1, ndim, ndim)
47 call disp%show("std")
48 call disp%show( std )
49 call disp%show("cov = getCovRand(1._TKG, std)")
50 cov = getCovRand(1._TKG, std)
51 call disp%show("cov")
52 call disp%show( cov )
53 call disp%show("cor = getCor(cov, uppDia)")
54 cor = getCor(cov, uppDia)
55 call disp%show("cor")
56 call disp%show( cor )
57 call disp%show("cor = getCor(cov, upp, stdinv = 1 / std)")
58 cor = getCor(cov, upp, stdinv = 1 / std)
59 call disp%show("cor")
60 call disp%show( cor )
61 call disp%show("cor = getCor(cov, lowDia)")
62 cor = getCor(cov, lowDia)
63 call disp%show("cor")
64 call disp%show( cor )
65 call disp%show("cor = getCor(cov, low, stdinv = 1 / std)")
66 cor = getCor(cov, low, stdinv = 1 / std)
67 call disp%show("cor")
68 call disp%show( cor )
69 call disp%show("getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.")
70 call disp%show( getCov(getCor(cov, lowDia), lowDia, std) , format = format )
71 call disp%show("getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.")
72 call disp%show( getCov(getCor(cov, uppDia), uppDia, std) , format = format )
73 call disp%skip()
74 end do
75 end block
76
77 call disp%skip()
78 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
79 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
80 call disp%show("!Compute the Pearson correlation matrix for a pair of time series.")
81 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
82 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
83 call disp%skip()
84
85 block
86 use pm_kind, only: TKG => RKS ! All other real types are also supported.
87 integer(IK) :: ndim, nsam, dim
88 real(TKG), allocatable :: sample(:,:), cor(:,:), mean(:)
89 format = getFormat(mold = [0._TKG], ed = SK_"es", signed = .true._LK)
90 call disp%show("ndim = 2; nsam = 10; dim = 2")
91 ndim = 2; nsam = 10; dim = 2
92 call disp%show("sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
93 sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
94 call disp%show("sample")
95 call disp%show( sample , format = format )
96 call disp%show("mean = getMean(sample, dim)")
97 mean = getMean(sample, dim)
98 call disp%show("mean")
99 call disp%show( mean , format = format )
100 call disp%show("cor = getCor(sample, dim) ! same result as above.")
101 cor = getCor(sample, dim)
102 call disp%show("cor")
103 call disp%show( cor , format = format )
104 call disp%skip()
105 call disp%show("Compute the sample correlation along the first dimension.", deliml = SK_'''')
106 call disp%skip()
107 call disp%show("dim = 1")
108 dim = 1
109 call disp%show("cor = getCor(transpose(sample), dim) ! same result as above.")
110 cor = getCor(transpose(sample), dim)
111 call disp%show("cor")
112 call disp%show( cor , format = format )
113 call disp%skip()
114 call disp%show("Compute the full sample correlation for a pair of time series.", deliml = SK_'''')
115 call disp%skip()
116 call disp%show("cor(1,1) = getCor(sample(1,:), sample(2,:)) ! same result as above.")
117 cor(1,1) = getCor(sample(1,:), sample(2,:))
118 call disp%show("cor(1,1)")
119 call disp%show( cor(1,1) , format = format )
120 call disp%skip()
121 end block
122
123 call disp%skip()
124 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
125 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
126 call disp%show("!Compute the Pearson correlation matrix for a weighted pair of time series.")
127 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
128 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
129 call disp%skip()
130
131 block
132 use pm_arrayVerbose, only: getVerbose
133 use pm_kind, only: TKG => RKS ! All other real types are also supported.
134 integer(IK), allocatable :: iweight(:)
135 real(TKG), allocatable :: rweight(:)
136 real(TKG) :: rweisum
137 integer(IK) :: iweisum
138 integer(IK) :: ndim, nsam, dim
139 real(TKG), allocatable :: sample(:,:), cor(:,:), mean(:)
140 format = getFormat(mold = [0._TKG], ed = SK_"es", signed = .true._LK)
141 call disp%show("ndim = 2; nsam = 10; dim = 2")
142 ndim = 2; nsam = 10; dim = 2
143 call disp%show("sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
144 sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
145 call disp%show("sample")
146 call disp%show( sample , format = format )
147 call disp%show("call setResized(mean, ndim)")
148 call setResized(mean, ndim)
149 call disp%show("iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.")
150 iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
151 call disp%show("iweight")
152 call disp%show( iweight )
153 call disp%show("call setMean(mean, sample, dim, iweight, iweisum)")
154 call setMean(mean, sample, dim, iweight, iweisum)
155 call disp%show("mean")
156 call disp%show( mean )
157 call disp%show("iweisum")
158 call disp%show( iweisum )
159 call disp%show("rweight = iweight ! or real-valued weights.")
160 rweight = iweight ! or real-valued weights.
161 call disp%show("iweight")
162 call disp%show( iweight )
163 call disp%show("call setMean(mean, sample, dim, rweight, rweisum)")
164 call setMean(mean, sample, dim, rweight, rweisum)
165 call disp%show("mean")
166 call disp%show( mean , format = format )
167 call disp%show("rweisum")
168 call disp%show( rweisum , format = format )
169
170 call disp%skip()
171 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
172 call disp%show("!Compute the correlation matrix with integer weights.")
173 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
174 call disp%skip()
175
176 call disp%show("cor = getCor(sample, dim, iweight) ! same result as above.")
177 cor = getCor(sample, dim, iweight)
178 call disp%show("cor")
179 call disp%show( cor , format = format )
180 call disp%skip()
181 call disp%show("Compute the sample correlation along the first dimension.", deliml = SK_'''')
182 call disp%skip()
183 call disp%show("dim = 1")
184 dim = 1
185 call disp%show("cor = getCor(transpose(sample), dim, iweight) ! same result as above.")
186 cor = getCor(transpose(sample), dim, iweight)
187 call disp%show("cor")
188 call disp%show( cor , format = format )
189 call disp%skip()
190 call disp%show("Compute the full sample correlation for a pair of time series.", deliml = SK_'''')
191 call disp%skip()
192 call disp%show("cor(1,1) = getCor(sample(1,:), sample(2,:), iweight) ! same result as above.")
193 cor(1,1) = getCor(sample(1,:), sample(2,:), iweight)
194 call disp%show("cor(1,1)")
195 call disp%show( cor(1,1) , format = format )
196 call disp%skip()
197
198 call disp%skip()
199 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
200 call disp%show("!Compute the correlation matrix with real weights.")
201 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
202 call disp%skip()
203
204 call disp%show("dim = 2")
205 dim = 2
206 call disp%show("cor = getCor(sample, dim, rweight) ! same result as above.")
207 cor = getCor(sample, dim, rweight)
208 call disp%show("cor")
209 call disp%show( cor , format = format )
210 call disp%skip()
211 call disp%show("Compute the sample correlation along the first dimension.", deliml = SK_'''')
212 call disp%skip()
213 call disp%show("dim = 1")
214 dim = 1
215 call disp%show("cor = getCor(transpose(sample), dim, rweight) ! same result as above.")
216 cor = getCor(transpose(sample), dim, rweight)
217 call disp%show("cor")
218 call disp%show( cor , format = format )
219 call disp%skip()
220 call disp%show("Compute the full sample correlation for a pair of time series.", deliml = SK_'''')
221 call disp%skip()
222 call disp%show("cor(1,1) = getCor(sample(1,:), sample(2,:), rweight) ! same result as above.")
223 cor(1,1) = getCor(sample(1,:), sample(2,:), rweight)
224 call disp%show("cor(1,1)")
225 call disp%show( cor(1,1) , format = format )
226 call disp%skip()
227 end block
228
229end 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 (optionally unbiased) covariance matrix of a pair of (potentially weighted) t...
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 properties related to the covariance ma...
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 covariance matrix.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7ndim = getUnifRand(0, 7)
8ndim
9+7
10std = getUnifRand(1, ndim, ndim)
11std
12+5.00000000, +1.00000000, +4.00000000, +6.00000000, +6.00000000, +7.00000000, +1.00000000
13cov = getCovRand(1._TKG, std)
14cov
15+25.0000000, +3.42286325, +14.8290634, +3.94746423, -20.8633156, -23.9637852, -2.16289282
16+3.42286325, +0.999999881, +1.67499459, -2.29277420, -4.99174309, -5.11780691, -0.634693623
17+14.8290634, +1.67499459, +15.9999981, +0.758326292, -15.0802774, -11.3248491, +0.506526232E-1
18+3.94746423, -2.29277420, +0.758326292, +36.0000000, +0.905534148, +8.03260994, -1.10761917
19-20.8633156, -4.99174309, -15.0802774, +0.905534148, +36.0000000, +20.2001801, +2.81178069
20-23.9637852, -5.11780691, -11.3248491, +8.03260994, +20.2001801, +49.0000000, +4.73817301
21-2.16289282, -0.634693623, +0.506526232E-1, -1.10761917, +2.81178069, +4.73817301, +1.00000000
22cor = getCor(cov, uppDia)
23cor
24+1.00000000, +0.684572756, +0.741453290, +0.131582141, -0.695443869, -0.684679627, -0.432578564
25+0.684572756, +1.00000000, +0.418748736, -0.382129073, -0.831957281, -0.731115401, -0.634693682
26+0.741453290, +0.418748736, +1.00000000, +0.315969326E-1, -0.628344953, -0.404458970, +0.126631577E-1
27+0.131582141, -0.382129073, +0.315969326E-1, +1.00000000, +0.251537282E-1, +0.191252634, -0.184603199
28-0.695443869, -0.831957281, -0.628344953, +0.251537282E-1, +1.00000000, +0.480956703, +0.468630135
29-0.684679627, -0.731115401, -0.404458970, +0.191252634, +0.480956703, +1.00000000, +0.676881909
30-0.432578564, -0.634693682, +0.126631577E-1, -0.184603199, +0.468630135, +0.676881909, +1.00000000
31cor = getCor(cov, upp, stdinv = 1 / std)
32cor
33+1.00000000, +0.684572637, +0.741453171, +0.131582141, -0.695443869, -0.684679627, -0.432578564
34+0.684572637, +1.00000000, +0.418748647, -0.382129043, -0.831957221, -0.731115282, -0.634693623
35+0.741453171, +0.418748647, +1.00000000, +0.315969288E-1, -0.628344893, -0.404458910, +0.126631558E-1
36+0.131582141, -0.382129043, +0.315969288E-1, +1.00000000, +0.251537282E-1, +0.191252634, -0.184603199
37-0.695443869, -0.831957221, -0.628344893, +0.251537282E-1, +1.00000000, +0.480956703, +0.468630135
38-0.684679627, -0.731115282, -0.404458910, +0.191252634, +0.480956703, +1.00000000, +0.676881909
39-0.432578564, -0.634693623, +0.126631558E-1, -0.184603199, +0.468630135, +0.676881909, +1.00000000
40cor = getCor(cov, lowDia)
41cor
42+1.00000000, +0.684572756, +0.741453290, +0.131582141, -0.695443869, -0.684679627, -0.432578564
43+0.684572756, +1.00000000, +0.418748736, -0.382129073, -0.831957281, -0.731115401, -0.634693682
44+0.741453290, +0.418748736, +1.00000000, +0.315969326E-1, -0.628344953, -0.404458970, +0.126631577E-1
45+0.131582141, -0.382129073, +0.315969326E-1, +1.00000000, +0.251537282E-1, +0.191252634, -0.184603199
46-0.695443869, -0.831957281, -0.628344953, +0.251537282E-1, +1.00000000, +0.480956703, +0.468630135
47-0.684679627, -0.731115401, -0.404458970, +0.191252634, +0.480956703, +1.00000000, +0.676881909
48-0.432578564, -0.634693682, +0.126631577E-1, -0.184603199, +0.468630135, +0.676881909, +1.00000000
49cor = getCor(cov, low, stdinv = 1 / std)
50cor
51+1.00000000, +0.684572637, +0.741453171, +0.131582141, -0.695443869, -0.684679627, -0.432578564
52+0.684572637, +1.00000000, +0.418748647, -0.382129043, -0.831957221, -0.731115282, -0.634693623
53+0.741453171, +0.418748647, +1.00000000, +0.315969288E-1, -0.628344893, -0.404458910, +0.126631558E-1
54+0.131582141, -0.382129043, +0.315969288E-1, +1.00000000, +0.251537282E-1, +0.191252634, -0.184603199
55-0.695443869, -0.831957221, -0.628344893, +0.251537282E-1, +1.00000000, +0.480956703, +0.468630135
56-0.684679627, -0.731115282, -0.404458910, +0.191252634, +0.480956703, +1.00000000, +0.676881909
57-0.432578564, -0.634693623, +0.126631558E-1, -0.184603199, +0.468630135, +0.676881909, +1.00000000
58getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
59+25.0000000, +3.42286372, +14.8290653, +3.94746423, -20.8633156, -23.9637871, -2.16289282
60+3.42286372, +1.00000000, +1.67499495, -2.29277444, -4.99174356, -5.11780787, -0.634693682
61+14.8290653, +1.67499495, +16.0000000, +0.758326411, -15.0802784, -11.3248510, +0.506526306E-1
62+3.94746423, -2.29277444, +0.758326411, +36.0000000, +0.905534208, +8.03261089, -1.10761917
63-20.8633156, -4.99174356, -15.0802784, +0.905534208, +36.0000000, +20.2001820, +2.81178093
64-23.9637871, -5.11780787, -11.3248510, +8.03261089, +20.2001820, +49.0000000, +4.73817348
65-2.16289282, -0.634693682, +0.506526306E-1, -1.10761917, +2.81178093, +4.73817348, +1.00000000
66getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
67+25.0000000, +3.42286372, +14.8290653, +3.94746423, -20.8633156, -23.9637871, -2.16289282
68+3.42286372, +1.00000000, +1.67499495, -2.29277444, -4.99174356, -5.11780787, -0.634693682
69+14.8290653, +1.67499495, +16.0000000, +0.758326411, -15.0802784, -11.3248510, +0.506526306E-1
70+3.94746423, -2.29277444, +0.758326411, +36.0000000, +0.905534208, +8.03261089, -1.10761917
71-20.8633156, -4.99174356, -15.0802784, +0.905534208, +36.0000000, +20.2001820, +2.81178093
72-23.9637871, -5.11780787, -11.3248510, +8.03261089, +20.2001820, +49.0000000, +4.73817348
73-2.16289282, -0.634693682, +0.506526306E-1, -1.10761917, +2.81178093, +4.73817348, +1.00000000
74
75
76ndim = getUnifRand(0, 7)
77ndim
78+2
79std = getUnifRand(1, ndim, ndim)
80std
81+2.00000000, +2.00000000
82cov = getCovRand(1._TKG, std)
83cov
84+4.00000000, -3.41578817
85-3.41578817, +4.00000048
86cor = getCor(cov, uppDia)
87cor
88+1.00000000, -0.853947043
89-0.853947043, +1.00000000
90cor = getCor(cov, upp, stdinv = 1 / std)
91cor
92+1.00000000, -0.853947043
93-0.853947043, +1.00000000
94cor = getCor(cov, lowDia)
95cor
96+1.00000000, -0.853947043
97-0.853947043, +1.00000000
98cor = getCor(cov, low, stdinv = 1 / std)
99cor
100+1.00000000, -0.853947043
101-0.853947043, +1.00000000
102getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
103+4.00000000, -3.41578817
104-3.41578817, +4.00000000
105getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
106+4.00000000, -3.41578817
107-3.41578817, +4.00000000
108
109
110ndim = getUnifRand(0, 7)
111ndim
112+2
113std = getUnifRand(1, ndim, ndim)
114std
115+2.00000000, +2.00000000
116cov = getCovRand(1._TKG, std)
117cov
118+4.00000000, -1.81410193
119-1.81410193, +3.99999976
120cor = getCor(cov, uppDia)
121cor
122+1.00000000, -0.453525543
123-0.453525543, +1.00000000
124cor = getCor(cov, upp, stdinv = 1 / std)
125cor
126+1.00000000, -0.453525484
127-0.453525484, +1.00000000
128cor = getCor(cov, lowDia)
129cor
130+1.00000000, -0.453525543
131-0.453525543, +1.00000000
132cor = getCor(cov, low, stdinv = 1 / std)
133cor
134+1.00000000, -0.453525484
135-0.453525484, +1.00000000
136getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
137+4.00000000, -1.81410217
138-1.81410217, +4.00000000
139getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
140+4.00000000, -1.81410217
141-1.81410217, +4.00000000
142
143
144ndim = getUnifRand(0, 7)
145ndim
146+7
147std = getUnifRand(1, ndim, ndim)
148std
149+1.00000000, +5.00000000, +3.00000000, +6.00000000, +6.00000000, +1.00000000, +1.00000000
150cov = getCovRand(1._TKG, std)
151cov
152+1.00000000, -4.86804104, +1.89701664, -1.37734020, +0.106371790, +0.503054142, -0.529540539
153-4.86804104, +25.0000019, -11.7432575, +11.8509970, +4.45921135, -2.91629577, +2.58521318
154+1.89701664, -11.7432575, +9.00000000, -11.7150049, -7.83632755, +2.07326889, -1.39819896
155-1.37734020, +11.8509970, -11.7150049, +36.0000076, +9.68124199, -1.39010274, +1.14766455
156+0.106371790, +4.45921135, -7.83632755, +9.68124199, +35.9999924, -2.36941385, -2.12545347
157+0.503054142, -2.91629577, +2.07326889, -1.39010274, -2.36941385, +1.00000000, -0.619132280
158-0.529540539, +2.58521318, -1.39819896, +1.14766455, -2.12545347, -0.619132280, +0.999999940
159cor = getCor(cov, uppDia)
160cor
161+1.00000000, -0.973608196, +0.632338881, -0.229556680, +0.177286342E-1, +0.503054142, -0.529540598
162-0.973608196, +1.00000000, -0.782883883, +0.395033211, +0.148640409, -0.583259165, +0.517042696
163+0.632338881, -0.782883883, +1.00000000, -0.650833547, -0.435351580, +0.691089630, -0.466066390
164-0.229556680, +0.395033211, -0.650833547, +1.00000000, +0.268923402, -0.231683776, +0.191277429
165+0.177286342E-1, +0.148640409, -0.435351580, +0.268923402, +1.00000000, -0.394902349, -0.354242325
166+0.503054142, -0.583259165, +0.691089630, -0.231683776, -0.394902349, +1.00000000, -0.619132340
167-0.529540598, +0.517042696, -0.466066390, +0.191277429, -0.354242325, -0.619132340, +1.00000000
168cor = getCor(cov, upp, stdinv = 1 / std)
169cor
170+1.00000000, -0.973608196, +0.632338881, -0.229556710, +0.177286323E-1, +0.503054142, -0.529540539
171-0.973608196, +1.00000000, -0.782883883, +0.395033240, +0.148640379, -0.583259165, +0.517042637
172+0.632338881, -0.782883883, +1.00000000, -0.650833666, -0.435351551, +0.691089630, -0.466066331
173-0.229556710, +0.395033240, -0.650833666, +1.00000000, +0.268923402, -0.231683791, +0.191277429
174+0.177286323E-1, +0.148640379, -0.435351551, +0.268923402, +1.00000000, -0.394902319, -0.354242265
175+0.503054142, -0.583259165, +0.691089630, -0.231683791, -0.394902319, +1.00000000, -0.619132280
176-0.529540539, +0.517042637, -0.466066331, +0.191277429, -0.354242265, -0.619132280, +1.00000000
177cor = getCor(cov, lowDia)
178cor
179+1.00000000, -0.973608196, +0.632338881, -0.229556680, +0.177286342E-1, +0.503054142, -0.529540598
180-0.973608196, +1.00000000, -0.782883883, +0.395033211, +0.148640409, -0.583259165, +0.517042696
181+0.632338881, -0.782883883, +1.00000000, -0.650833547, -0.435351580, +0.691089630, -0.466066390
182-0.229556680, +0.395033211, -0.650833547, +1.00000000, +0.268923402, -0.231683776, +0.191277429
183+0.177286342E-1, +0.148640409, -0.435351580, +0.268923402, +1.00000000, -0.394902349, -0.354242325
184+0.503054142, -0.583259165, +0.691089630, -0.231683776, -0.394902349, +1.00000000, -0.619132340
185-0.529540598, +0.517042696, -0.466066390, +0.191277429, -0.354242325, -0.619132340, +1.00000000
186cor = getCor(cov, low, stdinv = 1 / std)
187cor
188+1.00000000, -0.973608196, +0.632338881, -0.229556710, +0.177286323E-1, +0.503054142, -0.529540539
189-0.973608196, +1.00000000, -0.782883883, +0.395033240, +0.148640379, -0.583259165, +0.517042637
190+0.632338881, -0.782883883, +1.00000000, -0.650833666, -0.435351551, +0.691089630, -0.466066331
191-0.229556710, +0.395033240, -0.650833666, +1.00000000, +0.268923402, -0.231683791, +0.191277429
192+0.177286323E-1, +0.148640379, -0.435351551, +0.268923402, +1.00000000, -0.394902319, -0.354242265
193+0.503054142, -0.583259165, +0.691089630, -0.231683791, -0.394902319, +1.00000000, -0.619132280
194-0.529540539, +0.517042637, -0.466066331, +0.191277429, -0.354242265, -0.619132280, +1.00000000
195getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
196+1.00000000, -4.86804104, +1.89701664, -1.37734008, +0.106371805, +0.503054142, -0.529540598
197-4.86804104, +25.0000000, -11.7432585, +11.8509960, +4.45921230, -2.91629577, +2.58521342
198+1.89701664, -11.7432585, +9.00000000, -11.7150040, -7.83632851, +2.07326889, -1.39819920
199-1.37734008, +11.8509960, -11.7150040, +36.0000000, +9.68124199, -1.39010262, +1.14766455
200+0.106371805, +4.45921230, -7.83632851, +9.68124199, +36.0000000, -2.36941409, -2.12545395
201+0.503054142, -2.91629577, +2.07326889, -1.39010262, -2.36941409, +1.00000000, -0.619132340
202-0.529540598, +2.58521342, -1.39819920, +1.14766455, -2.12545395, -0.619132340, +1.00000000
203getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
204+1.00000000, -4.86804104, +1.89701664, -1.37734008, +0.106371805, +0.503054142, -0.529540598
205-4.86804104, +25.0000000, -11.7432585, +11.8509960, +4.45921230, -2.91629577, +2.58521342
206+1.89701664, -11.7432585, +9.00000000, -11.7150040, -7.83632851, +2.07326889, -1.39819920
207-1.37734008, +11.8509960, -11.7150040, +36.0000000, +9.68124199, -1.39010262, +1.14766455
208+0.106371805, +4.45921230, -7.83632851, +9.68124199, +36.0000000, -2.36941409, -2.12545395
209+0.503054142, -2.91629577, +2.07326889, -1.39010262, -2.36941409, +1.00000000, -0.619132340
210-0.529540598, +2.58521342, -1.39819920, +1.14766455, -2.12545395, -0.619132340, +1.00000000
211
212
213ndim = getUnifRand(0, 7)
214ndim
215+3
216std = getUnifRand(1, ndim, ndim)
217std
218+3.00000000, +2.00000000, +2.00000000
219cov = getCovRand(1._TKG, std)
220cov
221+9.00000000, -4.02719975, +0.983855188
222-4.02719975, +3.99999952, -2.94776940
223+0.983855188, -2.94776940, +3.99999976
224cor = getCor(cov, uppDia)
225cor
226+1.00000000, -0.671200037, +0.163975880
227-0.671200037, +1.00000000, -0.736942530
228+0.163975880, -0.736942530, +1.00000000
229cor = getCor(cov, upp, stdinv = 1 / std)
230cor
231+1.00000000, -0.671199977, +0.163975865
232-0.671199977, +1.00000000, -0.736942351
233+0.163975865, -0.736942351, +1.00000000
234cor = getCor(cov, lowDia)
235cor
236+1.00000000, -0.671200037, +0.163975880
237-0.671200037, +1.00000000, -0.736942530
238+0.163975880, -0.736942530, +1.00000000
239cor = getCor(cov, low, stdinv = 1 / std)
240cor
241+1.00000000, -0.671199977, +0.163975865
242-0.671199977, +1.00000000, -0.736942351
243+0.163975865, -0.736942351, +1.00000000
244getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
245+9.00000000, -4.02720022, +0.983855247
246-4.02720022, +4.00000000, -2.94777012
247+0.983855247, -2.94777012, +4.00000000
248getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
249+9.00000000, -4.02720022, +0.983855247
250-4.02720022, +4.00000000, -2.94777012
251+0.983855247, -2.94777012, +4.00000000
252
253
254ndim = getUnifRand(0, 7)
255ndim
256+2
257std = getUnifRand(1, ndim, ndim)
258std
259+1.00000000, +2.00000000
260cov = getCovRand(1._TKG, std)
261cov
262+1.00000000, +0.817710817
263+0.817710817, +4.00000000
264cor = getCor(cov, uppDia)
265cor
266+1.00000000, +0.408855408
267+0.408855408, +1.00000000
268cor = getCor(cov, upp, stdinv = 1 / std)
269cor
270+1.00000000, +0.408855408
271+0.408855408, +1.00000000
272cor = getCor(cov, lowDia)
273cor
274+1.00000000, +0.408855408
275+0.408855408, +1.00000000
276cor = getCor(cov, low, stdinv = 1 / std)
277cor
278+1.00000000, +0.408855408
279+0.408855408, +1.00000000
280getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
281+1.00000000, +0.817710817
282+0.817710817, +4.00000000
283getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
284+1.00000000, +0.817710817
285+0.817710817, +4.00000000
286
287
288ndim = getUnifRand(0, 7)
289ndim
290+3
291std = getUnifRand(1, ndim, ndim)
292std
293+1.00000000, +1.00000000, +1.00000000
294cov = getCovRand(1._TKG, std)
295cov
296+1.00000000, -0.839294791, +0.697833478
297-0.839294791, +1.00000000, -0.262053370
298+0.697833478, -0.262053370, +1.00000012
299cor = getCor(cov, uppDia)
300cor
301+1.00000000, -0.839294791, +0.697833478
302-0.839294791, +1.00000000, -0.262053370
303+0.697833478, -0.262053370, +1.00000000
304cor = getCor(cov, upp, stdinv = 1 / std)
305cor
306+1.00000000, -0.839294791, +0.697833478
307-0.839294791, +1.00000000, -0.262053370
308+0.697833478, -0.262053370, +1.00000000
309cor = getCor(cov, lowDia)
310cor
311+1.00000000, -0.839294791, +0.697833478
312-0.839294791, +1.00000000, -0.262053370
313+0.697833478, -0.262053370, +1.00000000
314cor = getCor(cov, low, stdinv = 1 / std)
315cor
316+1.00000000, -0.839294791, +0.697833478
317-0.839294791, +1.00000000, -0.262053370
318+0.697833478, -0.262053370, +1.00000000
319getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
320+1.00000000, -0.839294791, +0.697833478
321-0.839294791, +1.00000000, -0.262053370
322+0.697833478, -0.262053370, +1.00000000
323getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
324+1.00000000, -0.839294791, +0.697833478
325-0.839294791, +1.00000000, -0.262053370
326+0.697833478, -0.262053370, +1.00000000
327
328
329ndim = getUnifRand(0, 7)
330ndim
331+5
332std = getUnifRand(1, ndim, ndim)
333std
334+3.00000000, +5.00000000, +4.00000000, +1.00000000, +4.00000000
335cov = getCovRand(1._TKG, std)
336cov
337+9.00000000, +5.40125608, +8.68678856, -2.37290215, -7.43177319
338+5.40125608, +25.0000019, +0.494687259, -2.83118606, +0.835105181
339+8.68678856, +0.494687259, +16.0000019, -1.84764671, -10.5330791
340-2.37290215, -2.83118606, -1.84764671, +0.999999881, +2.36803794
341-7.43177319, +0.835105181, -10.5330791, +2.36803794, +16.0000000
342cor = getCor(cov, uppDia)
343cor
344+1.00000000, +0.360083759, +0.723899066, -0.790967464, -0.619314432
345+0.360083759, +1.00000000, +0.247343630E-1, -0.566237330, +0.417552590E-1
346+0.723899066, +0.247343630E-1, +1.00000000, -0.461911738, -0.658317447
347-0.790967464, -0.566237330, -0.461911738, +1.00000000, +0.592009544
348-0.619314432, +0.417552590E-1, -0.658317447, +0.592009544, +1.00000000
349cor = getCor(cov, upp, stdinv = 1 / std)
350cor
351+1.00000000, +0.360083759, +0.723899066, -0.790967405, -0.619314432
352+0.360083759, +1.00000000, +0.247343630E-1, -0.566237211, +0.417552590E-1
353+0.723899066, +0.247343630E-1, +1.00000000, -0.461911678, -0.658317447
354-0.790967405, -0.566237211, -0.461911678, +1.00000000, +0.592009485
355-0.619314432, +0.417552590E-1, -0.658317447, +0.592009485, +1.00000000
356cor = getCor(cov, lowDia)
357cor
358+1.00000000, +0.360083759, +0.723899066, -0.790967464, -0.619314432
359+0.360083759, +1.00000000, +0.247343630E-1, -0.566237330, +0.417552590E-1
360+0.723899066, +0.247343630E-1, +1.00000000, -0.461911738, -0.658317447
361-0.790967464, -0.566237330, -0.461911738, +1.00000000, +0.592009544
362-0.619314432, +0.417552590E-1, -0.658317447, +0.592009544, +1.00000000
363cor = getCor(cov, low, stdinv = 1 / std)
364cor
365+1.00000000, +0.360083759, +0.723899066, -0.790967405, -0.619314432
366+0.360083759, +1.00000000, +0.247343630E-1, -0.566237211, +0.417552590E-1
367+0.723899066, +0.247343630E-1, +1.00000000, -0.461911678, -0.658317447
368-0.790967405, -0.566237211, -0.461911678, +1.00000000, +0.592009485
369-0.619314432, +0.417552590E-1, -0.658317447, +0.592009485, +1.00000000
370getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
371+9.00000000, +5.40125656, +8.68678856, -2.37290239, -7.43177319
372+5.40125656, +25.0000000, +0.494687259, -2.83118677, +0.835105181
373+8.68678856, +0.494687259, +16.0000000, -1.84764695, -10.5330791
374-2.37290239, -2.83118677, -1.84764695, +1.00000000, +2.36803818
375-7.43177319, +0.835105181, -10.5330791, +2.36803818, +16.0000000
376getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
377+9.00000000, +5.40125656, +8.68678856, -2.37290239, -7.43177319
378+5.40125656, +25.0000000, +0.494687259, -2.83118677, +0.835105181
379+8.68678856, +0.494687259, +16.0000000, -1.84764695, -10.5330791
380-2.37290239, -2.83118677, -1.84764695, +1.00000000, +2.36803818
381-7.43177319, +0.835105181, -10.5330791, +2.36803818, +16.0000000
382
383
384ndim = getUnifRand(0, 7)
385ndim
386+3
387std = getUnifRand(1, ndim, ndim)
388std
389+2.00000000, +2.00000000, +1.00000000
390cov = getCovRand(1._TKG, std)
391cov
392+4.00000000, -3.08771920, +1.40654802
393-3.08771920, +3.99999976, -0.253569245
394+1.40654802, -0.253569245, +1.00000000
395cor = getCor(cov, uppDia)
396cor
397+1.00000000, -0.771929920, +0.703274012
398-0.771929920, +1.00000000, -0.126784638
399+0.703274012, -0.126784638, +1.00000000
400cor = getCor(cov, upp, stdinv = 1 / std)
401cor
402+1.00000000, -0.771929801, +0.703274012
403-0.771929801, +1.00000000, -0.126784623
404+0.703274012, -0.126784623, +1.00000000
405cor = getCor(cov, lowDia)
406cor
407+1.00000000, -0.771929920, +0.703274012
408-0.771929920, +1.00000000, -0.126784638
409+0.703274012, -0.126784638, +1.00000000
410cor = getCor(cov, low, stdinv = 1 / std)
411cor
412+1.00000000, -0.771929801, +0.703274012
413-0.771929801, +1.00000000, -0.126784623
414+0.703274012, -0.126784623, +1.00000000
415getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
416+4.00000000, -3.08771968, +1.40654802
417-3.08771968, +4.00000000, -0.253569275
418+1.40654802, -0.253569275, +1.00000000
419getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
420+4.00000000, -3.08771968, +1.40654802
421-3.08771968, +4.00000000, -0.253569275
422+1.40654802, -0.253569275, +1.00000000
423
424
425ndim = getUnifRand(0, 7)
426ndim
427+3
428std = getUnifRand(1, ndim, ndim)
429std
430+3.00000000, +3.00000000, +3.00000000
431cov = getCovRand(1._TKG, std)
432cov
433+9.00000000, -4.50277233, -1.55142355
434-4.50277233, +9.00000191, +5.46954679
435-1.55142355, +5.46954679, +9.00000191
436cor = getCor(cov, uppDia)
437cor
438+1.00000000, -0.500307977, -0.172380388
439-0.500307977, +1.00000000, +0.607727349
440-0.172380388, +0.607727349, +1.00000000
441cor = getCor(cov, upp, stdinv = 1 / std)
442cor
443+1.00000000, -0.500308096, -0.172380403
444-0.500308096, +1.00000000, +0.607727468
445-0.172380403, +0.607727468, +1.00000000
446cor = getCor(cov, lowDia)
447cor
448+1.00000000, -0.500307977, -0.172380388
449-0.500307977, +1.00000000, +0.607727349
450-0.172380388, +0.607727349, +1.00000000
451cor = getCor(cov, low, stdinv = 1 / std)
452cor
453+1.00000000, -0.500308096, -0.172380403
454-0.500308096, +1.00000000, +0.607727468
455-0.172380403, +0.607727468, +1.00000000
456getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
457+9.00000000, -4.50277185, -1.55142355
458-4.50277185, +9.00000000, +5.46954632
459-1.55142355, +5.46954632, +9.00000000
460getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
461+9.00000000, -4.50277185, -1.55142355
462-4.50277185, +9.00000000, +5.46954632
463-1.55142355, +5.46954632, +9.00000000
464
465
466!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
467!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
468!Compute the Pearson correlation matrix for a pair of time series.
469!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
470!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
471
472ndim = 2; nsam = 10; dim = 2
473sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
474sample
475+5.000000E+00, +2.000000E+01, +7.000000E+00, +1.400000E+01, +1.400000E+01, +2.000000E+00, +1.000000E+01, +1.700000E+01, +1.300000E+01, +1.600000E+01
476+6.000000E+00, +1.900000E+01, +8.000000E+00, +9.000000E+00, +3.000000E+00, +2.000000E+01, +1.200000E+01, +1.500000E+01, +1.900000E+01, +1.600000E+01
477mean = getMean(sample, dim)
478mean
479+1.180000E+01, +1.270000E+01
480cor = getCor(sample, dim) ! same result as above.
481cor
482+1.000000E+00, +1.737033E-01
483+1.737033E-01, +1.000000E+00
484
485'Compute the sample correlation along the first dimension.'
486
487dim = 1
488cor = getCor(transpose(sample), dim) ! same result as above.
489cor
490+1.000000E+00, +1.737033E-01
491+1.737033E-01, +1.000000E+00
492
493'Compute the full sample correlation for a pair of time series.'
494
495cor(1,1) = getCor(sample(1,:), sample(2,:)) ! same result as above.
496cor(1,1)
497+1.737033E-01
498
499
500!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
501!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
502!Compute the Pearson correlation matrix for a weighted pair of time series.
503!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
504!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
505
506ndim = 2; nsam = 10; dim = 2
507sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
508sample
509+1.100000E+01, +4.000000E+00, +1.400000E+01, +2.000000E+00, +2.000000E+01, +1.000000E+01, +1.400000E+01, +9.000000E+00, +1.800000E+01, +1.900000E+01
510+1.200000E+01, +9.000000E+00, +5.000000E+00, +2.000000E+01, +1.100000E+01, +1.000000E+01, +1.100000E+01, +4.000000E+00, +1.500000E+01, +1.200000E+01
511call setResized(mean, ndim)
512iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
513iweight
514+10, +10, +4, +6, +2, +1, +9, +10, +3, +5
515call setMean(mean, sample, dim, iweight, iweisum)
516mean
517+10.5500002, +10.4333344
518iweisum
519+60
520rweight = iweight ! or real-valued weights.
521iweight
522+10, +10, +4, +6, +2, +1, +9, +10, +3, +5
523call setMean(mean, sample, dim, rweight, rweisum)
524mean
525+1.055000E+01, +1.043333E+01
526rweisum
527+6.000000E+01
528
529!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
530!Compute the correlation matrix with integer weights.
531!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
532
533cor = getCor(sample, dim, iweight) ! same result as above.
534cor
535+1.000000E+00, -1.158815E-01
536-1.158815E-01, +1.000000E+00
537
538'Compute the sample correlation along the first dimension.'
539
540dim = 1
541cor = getCor(transpose(sample), dim, iweight) ! same result as above.
542cor
543+1.000000E+00, -1.158815E-01
544-1.158815E-01, +1.000000E+00
545
546'Compute the full sample correlation for a pair of time series.'
547
548cor(1,1) = getCor(sample(1,:), sample(2,:), iweight) ! same result as above.
549cor(1,1)
550-1.158815E-01
551
552
553!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
554!Compute the correlation matrix with real weights.
555!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
556
557dim = 2
558cor = getCor(sample, dim, rweight) ! same result as above.
559cor
560+1.000000E+00, -1.158815E-01
561-1.158815E-01, +1.000000E+00
562
563'Compute the sample correlation along the first dimension.'
564
565dim = 1
566cor = getCor(transpose(sample), dim, rweight) ! same result as above.
567cor
568+1.000000E+00, -1.158815E-01
569-1.158815E-01, +1.000000E+00
570
571'Compute the full sample correlation for a pair of time series.'
572
573cor(1,1) = getCor(sample(1,:), sample(2,:), rweight) ! same result as above.
574cor(1,1)
575-1.158815E-01
576
577
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

Definition at line 765 of file pm_sampleCor.F90.


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