Generate and return the Spearman rank correlation matrix of the input (weighted) sample of shape (ndim, nsam)
or (nsam, ndim)
or the Spearman rank correlation coefficient a pair of (weighted) time series x(1:nsam)
and y(1:nsam)
where ndim
is the number of data dimensions (the number of data attributes) and nsam
is the number of data points.
More...
Generate and return the Spearman rank correlation matrix of the input (weighted) sample of shape (ndim, nsam)
or (nsam, ndim)
or the Spearman rank correlation coefficient a pair of (weighted) time series x(1:nsam)
and y(1:nsam)
where ndim
is the number of data dimensions (the number of data attributes) and nsam
is the number of data points.
This generic interface performs one of the following computational tasks:
-
Compute the Spearman rank correlation coefficient corresponding to an input pair of time series
x
and y
of nsam
observations.
-
Compute the Spearman rank correlation matrix corresponding to an input multivariate sample of
nsam
observations each with ndim
attributes.
See the documentation of the parent module pm_sampleCor for algorithmic details and sample correlation matrix definition.
- Parameters
-
[in] | x | : The input contiguous vector of shape (nsam) of,
-
type
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU),
-
type
integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64),
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
-
type string container css_type,
-
type string PDT container css_pdt,
or a scalar of,
-
type
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) of arbitrary length type parameter,
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 rho and sample are missing.) |
[in] | y | : The input contiguous vector of shape (nsam) of the same type and kind as the input x , 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 rho and sample are missing.) |
[in] | sample | : The input contiguous array of shape (ndim, nsam) or (nsam, ndim) of,
-
type
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU),
-
type
integer of kind any supported by the processor (e.g., IK, IK8, IK16, IK32, or IK64),
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
-
type string container css_type,
-
type string PDT container css_pdt,
or a scalar of,
-
type
character of kind any supported by the processor (e.g., SK, SKA, SKD , or SKU) of arbitrary length type parameter,
containing the sample comprised of nsam observations each with ndim attributes.
If sample is a matrix, then the input argument dim dictates the direction along which the correlation matrix rho must be computed (i.e., the direction of individual observations).
(optional. It must be present if and only if the input arguments rho , 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.
-
If
dim = 1 , the input sample is assumed to have the shape (nsam, ndim) .
-
If
dim = 2 , the input sample is assumed to have the shape (ndim, nsam) .
(optional. It must be present if and only if the input argument sample is present and is of rank 2 .) |
[in] | weight | : The input contiguous vector of length nsam ,
-
type
integer of default kind IK,
-
type
real of default kind RK,
containing the corresponding weights of individual nsam observations in sample or the pair of vectors x and y .
Note that this default RK kind type parameter requirement on input weight of type real is unlike the other pm_sample* modules of the ParaMonte library.
This requirement is enforced by the default kind type parameter of the output of setRankFractional.
(optional. default = getFilled(1, nsam).) |
- Returns
rho
: The output positive semi-definite scalar or square matrix of shape (1 : ndim, 1 : ndim)
of,
-
type
real
of default kind RK,
containing the (Spearman rank) correlation coefficient or full correlation matrix corresponding to the input sample
or time series x
and y
, whichever is present.
-
If
x(:)
and y(:)
are present, then rho
shall be a scalara scalar of value \(r_{xy}\),
\begin{equation}
\ms{rho} =
\begin{bmatrix}
1 && r_{xy} \\
r_{yx} && 1 ~.
\end{bmatrix}
\end{equation}
-
If
sample
is present, then rho
shall be a square matrix of shape [size(sample, 3 - dim), size(sample, 3 - dim)]
.
Possible calling interfaces ⛓
rho
= getRho(x(
1:nsam), y(
1:nsam) )
rho
= getRho(x(
1:nsam), y(
1:nsam), weight(
1:nsam))
rho(
1:ndim,
1:ndim)
= getRho(sample(:,:), dim)
rho(
1:ndim,
1:ndim)
= getRho(sample(:,:), dim, weight(
1:nsam))
Generate and return the Spearman rank correlation matrix of the input (weighted) sample of shape (ndi...
This module contains classes and procedures for computing properties related to the correlation matri...
- Warning
- All conditions that must hold for the generic interface setRho must equally hold for this generic interface.
These conditions are verified only if the library is built with the preprocessor macro CHECK_ENABLED=1
.
-
The
pure
procedure(s) documented herein become impure
when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1
.
By default, these procedures are pure
in release
build and impure
in debug
and testing
builds.
- See also
- getCor
setCor
getRho
setRho
getCov
setCov
setECDF
getMean
setMean
getShifted
setShifted
getVar
setVar
Example usage ⛓
23 type(display_type) :: disp
24 integer(IK) :: itry, ntry
= 10
25 character(:),
allocatable :: format
26 real(RK),
allocatable :: rho(:,:), rweight(:)
30 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
31 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
32 call disp%show(
"!Compute the Spearman correlation matrix for a pair of character-sequence time series.")
33 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
34 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
41 character(:,SKG),
allocatable :: x, y
44 call disp%show(
"x = getUnifRand(repeat('A', nsam), repeat('Z', nsam))")
45 x
= getUnifRand(
repeat(
'A', nsam),
repeat(
'Z', nsam))
47 call disp%show( x , deliml
= SK_
'''' )
48 call disp%show(
"y = getUnifRand(repeat('A', nsam), repeat('Z', nsam))")
49 y
= getUnifRand(
repeat(
'A', nsam),
repeat(
'Z', nsam))
51 call disp%show( y , deliml
= SK_
'''' )
60 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
61 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
62 call disp%show(
"!Compute the Spearman correlation matrix for a pair of character-valued time series.")
63 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
64 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
69 integer(IK) :: ndim, nsam
70 character(
2,SKG),
allocatable :: sample(:,:)
71 call disp%show(
"ndim = 2; nsam = 10")
73 call disp%show(
"call setResized(sample, [ndim, nsam])")
75 call disp%show(
"call setUnifRand(sample, SKG_'AA', SKG_'ZZ')")
78 call disp%show( sample , deliml
= SK_
'''' )
79 call disp%show(
"rho = getFilled(0., ndim, ndim)")
81 call disp%show(
"rho = getRho(sample, dim = 2_IK)")
82 rho
= getRho(sample,
dim = 2_IK)
86 call disp%show(
"Compute the sample correlation along the first dimension.", deliml
= SK_
'''')
88 call disp%show(
"rho = getFilled(0., ndim, ndim)")
90 call disp%show(
"rho = getRho(transpose(sample), dim = 1_IK)")
91 rho
= getRho(
transpose(sample),
dim = 1_IK)
95 call disp%show(
"Compute the full sample correlation for a pair of time series.", deliml
= SK_
'''')
97 call disp%show(
"rho(1,1) = getRho(sample(1,:), sample(2,:))")
98 rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:))
105 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
106 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
107 call disp%show(
"!Compute the Spearman correlation matrix for a pair of integer-valued time series.")
108 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
109 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
114 integer(IK) :: ndim, nsam
115 integer(IKG),
allocatable :: sample(:,:)
116 call disp%show(
"ndim = 2; nsam = 10")
118 call disp%show(
"sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
119 sample
= reshape(
getUnifRand(
1,
20, ndim
* nsam), shape
= [ndim, nsam], order
= [
2,
1])
122 call disp%show(
"rho = getFilled(0., ndim, ndim)")
124 call disp%show(
"rho = getRho(sample, dim = 2_IK)")
125 rho
= getRho(sample,
dim = 2_IK)
129 call disp%show(
"Compute the sample correlation along the first dimension.", deliml
= SK_
'''')
131 call disp%show(
"rho = getFilled(0., ndim, ndim)")
133 call disp%show(
"rho = getRho(transpose(sample), dim = 1_IK)")
134 rho
= getRho(
transpose(sample),
dim = 1_IK)
138 call disp%show(
"Compute the full sample correlation for a pair of time series.", deliml
= SK_
'''')
140 call disp%show(
"rho(1,1) = getRho(sample(1,:), sample(2,:))")
141 rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:))
148 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
149 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
150 call disp%show(
"!Compute the Spearman correlation matrix for a pair of real-valued time series.")
151 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
152 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
157 integer(IK) :: ndim, nsam
158 real(RKG),
allocatable :: sample(:,:)
159 format = getFormat(mold
= [
0._RKG], ed
= SK_
"es", signed
= .true._LK)
160 call disp%show(
"ndim = 2; nsam = 10")
162 call disp%show(
"sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
163 sample
= reshape(
getUnifRand(
1,
20, ndim
* nsam), shape
= [ndim, nsam], order
= [
2,
1])
165 call disp%show( sample ,
format = format )
166 call disp%show(
"rho = getFilled(0., ndim, ndim)")
168 call disp%show(
"rho = getRho(sample, dim = 2_IK)")
169 rho
= getRho(sample,
dim = 2_IK)
171 call disp%show( rho ,
format = format )
173 call disp%show(
"Compute the sample correlation along the first dimension.", deliml
= SK_
'''')
175 call disp%show(
"rho = getFilled(0., ndim, ndim)")
177 call disp%show(
"rho = getRho(transpose(sample), dim = 1_IK)")
178 rho
= getRho(
transpose(sample),
dim = 1_IK)
180 call disp%show( rho ,
format = format )
182 call disp%show(
"Compute the full sample correlation for a pair of time series.", deliml
= SK_
'''')
184 call disp%show(
"rho(1,1) = getRho(sample(1,:), sample(2,:))")
185 rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:))
187 call disp%show( rho(
1,
1) ,
format = format )
192 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
193 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
194 call disp%show(
"!Compute the Spearman correlation matrix for a weighted pair of time series.")
195 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
196 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
202 integer(IK) :: ndim, nsam
203 integer(IK),
allocatable :: iweight(:)
204 real(RKG),
allocatable :: sample(:,:)
205 format = getFormat(mold
= [
0._RKG], ed
= SK_
"es", signed
= .true._LK)
206 call disp%show(
"ndim = 2; nsam = 10")
208 call disp%show(
"sample = reshape(getUnifRand(1, 20, ndim * nsam), shape = [ndim, nsam], order = [2, 1])")
209 sample
= reshape(
getUnifRand(
1,
20, ndim
* nsam), shape
= [ndim, nsam], order
= [
2,
1])
211 call disp%show( sample ,
format = format )
212 call disp%show(
"iweight = getUnifRand(1, 10, nsam) ! integer-valued weights.")
216 call disp%show(
"rweight = iweight ! or real-valued weights.")
222 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
223 call disp%show(
"!Compute the correlation matrix with integer weights.")
224 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
227 call disp%show(
"rho = getFilled(0., ndim, ndim)")
229 call disp%show(
"rho = getRho(sample, 2_IK, iweight)")
230 rho
= getRho(sample,
2_IK, iweight)
232 call disp%show( rho ,
format = format )
234 call disp%show(
"Compute the sample correlation along the first dimension.", deliml
= SK_
'''')
236 call disp%show(
"rho = getFilled(0., ndim, ndim)")
238 call disp%show(
"rho = getRho(transpose(sample), 1_IK, iweight)")
239 rho
= getRho(
transpose(sample),
1_IK, iweight)
241 call disp%show( rho ,
format = format )
243 call disp%show(
"rho = getFilled(0., ndim, ndim)")
246 call disp%show(
"Compute the full sample correlation for a pair of time series.", deliml
= SK_
'''')
248 call disp%show(
"rho(1,1) = getRho(sample(1,:), sample(2,:), iweight)")
249 rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:), iweight)
251 call disp%show( rho(
1,
1) ,
format = format )
255 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
256 call disp%show(
"!Compute the correlation matrix with real weights.")
257 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
260 call disp%show(
"rho = getFilled(0., ndim, ndim)")
262 call disp%show(
"rho = getRho(sample, 2_IK, rweight)")
263 rho
= getRho(sample,
2_IK, rweight)
265 call disp%show( rho ,
format = format )
267 call disp%show(
"Compute the sample correlation along the first dimension.", deliml
= SK_
'''')
269 call disp%show(
"rho = getFilled(0., ndim, ndim)")
271 call disp%show(
"rho = getRho(transpose(sample), 1_IK, rweight)")
272 rho
= getRho(
transpose(sample),
1_IK, rweight)
274 call disp%show( rho ,
format = format )
276 call disp%show(
"rho = getFilled(0., ndim, ndim)")
279 call disp%show(
"Compute the full sample correlation for a pair of time series.", deliml
= SK_
'''')
281 call disp%show(
"rho(1,1) = getRho(sample(1,:), sample(2,:), rweight)")
282 rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:), rweight)
284 call disp%show( rho(
1,
1) ,
format = format )
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 and return an output array whose elements are the reversed-order elements of the input array...
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...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Return a uniform random scalar or contiguous array of arbitrary rank of randomly uniformly distribute...
This is a generic method of the derived type display_type with pass attribute.
This is a generic method of the derived type display_type with pass attribute.
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 reversing the order of elements in arrays ...
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...
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...
type(display_type) disp
This is a scalar module variable an object of type display_type for general display.
This module defines the relevant Fortran kind type-parameters frequently used in the ParaMonte librar...
integer, parameter RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
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.
Example Unix compile command via Intel ifort
compiler ⛓
3ifort -fpp -standard-semantics -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example Windows Batch compile command via Intel ifort
compiler ⛓
2set PATH=..\..\..\lib;%PATH%
3ifort /fpp /standard-semantics /O3 /I:..\..\..\include main.F90 ..\..\..\lib\libparamonte*.lib /exe:main.exe
Example Unix / MinGW compile command via GNU gfortran
compiler ⛓
3gfortran -cpp -ffree-line-length-none -O3 -Wl,-rpath,../../../lib -I../../../inc main.F90 ../../../lib/libparamonte* -o main.exe
Example output ⛓
12y
= getUnifRand(
repeat(
'A', nsam),
repeat(
'Z', nsam))
30'HK',
'GO',
'FK',
'ND',
'LM',
'AH',
'FJ',
'FU',
'IY',
'WW'
31'IN',
'GK',
'IY',
'UA',
'FK',
'WP',
'GG',
'MA',
'OL',
'EQ'
33rho
= getRho(sample,
dim = 2_IK)
35+1.0000000000000000,
-0.33333333333333337
36-0.33333333333333337,
+1.0000000000000000
38'Compute the sample correlation along the first dimension.'
41rho
= getRho(
transpose(sample),
dim = 1_IK)
43+1.0000000000000000,
-0.33333333333333337
44-0.33333333333333337,
+1.0000000000000000
46'Compute the full sample correlation for a pair of time series.'
48rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:))
60sample
= reshape(
getUnifRand(
1,
20, ndim
* nsam), shape
= [ndim, nsam], order
= [
2,
1])
62+9,
+1,
+15,
+14,
+1,
+18,
+16,
+19,
+10,
+16
63+5,
+7,
+17,
+17,
+4,
+19,
+19,
+11,
+9,
+17
65rho
= getRho(sample,
dim = 2_IK)
67+1.0000000000000000,
+0.76473886837914540
68+0.76473886837914540,
+1.0000000000000000
70'Compute the sample correlation along the first dimension.'
73rho
= getRho(
transpose(sample),
dim = 1_IK)
75+1.0000000000000000,
+0.76473886837914540
76+0.76473886837914540,
+1.0000000000000000
78'Compute the full sample correlation for a pair of time series.'
80rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:))
92sample
= reshape(
getUnifRand(
1,
20, ndim
* nsam), shape
= [ndim, nsam], order
= [
2,
1])
94+2.000000E+01,
+1.800000E+01,
+9.000000E+00,
+1.100000E+01,
+1.000000E+01,
+1.100000E+01,
+1.900000E+01,
+1.200000E+01,
+1.300000E+01,
+7.000000E+00
95+2.000000E+00,
+6.000000E+00,
+1.700000E+01,
+1.300000E+01,
+1.100000E+01,
+1.900000E+01,
+1.900000E+01,
+6.000000E+00,
+1.700000E+01,
+3.000000E+00
97rho
= getRho(sample,
dim = 2_IK)
99+1.000000E+00,
-4.601314E-02
100-4.601314E-02,
+1.000000E+00
102'Compute the sample correlation along the first dimension.'
105rho
= getRho(
transpose(sample),
dim = 1_IK)
107+1.000000E+00,
-4.601314E-02
108-4.601314E-02,
+1.000000E+00
110'Compute the full sample correlation for a pair of time series.'
112rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:))
124sample
= reshape(
getUnifRand(
1,
20, ndim
* nsam), shape
= [ndim, nsam], order
= [
2,
1])
126+1.400000E+01,
+8.000000E+00,
+1.500000E+01,
+1.400000E+01,
+1.500000E+01,
+5.000000E+00,
+4.000000E+00,
+5.000000E+00,
+1.500000E+01,
+6.000000E+00
127+9.000000E+00,
+5.000000E+00,
+9.000000E+00,
+1.500000E+01,
+1.200000E+01,
+1.300000E+01,
+1.200000E+01,
+1.800000E+01,
+5.000000E+00,
+1.000000E+01
130+6,
+9,
+8,
+2,
+9,
+5,
+5,
+7,
+2,
+5
133+6,
+9,
+8,
+2,
+9,
+5,
+5,
+7,
+2,
+5
140rho
= getRho(sample,
2_IK, iweight)
142+1.000000E+00,
-4.184408E-01
143-4.184408E-01,
+1.000000E+00
145'Compute the sample correlation along the first dimension.'
148rho
= getRho(
transpose(sample),
1_IK, iweight)
150+1.000000E+00,
-4.184408E-01
151-4.184408E-01,
+1.000000E+00
155'Compute the full sample correlation for a pair of time series.'
157rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:), iweight)
167rho
= getRho(sample,
2_IK, rweight)
169+1.000000E+00,
-4.184408E-01
170-4.184408E-01,
+1.000000E+00
172'Compute the sample correlation along the first dimension.'
175rho
= getRho(
transpose(sample),
1_IK, rweight)
177+1.000000E+00,
-4.184408E-01
178-4.184408E-01,
+1.000000E+00
182'Compute the full sample correlation for a pair of time series.'
184rho(
1,
1)
= getRho(sample(
1,:), sample(
2,:), rweight)
- 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.
-
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.
-
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.
- Copyright
- Computational Data Science Lab
- Author:
- Fatemeh Bagheri, Monday 02:15 AM, September 27, 2021, Dallas, TX
Amir Shahmoradi, Monday March 6, 2017, 3:22 pm, Institute for Computational Engineering and Sciences (ICES), The University of Texas at Austin.
Definition at line 7643 of file pm_sampleCor.F90.