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

Return the (weighted) sample variance and mean of an input time series of nsam observations, or of an input sample of nsam observations with ndim attributes optionally weighted by the input weight, optionally also sum(weight).
More...

Detailed Description

Return the (weighted) sample variance and mean of an input time series of nsam observations, or of an input sample of nsam observations with ndim attributes optionally weighted by the input weight, optionally also sum(weight).

This generic interface is developed to specifically improve the performance of other procedures where sum(weight) are also needed in addition to the mean of the weighted sample.
It uses a one-pass approach to simultaneously compute the sample mean and variance.

Parameters
[out]var: The output object of rank rank(sample) - 1 of type real of the same kind as the input sample representing its variance:
  1. If sample is a vector, the output var must be a scalar.
  2. If sample is a matrix, the output var must be a vector of size ndim = size(sample, 3 - dim).
[out]mean: The output object of the same type, kind, and rank as the output var containing the sample mean.
[in]sample: The input contiguous array of shape (nsam), (ndim, nsam) or (nsam, ndim) of the same type and kind as the output var, containing the sample whose mean is to be computed.
[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.
  1. If dim = 1, the input sample of rank 2 is assumed to have the shape (nsam, ndim).
  2. 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 input contiguous vector of length nsam of,
  1. type integer of default kind IK, or
  2. type real of the same kind as the input sample,
containing the corresponding weights of the data points in the input sample.
[out]weisum: The output scalar of the same type and kind as the input weight, containing sum(weight).
(optional. It must be present if and only if the input argument weight is also present.)
[in]meang: The input object of the same type, kind, and rank as the output var containing the best guess for the sample mean.
If no good guess is known a priori, meang can be set to any observation in sample.
See the example usage below.


Possible calling interfaces

! 1D sample
call setVarMean(var, mean, sample(1 : nsam), meang)
call setVarMean(var, mean, sample(1 : nsam), weight(1 : nsam), weisum, meang)
call setVarMean(var, mean, sample(1 : nsam), dim, meang)
call setVarMean(var, mean, sample(1 : nsam), dim, weight(1 : nsam), weisum, meang)
! 2D sample
call setVarMean(var, mean(1 : ndim), sample(:,:), meang(1 : ndim))
call setVarMean(var, mean(1 : ndim), sample(:,:), weight(1 : size(sample)), weisum, meang(1 : ndim))
call setVarMean(var, mean(1 : ndim), sample(:,:), dim, meang(1 : ndim))
call setVarMean(var, mean(1 : ndim), sample(:,:), dim, weight(1 : size(sample, dim)), weisum, meang(1 : ndim))
Return the (weighted) sample variance and mean of an input time series of nsam observations,...
This module contains classes and procedures for computing the properties related to the covariance ma...
Warning
The condition all(0. <= weight) must hold for the corresponding input arguments.
The condition 1 <= dim .and. dim <= rank(sample) must hold for the corresponding input arguments.
The condition size(sample, dim) == size(weight, 1) must hold for the corresponding input arguments.
The condition size(sample, 3 - dim) == size(var, 1) must hold for the corresponding input arguments.
The condition size(sample, 3 - dim) == size(mean, 1) must hold for the corresponding input arguments.
The condition size(sample, 3 - dim) == size(meang, 1) must hold for the corresponding input arguments.
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.
Note
If the input sample is to be an array of type integer, simply convert the sample to an array of type real of the desired kind for the output real mean of the sample.
There is no point in accepting an input sample of type integer since it will be inevitably converted to an array of type real within the procedure to avoid potential integer overflow.
Furthermore, an sample of type integer creates ambiguity about the kind of the real-valued returned mean by the procedure.
See the notes in the description of the pm_sampleMean.
Note that the mean of any one or two-dimensional sample can be simply computed via the Fortran intrinsic routine sum():
integer(IK) :: i
integer(IK) , parameter :: NDIM = 3_IK
integer(IK) , parameter :: NSAM = 1000_IK
real(TKG) , parameter :: sample(NDIM,NSAM) = reshape([( real(i,RK), i = 1, NSAM )], shape = shape(sample))
real(TKG) , allocatable :: mean(:)
mean = sum(sample, dim = 1) / size(transpose(sample), dim = 1) ! assuming the first dimension represents the observations
mean = sum(sample, dim = 2) / size(sample, dim = 2) ! assuming the second dimension represents the observations
The mean of a whole multidimensional array can be obtained by either,
  1. reshaping the array to a vector form and passing it to this procedure, or
  2. mapping the array to a 1-dimensional pointer array of the same size as the ndim dimensional array.
See the examples below.
See also
getVar
setVar
getCov
setCov
getMean
setMean
setCovMean


Example usage

1program example
2
3 use pm_kind, only: SK, IK, LK, RK
4 use pm_kind, only: TKG => RKS ! All other real types are also supported.
5 use pm_io, only: display_type
6 use pm_sampleVar, only: getVar
7 use pm_sampleMean, only: getMean
8 use pm_sampleVar, only: setVarMean
10 use pm_sampleCov, only: fweight, rweight
11 use pm_sampleShift, only: getShifted
12 use pm_arrayResize, only: setResized
13 use pm_distUnif, only: getUnifRand
14
15 implicit none
16
17 integer(IK) :: dim
18 type(display_type) :: disp
19 disp = display_type(file = "main.out.F90")
20
21 call disp%skip()
22 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
23 call disp%show("!Compute the variance of a 1-D sample.")
24 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
25 call disp%skip()
26
27 block
28 real(TKG) :: var, mean
29 real(TKG), allocatable :: sample(:)
30 call disp%show("sample = getLinSpace(1._TKG, 9._TKG, 5_IK)")
31 sample = getLinSpace(1._TKG, 9._TKG, 5_IK)
32 call disp%show("sample")
33 call disp%show( sample )
34 call disp%show("call setVarMean(var, mean, sample, meang = sample(1))")
35 call setVarMean(var, mean, sample, meang = sample(1))
36 call disp%show("mean")
37 call disp%show( mean )
38 call disp%show("getMean(sample) ! for comparison.")
39 call disp%show( getMean(sample) )
40 call disp%show("var")
41 call disp%show( var )
42 call disp%show("getVar(sample) ! for comparison.")
43 call disp%show( getVar(sample) )
44 call disp%show("call setVarMean(var, mean, sample, dim = 1_IK, meang = sample(1))")
45 call setVarMean(var, mean, sample, dim = 1_IK, meang = sample(1))
46 call disp%show("mean")
47 call disp%show( mean )
48 call disp%show("getMean(sample, dim = 1_IK) ! for comparison.")
49 call disp%show( getMean(sample, dim = 1_IK) )
50 call disp%show("var")
51 call disp%show( var )
52 call disp%show("getVar(sample, dim = 1_IK) ! for comparison.")
53 call disp%show( getVar(sample, dim = 1_IK) )
54 call disp%skip()
55 end block
56
57 block
58 real(TKG) :: var
59 complex(TKG) :: mean
60 complex(TKG), allocatable :: sample(:)
61 call disp%show("sample = cmplx(getLinSpace(1., 9., 5_IK), -getLinSpace(1., 9., 5_IK), TKG)")
62 sample = cmplx(getLinSpace(1., 9., 5_IK), -getLinSpace(1., 9., 5_IK), TKG)
63 call disp%show("sample")
64 call disp%show( sample )
65 call disp%show("call setVarMean(var, mean, sample, meang = sample(1))")
66 call setVarMean(var, mean, sample, meang = sample(1))
67 call disp%show("mean")
68 call disp%show( mean )
69 call disp%show("getMean(sample) ! for comparison.")
70 call disp%show( getMean(sample) )
71 call disp%show("var")
72 call disp%show( var )
73 call disp%show("getVar(sample) ! for comparison.")
74 call disp%show( getVar(sample) )
75 call disp%show("call setVarMean(var, mean, sample, dim = 1_IK, meang = sample(1))")
76 call setVarMean(var, mean, sample, dim = 1_IK, meang = sample(1))
77 call disp%show("mean")
78 call disp%show( mean )
79 call disp%show("getMean(sample, dim = 1_IK) ! for comparison.")
80 call disp%show( getMean(sample, dim = 1_IK) )
81 call disp%show("var")
82 call disp%show( var )
83 call disp%show("getVar(sample, dim = 1_IK) ! for comparison.")
84 call disp%show( getVar(sample, dim = 1_IK) )
85 call disp%skip()
86 end block
87
88 call disp%skip()
89 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
90 call disp%show("!Compute the variance of a 1-D weighted sample.")
91 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
92 call disp%skip()
93
94 block
95 real(TKG) :: var, weisum
96 real(TKG), allocatable :: weight(:)
97 real(TKG), allocatable :: sample(:)
98 real(TKG) :: mean
99 call disp%show("sample = getLinSpace(1._TKG, 9._TKG, 5_IK)")
100 sample = getLinSpace(1._TKG, 9._TKG, 5_IK)
101 call disp%show("sample")
102 call disp%show( sample )
103 call disp%show("weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))")
104 weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))
105 call disp%show("weight")
106 call disp%show( weight )
107 call disp%show("call setVarMean(var, mean, sample, weight, weisum, meang = sample(1))")
108 call setVarMean(var, mean, sample, weight, weisum, meang = sample(1))
109 call disp%show("mean")
110 call disp%show( mean )
111 call disp%show("getMean(sample, weight) ! for comparison.")
112 call disp%show( getMean(sample, weight) )
113 call disp%show("var")
114 call disp%show( var )
115 call disp%show("getVar(sample, weight) ! for comparison.")
116 call disp%show( getVar(sample, weight) )
117 call disp%show("call setVarMean(var, mean, sample, 1_IK, weight, weisum, meang = sample(1))")
118 call setVarMean(var, mean, sample, 1_IK, weight, weisum, meang = sample(1))
119 call disp%show("mean")
120 call disp%show( mean )
121 call disp%show("getMean(sample, 1_IK, weight) ! for comparison.")
122 call disp%show( getMean(sample, 1_IK, weight) )
123 call disp%show("var")
124 call disp%show( var )
125 call disp%show("getVar(sample, 1_IK, weight) ! for comparison.")
126 call disp%show( getVar(sample, 1_IK, weight) )
127 call disp%skip()
128 end block
129
130 block
131 real(TKG) :: var, weisum
132 real(TKG), allocatable :: weight(:)
133 complex(TKG), allocatable :: sample(:)
134 complex(TKG) :: mean
135 call disp%show("sample = cmplx(getLinSpace(1., 9., 5_IK), -getLinSpace(1., 9., 5_IK), TKG)")
136 sample = cmplx(getLinSpace(1., 9., 5_IK), -getLinSpace(1., 9., 5_IK), TKG)
137 call disp%show("sample")
138 call disp%show( sample )
139 call disp%show("weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))")
140 weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))
141 call disp%show("weight")
142 call disp%show( weight )
143 call disp%show("call setVarMean(var, mean, sample, weight, weisum, meang = sample(1))")
144 call setVarMean(var, mean, sample, weight, weisum, meang = sample(1))
145 call disp%show("mean")
146 call disp%show( mean )
147 call disp%show("getMean(sample, weight) ! for comparison.")
148 call disp%show( getMean(sample, weight) )
149 call disp%show("var")
150 call disp%show( var )
151 call disp%show("getVar(sample, weight) ! for comparison.")
152 call disp%show( getVar(sample, weight) )
153 call disp%show("call setVarMean(var, mean, sample, 1_IK, weight, weisum, meang = sample(1))")
154 call setVarMean(var, mean, sample, 1_IK, weight, weisum, meang = sample(1))
155 call disp%show("mean")
156 call disp%show( mean )
157 call disp%show("getMean(sample, 1_IK, weight) ! for comparison.")
158 call disp%show( getMean(sample, 1_IK, weight) )
159 call disp%show("var")
160 call disp%show( var )
161 call disp%show("getVar(sample, 1_IK, weight) ! for comparison.")
162 call disp%show( getVar(sample, 1_IK, weight) )
163 call disp%skip()
164 end block
165
166 call disp%skip()
167 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
168 call disp%show("!Compute the variance of a 2-D array.")
169 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
170 call disp%skip()
171
172 block
173 real(TKG), allocatable :: var(:)
174 real(TKG), allocatable :: meang(:), mean(:), sample(:,:)
175 call disp%skip()
176 call disp%show("sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
177 sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
178 call disp%show("sample")
179 call disp%show( sample )
180 call disp%show("call setResized(var, 1_IK)")
181 call setResized(var, 1_IK)
182 call disp%show("call setResized(mean, 1_IK)")
183 call setResized(mean, 1_IK)
184 call disp%show("call setVarMean(var(1), mean(1), sample, meang = sample(1,1))")
185 call setVarMean(var(1), mean(1), sample, meang = sample(1,1))
186 call disp%show("mean")
187 call disp%show( mean )
188 call disp%show("getMean(sample) ! for comparison.")
189 call disp%show( getMean(sample) )
190 call disp%show("var")
191 call disp%show( var )
192 call disp%show("getVar(sample) ! for comparison.")
193 call disp%show( getVar(sample) )
194 call disp%skip()
195 do dim = 1, 2
196 call disp%show("dim ! The observations axis.")
197 call disp%show( dim )
198 call disp%show("sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
199 sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
200 call disp%show("sample")
201 call disp%show( sample )
202 call disp%show("call setResized(var, size(sample, 3 - dim, IK))")
203 call setResized(var, size(sample, 3 - dim, IK))
204 call disp%show("call setResized(mean, size(sample, 3 - dim, IK))")
205 call setResized(mean, size(sample, 3 - dim, IK))
206 call disp%show("if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if")
207 if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
208 call disp%show("call setVarMean(var, mean, sample, dim, meang)")
209 call setVarMean(var, mean, sample, dim, meang)
210 call disp%show("mean")
211 call disp%show( mean )
212 call disp%show("getMean(sample, dim) ! for comparison.")
213 call disp%show( getMean(sample, dim) )
214 call disp%show("var")
215 call disp%show( var )
216 call disp%show("getVar(sample, dim) ! for comparison.")
217 call disp%show( getVar(sample, dim) )
218 call disp%skip()
219 end do
220 end block
221
222 block
223 real(TKG), allocatable :: var(:)
224 complex(TKG), allocatable :: meang(:), mean(:), sample(:,:)
225 call disp%skip()
226 call disp%show("sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)")
227 sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)
228 call disp%show("sample")
229 call disp%show( sample )
230 call disp%show("call setResized(var, 1_IK)")
231 call setResized(var, 1_IK)
232 call disp%show("call setResized(mean, 1_IK)")
233 call setResized(mean, 1_IK)
234 call disp%show("call setVarMean(var(1), mean(1), sample, meang = sample(1,1))")
235 call setVarMean(var(1), mean(1), sample, meang = sample(1,1))
236 call disp%show("mean")
237 call disp%show( mean )
238 call disp%show("getMean(sample) ! for comparison.")
239 call disp%show( getMean(sample) )
240 call disp%show("var")
241 call disp%show( var )
242 call disp%show("getVar(sample) ! for comparison.")
243 call disp%show( getVar(sample) )
244 call disp%skip()
245 do dim = 1, 2
246 call disp%show("dim ! The observations axis.")
247 call disp%show( dim )
248 call disp%show("sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)")
249 sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)
250 call disp%show("sample")
251 call disp%show( sample )
252 call disp%show("call setResized(var, size(sample, 3 - dim, IK))")
253 call setResized(var, size(sample, 3 - dim, IK))
254 call disp%show("call setResized(mean, size(sample, 3 - dim, IK))")
255 call setResized(mean, size(sample, 3 - dim, IK))
256 call disp%show("if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if")
257 if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
258 call disp%show("call setVarMean(var, mean, sample, dim, meang)")
259 call setVarMean(var, mean, sample, dim, meang)
260 call disp%show("mean")
261 call disp%show( mean )
262 call disp%show("getMean(sample, dim) ! for comparison.")
263 call disp%show( getMean(sample, dim) )
264 call disp%show("var")
265 call disp%show( var )
266 call disp%show("getVar(sample, dim) ! for comparison.")
267 call disp%show( getVar(sample, dim) )
268 call disp%skip()
269 end do
270 end block
271
272 call disp%skip()
273 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
274 call disp%show("!Compute the variance of a 2-D weighted sample.")
275 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
276 call disp%skip()
277
278 block
279 real(TKG):: rweisum
280 real(TKG), allocatable :: var(:)
281 real(TKG), allocatable :: rweight(:)
282 real(TKG), allocatable :: meang(:), mean(:), sample(:,:)
283 integer(IK), allocatable :: iweight(:)
284 integer(IK) :: iweisum
285 call disp%skip()
286 call disp%show("sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
287 sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
288 call disp%show("sample")
289 call disp%show( sample )
290 call disp%show("iweight = getUnifRand(1, 9, size(sample, kind = IK))")
291 iweight = getUnifRand(1, 9, size(sample, kind = IK))
292 call disp%show("iweight")
293 call disp%show( iweight )
294 call disp%show("rweight = iweight")
295 rweight = iweight
296 call disp%show("call setResized(var, 1_IK)")
297 call setResized(var, 1_IK)
298 call disp%show("call setResized(mean, 1_IK)")
299 call setResized(mean, 1_IK)
300 call disp%show("call setVarMean(var(1), mean(1), sample, iweight, iweisum, meang = sample(1,1))")
301 call setVarMean(var(1), mean(1), sample, iweight, iweisum, meang = sample(1,1))
302 call disp%show("mean")
303 call disp%show( mean )
304 call disp%show("getMean(sample, iweight) ! for comparison.")
305 call disp%show( getMean(sample, iweight) )
306 call disp%show("var")
307 call disp%show( var )
308 call disp%show("getVar(sample, iweight) ! for comparison.")
309 call disp%show( getVar(sample, iweight) )
310 call disp%show("[iweisum, sum(iweight)] ! for comparison.")
311 call disp%show( [iweisum, sum(iweight)] )
312 call disp%show("call setVarMean(var(1), mean(1), sample, rweight, rweisum, meang = sample(1,1))")
313 call setVarMean(var(1), mean(1), sample, rweight, rweisum, meang = sample(1,1))
314 call disp%show("mean")
315 call disp%show( mean )
316 call disp%show("getMean(sample, rweight) ! for comparison.")
317 call disp%show( getMean(sample, rweight) )
318 call disp%show("var")
319 call disp%show( var )
320 call disp%show("getVar(sample, rweight) ! for comparison.")
321 call disp%show( getVar(sample, rweight) )
322 call disp%show("[rweisum, sum(rweight)] ! for comparison.")
323 call disp%show( [rweisum, sum(rweight)] )
324 call disp%skip()
325 do dim = 1, 2
326 call disp%show("dim ! The observations axis.")
327 call disp%show( dim )
328 call disp%show("sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)")
329 sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
330 call disp%show("sample")
331 call disp%show( sample )
332 call disp%show("iweight = getUnifRand(1, 9, size(sample, dim, IK))")
333 iweight = getUnifRand(1, 9, size(sample, dim, IK))
334 call disp%show("iweight")
335 call disp%show( iweight )
336 call disp%show("call setResized(var, size(sample, 3 - dim, IK))")
337 call setResized(var, size(sample, 3 - dim, IK))
338 call disp%show("call setResized(mean, size(sample, 3 - dim, IK))")
339 call setResized(mean, size(sample, 3 - dim, IK))
340 call disp%show("if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if")
341 if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
342 call disp%show("call setVarMean(var, mean, sample, dim, iweight, iweisum, meang)")
343 call setVarMean(var, mean, sample, dim, iweight, iweisum, meang)
344 call disp%show("mean")
345 call disp%show( mean )
346 call disp%show("getMean(sample, dim, iweight) ! for comparison.")
347 call disp%show( getMean(sample, dim, iweight) )
348 call disp%show("var")
349 call disp%show( var )
350 call disp%show("getVar(sample, dim, iweight) ! for comparison.")
351 call disp%show( getVar(sample, dim, iweight) )
352 call disp%show("[iweisum, sum(iweight)] ! for comparison.")
353 call disp%show( [iweisum, sum(iweight)] )
354 call disp%skip()
355 end do
356 end block
357
358end program example
Allocate or resize (shrink or expand) an input allocatable scalar string or array of rank 1....
Generate count evenly spaced points over the interval [x1, x2] if x1 < x2, or [x2,...
Generate and return a 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.
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 (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...
Generate and return the variance of the input sample of type complex or real of shape (nsam) or (ndim...
This module contains procedures and generic interfaces for resizing allocatable arrays of various typ...
This module contains procedures and generic interfaces for generating arrays with linear or logarithm...
This module contains classes and procedures for 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 RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
Definition: pm_kind.F90:543
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 variance of a 1-D sample.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6sample = getLinSpace(1._TKG, 9._TKG, 5_IK)
7sample
8+1.00000000, +3.00000000, +5.00000000, +7.00000000, +9.00000000
9call setVarMean(var, mean, sample, meang = sample(1))
10mean
11+5.00000000
12getMean(sample) ! for comparison.
13+5.00000000
14var
15+8.00000000
16getVar(sample) ! for comparison.
17+8.00000000
18call setVarMean(var, mean, sample, dim = 1_IK, meang = sample(1))
19mean
20+5.00000000
21getMean(sample, dim = 1_IK) ! for comparison.
22+5.00000000
23var
24+8.00000000
25getVar(sample, dim = 1_IK) ! for comparison.
26+8.00000000
27
28sample = cmplx(getLinSpace(1., 9., 5_IK), -getLinSpace(1., 9., 5_IK), TKG)
29sample
30(+1.00000000, -1.00000000), (+3.00000000, -3.00000000), (+5.00000000, -5.00000000), (+7.00000000, -7.00000000), (+9.00000000, -9.00000000)
31call setVarMean(var, mean, sample, meang = sample(1))
32mean
33(+5.00000000, -5.00000000)
34getMean(sample) ! for comparison.
35(+5.00000000, -5.00000000)
36var
37+16.0000000
38getVar(sample) ! for comparison.
39+16.0000000
40call setVarMean(var, mean, sample, dim = 1_IK, meang = sample(1))
41mean
42(+5.00000000, -5.00000000)
43getMean(sample, dim = 1_IK) ! for comparison.
44(+5.00000000, -5.00000000)
45var
46+16.0000000
47getVar(sample, dim = 1_IK) ! for comparison.
48+16.0000000
49
50
51!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
52!Compute the variance of a 1-D weighted sample.
53!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
54
55sample = getLinSpace(1._TKG, 9._TKG, 5_IK)
56sample
57+1.00000000, +3.00000000, +5.00000000, +7.00000000, +9.00000000
58weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))
59weight
60+1.00000000, +3.00000000, +5.00000000, +7.00000000, +9.00000000
61call setVarMean(var, mean, sample, weight, weisum, meang = sample(1))
62mean
63+6.59999990
64getMean(sample, weight) ! for comparison.
65+6.59999990
66var
67+5.44000244
68getVar(sample, weight) ! for comparison.
69+5.44000006
70call setVarMean(var, mean, sample, 1_IK, weight, weisum, meang = sample(1))
71mean
72+6.59999990
73getMean(sample, 1_IK, weight) ! for comparison.
74+6.59999990
75var
76+5.44000244
77getVar(sample, 1_IK, weight) ! for comparison.
78+5.44000006
79
80sample = cmplx(getLinSpace(1., 9., 5_IK), -getLinSpace(1., 9., 5_IK), TKG)
81sample
82(+1.00000000, -1.00000000), (+3.00000000, -3.00000000), (+5.00000000, -5.00000000), (+7.00000000, -7.00000000), (+9.00000000, -9.00000000)
83weight = getLinSpace(1._TKG, 9._TKG, size(sample, kind = IK))
84weight
85+1.00000000, +3.00000000, +5.00000000, +7.00000000, +9.00000000
86call setVarMean(var, mean, sample, weight, weisum, meang = sample(1))
87mean
88(+6.59999990, -6.59999990)
89getMean(sample, weight) ! for comparison.
90(+6.59999990, -6.59999990)
91var
92+10.8800049
93getVar(sample, weight) ! for comparison.
94+10.8800001
95call setVarMean(var, mean, sample, 1_IK, weight, weisum, meang = sample(1))
96mean
97(+6.59999990, -6.59999990)
98getMean(sample, 1_IK, weight) ! for comparison.
99(+6.59999990, -6.59999990)
100var
101+10.8800049
102getVar(sample, 1_IK, weight) ! for comparison.
103+10.8800001
104
105
106!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
107!Compute the variance of a 2-D array.
108!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109
110
111sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
112sample
113+2.71549416, +8.70144176, +5.82819033, +1.19167852, +8.32015133
114+7.69862223, +5.48890543, +3.90452385, +5.66230059, +3.40289307
115+2.02779531, +1.91145182, +4.23033762, +8.16210365, +1.42035627
116+4.06375647, +2.13256884, +2.12689304, +6.76190281, +5.46591902
117call setResized(var, 1_IK)
118call setResized(mean, 1_IK)
119call setVarMean(var(1), mean(1), sample, meang = sample(1,1))
120mean
121+4.56086397
122getMean(sample) ! for comparison.
123+4.56086445
124var
125+5.77198172
126getVar(sample) ! for comparison.
127+5.77198029
128
129dim ! The observations axis.
130+1
131sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
132sample
133+7.82783556, +6.55134678, +8.08672714, +1.38690710, +4.06228542
134+1.49205065, +1.23432016, +3.60968208, +7.44220686, +2.50107288
135+3.85657740, +6.44197893, +3.09694004, +8.52971458, +6.62574482
136+1.98181057, +3.74114513, +5.61530495, +1.97345972, +3.08116484
137call setResized(var, size(sample, 3 - dim, IK))
138call setResized(mean, size(sample, 3 - dim, IK))
139if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
140call setVarMean(var, mean, sample, dim, meang)
141mean
142+3.78956842, +4.49219799, +5.10216331, +4.83307171, +4.06756687
143getMean(sample, dim) ! for comparison.
144+3.78956866, +4.49219799, +5.10216331, +4.83307219, +4.06756687
145var
146+6.21466446, +4.80489683, +3.85484028, +10.1315489, +2.49279881
147getVar(sample, dim) ! for comparison.
148+6.21466780, +4.80489683, +3.85483909, +10.1315460, +2.49279857
149
150dim ! The observations axis.
151+2
152sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
153sample
154+4.86519575, +5.17080021, +6.07770634, +4.70440483, +3.72420311
155+8.93040085, +7.03005695, +4.77905560, +5.46933937, +7.68396139
156+2.40399837, +8.78949642, +8.46528149, +3.81796646, +2.42425966
157+2.52637243, +6.82366037, +6.51675224, +2.33033037, +2.08296776
158call setResized(var, size(sample, 3 - dim, IK))
159call setResized(mean, size(sample, 3 - dim, IK))
160if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
161call setVarMean(var, mean, sample, dim, meang)
162mean
163+4.90846205, +6.77856255, +5.18020058, +4.05601692
164getMean(sample, dim) ! for comparison.
165+4.90846205, +6.77856302, +5.18020010, +4.05601645
166var
167+0.576386869, +2.24509978, +8.19539452, +4.58516026
168getVar(sample, dim) ! for comparison.
169+0.576386869, +2.24509954, +8.19539261, +4.58516026
170
171
172sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)
173sample
174(+8.19689751, -6.26212835), (+6.32065010, -7.55263567), (+4.96137047, -2.30146790), (+2.67209959, -7.61132574), (+5.45345354, -5.85012150)
175(+2.02113676, -1.36360502), (+1.95255613, -4.96722412), (+5.40246725, -3.70249081), (+3.66531849, -1.05352640), (+2.48270464, -5.66561127)
176(+6.91553974, -8.15740585), (+5.44673347, -2.84762049), (+2.15664864, -5.02132463), (+4.30838060, -3.00564957), (+8.53755379, -4.40673399)
177(+7.20681238, -3.82577229), (+2.12556839, -1.73339558), (+4.46287584, -8.44206905), (+4.57033920, -3.33061457), (+3.00900078, -7.71163940)
178call setResized(var, 1_IK)
179call setResized(mean, 1_IK)
180call setVarMean(var(1), mean(1), sample, meang = sample(1,1))
181mean
182(+4.59340572, -4.74061823)
183getMean(sample) ! for comparison.
184(+4.59340525, -4.74061775)
185var
186+9.43611431
187getVar(sample) ! for comparison.
188+9.43611622
189
190dim ! The observations axis.
191+1
192sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)
193sample
194(+8.24742699, -4.64505148), (+1.69802570, -5.62454796), (+4.95255470, -8.15800667), (+3.25329256, -1.54981565), (+6.34921980, -7.00323200)
195(+4.49966097, -1.22723103), (+6.68869209, -3.12893295), (+5.38118982, -6.65558290), (+1.00942326, -7.64161825), (+6.76740026, -8.24177361)
196(+6.53018951, -5.07036495), (+3.41584492, -8.42527962), (+1.08111382, -7.82675076), (+1.53410196, -4.53349733), (+5.77184248, -2.81533957)
197(+8.77171135, -1.20631790), (+1.57522583, -1.31798840), (+6.65175056, -4.27846575), (+3.71951199, -1.02973318), (+2.39030981, -2.11225986)
198call setResized(var, size(sample, 3 - dim, IK))
199call setResized(mean, size(sample, 3 - dim, IK))
200if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
201call setVarMean(var, mean, sample, dim, meang)
202mean
203(+7.01224709, -3.03724146), (+3.34444714, -4.62418747), (+4.51665211, -6.72970152), (+2.37908244, -3.68866611), (+5.31969309, -5.04315138)
204getMean(sample, dim) ! for comparison.
205(+7.01224709, -3.03724146), (+3.34444714, -4.62418747), (+4.51665211, -6.72970104), (+2.37908244, -3.68866634), (+5.31969309, -5.04315090)
206var
207+6.12847805, +11.4114208, +6.63915730, +8.28378201, +9.89196777
208getVar(sample, dim) ! for comparison.
209+6.12847900, +11.4114199, +6.63915682, +8.28378201, +9.89196968
210
211dim ! The observations axis.
212+2
213sample = cmplx(getUnifRand(1., 9., 4_IK, 5_IK), -getUnifRand(1., 9., 4_IK, 5_IK), TKG)
214sample
215(+8.91586399, -3.28891897), (+5.00160074, -7.15572786), (+1.14957523, -2.88127565), (+4.83427382, -2.08714008), (+3.41051817, -3.86440372)
216(+1.30766630, -7.21664619), (+7.59619904, -3.38078403), (+7.66982651, -7.94244242), (+3.81181669, -5.91459703), (+6.09494305, -3.96149969)
217(+4.46235991, -5.73899221), (+4.33346128, -7.77805328), (+8.63883686, -7.79576492), (+6.03275871, -2.11102772), (+1.02897978, -4.26668024)
218(+2.71071386, -3.30129671), (+1.61056662, -1.09335375), (+7.74437094, -8.18033218), (+6.72093534, -8.79166794), (+1.56245613, -8.96508217)
219call setResized(var, size(sample, 3 - dim, IK))
220call setResized(mean, size(sample, 3 - dim, IK))
221if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
222call setVarMean(var, mean, sample, dim, meang)
223mean
224(+4.66236639, -3.85549331), (+5.29609013, -5.68319368), (+4.89927912, -5.53810358), (+4.06980848, -6.06634617)
225getMean(sample, dim) ! for comparison.
226(+4.66236639, -3.85549331), (+5.29609060, -5.68319416), (+4.89927959, -5.53810406), (+4.06980848, -6.06634665)
227var
228+9.48650265, +9.08968544, +10.8550940, +17.4775696
229getVar(sample, dim) ! for comparison.
230+9.48650169, +9.08968639, +10.8550940, +17.4775696
231
232
233!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
234!Compute the variance of a 2-D weighted sample.
235!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
236
237
238sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
239sample
240+4.51960468, +8.71387768, +1.69439554, +3.91203833, +6.96232891
241+2.96419239, +5.04764128, +6.34306145, +4.60577965, +5.04803562
242+7.45471573, +1.78127575, +6.18107462, +7.87240362, +8.38695335
243+4.97915459, +8.43780994, +5.47503328, +2.24818325, +2.47297287
244iweight = getUnifRand(1, 9, size(sample, kind = IK))
245iweight
246+7, +8, +1, +4, +2, +3, +7, +7, +1, +8, +3, +2, +5, +3, +8, +4, +4, +5, +4, +7
247rweight = iweight
248call setResized(var, 1_IK)
249call setResized(mean, 1_IK)
250call setVarMean(var(1), mean(1), sample, iweight, iweisum, meang = sample(1,1))
251mean
252+5.14038897
253getMean(sample, iweight) ! for comparison.
254+5.14038897
255var
256+4.95144939
257getVar(sample, iweight) ! for comparison.
258+4.95144939
259[iweisum, sum(iweight)] ! for comparison.
260+93, +93
261call setVarMean(var(1), mean(1), sample, rweight, rweisum, meang = sample(1,1))
262mean
263+5.14038897
264getMean(sample, rweight) ! for comparison.
265+5.14038897
266var
267+4.95144939
268getVar(sample, rweight) ! for comparison.
269+4.95144939
270[rweisum, sum(rweight)] ! for comparison.
271+93.0000000, +93.0000000
272
273dim ! The observations axis.
274+1
275sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
276sample
277+1.88925505, +5.69060040, +4.45722532, +8.87386036, +3.64197922
278+5.64725828, +3.45551586, +5.88054228, +3.79754496, +5.27457094
279+2.92074537, +5.57918119, +5.13014936, +5.57522345, +4.37556982
280+6.04048061, +7.83769131, +2.45026779, +5.34664249, +3.24765253
281iweight = getUnifRand(1, 9, size(sample, dim, IK))
282iweight
283+3, +2, +6, +3
284call setResized(var, size(sample, 3 - dim, IK))
285call setResized(mean, size(sample, 3 - dim, IK))
286if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
287call setVarMean(var, mean, sample, dim, iweight, iweisum, meang)
288mean
289+3.75772834, +5.78364229, +4.51888990, +5.97913837, +4.10510397
290getMean(sample, dim, iweight) ! for comparison.
291+3.75772858, +5.78364277, +4.51888990, +5.97913933, +4.10510397
292var
293+2.67502379, +1.69817841, +1.34278738, +2.63114238, +0.430238515
294getVar(sample, dim, iweight) ! for comparison.
295+2.67502451, +1.69817817, +1.34278738, +2.63114166, +0.430238545
296[iweisum, sum(iweight)] ! for comparison.
297+14, +14
298
299dim ! The observations axis.
300+2
301sample = getUnifRand(1._TKG, 9._TKG, 4_IK, 5_IK)
302sample
303+1.78249931, +2.72149134, +1.68096256, +3.83078337, +6.70195627
304+3.64154243, +4.76133633, +2.12173319, +1.54790163, +4.58819962
305+1.90334749, +3.05829954, +2.15340805, +3.08271360, +3.97181892
306+7.09604263, +7.57688141, +2.75400877, +5.64762306, +5.78451347
307iweight = getUnifRand(1, 9, size(sample, dim, IK))
308iweight
309+2, +1, +1, +4, +5
310call setResized(var, size(sample, 3 - dim, IK))
311call setResized(mean, size(sample, 3 - dim, IK))
312if (dim == 1) then; meang = sample(1,:); else; meang = sample(:,1); end if
313call setVarMean(var, mean, sample, dim, iweight, iweisum, meang)
314mean
315+4.36925888, +3.33067369, +3.16987324, +5.84892559
316getMean(sample, dim, iweight) ! for comparison.
317+4.36925936, +3.33067393, +3.16987324, +5.84892607
318var
319+3.97630215, +1.87089038, +0.576907396, +1.21982884
320getVar(sample, dim, iweight) ! for comparison.
321+3.97630262, +1.87089038, +0.576907456, +1.21982920
322[iweisum, sum(iweight)] ! for comparison.
323+13, +13
324
325
Test:
test_pm_sampleMean


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:
Amir Shahmoradi, April 21, 2017, 1:54 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin

Definition at line 6304 of file pm_sampleVar.F90.


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