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

Generate and return the (weighted) mean of an input sample of nsam observations with ndim = 1 or 2 attributes, optionally weighted by the input weight.
More...

Detailed Description

Generate and return the (weighted) mean of an input sample of nsam observations with ndim = 1 or 2 attributes, optionally weighted by the input weight.

Parameters
[in]sample: The contiguous array of shape (nsam), (ndim, nsam) or (nsam, 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),
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 mean of the whole input sample is computed.)
[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 that of the output mean,
containing the corresponding weights of the data points in the input sample.
(optional, default = a vector of ones.)
Returns
mean : The output object of the same type and kind as the input sample.
  1. When the input sample has the shape (nsam), or (:,:) and dim is missing, the output mean is a scalar.
  2. When the input sample has the shape (nsam, ndim) or (ndim, nsam), the output mean is a vector of length ndim.


Possible calling interfaces

! 1D sample.
mean = getMean(sample(1 : nsam))
mean = getMean(sample(1 : nsam), weight(1 : nsam))
mean = getMean(sample(1 : nsam), dim)
mean = getMean(sample(1 : nsam), dim, weight(1 : nsam))
! 2D sample.
mean(1 : ndim) = getMean(sample(:,:))
mean(1 : ndim) = getMean(sample(:,:), weight(1 : nsam))
mean(1 : ndim) = getMean(sample(:,:), dim)
mean(1 : ndim) = getMean(sample(:,:), dim, weight(1 : nsam))
Generate and return the (weighted) mean of an input sample of nsam observations with ndim = 1 or 2 at...
This module contains classes and procedures for computing the first moment (i.e., the statistical mea...
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.
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.
It is therefore sensible to offer interfaces for only weights of type integer of default kind and real of the same kind as the sample.
Note that the mean of any one or two-dimensional sample can be simply computed via the Fortran intrinsic routine sum():
integer :: i
integer , parameter :: NDIM = 3_IK
integer , parameter :: NSAM = 1000_IK
real , parameter :: sample(NDIM,NSAM) = reshape([( real(i,RK), i = 1, NSAM )], shape = shape(sample))
real , 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.
Developer Remark:
An XY input sample interface is impossible due to ambiguity with existing interfaces.
See also
getVar
setMean
setVarMean


Example usage

1program example
2
3 use pm_kind, only: SK, IK
4 use pm_kind, only: TKG => RKS ! All other real types are also supported.
5 use pm_sampleMean, only: getMean
6 use pm_arrayRange, only: getRange
7 use pm_distUnif, only: getUnifRand
8 use pm_io, only: display_type
9
10 implicit none
11
12 integer(IK) :: idim, ndim, nsam
13 type(display_type) :: disp
14 disp = display_type(file = "main.out.F90")
15
16 call disp%skip()
17 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
18 call disp%show("!Compute the mean of a 1-D array.")
19 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
20 call disp%skip()
21
22 block
23 real(TKG) :: mean
24 real(TKG), allocatable :: sample(:)
25 call disp%skip()
26 call disp%show("nsam = getUnifRand(1, 5)")
27 nsam = getUnifRand(1, 5)
28 call disp%show("sample = getUnifRand(0., 1., nsam)")
29 sample = getUnifRand(0., 1., nsam)
30 call disp%show("sample")
31 call disp%show( sample )
32 call disp%show("getMean(sample)")
33 call disp%show( getMean(sample) )
34 call disp%show("mean = getMean(sample, dim = 1_IK)")
35 mean = getMean(sample, dim = 1_IK)
36 call disp%show("mean")
37 call disp%show( mean )
38 call disp%skip()
39 end block
40
41 call disp%skip()
42 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
43 call disp%show("!Compute the mean of a 2-D array along a specific dimension.")
44 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
45 call disp%skip()
46
47 block
48 real(TKG), allocatable :: mean(:)
49 real(TKG), allocatable :: sample(:,:)
50 call disp%skip()
51 call disp%show("ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)")
52 ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)
53 call disp%show("sample = getUnifRand(0., 1., ndim, nsam)")
54 sample = getUnifRand(0., 1., ndim, nsam)
55 call disp%show("sample")
56 call disp%show( sample )
57 call disp%show("getMean(sample)")
58 call disp%show( getMean(sample) )
59 call disp%show("mean = getMean(sample, dim = 2_IK)")
60 mean = getMean(sample, dim = 2_IK)
61 call disp%show("mean")
62 call disp%show( mean )
63 call disp%show("mean = getMean(transpose(sample), dim = 1_IK)")
64 mean = getMean(transpose(sample), dim = 1_IK)
65 call disp%show("mean")
66 call disp%show( mean )
67 call disp%skip()
68 end block
69
70 call disp%skip()
71 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
72 call disp%show("!Compute the mean of a 1-D weighted array.")
73 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
74 call disp%skip()
75
76 block
77 real(TKG) :: mean
78 real(TKG), allocatable :: sample(:)
79 integer(IK), allocatable :: weight(:)
80 call disp%skip()
81 call disp%show("nsam = getUnifRand(1, 5)")
82 nsam = getUnifRand(1, 5)
83 call disp%show("sample = getUnifRand(0., 1., nsam)")
84 sample = getUnifRand(0., 1., nsam)
85 call disp%show("sample")
86 call disp%show( sample )
87 call disp%show("weight = getUnifRand(1, 9, nsam)")
88 weight = getUnifRand(1, 9, nsam)
89 call disp%show("weight")
90 call disp%show( weight )
91 call disp%show("mean = getMean(sample)")
92 mean = getMean(sample)
93 call disp%show("mean")
94 call disp%show( mean )
95 call disp%show("mean = getMean(sample, weight)")
96 mean = getMean(sample, weight)
97 call disp%show("mean")
98 call disp%show( mean )
99 call disp%show("mean = getMean(sample, 1_IK, weight)")
100 mean = getMean(sample, 1_IK, weight)
101 call disp%show("mean")
102 call disp%show( mean )
103 call disp%skip()
104 end block
105
106 call disp%skip()
107 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
108 call disp%show("!Compute the mean of a 1-D weighted array.")
109 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
110 call disp%skip()
111
112 block
113 real(TKG) :: mean
114 real(TKG), allocatable :: sample(:)
115 real(TKG), allocatable :: weight(:)
116 call disp%skip()
117 call disp%show("nsam = getUnifRand(1, 5)")
118 nsam = getUnifRand(1, 5)
119 call disp%show("sample = getUnifRand(0., 1., nsam)")
120 sample = getUnifRand(0., 1., nsam)
121 call disp%show("sample")
122 call disp%show( sample )
123 call disp%show("weight = getUnifRand(0., 1., nsam)")
124 weight = getUnifRand(0., 1., nsam)
125 call disp%show("weight")
126 call disp%show( weight )
127 call disp%show("mean = getMean(sample)")
128 mean = getMean(sample)
129 call disp%show("mean")
130 call disp%show( mean )
131 call disp%show("mean = getMean(sample, weight)")
132 mean = getMean(sample, weight)
133 call disp%show("mean")
134 call disp%show( mean )
135 call disp%show("mean = getMean(sample, 1_IK, weight)")
136 mean = getMean(sample, 1_IK, weight)
137 call disp%show("mean")
138 call disp%show( mean )
139 call disp%skip()
140 end block
141
142 call disp%skip()
143 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
144 call disp%show("!Compute the mean of a 2-D weighted array along a specific dimension.")
145 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
146 call disp%skip()
147
148 block
149 real(TKG), allocatable :: mean(:)
150 real(TKG), allocatable :: sample(:,:)
151 real(TKG), allocatable :: weight(:)
152 call disp%skip()
153 call disp%show("ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)")
154 ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)
155 call disp%show("sample = getUnifRand(0., 1., ndim, nsam)")
156 sample = getUnifRand(0., 1., ndim, nsam)
157 call disp%show("sample")
158 call disp%show( sample )
159 call disp%show("weight = getUnifRand(0., 1., nsam)")
160 weight = getUnifRand(0., 1., nsam)
161 call disp%show("weight")
162 call disp%show( weight )
163 call disp%show("getMean(sample, [(weight, idim = 1, ndim)])")
164 call disp%show( getMean(sample, [(weight, idim = 1, ndim)]) )
165 call disp%show("mean = getMean(sample, 2_IK, weight)")
166 mean = getMean(sample, 2_IK, weight)
167 call disp%show("mean")
168 call disp%show( mean )
169 call disp%show("mean = getMean(transpose(sample), 1_IK, weight)")
170 mean = getMean(transpose(sample), 1_IK, weight)
171 call disp%show("mean")
172 call disp%show( mean )
173 call disp%skip()
174 end block
175
176 call disp%skip()
177 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
178 call disp%show("!Compute the mean of a multidimensional array by associating it with a 1D pointer.")
179 call disp%show("!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
180 call disp%skip()
181
182 block
183 integer(IK) :: nslice
184 real(TKG), allocatable :: mean
185 real(TKG), allocatable, target :: sample(:,:,:)
186 real(TKG), pointer :: samptr(:)
187 call disp%skip()
188 call disp%show("ndim = getUnifRand(1, 2); nsam = getUnifRand(1, 3); nslice = getUnifRand(1, 4)")
189 ndim = getUnifRand(1, 2); nsam = getUnifRand(1, 3); nslice = getUnifRand(1, 4)
190 call disp%show("sample = getUnifRand(0., 1., ndim, nsam, nslice)")
191 sample = getUnifRand(0., 1., ndim, nsam, nslice)
192 call disp%show("sample")
193 call disp%show( sample )
194 call disp%show("shape(sample)")
195 call disp%show( shape(sample) )
196 call disp%show("samptr(1:product(shape(sample))) => sample")
197 samptr(1:product(shape(sample))) => sample
198 call disp%show("mean = getMean(samptr)")
199 mean = getMean(samptr)
200 call disp%show("nullify(samptr)")
201 nullify(samptr)
202 call disp%show("mean")
203 call disp%show( mean )
204 call disp%show("mean = getMean(reshape(sample, [product(shape(sample))]))")
205 mean = getMean(reshape(sample, [product(shape(sample))]))
206 call disp%show("mean")
207 call disp%show( mean )
208 call disp%skip()
209 end block
210
211end program example
Generate minimally-spaced character, integer, real sequences or sequences at fixed intervals of size ...
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
This module contains procedures and generic interfaces for generating ranges of discrete character,...
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 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
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 mean of a 1-D array.
4!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
5
6
7nsam = getUnifRand(1, 5)
8sample = getUnifRand(0., 1., nsam)
9sample
10+0.595203698, +0.475735009, +0.863631487, +0.700949430E-1
11getMean(sample)
12+0.501166284
13mean = getMean(sample, dim = 1_IK)
14mean
15+0.501166284
16
17
18!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
19!Compute the mean of a 2-D array along a specific dimension.
20!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
21
22
23ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)
24sample = getUnifRand(0., 1., ndim, nsam)
25sample
26+0.398384452
27+0.523957610E-1
28getMean(sample)
29+0.225390106
30mean = getMean(sample, dim = 2_IK)
31mean
32+0.398384452, +0.523957610E-1
33mean = getMean(transpose(sample), dim = 1_IK)
34mean
35+0.398384452, +0.523957610E-1
36
37
38!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
39!Compute the mean of a 1-D weighted array.
40!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
41
42
43nsam = getUnifRand(1, 5)
44sample = getUnifRand(0., 1., nsam)
45sample
46+0.315453410
47weight = getUnifRand(1, 9, nsam)
48weight
49+6
50mean = getMean(sample)
51mean
52+0.315453410
53mean = getMean(sample, weight)
54mean
55+0.315453410
56mean = getMean(sample, 1_IK, weight)
57mean
58+0.315453410
59
60
61!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
62!Compute the mean of a 1-D weighted array.
63!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
64
65
66nsam = getUnifRand(1, 5)
67sample = getUnifRand(0., 1., nsam)
68sample
69+0.135345101
70weight = getUnifRand(0., 1., nsam)
71weight
72+0.745889008
73mean = getMean(sample)
74mean
75+0.135345101
76mean = getMean(sample, weight)
77mean
78+0.135345101
79mean = getMean(sample, 1_IK, weight)
80mean
81+0.135345101
82
83
84!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
85!Compute the mean of a 2-D weighted array along a specific dimension.
86!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
87
88
89ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)
90sample = getUnifRand(0., 1., ndim, nsam)
91sample
92+0.911404252, +0.139900982, +0.618470728, +0.399185598, +0.119258225
93weight = getUnifRand(0., 1., nsam)
94weight
95+0.521258175, +0.984511793, +0.792486846, +0.788852513, +0.254609108
96getMean(sample, [(weight, idim = 1, ndim)])
97+0.433370978
98mean = getMean(sample, 2_IK, weight)
99mean
100+0.433370978
101mean = getMean(transpose(sample), 1_IK, weight)
102mean
103+0.433370978
104
105
106!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
107!Compute the mean of a multidimensional array by associating it with a 1D pointer.
108!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%
109
110
111ndim = getUnifRand(1, 2); nsam = getUnifRand(1, 3); nslice = getUnifRand(1, 4)
112sample = getUnifRand(0., 1., ndim, nsam, nslice)
113sample
114slice(:,:,1) =
115+0.758918405
116slice(:,:,2) =
117+0.409432650
118slice(:,:,3) =
119+0.542795658
120shape(sample)
121+1, +1, +3
122samptr(1:product(shape(sample))) => sample
123mean = getMean(samptr)
124nullify(samptr)
125mean
126+0.570382237
127mean = getMean(reshape(sample, [product(shape(sample))]))
128mean
129+0.570382237
130
131
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 at Austin

Definition at line 229 of file pm_sampleMean.F90.


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