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+1
10std = getUnifRand(1, ndim, ndim)
11std
12+1.00000000
13cov = getCovRand(1._TKG, std)
14cov
15+1.00000000
16cor = getCor(cov, uppDia)
17cor
18+1.00000000
19cor = getCor(cov, upp, stdinv = 1 / std)
20cor
21+1.00000000
22cor = getCor(cov, lowDia)
23cor
24+1.00000000
25cor = getCor(cov, low, stdinv = 1 / std)
26cor
27+1.00000000
28getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
29+1.00000000
30getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
31+1.00000000
32
33
34ndim = getUnifRand(0, 7)
35ndim
36+5
37std = getUnifRand(1, ndim, ndim)
38std
39+4.00000000, +4.00000000, +3.00000000, +5.00000000, +2.00000000
40cov = getCovRand(1._TKG, std)
41cov
42+16.0000000, -3.50091743, -9.06825066, +7.43039036, +5.47583580
43-3.50091743, +16.0000000, -4.15841007, -8.57373238, -4.47891235
44-9.06825066, -4.15841007, +9.00000191, -5.13605833, -1.18745136
45+7.43039036, -8.57373238, -5.13605833, +25.0000000, +5.46598053
46+5.47583580, -4.47891235, -1.18745136, +5.46598053, +3.99999976
47cor = getCor(cov, uppDia)
48cor
49+1.00000000, -0.218807340, -0.755687535, +0.371519536, +0.684479535
50-0.218807340, +1.00000000, -0.346534163, -0.428686619, -0.559864104
51-0.755687535, -0.346534163, +1.00000000, -0.342403859, -0.197908565
52+0.371519536, -0.428686619, -0.342403859, +1.00000000, +0.546598136
53+0.684479535, -0.559864104, -0.197908565, +0.546598136, +1.00000000
54cor = getCor(cov, upp, stdinv = 1 / std)
55cor
56+1.00000000, -0.218807340, -0.755687594, +0.371519536, +0.684479475
57-0.218807340, +1.00000000, -0.346534193, -0.428686619, -0.559864044
58-0.755687594, -0.346534193, +1.00000000, -0.342403919, -0.197908565
59+0.371519536, -0.428686619, -0.342403919, +1.00000000, +0.546598077
60+0.684479475, -0.559864044, -0.197908565, +0.546598077, +1.00000000
61cor = getCor(cov, lowDia)
62cor
63+1.00000000, -0.218807340, -0.755687535, +0.371519536, +0.684479535
64-0.218807340, +1.00000000, -0.346534163, -0.428686619, -0.559864104
65-0.755687535, -0.346534163, +1.00000000, -0.342403859, -0.197908565
66+0.371519536, -0.428686619, -0.342403859, +1.00000000, +0.546598136
67+0.684479535, -0.559864104, -0.197908565, +0.546598136, +1.00000000
68cor = getCor(cov, low, stdinv = 1 / std)
69cor
70+1.00000000, -0.218807340, -0.755687594, +0.371519536, +0.684479475
71-0.218807340, +1.00000000, -0.346534193, -0.428686619, -0.559864044
72-0.755687594, -0.346534193, +1.00000000, -0.342403919, -0.197908565
73+0.371519536, -0.428686619, -0.342403919, +1.00000000, +0.546598077
74+0.684479475, -0.559864044, -0.197908565, +0.546598077, +1.00000000
75getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
76+16.0000000, -3.50091743, -9.06825066, +7.43039083, +5.47583628
77-3.50091743, +16.0000000, -4.15841007, -8.57373238, -4.47891283
78-9.06825066, -4.15841007, +9.00000000, -5.13605785, -1.18745136
79+7.43039083, -8.57373238, -5.13605785, +25.0000000, +5.46598148
80+5.47583628, -4.47891283, -1.18745136, +5.46598148, +4.00000000
81getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
82+16.0000000, -3.50091743, -9.06825066, +7.43039083, +5.47583628
83-3.50091743, +16.0000000, -4.15841007, -8.57373238, -4.47891283
84-9.06825066, -4.15841007, +9.00000000, -5.13605785, -1.18745136
85+7.43039083, -8.57373238, -5.13605785, +25.0000000, +5.46598148
86+5.47583628, -4.47891283, -1.18745136, +5.46598148, +4.00000000
87
88
89ndim = getUnifRand(0, 7)
90ndim
91+5
92std = getUnifRand(1, ndim, ndim)
93std
94+3.00000000, +3.00000000, +3.00000000, +5.00000000, +3.00000000
95cov = getCovRand(1._TKG, std)
96cov
97+9.00000000, +6.81052399, +6.75230598, +11.4625969, +0.963286519
98+6.81052399, +8.99999905, +8.17658234, +6.70068026, -2.21397209
99+6.75230598, +8.17658234, +9.00000000, +6.10922766, -3.61728859
100+11.4625969, +6.70068026, +6.10922766, +24.9999981, -1.81004620
101+0.963286519, -2.21397209, -3.61728859, -1.81004620, +9.00000191
102cor = getCor(cov, uppDia)
103cor
104+1.00000000, +0.756725013, +0.750256300, +0.764173150, +0.107031830
105+0.756725013, +1.00000000, +0.908509254, +0.446712077, -0.245996922
106+0.750256300, +0.908509254, +1.00000000, +0.407281876, -0.401920944
107+0.764173150, +0.446712077, +0.407281876, +1.00000000, -0.120669737
108+0.107031830, -0.245996922, -0.401920944, -0.120669737, +1.00000000
109cor = getCor(cov, upp, stdinv = 1 / std)
110cor
111+1.00000000, +0.756724954, +0.750256300, +0.764173150, +0.107031845
112+0.756724954, +1.00000000, +0.908509195, +0.446712047, -0.245996922
113+0.750256300, +0.908509195, +1.00000000, +0.407281876, -0.401920974
114+0.764173150, +0.446712047, +0.407281876, +1.00000000, -0.120669752
115+0.107031845, -0.245996922, -0.401920974, -0.120669752, +1.00000000
116cor = getCor(cov, lowDia)
117cor
118+1.00000000, +0.756725013, +0.750256300, +0.764173150, +0.107031830
119+0.756725013, +1.00000000, +0.908509254, +0.446712077, -0.245996922
120+0.750256300, +0.908509254, +1.00000000, +0.407281876, -0.401920944
121+0.764173150, +0.446712077, +0.407281876, +1.00000000, -0.120669737
122+0.107031830, -0.245996922, -0.401920944, -0.120669737, +1.00000000
123cor = getCor(cov, low, stdinv = 1 / std)
124cor
125+1.00000000, +0.756724954, +0.750256300, +0.764173150, +0.107031845
126+0.756724954, +1.00000000, +0.908509195, +0.446712047, -0.245996922
127+0.750256300, +0.908509195, +1.00000000, +0.407281876, -0.401920974
128+0.764173150, +0.446712047, +0.407281876, +1.00000000, -0.120669752
129+0.107031845, -0.245996922, -0.401920974, -0.120669752, +1.00000000
130getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
131+9.00000000, +6.81052494, +6.75230694, +11.4625969, +0.963286459
132+6.81052494, +9.00000000, +8.17658329, +6.70068121, -2.21397233
133+6.75230694, +8.17658329, +9.00000000, +6.10922813, -3.61728859
134+11.4625969, +6.70068121, +6.10922813, +25.0000000, -1.81004608
135+0.963286459, -2.21397233, -3.61728859, -1.81004608, +9.00000000
136getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
137+9.00000000, +6.81052494, +6.75230694, +11.4625969, +0.963286459
138+6.81052494, +9.00000000, +8.17658329, +6.70068121, -2.21397233
139+6.75230694, +8.17658329, +9.00000000, +6.10922813, -3.61728859
140+11.4625969, +6.70068121, +6.10922813, +25.0000000, -1.81004608
141+0.963286459, -2.21397233, -3.61728859, -1.81004608, +9.00000000
142
143
144ndim = getUnifRand(0, 7)
145ndim
146+3
147std = getUnifRand(1, ndim, ndim)
148std
149+1.00000000, +2.00000000, +1.00000000
150cov = getCovRand(1._TKG, std)
151cov
152+1.00000000, +1.36172950, +0.170969695
153+1.36172950, +4.00000000, -0.695341766
154+0.170969695, -0.695341766, +0.999999881
155cor = getCor(cov, uppDia)
156cor
157+1.00000000, +0.680864751, +0.170969710
158+0.680864751, +1.00000000, -0.347670913
159+0.170969710, -0.347670913, +1.00000000
160cor = getCor(cov, upp, stdinv = 1 / std)
161cor
162+1.00000000, +0.680864751, +0.170969695
163+0.680864751, +1.00000000, -0.347670883
164+0.170969695, -0.347670883, +1.00000000
165cor = getCor(cov, lowDia)
166cor
167+1.00000000, +0.680864751, +0.170969710
168+0.680864751, +1.00000000, -0.347670913
169+0.170969710, -0.347670913, +1.00000000
170cor = getCor(cov, low, stdinv = 1 / std)
171cor
172+1.00000000, +0.680864751, +0.170969695
173+0.680864751, +1.00000000, -0.347670883
174+0.170969695, -0.347670883, +1.00000000
175getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
176+1.00000000, +1.36172950, +0.170969710
177+1.36172950, +4.00000000, -0.695341825
178+0.170969710, -0.695341825, +1.00000000
179getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
180+1.00000000, +1.36172950, +0.170969710
181+1.36172950, +4.00000000, -0.695341825
182+0.170969710, -0.695341825, +1.00000000
183
184
185ndim = getUnifRand(0, 7)
186ndim
187+3
188std = getUnifRand(1, ndim, ndim)
189std
190+1.00000000, +2.00000000, +2.00000000
191cov = getCovRand(1._TKG, std)
192cov
193+1.00000000, +1.75050652, +1.39740634
194+1.75050652, +4.00000000, +3.81495476
195+1.39740634, +3.81495476, +4.00000048
196cor = getCor(cov, uppDia)
197cor
198+1.00000000, +0.875253260, +0.698703170
199+0.875253260, +1.00000000, +0.953738689
200+0.698703170, +0.953738689, +1.00000000
201cor = getCor(cov, upp, stdinv = 1 / std)
202cor
203+1.00000000, +0.875253260, +0.698703170
204+0.875253260, +1.00000000, +0.953738689
205+0.698703170, +0.953738689, +1.00000000
206cor = getCor(cov, lowDia)
207cor
208+1.00000000, +0.875253260, +0.698703170
209+0.875253260, +1.00000000, +0.953738689
210+0.698703170, +0.953738689, +1.00000000
211cor = getCor(cov, low, stdinv = 1 / std)
212cor
213+1.00000000, +0.875253260, +0.698703170
214+0.875253260, +1.00000000, +0.953738689
215+0.698703170, +0.953738689, +1.00000000
216getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
217+1.00000000, +1.75050652, +1.39740634
218+1.75050652, +4.00000000, +3.81495476
219+1.39740634, +3.81495476, +4.00000000
220getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
221+1.00000000, +1.75050652, +1.39740634
222+1.75050652, +4.00000000, +3.81495476
223+1.39740634, +3.81495476, +4.00000000
224
225
226ndim = getUnifRand(0, 7)
227ndim
228+0
229std = getUnifRand(1, ndim, ndim)
230std
231
232cov = getCovRand(1._TKG, std)
233cov
234cor = getCor(cov, uppDia)
235cor
236cor = getCor(cov, upp, stdinv = 1 / std)
237cor
238cor = getCor(cov, lowDia)
239cor
240cor = getCor(cov, low, stdinv = 1 / std)
241cor
242getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
243getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
244
245
246ndim = getUnifRand(0, 7)
247ndim
248+6
249std = getUnifRand(1, ndim, ndim)
250std
251+2.00000000, +5.00000000, +6.00000000, +2.00000000, +2.00000000, +1.00000000
252cov = getCovRand(1._TKG, std)
253cov
254+3.99999952, +7.70508003, -6.53893089, +1.20600760, -1.48390090, +0.696430922
255+7.70508003, +24.9999962, -7.74564838, +3.88635874, -6.25388718, -0.214107186
256-6.53893089, -7.74564838, +36.0000000, +4.85372639, -4.11909485, +0.917261839E-2
257+1.20600760, +3.88635874, +4.85372639, +4.00000000, -3.05074596, +0.534889698
258-1.48390090, -6.25388718, -4.11909485, -3.05074596, +3.99999952, -0.616045356
259+0.696430922, -0.214107186, +0.917261839E-2, +0.534889698, -0.616045356, +1.00000000
260cor = getCor(cov, uppDia)
261cor
262+1.00000000, +0.770508170, -0.544910967, +0.301501930, -0.370975316, +0.348215491
263+0.770508170, +1.00000000, -0.258188307, +0.388635904, -0.625388861, -0.428214408E-1
264-0.544910967, -0.258188307, +1.00000000, +0.404477209, -0.343257934, +0.152876973E-2
265+0.301501930, +0.388635904, +0.404477209, +1.00000000, -0.762686610, +0.267444849
266-0.370975316, -0.625388861, -0.343257934, -0.762686610, +1.00000000, -0.308022708
267+0.348215491, -0.428214408E-1, +0.152876973E-2, +0.267444849, -0.308022708, +1.00000000
268cor = getCor(cov, upp, stdinv = 1 / std)
269cor
270+1.00000000, +0.770507991, -0.544910908, +0.301501900, -0.370975226, +0.348215461
271+0.770507991, +1.00000000, -0.258188307, +0.388635874, -0.625388741, -0.428214371E-1
272-0.544910908, -0.258188307, +1.00000000, +0.404477209, -0.343257904, +0.152876973E-2
273+0.301501900, +0.388635874, +0.404477209, +1.00000000, -0.762686491, +0.267444849
274-0.370975226, -0.625388741, -0.343257904, -0.762686491, +1.00000000, -0.308022678
275+0.348215461, -0.428214371E-1, +0.152876973E-2, +0.267444849, -0.308022678, +1.00000000
276cor = getCor(cov, lowDia)
277cor
278+1.00000000, +0.770508170, -0.544910967, +0.301501930, -0.370975316, +0.348215491
279+0.770508170, +1.00000000, -0.258188307, +0.388635904, -0.625388861, -0.428214408E-1
280-0.544910967, -0.258188307, +1.00000000, +0.404477209, -0.343257934, +0.152876973E-2
281+0.301501930, +0.388635904, +0.404477209, +1.00000000, -0.762686610, +0.267444849
282-0.370975316, -0.625388861, -0.343257934, -0.762686610, +1.00000000, -0.308022708
283+0.348215491, -0.428214408E-1, +0.152876973E-2, +0.267444849, -0.308022708, +1.00000000
284cor = getCor(cov, low, stdinv = 1 / std)
285cor
286+1.00000000, +0.770507991, -0.544910908, +0.301501900, -0.370975226, +0.348215461
287+0.770507991, +1.00000000, -0.258188307, +0.388635874, -0.625388741, -0.428214371E-1
288-0.544910908, -0.258188307, +1.00000000, +0.404477209, -0.343257904, +0.152876973E-2
289+0.301501900, +0.388635874, +0.404477209, +1.00000000, -0.762686491, +0.267444849
290-0.370975226, -0.625388741, -0.343257904, -0.762686491, +1.00000000, -0.308022678
291+0.348215461, -0.428214371E-1, +0.152876973E-2, +0.267444849, -0.308022678, +1.00000000
292getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
293+4.00000000, +7.70508194, -6.53893185, +1.20600772, -1.48390126, +0.696430981
294+7.70508194, +25.0000000, -7.74564934, +3.88635898, -6.25388861, -0.214107201
295-6.53893185, -7.74564934, +36.0000000, +4.85372639, -4.11909533, +0.917261839E-2
296+1.20600772, +3.88635898, +4.85372639, +4.00000000, -3.05074644, +0.534889698
297-1.48390126, -6.25388861, -4.11909533, -3.05074644, +4.00000000, -0.616045415
298+0.696430981, -0.214107201, +0.917261839E-2, +0.534889698, -0.616045415, +1.00000000
299getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
300+4.00000000, +7.70508194, -6.53893185, +1.20600772, -1.48390126, +0.696430981
301+7.70508194, +25.0000000, -7.74564934, +3.88635898, -6.25388861, -0.214107201
302-6.53893185, -7.74564934, +36.0000000, +4.85372639, -4.11909533, +0.917261839E-2
303+1.20600772, +3.88635898, +4.85372639, +4.00000000, -3.05074644, +0.534889698
304-1.48390126, -6.25388861, -4.11909533, -3.05074644, +4.00000000, -0.616045415
305+0.696430981, -0.214107201, +0.917261839E-2, +0.534889698, -0.616045415, +1.00000000
306
307
308ndim = getUnifRand(0, 7)
309ndim
310+1
311std = getUnifRand(1, ndim, ndim)
312std
313+1.00000000
314cov = getCovRand(1._TKG, std)
315cov
316+1.00000000
317cor = getCor(cov, uppDia)
318cor
319+1.00000000
320cor = getCor(cov, upp, stdinv = 1 / std)
321cor
322+1.00000000
323cor = getCor(cov, lowDia)
324cor
325+1.00000000
326cor = getCor(cov, low, stdinv = 1 / std)
327cor
328+1.00000000
329getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
330+1.00000000
331getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
332+1.00000000
333
334
335ndim = getUnifRand(0, 7)
336ndim
337+4
338std = getUnifRand(1, ndim, ndim)
339std
340+3.00000000, +1.00000000, +4.00000000, +4.00000000
341cov = getCovRand(1._TKG, std)
342cov
343+9.00000000, -2.91477394, +9.95705700, +3.87634468
344-2.91477394, +1.00000000, -2.91490459, -1.73088920
345+9.95705700, -2.91490459, +16.0000000, +4.95686865
346+3.87634468, -1.73088920, +4.95686865, +16.0000000
347cor = getCor(cov, uppDia)
348cor
349+1.00000000, -0.971591353, +0.829754770, +0.323028743
350-0.971591353, +1.00000000, -0.728726149, -0.432722300
351+0.829754770, -0.728726149, +1.00000000, +0.309804291
352+0.323028743, -0.432722300, +0.309804291, +1.00000000
353cor = getCor(cov, upp, stdinv = 1 / std)
354cor
355+1.00000000, -0.971591353, +0.829754770, +0.323028743
356-0.971591353, +1.00000000, -0.728726149, -0.432722300
357+0.829754770, -0.728726149, +1.00000000, +0.309804291
358+0.323028743, -0.432722300, +0.309804291, +1.00000000
359cor = getCor(cov, lowDia)
360cor
361+1.00000000, -0.971591353, +0.829754770, +0.323028743
362-0.971591353, +1.00000000, -0.728726149, -0.432722300
363+0.829754770, -0.728726149, +1.00000000, +0.309804291
364+0.323028743, -0.432722300, +0.309804291, +1.00000000
365cor = getCor(cov, low, stdinv = 1 / std)
366cor
367+1.00000000, -0.971591353, +0.829754770, +0.323028743
368-0.971591353, +1.00000000, -0.728726149, -0.432722300
369+0.829754770, -0.728726149, +1.00000000, +0.309804291
370+0.323028743, -0.432722300, +0.309804291, +1.00000000
371getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
372+9.00000000, -2.91477394, +9.95705700, +3.87634492
373-2.91477394, +1.00000000, -2.91490459, -1.73088920
374+9.95705700, -2.91490459, +16.0000000, +4.95686865
375+3.87634492, -1.73088920, +4.95686865, +16.0000000
376getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
377+9.00000000, -2.91477394, +9.95705700, +3.87634492
378-2.91477394, +1.00000000, -2.91490459, -1.73088920
379+9.95705700, -2.91490459, +16.0000000, +4.95686865
380+3.87634492, -1.73088920, +4.95686865, +16.0000000
381
382
383ndim = getUnifRand(0, 7)
384ndim
385+3
386std = getUnifRand(1, ndim, ndim)
387std
388+1.00000000, +3.00000000, +2.00000000
389cov = getCovRand(1._TKG, std)
390cov
391+1.00000000, +2.09673786, +1.46986127
392+2.09673786, +9.00000191, +0.369257390
393+1.46986127, +0.369257390, +4.00000048
394cor = getCor(cov, uppDia)
395cor
396+1.00000000, +0.698912561, +0.734930634
397+0.698912561, +1.00000000, +0.615428947E-1
398+0.734930634, +0.615428947E-1, +1.00000000
399cor = getCor(cov, upp, stdinv = 1 / std)
400cor
401+1.00000000, +0.698912621, +0.734930634
402+0.698912621, +1.00000000, +0.615428984E-1
403+0.734930634, +0.615428984E-1, +1.00000000
404cor = getCor(cov, lowDia)
405cor
406+1.00000000, +0.698912561, +0.734930634
407+0.698912561, +1.00000000, +0.615428947E-1
408+0.734930634, +0.615428947E-1, +1.00000000
409cor = getCor(cov, low, stdinv = 1 / std)
410cor
411+1.00000000, +0.698912621, +0.734930634
412+0.698912621, +1.00000000, +0.615428984E-1
413+0.734930634, +0.615428984E-1, +1.00000000
414getCov(getCor(cov, lowDia), lowDia, std) ! reconstruct the original covariance matrix.
415+1.00000000, +2.09673762, +1.46986127
416+2.09673762, +9.00000000, +0.369257361
417+1.46986127, +0.369257361, +4.00000000
418getCov(getCor(cov, uppDia), uppDia, std) ! reconstruct the original covariance matrix.
419+1.00000000, +2.09673762, +1.46986127
420+2.09673762, +9.00000000, +0.369257361
421+1.46986127, +0.369257361, +4.00000000
422
423
424!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
425!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
426!Compute the Pearson correlation matrix for a pair of time series.
427!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
428!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
429
430ndim = 2; nsam = 10; dim = 2
431sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
432sample
433+1.000000E+00, +1.600000E+01, +1.500000E+01, +2.000000E+00, +2.000000E+01, +1.900000E+01, +1.000000E+01, +1.400000E+01, +4.000000E+00, +1.800000E+01
434+1.600000E+01, +1.700000E+01, +1.200000E+01, +1.700000E+01, +1.000000E+00, +1.700000E+01, +4.000000E+00, +3.000000E+00, +1.700000E+01, +1.000000E+01
435mean = getMean(sample, dim)
436mean
437+1.190000E+01, +1.140000E+01
438cor = getCor(sample, dim) ! same result as above.
439cor
440+1.000000E+00, -4.297788E-01
441-4.297788E-01, +1.000000E+00
442
443'Compute the sample correlation along the first dimension.'
444
445dim = 1
446cor = getCor(transpose(sample), dim) ! same result as above.
447cor
448+1.000000E+00, -4.297788E-01
449-4.297788E-01, +1.000000E+00
450
451'Compute the full sample correlation for a pair of time series.'
452
453cor(1,1) = getCor(sample(1,:), sample(2,:)) ! same result as above.
454cor(1,1)
455-4.297788E-01
456
457
458!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
459!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
460!Compute the Pearson correlation matrix for a weighted pair of time series.
461!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
462!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
463
464ndim = 2; nsam = 10; dim = 2
465sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])
466sample
467+4.000000E+00, +5.000000E+00, +1.600000E+01, +1.400000E+01, +1.000000E+01, +2.000000E+00, +4.000000E+00, +1.900000E+01, +4.000000E+00, +2.000000E+00
468+1.100000E+01, +1.600000E+01, +3.000000E+00, +7.000000E+00, +1.100000E+01, +1.300000E+01, +7.000000E+00, +2.000000E+00, +4.000000E+00, +1.700000E+01
469call setResized(mean, ndim)
470iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.
471iweight
472+10, +10, +8, +2, +1, +8, +10, +2, +4, +9
473call setMean(mean, sample, dim, iweight, iweisum)
474mean
475+6.00000000, +10.4062500
476iweisum
477+64
478rweight = iweight ! or real-valued weights.
479iweight
480+10, +10, +8, +2, +1, +8, +10, +2, +4, +9
481call setMean(mean, sample, dim, rweight, rweisum)
482mean
483+6.000000E+00, +1.040625E+01
484rweisum
485+6.400000E+01
486
487!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
488!Compute the correlation matrix with integer weights.
489!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
490
491cor = getCor(sample, dim, iweight) ! same result as above.
492cor
493+1.000000E+00, -6.849387E-01
494-6.849387E-01, +1.000000E+00
495
496'Compute the sample correlation along the first dimension.'
497
498dim = 1
499cor = getCor(transpose(sample), dim, iweight) ! same result as above.
500cor
501+1.000000E+00, -6.849387E-01
502-6.849387E-01, +1.000000E+00
503
504'Compute the full sample correlation for a pair of time series.'
505
506cor(1,1) = getCor(sample(1,:), sample(2,:), iweight) ! same result as above.
507cor(1,1)
508-6.849387E-01
509
510
511!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
512!Compute the correlation matrix with real weights.
513!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
514
515dim = 2
516cor = getCor(sample, dim, rweight) ! same result as above.
517cor
518+1.000000E+00, -6.849387E-01
519-6.849387E-01, +1.000000E+00
520
521'Compute the sample correlation along the first dimension.'
522
523dim = 1
524cor = getCor(transpose(sample), dim, rweight) ! same result as above.
525cor
526+1.000000E+00, -6.849387E-01
527-6.849387E-01, +1.000000E+00
528
529'Compute the full sample correlation for a pair of time series.'
530
531cor(1,1) = getCor(sample(1,:), sample(2,:), rweight) ! same result as above.
532cor(1,1)
533-6.849387E-01
534
535
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: