Generate and return the bias correction factor for the computation of the variance of a (weighted) sample.
More...
Generate and return the bias correction factor for the computation of the variance of a (weighted) sample.
When the population mean is replaced with sample mean, multiplying the output of this generic interface with the output (co)variance from setVar or setCov yields an unbiased estimate of the sample variance.
See the documentation of the parent module pm_sampleVar for algorithmic details and variance definitions.
- Parameters
-
[in] | weisum | : The input scalar of,
-
type
real of kind any supported by the processor (e.g., RK, RK32, RK64, or RK128),
containing either,
-
the size of an unweighted sample, only if the input argument
correction is missing or is not rweight_type,
-
the quantity sum(weight) of a weighted sample, only if the input argument correction is missing or is not rweight_type,
based upon which the bias correction is computed.
The vector weight(1:nsam) refers to the weights of a sample of size nsam .
|
[in] | weisqs | : The input scalar of the same type and kind as weisum , containing the sum of squared sample weight: sum(weight**2) .
(optional. If missing, the output corresponds to the Bessel correction (i.e., frequency weight is assumed).
If present, the correction corresponding to reliability weights is returned.) |
- Returns
normfac
: The output scalar of the same type and kind as weisum
containing the bias correction factor of variance equation.
Multiplying normfac
with the sum of the sample weight (or its size, if unweighted) yields in the Bias Correction factor.
Possible calling interfaces ⛓
!
Generate and return the bias correction factor for the computation of the variance of a (weighted) sa...
This module contains classes and procedures for computing the properties related to the covariance ma...
- Warning
- The condition
0 < weisum
must hold for the corresponding input arguments.
The condition 0 < weisqs
must hold for the corresponding input arguments.
This condition is 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
- Note the effects of bias-correction in computing the variance become noticeable only for very sample sample sizes (i.e., when
nsam
is small).
- See also
- getVar
setVar
getCov
setCov
Example usage ⛓
10 integer(IK) :: itry, ntry
= 10
11 type(display_type) :: disp
16 integer(RKG),
allocatable :: weight(:)
17 real(RKG) :: correction
20 call disp%show(
"weight = getUnifRand(1, 9, getUnifRand(2_IK, 20_IK))")
24 call disp%show(
"correction = getVarCorrection(real(sum(weight), RKG))")
28 call disp%show(
"correction = getVarCorrection(real(sum(weight), RKG), real(sum(weight**2), RKG))")
29 correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
44 real(RKG),
allocatable :: sampleSize(:)
45 sampleSize
= getLogSpace(
log(
1.1_RKG),
log(
1001._RKG),
1000_IK)
46 if (
0 /= getErrTableWrite(
"getVarCorrection.RK.txt",
reshape([sampleSize,
getVarCorrection(sampleSize)], [
size(sampleSize),
2]), header
= SK_
"Sample Size,Bessel Correction"))
error stop 'table write failed.'
Generate count evenly-logarithmically-spaced points over the interval [base**logx1,...
Generate and return a scalar or a contiguous array of rank 1 of length s1 of randomly uniformly distr...
Generate and return the iostat code resulting from writing the input table of rank 1 or 2 to the spec...
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 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 RK
The default real kind in the ParaMonte library: real64 in Fortran, c_double in C-Fortran Interoperati...
integer, parameter LK
The default logical kind in the ParaMonte library: kind(.true.) in Fortran, kind(....
integer, parameter IK
The default integer kind in the ParaMonte library: int32 in Fortran, c_int32_t in C-Fortran Interoper...
integer, parameter SK
The default character kind in the ParaMonte library: kind("a") in Fortran, c_char in C-Fortran Intero...
integer, parameter RKS
The single-precision real kind in Fortran mode. On most platforms, this is an 32-bit real kind.
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 ⛓
4+1,
+7,
+6,
+3,
+2,
+6,
+8,
+6,
+1,
+1,
+8,
+2,
+8,
+4,
+5
8correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
15+2,
+2,
+6,
+6,
+7,
+9,
+2,
+9,
+9,
+8,
+1,
+9,
+6,
+5,
+2,
+1,
+1,
+3,
+7,
+5
19correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
30correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
37+5,
+1,
+5,
+1,
+7,
+7,
+5,
+3,
+2,
+2,
+9,
+4,
+7,
+3,
+7,
+7,
+1,
+6,
+5,
+3
41correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
48+8,
+7,
+1,
+7,
+1,
+3,
+8
52correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
59+1,
+3,
+3,
+1,
+5,
+3,
+9,
+9,
+7,
+3,
+2,
+2,
+7,
+5
63correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
70+1,
+8,
+5,
+5,
+2,
+6,
+6
74correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
81+4,
+8,
+2,
+6,
+4,
+4,
+4,
+2,
+5,
+3,
+1,
+2,
+1,
+4,
+4,
+2,
+4,
+3,
+9,
+9
85correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
92+9,
+1,
+6,
+9,
+4,
+5,
+2,
+8,
+1,
+1,
+3,
+2,
+6,
+5,
+7,
+2,
+9,
+4
96correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
103+9,
+3,
+2,
+5,
+9,
+5,
+5,
+2,
+7,
+2,
+7,
+7,
+7,
+5,
+1,
+4,
+6,
+5,
+6
107correction
= getVarCorrection(
real(
sum(weight), RKG),
real(
sum(weight
**2), RKG))
Postprocessing of the example output ⛓
3import matplotlib.pyplot
as plt
13fileList = glob.glob(pattern)
16 df = pd.read_csv(fileList[0], delimiter =
",")
18 fig = plt.figure(figsize = 1.25 * np.array([6.4, 4.8]), dpi = 200)
21 plt.plot( df.values[:, 0]
23 , linewidth = linewidth
26 ax.legend ( list(df.columns.values[1:])
30 plt.xticks(fontsize = fontsize - 2)
31 plt.yticks(fontsize = fontsize - 2)
34 ax.set_xlabel(df.columns.values[0], fontsize = 17)
35 ax.set_ylabel(
"Correction", fontsize = 17)
37 plt.grid(visible =
True, which =
"both", axis =
"both", color =
"0.85", linestyle =
"-")
38 ax.tick_params(axis =
"y", which =
"minor")
39 ax.tick_params(axis =
"x", which =
"minor")
41 plt.savefig(fileList[0].replace(
".txt",
".png"))
43elif len(fileList) > 1:
45 sys.exit(
"Ambiguous file list exists.")
Visualization of the example output ⛓
- Test:
- test_pm_sampleVar
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 481 of file pm_sampleVar.F90.