Generate and return the variance of the input sample of type complex
or real
of shape (nsam)
or (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.
See the documentation of the parent module pm_sampleVar for algorithmic details and variance definitions.
- Parameters
-
[in] | sample | : The input contiguous array of shape (nsam) , (ndim, nsam) , or (nsam, ndim) of,
-
type
complex of kind any supported by the processor (e.g., CK, CK32, CK64, or CK128),
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing the sample comprised of nsam observations each with ndim attributes.
If sample is a matrix, then the direction along which variance is computed is dictated by the input argument dim .
|
[in] | dim | : The input scalar integer of default kind IK representing the dimension (1 or 2 ) of the input sample along which the mean must be computed.
-
If
dim = 1 , the input sample of rank 2 is assumed to have the shape (nsam, ndim) .
-
If
dim = 2 , the input sample of rank 2 is assumed to have the shape (ndim, nsam) .
The input dim must always be 1 or missing for an input sample of rank 1 .
(optional. If missing, the variance of the whole input sample is computed.) |
[in] | weight | : The contiguous vector of length nsam of,
-
type
integer of default kind IK, or
-
type
real of the same kind as the input sample ,
containing the corresponding weights of individual nsam observations in sample .
(optional. default = getFilled(1, size(sample, dim)) if dim is present, or getFilled(1, size(sample)) if dim is missing.) |
[in] | correction | : The input scalar object that can be the following:
-
The constant fweight or an object of type fweight implying a bias correction based on the assumption of frequency weights for the sample observations, even if the
weight argument is missing.
This is the most popular type of correction, also known as the Bessel correction.
-
The constant rweight or an object of type rweight_type implying a bias correction based on the assumption of reliability weights for the sample observations.
(optional. If missing, no bias-correction will be applied to the output var .) |
- Returns
var
: The output object of type real
of the same kind as the input sample
of rank rank(sample) - 1
, representing its variance:
-
If
sample
is a vector, the output var
is a scalar.
-
If
sample
is a matrix, the output var
is an allocatable
vector of size ndim = size(sample, dim)
.
Possible calling interfaces ⛓
var
= getVar(sample(:), weight
= weight(
1 :
size(sample)), correction
= correction)
var
= getVar(sample(:,:), weight
= weight(
1 :
size(sample)), correction
= correction)
var
= getVar(sample(:), weight
= weight(
1 :
size(sample, dim)), correction
= correction)
var(:)
= getVar(sample(:,:), dim, weight
= weight(
1 :
size(sample, dim)), correction
= correction)
Generate and return the variance of the input sample of type complex or real of shape (nsam) or (ndim...
This module contains classes and procedures for computing the properties related to the covariance ma...
- Warning
- This generic interface is merely a flexible wrapper around the generic
subroutine
interface setVar.
As such, all conditions and warnings associated with setVar equally hold for this generic interface.
-
The
pure
procedure(s) documented herein become impure
when the ParaMonte library is compiled with preprocessor macro CHECK_ENABLED=1
.
By default, these procedures are pure
in release
build and impure
in debug
and testing
builds.
- Note
- Note the effects of bias-correction in computing the variance become noticeable only for very sample sample sizes (i.e., when
nsam
is small).
-
For a one-dimensional
sample
, one can also use the concise Fortran syntax to achieve the same goal as with the interface var = getVar(sample(:), mean, correction)
with integer weight
:
mean = sum(sample) / size(sample)
var = sum((sample-mean)**2) / (size(sample) - 1)
But the above concise version will be likely slightly slower as it potentially involves more loops.
-
For a two or higher-dimensional
sample
, if the variance is to be computed for the entire sample
(as opposed to computing it along a particular dimension), simply pass reshape(sample, shape = size(sample))
to the appropriate getVar interface.
Alternatively, a 1D pointer of the same size as the multidimensional sample can be passed to the procedure.
- See also
- getMean
setMean
getShifted
setShifted
getVar
setVar
Example usage ⛓
18 type(display_type) :: disp
22 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%show(
"!Compute the variance of a 1-D array.")
24 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
29 real(TKG),
allocatable :: sample(:)
30 call disp%show(
"sample = getLinSpace(1._TKG, 9._TKG, 5_IK)")
34 call disp%show(
"var = getVar(sample)")
38 call disp%show(
"var = getVar(sample, correction = fweight)")
39 var
= getVar(sample, correction
= fweight)
42 call disp%show(
"var = getVar(sample, correction = rweight)")
43 var
= getVar(sample, correction
= rweight)
46 call disp%show(
"var = getVar(sample, dim = 1_IK)")
47 var
= getVar(sample,
dim = 1_IK)
50 call disp%show(
"var = getVar(sample, dim = 1_IK, correction = fweight)")
51 var
= getVar(sample,
dim = 1_IK, correction
= fweight)
54 call disp%show(
"var = getVar(sample, dim = 1_IK, correction = rweight)")
55 var
= getVar(sample,
dim = 1_IK, correction
= rweight)
63 complex(TKG) ,
allocatable :: sample(:)
64 call disp%show(
"sample = cmplx(getShuffled(getLinSpace(1., 9., 5_IK)), getShuffled(getLinSpace(1., 9., 5_IK)), TKG)")
68 call disp%show(
"var = getVar(sample)")
72 call disp%show(
"var = getVar(sample, correction = fweight)")
73 var
= getVar(sample, correction
= fweight)
76 call disp%show(
"var = getVar(sample, correction = rweight)")
77 var
= getVar(sample, correction
= rweight)
80 call disp%show(
"var = getVar(sample, dim = 1_IK)")
81 var
= getVar(sample,
dim = 1_IK)
84 call disp%show(
"var = getVar(sample, dim = 1_IK, correction = fweight)")
85 var
= getVar(sample,
dim = 1_IK, correction
= fweight)
88 call disp%show(
"var = getVar(sample, dim = 1_IK, correction = rweight)")
89 var
= getVar(sample,
dim = 1_IK, correction
= rweight)
96 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
97 call disp%show(
"!Compute the variance of a 2-D array.")
98 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
102 real(TKG),
allocatable :: var(:), sample(:,:)
104 call disp%show(
"sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
108 call disp%show(
"call setResized(var, 1_IK)")
110 call disp%show(
"var(1) = getVar(sample)")
114 call disp%show(
"var(1) = getVar(sample, correction = fweight)")
115 var(
1)
= getVar(sample, correction
= fweight)
118 call disp%show(
"var(1) = getVar(sample, correction = rweight)")
119 var(
1)
= getVar(sample, correction
= rweight)
124 call disp%show(
"dim ! The observations axis.")
126 call disp%show(
"sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
128 call disp%show(
"var = getVar(sample, dim, correction = fweight)")
129 var
= getVar(sample, dim, correction
= fweight)
132 call disp%show(
"var = getVar(sample, dim, correction = rweight)")
133 var
= getVar(sample, dim, correction
= rweight)
141 real(TKG),
allocatable :: var(:)
142 complex(TKG),
allocatable :: sample(:,:)
144 call disp%show(
"sample = cmplx(getUnifRand(-9, 9, 4_IK, 5_IK), -getUnifRand(-9, 9, 4_IK, 5_IK), TKG)")
148 call disp%show(
"call setResized(var, 1_IK)")
150 call disp%show(
"var(1) = getVar(sample)")
154 call disp%show(
"var(1) = getVar(sample, correction = fweight)")
155 var(
1)
= getVar(sample, correction
= fweight)
158 call disp%show(
"var(1) = getVar(sample, correction = rweight)")
159 var(
1)
= getVar(sample, correction
= rweight)
164 call disp%show(
"dim ! The observations axis.")
166 call disp%show(
"sample = cmplx(getUnifRand(-9, 9, 4_IK, 5_IK), -getUnifRand(-9, 9, 4_IK, 5_IK), TKG)")
170 call disp%show(
"var = getVar(sample, dim)")
174 call disp%show(
"var = getVar(sample, dim, correction = fweight)")
175 var
= getVar(sample, dim, correction
= fweight)
178 call disp%show(
"var = getVar(sample, dim, correction = rweight)")
179 var
= getVar(sample, dim, correction
= rweight)
187 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
188 call disp%show(
"!Compute the variance of a 1-D weighted array.")
189 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
194 real(TKG),
allocatable :: sample(:)
195 real(TKG) ,
allocatable :: weight(:)
196 call disp%show(
"sample = getLinSpace(1._TKG, 9._TKG, 5_IK)")
200 call disp%show(
"weight = getLinSpace(1._TKG, 9._TKG, size(sample, 1, IK))")
204 call disp%show(
"var = getVar(sample, weight = weight)")
205 var
= getVar(sample, weight
= weight)
208 call disp%show(
"var = getVar(sample, weight = weight, correction = fweight)")
209 var
= getVar(sample, weight
= weight, correction
= fweight)
212 call disp%show(
"var = getVar(sample, weight = weight, correction = rweight)")
213 var
= getVar(sample, weight
= weight, correction
= rweight)
216 call disp%show(
"var = getVar(sample, dim = 1_IK, weight = weight)")
217 var
= getVar(sample,
dim = 1_IK, weight
= weight)
220 call disp%show(
"var = getVar(sample, dim = 1_IK, weight = weight, correction = fweight)")
221 var
= getVar(sample,
dim = 1_IK, weight
= weight, correction
= fweight)
224 call disp%show(
"var = getVar(sample, dim = 1_IK, weight = weight, correction = rweight)")
225 var
= getVar(sample,
dim = 1_IK, weight
= weight, correction
= rweight)
232 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
233 call disp%show(
"!Compute the variance of a 2-D weighted array.")
234 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
238 real(TKG),
allocatable :: var(:), sample(:,:)
239 real(TKG),
allocatable :: weight(:)
241 call disp%show(
"sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
245 call disp%show(
"weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))")
246 weight
= getLinSpace(
1._TKG,
9._TKG,
size(sample,
kind = IK))
249 call disp%show(
"call setResized(var, 1_IK)")
251 call disp%show(
"var(1) = getVar(sample, weight = weight)")
252 var(
1)
= getVar(sample, weight
= weight)
255 call disp%show(
"var(1) = getVar(sample, weight = weight, correction = fweight)")
256 var(
1)
= getVar(sample, weight
= weight, correction
= fweight)
259 call disp%show(
"var(1) = getVar(sample, weight = weight, correction = rweight)")
260 var(
1)
= getVar(sample, weight
= weight, correction
= rweight)
265 call disp%show(
"dim ! The observations axis.")
267 call disp%show(
"sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
271 call disp%show(
"weight = getLinSpace(1._TKG, 9._TKG, size(sample, dim, IK))")
275 call disp%show(
"var = getVar(sample, dim, weight = weight)")
276 var
= getVar(sample, dim, weight
= weight)
279 call disp%show(
"var = getVar(sample, dim, weight = weight, correction = fweight)")
280 var
= getVar(sample, dim, weight
= weight, correction
= fweight)
283 call disp%show(
"var = getVar(sample, dim, weight = weight, correction = rweight)")
284 var
= getVar(sample, dim, weight
= weight, correction
= rweight)
292 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
293 call disp%show(
"!Compute the variance of a multidimensional array by reshaping the array.")
294 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
298 real(TKG),
allocatable :: var
299 real(TKG),
allocatable,
target :: sample(:,:,:)
300 real(TKG),
pointer :: samptr(:)
301 call disp%show(
"sample = getUnifRand(0., 1., 3_IK, 3_IK, 3_IK)")
305 call disp%show(
"samptr(1:product(shape(sample))) => sample")
306 samptr(
1:
product(
shape(sample)))
=> sample
307 call disp%show(
"var = getVar(samptr)")
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Perform an unbiased random shuffling of the input array, known as the Knuth or Fisher-Yates shuffle,...
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
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...
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 resizing allocatable arrays of various typ...
This module contains procedures and generic interfaces for shuffling arrays of various types.
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
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 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 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.
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 ⛓
8+1.00000000,
+3.00000000,
+5.00000000,
+7.00000000,
+9.00000000
12var
= getVar(sample, correction
= fweight)
15var
= getVar(sample, correction
= rweight)
18var
= getVar(sample,
dim = 1_IK)
21var
= getVar(sample,
dim = 1_IK, correction
= fweight)
24var
= getVar(sample,
dim = 1_IK, correction
= rweight)
30(
+5.00000000,
+1.00000000), (
+9.00000000,
+7.00000000), (
+3.00000000,
+5.00000000), (
+1.00000000,
+3.00000000), (
+7.00000000,
+9.00000000)
34var
= getVar(sample, correction
= fweight)
37var
= getVar(sample, correction
= rweight)
40var
= getVar(sample,
dim = 1_IK)
43var
= getVar(sample,
dim = 1_IK, correction
= fweight)
46var
= getVar(sample,
dim = 1_IK, correction
= rweight)
58+6.52446842,
+3.83352661,
+1.77059317,
+8.38501263,
+8.95649338
59+5.00342083,
+7.77162886,
+4.97027016,
+1.27960253,
+7.32357645
60+3.90767574,
+4.08084774,
+5.77870083,
+1.82107782,
+5.74156380
61+8.22197819,
+7.56458569,
+3.83126736,
+6.85559607,
+3.62764645
66var(
1)
= getVar(sample, correction
= fweight)
69var(
1)
= getVar(sample, correction
= rweight)
76var
= getVar(sample, dim, correction
= fweight)
78+4.97601366,
+4.45858145,
+4.38417006,
+3.08473539,
+3.36357975
79var
= getVar(sample, dim, correction
= rweight)
81+4.97601366,
+4.45858145,
+4.38417006,
+3.08473539,
+3.36357975
86var
= getVar(sample, dim, correction
= fweight)
88+8.88743401,
+6.12243223,
+7.31469440,
+0.990785539
89var
= getVar(sample, dim, correction
= rweight)
91+8.88743401,
+6.12243223,
+7.31469440,
+0.990785539
96(
-8.00000000,
+6.00000000), (
-7.00000000,
+2.00000000), (
+4.00000000,
-1.00000000), (
+0.00000000,
+2.00000000), (
+2.00000000,
-3.00000000)
97(
-2.00000000,
+0.00000000), (
+9.00000000,
-4.00000000), (
+7.00000000,
+7.00000000), (
+0.00000000,
-7.00000000), (
+1.00000000,
+6.00000000)
98(
+3.00000000,
-8.00000000), (
-5.00000000,
+6.00000000), (
+8.00000000,
+4.00000000), (
+8.00000000,
-8.00000000), (
+4.00000000,
+7.00000000)
99(
-2.00000000,
+2.00000000), (
-1.00000000,
+5.00000000), (
-4.00000000,
+3.00000000), (
+1.00000000,
+6.00000000), (
+5.00000000,
+1.00000000)
104var(
1)
= getVar(sample, correction
= fweight)
107var(
1)
= getVar(sample, correction
= rweight)
115(
-1.00000000,
+3.00000000), (
+3.00000000,
+1.00000000), (
-2.00000000,
+6.00000000), (
+5.00000000,
+1.00000000), (
-8.00000000,
-7.00000000)
116(
+8.00000000,
-2.00000000), (
-3.00000000,
+8.00000000), (
+1.00000000,
+3.00000000), (
+5.00000000,
+3.00000000), (
-5.00000000,
-2.00000000)
117(
+7.00000000,
-3.00000000), (
-5.00000000,
+8.00000000), (
-2.00000000,
-9.00000000), (
-6.00000000,
-7.00000000), (
+9.00000000,
+4.00000000)
118(
-8.00000000,
-5.00000000), (
+3.00000000,
+8.00000000), (
-3.00000000,
+0.00000000), (
+6.00000000,
-3.00000000), (
+2.00000000,
+3.00000000)
121+50.9375000,
+21.9375000,
+33.7500000,
+39.0000000,
+62.5000000
122var
= getVar(sample, dim, correction
= fweight)
124+53.6184235,
+23.0921059,
+35.5263176,
+41.0526314,
+65.7894745
125var
= getVar(sample, dim, correction
= rweight)
127+53.6184235,
+23.0921059,
+35.5263176,
+41.0526314,
+65.7894745
133(
+4.00000000,
+6.00000000), (
-9.00000000,
+5.00000000), (
-1.00000000,
+0.00000000), (
-2.00000000,
-1.00000000), (
+5.00000000,
+3.00000000)
134(
+9.00000000,
+8.00000000), (
-8.00000000,
-2.00000000), (
-1.00000000,
-4.00000000), (
+8.00000000,
-9.00000000), (
+7.00000000,
+3.00000000)
135(
-4.00000000,
+0.00000000), (
+2.00000000,
+8.00000000), (
-7.00000000,
-8.00000000), (
-6.00000000,
+3.00000000), (
-5.00000000,
+7.00000000)
136(
-2.00000000,
+0.00000000), (
+3.00000000,
-3.00000000), (
+9.00000000,
-4.00000000), (
-7.00000000,
+7.00000000), (
-3.00000000,
+2.00000000)
139+32.4799995,
+76.9599991,
+43.2000008,
+45.8400002
140var
= getVar(sample, dim, correction
= fweight)
142+34.1894760,
+81.0105286,
+45.4736862,
+48.2526321
143var
= getVar(sample, dim, correction
= rweight)
145+34.1894760,
+81.0105286,
+45.4736862,
+48.2526321
154+1.00000000,
+3.00000000,
+5.00000000,
+7.00000000,
+9.00000000
157+1.00000000,
+3.00000000,
+5.00000000,
+7.00000000,
+9.00000000
158var
= getVar(sample, weight
= weight)
161var
= getVar(sample, weight
= weight, correction
= fweight)
164var
= getVar(sample, weight
= weight, correction
= rweight)
167var
= getVar(sample,
dim = 1_IK, weight
= weight)
170var
= getVar(sample,
dim = 1_IK, weight
= weight, correction
= fweight)
173var
= getVar(sample,
dim = 1_IK, weight
= weight, correction
= rweight)
185+4.92118597,
+3.93205404,
+2.04195070,
+1.08013439,
+8.42021465
186+6.33719110,
+4.58371544,
+3.42720079,
+5.51021671,
+4.78841782
187+6.35301495,
+5.05011415,
+7.05587339,
+5.93136883,
+3.38336372
188+7.58829355,
+5.07982922,
+3.44624376,
+7.93914604,
+1.88381815
191+1.00000000,
+1.42105269,
+1.84210527,
+2.26315784,
+2.68421054,
+3.10526323,
+3.52631569,
+3.94736838,
+4.36842108,
+4.78947353,
+5.21052647,
+5.63157892,
+6.05263138,
+6.47368431,
+6.89473677,
+7.31578970,
+7.73684216,
+8.15789413,
+8.57894707,
+9.00000000
193var(
1)
= getVar(sample, weight
= weight)
196var(
1)
= getVar(sample, weight
= weight, correction
= fweight)
199var(
1)
= getVar(sample, weight
= weight, correction
= rweight)
207+4.29987097,
+4.51889372,
+2.49617767,
+4.28629637,
+4.84596252
208+1.13063669,
+6.40167713,
+1.41440201,
+4.65148449,
+1.05075884
209+2.81298780,
+5.53477764,
+3.52329731,
+1.20692921,
+8.21473789
210+8.45940113,
+8.63117886,
+5.85738850,
+7.98474264,
+4.74142790
213+1.00000000,
+3.66666675,
+6.33333349,
+9.00000000
214var
= getVar(sample, dim, weight
= weight)
216+9.65500069,
+2.24928784,
+2.94472671,
+8.61085415,
+6.13436651
217var
= getVar(sample, dim, weight
= weight, correction
= fweight)
219+10.1631594,
+2.36767149,
+3.09971237,
+9.06405735,
+6.45722818
220var
= getVar(sample, dim, weight
= weight, correction
= rweight)
222+14.6042023,
+3.40228391,
+4.45420837,
+13.0248203,
+9.27887344
228+5.95939255,
+8.73706532,
+5.65549803,
+2.86419582,
+1.95827341
229+6.32004738,
+4.42536259,
+2.82742453,
+8.56602859,
+3.90013361
230+6.40244722,
+6.10307646,
+3.61640501,
+8.00912476,
+6.62316418
231+1.59151077,
+2.09931469,
+5.49416971,
+3.51770878,
+7.18769121
234+1.00000000,
+3.00000000,
+5.00000000,
+7.00000000,
+9.00000000
235var
= getVar(sample, dim, weight
= weight)
237+5.25076151,
+5.02641058,
+2.29942298,
+3.86122894
238var
= getVar(sample, dim, weight
= weight, correction
= fweight)
240+5.46954298,
+5.23584414,
+2.39523220,
+4.02211332
241var
= getVar(sample, dim, weight
= weight, correction
= rweight)
243+7.13418674,
+6.82936192,
+3.12421584,
+5.24623489
253+0.670200706,
+0.245780528,
+0.691632628E-1
254+0.741564631,
+0.370370686,
+0.645985901
255+0.567170262,
+0.154305458,
+0.208736062
257+0.954846084,
+0.512954772,
+0.562396646E-1
258+0.101873875,
+0.967078269,
+0.719123542
259+0.322870612E-1,
+0.252239347,
+0.973619521
261+0.739393175,
+0.216321886,
+0.134367526
262+0.837722719,
+0.347896039,
+0.566112757
263+0.416002929,
+0.695480704,
+0.562496781E-1
264samptr(
1:
product(
shape(sample)))
=> sample
- Test:
- test_pm_sampleVar
- Todo:
- Very Low Priority: The functionality of this interface can be expanded in the future to include the computation of the variance of higher dimensional input
sample
and whole sample
input arrays of arbitrary shape.
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
Definition at line 715 of file pm_sampleVar.F90.