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,
-
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 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.
-
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 mean 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 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
.
-
When the input
sample
has the shape (nsam)
, or (:,:)
and dim
is missing, the output mean is a scalar.
-
When the input
sample
has the shape (nsam, ndim)
or (ndim, nsam)
, the output mean is a vector of length ndim
.
Possible calling interfaces ⛓
mean
= getMean(sample(
1 : nsam), weight(
1 : nsam))
mean
= getMean(sample(
1 : nsam), dim)
mean
= getMean(sample(
1 : nsam), dim, weight(
1 : nsam))
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)
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
-
reshaping the array to a vector form and passing it to this procedure, or
-
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 ⛓
12 integer(IK) :: idim, ndim, nsam
13 type(display_type) :: disp
17 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
18 call disp%show(
"!Compute the mean of a 1-D array.")
19 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
24 real(TKG),
allocatable :: sample(:)
26 call disp%show(
"nsam = getUnifRand(1, 5)")
28 call disp%show(
"sample = getUnifRand(0., 1., nsam)")
34 call disp%show(
"mean = getMean(sample, dim = 1_IK)")
35 mean
= getMean(sample,
dim = 1_IK)
42 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
43 call disp%show(
"!Compute the mean of a 2-D array along a specific dimension.")
44 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
48 real(TKG),
allocatable :: mean(:)
49 real(TKG),
allocatable :: sample(:,:)
51 call disp%show(
"ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)")
53 call disp%show(
"sample = getUnifRand(0., 1., ndim, nsam)")
59 call disp%show(
"mean = getMean(sample, dim = 2_IK)")
60 mean
= getMean(sample,
dim = 2_IK)
63 call disp%show(
"mean = getMean(transpose(sample), dim = 1_IK)")
64 mean
= getMean(
transpose(sample),
dim = 1_IK)
71 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
72 call disp%show(
"!Compute the mean of a 1-D weighted array.")
73 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
78 real(TKG),
allocatable :: sample(:)
79 integer(IK),
allocatable :: weight(:)
81 call disp%show(
"nsam = getUnifRand(1, 5)")
83 call disp%show(
"sample = getUnifRand(0., 1., nsam)")
87 call disp%show(
"weight = getUnifRand(1, 9, nsam)")
91 call disp%show(
"mean = getMean(sample)")
95 call disp%show(
"mean = getMean(sample, weight)")
99 call disp%show(
"mean = getMean(sample, 1_IK, weight)")
100 mean
= getMean(sample,
1_IK, weight)
107 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
108 call disp%show(
"!Compute the mean of a 1-D weighted array.")
109 call disp%show(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
114 real(TKG),
allocatable :: sample(:)
115 real(TKG),
allocatable :: weight(:)
117 call disp%show(
"nsam = getUnifRand(1, 5)")
119 call disp%show(
"sample = getUnifRand(0., 1., nsam)")
123 call disp%show(
"weight = getUnifRand(0., 1., nsam)")
127 call disp%show(
"mean = getMean(sample)")
131 call disp%show(
"mean = getMean(sample, weight)")
135 call disp%show(
"mean = getMean(sample, 1_IK, weight)")
136 mean
= getMean(sample,
1_IK, weight)
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
149 real(TKG),
allocatable :: mean(:)
150 real(TKG),
allocatable :: sample(:,:)
151 real(TKG),
allocatable :: weight(:)
153 call disp%show(
"ndim = getUnifRand(1, 3); nsam = getUnifRand(1, 5)")
155 call disp%show(
"sample = getUnifRand(0., 1., ndim, nsam)")
159 call disp%show(
"weight = getUnifRand(0., 1., nsam)")
163 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)
169 call disp%show(
"mean = getMean(transpose(sample), 1_IK, weight)")
170 mean
= getMean(
transpose(sample),
1_IK, weight)
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(
"!%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%%")
183 integer(IK) :: nslice
184 real(TKG),
allocatable :: mean
185 real(TKG),
allocatable,
target :: sample(:,:,:)
186 real(TKG),
pointer :: samptr(:)
188 call disp%show(
"ndim = getUnifRand(1, 2); nsam = getUnifRand(1, 3); nslice = getUnifRand(1, 4)")
190 call disp%show(
"sample = getUnifRand(0., 1., ndim, nsam, nslice)")
196 call disp%show(
"samptr(1:product(shape(sample))) => sample")
197 samptr(
1:
product(
shape(sample)))
=> sample
198 call disp%show(
"mean = getMean(samptr)")
204 call disp%show(
"mean = getMean(reshape(sample, [product(shape(sample))]))")
205 mean
= getMean(
reshape(sample, [
product(
shape(sample))]))
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.
This is a generic method of the derived type display_type with pass attribute.
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...
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 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.
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 ⛓
10+0.859427989,
+0.353745699,
+0.293616951,
+0.262036860
13mean
= getMean(sample,
dim = 1_IK)
26+0.336883903,
+0.978171170,
+0.686017036
27+0.551356018,
+0.250263155,
+0.141360044
30mean
= getMean(sample,
dim = 2_IK)
32+0.667024076,
+0.314326406
33mean
= getMean(
transpose(sample),
dim = 1_IK)
35+0.667024076,
+0.314326406
46+0.307034671,
+0.588855505,
+0.269147933,
+0.759839654,
+0.727329195
56mean
= getMean(sample,
1_IK, weight)
79mean
= getMean(sample,
1_IK, weight)
92+0.613177538,
+0.954005778,
+0.143002987,
+0.308543265,
+0.799150944
95+0.558582485,
+0.612668514,
+0.499010265,
+0.254969299,
+0.483779669
96getMean(sample, [(weight, idim
= 1, ndim)])
98mean
= getMean(sample,
2_IK, weight)
101mean
= getMean(
transpose(sample),
1_IK, weight)
118samptr(
1:
product(
shape(sample)))
=> sample
123mean
= getMean(
reshape(sample, [
product(
shape(sample))]))
- 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.
-
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:
- Amir Shahmoradi, April 21, 2017, 1:54 AM, Institute for Computational Engineering and Sciences (ICES), The University of Texas Austin
Definition at line 229 of file pm_sampleMean.F90.