STK++ 0.9.13
The descriptive statistics sub-project.

In this sub-project, we compute the usual descriptive statistics of variables. More...

Collaboration diagram for The descriptive statistics sub-project.:

Classes

class  STK::Stat::Bivariate< XTYPE, YTYPE, TContainer1D >
 Computation of the bivariate Statistics of two Variables. More...
 
class  STK::Stat::ConfusionMatrix< TrueArray, PredArray >
 Class computing the confusion matrix given two vectors of integer. More...
 
class  STK::Stat::Factor< Array >
 Computation of the Factors of a 1D Container. More...
 
class  STK::Stat::MultiFactor< Array >
 Computation of the MultiFactors of a 2D Container. More...
 
class  STK::Stat::Multivariate< Array, WColVector, Type >
 Computation of the multivariate statistics of a Variable. More...
 
class  STK::Stat::Multivariate< Array, Real >
 Computation of the Multivariate Statistics of a 2D Container of Real. More...
 
struct  STK::Stat::Online< Array, Real >
 Computation online of the statistics of a Real set of variable. More...
 
struct  STK::Stat::Online< Real, Real >
 Computation online of the statistics of a Real set of variable. More...
 
class  STK::Stat::Univariate< TContainer1D, Real >
 Computation of the univariate Statistics of a Real Variable. More...
 

Functions

template<class TrueArray , class PredArray >
CSquareXi STK::Stat::confusionMatrix (TrueArray const &trueClass, PredArray const &predictClass)
 Utility function computing the confusion matrix given two vectors of integer.
 
template<class XArray , class YArray >
Real STK::Stat::covariance (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, bool unbiased=false)
 Compute the covariance of the variables X and Y.
 
template<class XArray , class YArray >
Real STK::Stat::covarianceSafe (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, bool unbiased=false)
 Compute the covariance of the variables X and Y.
 
template<class XArray , class YArray , class Weights >
Real STK::Stat::covariance (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, ExprBase< Weights > const &W, bool unbiased=false)
 The covariance between X and Y will be computed using the formula.
 
template<class XArray , class YArray , class Weights >
Real STK::Stat::covarianceSafe (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, ExprBase< Weights > const &W, bool unbiased=false)
 The covariance between X and Y will be computed using the formula.
 
template<class XArray , class YArray >
Real STK::Stat::covarianceWithFixedMean (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, typename hidden::Traits< XArray >::Type const &xMu, typename hidden::Traits< YArray >::Type const &yMu, bool unbiased=false)
 Compute the covariance of the variables X and Y with fixed means.
 
template<class XArray , class YArray >
Real STK::Stat::covarianceWithFixedMeanSafe (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, typename hidden::Traits< XArray >::Type const &xMu, typename hidden::Traits< YArray >::Type const &yMu, bool unbiased=false)
 Compute the covariance of the variables X and Y with fixed means.
 
template<class XArray , class YArray , class Weights >
Real STK::Stat::covarianceWithFixedMean (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, ExprBase< Weights > const &W, typename hidden::Traits< XArray >::Type const &xMean, typename hidden::Traits< YArray >::Type const &yMean, bool unbiased=false)
 The covariance between X and Y will be computed using the formula.
 
template<class XArray , class YArray , class Weights >
Real STK::Stat::covarianceWithFixedMeanSafe (ExprBase< XArray > const &X, ExprBase< YArray > const &Y, ExprBase< Weights > const &W, typename hidden::Traits< XArray >::Type const &xMu, typename hidden::Traits< YArray >::Type const &yMu, bool unbiased=false)
 The covariance between X and Y will be computed using the formula.
 
template<class Array >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covariance (ExprBase< Array > const &V, bool unbiased=false)
 Compute the covariance matrix using the column of the data set V.
 
template<class Array >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceByRow (ExprBase< Array > const &V, bool unbiased=false)
 Compute the covariance matrix using the rows of the data set V.
 
template<class Array , class Weights >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covariance (ExprBase< Array > const &V, Weights const &W, bool unbiased=false)
 Compute the weighted covariance matrix using the column of the data set V.
 
template<class Array , class Weights >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceByRow (ExprBase< Array > const &V, Weights const &W, bool unbiased=false)
 Compute the weighted covariance matrix using the rows of the data set V.
 
template<class Array , class Mean >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covarianceWithFixedMean (ExprBase< Array > const &V, ExprBase< Mean > const &mean, bool unbiased=false)
 Compute the covariance matrix using the columns of the data set V.
 
template<class Array , class Mean >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceWithFixedMeanByRow (ExprBase< Array > const &V, ExprBase< Mean > const &mean, bool unbiased=false)
 Compute the covariance matrix using the rows of the data set V.
 
template<class Array , class Weights , class Mean >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covarianceWithFixedMean (ExprBase< Array > const &V, Weights const &W, ExprBase< Mean > const &mean, bool unbiased=false)
 Compute the weighted covariance matrix using the column of the data set V.
 
template<class Array , class Weights , class Mean >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceWithFixedMeanByRow (ExprBase< Array > const &V, Weights const &W, ExprBase< Mean > const &mean, bool unbiased=false)
 Compute the weighted covariance matrix using the rows of the data set V.
 
template<class Derived >
hidden::FunctorTraits< Derived, MinOp >::Row STK::Stat::min (Derived const &A)
 Compute the minimal(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, MinSafeOp >::Row STK::Stat::minSafe (Derived const &A)
 Compute safely the minimal(s) [weighted] value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, MaxOp >::Row STK::Stat::max (Derived const &A)
 Compute the maximal(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, MaxSafeOp >::Row STK::Stat::maxSafe (Derived const &A)
 Compute safely the maximal(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, SumOp >::Row STK::Stat::sum (Derived const &A)
 Compute the sum of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, SumSafeOp >::Row STK::Stat::sumSafe (Derived const &A)
 Compute safely the mean(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, MeanOp >::Row STK::Stat::mean (Derived const &A)
 Compute the mean(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, MeanSafeOp >::Row STK::Stat::meanSafe (Derived const &A)
 Compute safely the mean(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, VarianceOp >::Row STK::Stat::variance (Derived const &A, bool unbiased=false)
 Compute the variance(s) value(s) of A.
 
template<class Derived >
hidden::FunctorTraits< Derived, VarianceSafeOp >::Row STK::Stat::varianceSafe (Derived const &A, bool unbiased=false)
 Compute safely the variance(s) value(s) of A.
 
template<class Derived , class MeanType >
hidden::FunctorTraits< Derived, VarianceWithFixedMeanOp >::Row STK::Stat::varianceWithFixedMean (Derived const &A, MeanType const &mean, bool unbiased)
 Compute the VarianceWithFixedMean(s) value(s) of A.
 
template<class Derived , class MeanType >
hidden::FunctorTraits< Derived, VarianceWithFixedMeanSafeOp >::Row STK::Stat::varianceWithFixedMeanSafe (Derived const &A, MeanType const &mean, bool unbiased)
 Compute safely the VarianceWithFixedMean(s) value(s) of A.
 
template<class Array , class RowVector >
void STK::Stat::center (Array &m, RowVector &mean)
 Compute the mean by column of the variables in the container m and center it.
 
template<class Array , class ColVector >
void STK::Stat::centerByRow (Array &m, ColVector &mean)
 Compute the mean by row of the variables in the container V and center V.
 
template<class Array , class Weights , class RowVector >
void STK::Stat::center (Array &m, Weights const &W, RowVector &mean)
 Compute the weighted mean by column of the variables and center m.
 
template<class Array , class Weights , class ColVector >
void STK::Stat::centerByRow (Array &m, Weights const &W, ColVector &mean)
 Compute the weighted mean of the variables by rows and center m.
 
template<class Array , class RowVector >
void STK::Stat::standardize (Array &m, RowVector &mean, RowVector &std, bool unbiased=false)
 Compute the mean and the standard deviation by columns of the variable m and standardize it.
 
template<class Array , class ColVector >
void STK::Stat::standardizeByRow (Array &m, ColVector &mean, ColVector &std, bool unbiased=false)
 Compute the mean and the standard deviation by rows of the variable m and standardize it.
 
template<class Array , class Weights , class RowVector >
void STK::Stat::standardize (Array &m, Weights const &W, RowVector &mean, RowVector &std, bool unbiased=false)
 Compute the weighted means and standard deviations by columns of the variable m and standardize it.
 
template<class Array , class Weights , class ColVector >
void STK::Stat::standardizeByRow (Array &m, Weights const &W, ColVector &mean, ColVector &std, bool unbiased=false)
 Compute the weighted means and standard deviations by rows of the variable m and standardize it.
 
template<class Array , class RowVector >
void STK::Stat::unstandardize (Array &m, RowVector const &std)
 undo the standardization by columns of the standardized variable m.
 
template<class Array , class ColVector >
void STK::Stat::unstandardizeByRow (Array &m, ColVector const &std)
 undo the standardization by rows of the standardized variable m.
 
template<class Array , class RowVector >
void STK::Stat::unstandardize (Array &m, RowVector const &mean, RowVector const &std)
 undo the standardization by columns of the standardized variable m
 
template<class Array , class ColVector >
void STK::Stat::unstandardizeByRow (Array &m, ColVector const &mean, ColVector const &std)
 undo the standardization by rows of the standardized variable m
 

Detailed Description

In this sub-project, we compute the usual descriptive statistics of variables.

Function Documentation

◆ center() [1/2]

template<class Array , class RowVector >
void STK::Stat::center ( Array &  m,
RowVector &  mean 
)

Compute the mean by column of the variables in the container m and center it.

Parameters
m,meandata and means arrays

Definition at line 52 of file STK_Stat_Transform.h.

53{
54 typedef typename Array::Type Type;
55 enum
56 {
57 sizeRows_ = Array::sizeRows_,
58 sizeCols_ = Array::sizeCols_
59 };
60 m -= Const::Vector<Type, sizeRows_>(m.rows())
61 * (mean = Stat::meanByCol(m));
62}
hidden::SliceVisitorSelector< Derived, hidden::MeanVisitor, Arrays::by_col_ >::type_result mean(Derived const &A)
If A is a row-vector or a column-vector then the function will return the usual mean value of the vec...

References STK::Stat::mean(), and STK::Stat::meanByCol().

Referenced by STK::IAAModel< Array >::center(), STK::IAAModel< Array >::center(), STK::Stat::centerByCol(), and STK::Stat::centerByCol().

◆ center() [2/2]

template<class Array , class Weights , class RowVector >
void STK::Stat::center ( Array &  m,
Weights const W,
RowVector &  mean 
)

Compute the weighted mean by column of the variables and center m.

Parameters
m,W,meancontainer of the data, weights and means

Definition at line 89 of file STK_Stat_Transform.h.

90{
91 typedef typename Array::Type Type;
92 enum
93 {
94 sizeRows_ = Array::sizeRows_,
95 sizeCols_ = Array::sizeCols_
96 };
97 m -= Const::Vector<Type, sizeRows_>(m.rows())
98 * (mean = Stat::meanByCol(m, W));
99}

References STK::Stat::mean(), and STK::Stat::meanByCol().

◆ centerByRow() [1/2]

template<class Array , class ColVector >
void STK::Stat::centerByRow ( Array &  m,
ColVector &  mean 
)

Compute the mean by row of the variables in the container V and center V.

Parameters
m,meandata and means arrays

Definition at line 72 of file STK_Stat_Transform.h.

73{
74 typedef typename Array::Type Type;
75 enum
76 {
77 sizeRows_ = Array::sizeRows_,
78 sizeCols_ = Array::sizeCols_
79 };
80 m -= (mean = Stat::meanByRow(m))
81 * Const::Point<Type, sizeCols_>(m.cols());
82}

References STK::Stat::mean(), and STK::Stat::meanByRow().

Referenced by STK::Stat::standardizeByRow(), and STK::Stat::standardizeByRow().

◆ centerByRow() [2/2]

template<class Array , class Weights , class ColVector >
void STK::Stat::centerByRow ( Array &  m,
Weights const W,
ColVector &  mean 
)

Compute the weighted mean of the variables by rows and center m.

Parameters
m,W,meancontainer of the data, weights and means

Definition at line 109 of file STK_Stat_Transform.h.

110{
111 typedef typename Array::Type Type;
112 enum
113 {
114 sizeRows_ = Array::sizeRows_,
115 sizeCols_ = Array::sizeCols_
116 };
117 m -= (mean = Stat::meanByRow(m, W))
118 * Const::Point<Type, sizeCols_>(m.cols());
119}

References STK::Stat::mean(), and STK::Stat::meanByRow().

◆ confusionMatrix()

CSquareXi STK::Stat::confusionMatrix ( TrueArray const trueClass,
PredArray const predictClass 
)

Utility function computing the confusion matrix given two vectors of integer.

Definition at line 87 of file STK_Stat_ConfusionMatrix.h.

88{ return ConfusionMatrix<TrueArray, PredArray>(trueClass,predictClass)();}

Referenced by STK::Stat::ConfusionMatrix< TrueArray, PredArray >::ConfusionMatrix().

◆ covariance() [1/4]

template<class Array >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covariance ( ExprBase< Array > const V,
bool  unbiased = false 
)

Compute the covariance matrix using the column of the data set V.

\[ \hat{\Sigma} = \frac{1}{n}
                     \sum_{i=1}^n (V_i-\hat{\mu}) (V_i-\hat{\mu})^T.
\]

Parameters
Vvariable
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 484 of file STK_Stat_Covariance.h.

485{
486 CArraySquare<typename Array::Type, Array::sizeCols_, Array::orient_> cov_(V.cols());
487 // typename Array::Row mean;
488 typename hidden::FunctorTraits<Array, MeanOp>::Row mean;
489 // compute the mean
490 mean.move(Stat::mean(V.asDerived()));
491 for (int j= cov_.begin(); j< cov_.end(); j++)
492 {
493 cov_(j, j) = varianceWithFixedMean(V.col(j), mean[j], unbiased);
494 for (int i= cov_.begin(); i<j; i++)
495 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.col(i), V.col(j), mean[i], mean[j], unbiased));}
496 }
497 return cov_;
498}
hidden::FunctorTraits< Derived, VarianceWithFixedMeanOp >::Row varianceWithFixedMean(Derived const &A, MeanType const &mean, bool unbiased)
Compute the VarianceWithFixedMean(s) value(s) of A.
Real covarianceWithFixedMean(ExprBase< XArray > const &X, ExprBase< YArray > const &Y, typename hidden::Traits< XArray >::Type const &xMu, typename hidden::Traits< YArray >::Type const &yMu, bool unbiased=false)
Compute the covariance of the variables X and Y with fixed means.

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), and STK::Stat::varianceWithFixedMean().

◆ covariance() [2/4]

template<class Array , class Weights >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covariance ( ExprBase< Array > const V,
Weights const W,
bool  unbiased = false 
)

Compute the weighted covariance matrix using the column of the data set V.

\[ \hat{\sigma}^2 = \frac{1}{n}
                     \sum_{i=1}^n  w_{i} w_j (V_i-\hat{\mu}) (V_i-\hat{\mu})^T.
\]

Parameters
V,Wthe variables and the weights
unbiasedtrue if we want an unbiased estimate of the covariance, false otherwise (default is false)

Definition at line 537 of file STK_Stat_Covariance.h.

538{
539 CArraySquare<typename Array::Type, Array::sizeCols_, Array::orient_> cov_(V.cols());
540 typename hidden::FunctorTraits<Array, MeanOp>::Row mean;
541 // compute the mean
542 mean.move(Stat::mean(V.asDerived(), W.asDerived()));
543 for (int j= cov_.begin(); j< cov_.end(); j++)
544 {
545 cov_(j, j) = varianceWithFixedMean(V.col(j), W, mean[j], unbiased);
546 for (int i= cov_.begin(); i<j; i++)
547 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.col(i), V.col(j), W, mean[i], mean[j], unbiased));}
548 }
549 return cov_;
550}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), and STK::Stat::varianceWithFixedMean().

◆ covariance() [3/4]

template<class XArray , class YArray >
Real STK::Stat::covariance ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
bool  unbiased = false 
)

Compute the covariance of the variables X and Y.

\[ \hat{cov}(X,Y) = \frac{1}{n}
                \sum_{i=1}^n (X_i - \bar{X})(Y_i - \bar{Y}).
\]

Parameters
X,Yvariables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 55 of file STK_Stat_Covariance.h.

59{
62 typedef typename hidden::Traits<XArray>::Type XType;
63 typedef typename hidden::Traits<YArray>::Type YType;
64 typedef typename hidden::Promote<XType, YType>::result_type Type;
65
66#ifdef STK_BOUNDS_CHECK
67 if (X.range() != Y.range())
68 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covariance(X,Y,unbiased),ranges are not the sames);
69#endif
70 Type xMean = mean(X.asDerived()), yMean = mean(Y.asDerived());
71 int nobs = X.size();
72 Type xsum = 0.0, ysum = 0.0, cov = 0.0, xdev, ydev;
73 for (int i=X.begin(); i<X.end(); i++)
74 {
75 xsum += (xdev = X[i] - xMean);
76 ysum += (ydev = Y[i] - yMean);
77 cov += (xdev*ydev);
78 }
79 // compute the covariance with corrected bias
80 return (unbiased) ?
81 (nobs > 1) ? ((cov - (xsum*ysum)/Type(nobs))/Type(nobs -1)) : 0.
82
83 :
84 (nobs > 0) ? (cov - (xsum*ysum)/Type(nobs))/Type(nobs) : 0.;
85}
#define STKRUNTIME_ERROR_NO_ARG(Where, Error)
Definition STK_Macros.h:138
#define STK_STATIC_ASSERT_ONE_DIMENSION_ONLY(EXPR)

References STK::Stat::covariance(), STK::Stat::mean(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

Referenced by STK::GaussianModel< Array >::compCovariance(), STK::GaussianAAModel< Array >::computeProjectedCovariance(), STK::GaussianAAModel< Array >::computeResidualCovariance(), STK::GaussianModel< Array >::compWeightedCovariance(), STK::Stat::covariance(), main(), STK::ProjectedVariance< Array >::maximizeStep(), and STK::ProjectedVariance< Array >::maximizeStep().

◆ covariance() [4/4]

Real STK::Stat::covariance ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
ExprBase< Weights > const W,
bool  unbiased = false 
)

The covariance between X and Y will be computed using the formula.

See also
Mark Galassi, Jim Davies, James Theiler, Brian Gough, Gerard Jungman, Michael Booth, and Fabrice Rossi. http://www.gnu.org/software/gsl/manual GNU Scientific Library - Reference manual, Version 1.15, 2011. http://www.gnu.org/software/gsl/manual/html_node/Weighted-Samples.html Sec. 21.7 Weighted Samples

$ cov(X,Y)=\frac{\sum_{i=1}^{N}w_i}{\left(\sum_{i=1}^{N}w_i\right)^2-\sum_{i=1}^{N}w_i^2}
\sum_{i=1}^N w_i \left(  x_{i}-\bar{x} \right)  \left( y_{i}-\bar{y} \right). $

Parameters
X,Y,Wvariables and weights
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 153 of file STK_Stat_Covariance.h.

158{
162 typedef typename hidden::Traits<XArray>::Type XType;
163 typedef typename hidden::Traits<YArray>::Type YType;
164 typedef typename hidden::Promote<XType, YType>::result_type Type;
165
166#ifdef STK_BOUNDS_CHECK
167 if (X.range() != Y.range())
168 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceWithFixedMeanSafe(X,Y,W,xMu,yMu,unbiased),ranges are not the sames);
169#endif
170
171 Type xMean = mean(X.asDerived(),W.asDerived()), yMean = mean(Y.asDerived(),W.asDerived());
172 Real xsum = 0.0, ysum = 0.0, xdev, ydev, sumWeights = 0.0, sum2Weights = 0.0, cov = 0.0;
173 for (int i=X.begin(); i<X.end(); i++)
174 {
175 xsum += (xdev = X[i] - xMean);
176 ysum += (ydev = Y[i] - yMean);
177 Real Wi = std::abs((Real)W[i]);
178 cov += Wi * (xdev*ydev);
179 sumWeights += Wi;
180 sum2Weights += Wi * Wi;
181 }
182 // compute the variance
183 return (unbiased) ?
184 (sumWeights*sumWeights > sum2Weights) ?
185 (cov - xsum*ysum)/(sumWeights - sum2Weights/sumWeights) : 0.
186// (cov - xsum*ysum/sumWeights)/(sumWeights - sum2Weights/sumWeights) : 0.
187 : (sumWeights) ? (cov - xsum*ysum)/(sumWeights) : 0.;
188}
double Real
STK fundamental type of Real values.

References STK::Stat::covarianceWithFixedMeanSafe(), STK::Stat::mean(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

◆ covarianceByRow() [1/2]

template<class Array >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceByRow ( ExprBase< Array > const V,
bool  unbiased = false 
)

Compute the covariance matrix using the rows of the data set V.

\[ \hat{\Sigma} = \frac{1}{p}
                     \sum_{j=1}^p (V_j-\hat{\mu})^T (V_j-\hat{\mu}).
\]

Parameters
Vvariable
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 511 of file STK_Stat_Covariance.h.

512{
513 CArraySquare<typename Array::Type, Array::sizeRows_, Array::orient_> cov_(V.rows());
514 typename hidden::FunctorTraits<Array, MeanOp>::Col mean;
515 // compute the mean
516 mean.move(Stat::meanByRow(V.asDerived()));
517 for (int j= cov_.begin(); j< cov_.end(); j++)
518 {
519 cov_(j, j) = varianceWithFixedMean(V.row(j), mean[j], unbiased);
520 for (int i= cov_.begin(); i<j; i++)
521 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.row(i), V.row(j), mean[i], mean[j], unbiased));}
522 }
523 return cov_;
524}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), STK::Stat::meanByRow(), and STK::Stat::varianceWithFixedMean().

◆ covarianceByRow() [2/2]

template<class Array , class Weights >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceByRow ( ExprBase< Array > const V,
Weights const W,
bool  unbiased = false 
)

Compute the weighted covariance matrix using the rows of the data set V.

\[ \hat{\sigma}^2 = \frac{1}{n}
                     \sum_{i=1}^n  w_{i} w_j (V_i-\hat{\mu}) (V_i-\hat{\mu})^T.
\]

Parameters
V,Wthe variables and the weights
unbiasedtrue if we want an unbiased estimate of the covariance, false otherwise (default is false)

Definition at line 563 of file STK_Stat_Covariance.h.

564{
565 CArraySquare<typename Array::Type, Array::sizeRows_, Array::orient_> cov_(V.rows());
566 typename hidden::FunctorTraits<Array, MeanOp>::Col mean;
567 mean.move(Stat::meanByRow(V.asDerived(), W.asDerived()));
568 for (int j= cov_.begin(); j< cov_.end(); j++)
569 {
570 cov_(j, j) = varianceWithFixedMean(V.row(j), W, mean[j], unbiased);
571 for (int i= cov_.begin(); i<j; i++)
572 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.row(i), V.row(j), W, mean[i], mean[j], unbiased));}
573 }
574 return cov_;
575}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), STK::Stat::meanByRow(), and STK::Stat::varianceWithFixedMean().

◆ covarianceSafe() [1/2]

template<class XArray , class YArray >
Real STK::Stat::covarianceSafe ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
bool  unbiased = false 
)

Compute the covariance of the variables X and Y.

\[ \hat{cov}(X,Y) = \frac{1}{n}
                \sum_{i=1}^n (X_i - \hat{X})(Y_i - \hat{Y}).
\]

Parameters
X,Yvariables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 98 of file STK_Stat_Covariance.h.

102{
105 typedef typename hidden::Traits<XArray>::Type XType;
106 typedef typename hidden::Traits<YArray>::Type YType;
107 typedef typename hidden::Promote<XType, YType>::result_type Type;
108
109 #ifdef STK_BOUNDS_CHECK
110 if (X.range() != Y.range())
111 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceSafe(X,Y,unbiased),ranges are not the sames);
112 #endif
113
114 int nobs = X.size();
115 Type xMean = mean(X.asDerived()), yMean = mean(Y.asDerived());
116 Type xsum = 0.0, ysum = 0.0, cov = 0.0, xdev, ydev;
117 for (int i=X.begin(); i<X.end(); i++)
118 {
119 if (Arithmetic<Type>::isFinite(X[i]) && Arithmetic<Type>::isFinite(Y[i]))
120 {
121 xsum += (xdev = X[i] - xMean);
122 ysum += (ydev = Y[i] - yMean);
123 cov += (xdev*ydev);
124 }
125 else nobs--;
126 }
127 // compute the covariance with corrected bias
128 return (unbiased) ?
129 (nobs > 1) ? ((cov - (xsum*ysum)/Type(nobs))/Type(nobs -1)) : 0.
130
131 :
132 (nobs > 0) ? (cov - (xsum*ysum)/Type(nobs))/Type(nobs) : 0.;
133}

References STK::Stat::covarianceSafe(), STK::Arithmetic< Type >::isFinite(), STK::Stat::mean(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

Referenced by STK::Stat::covarianceSafe().

◆ covarianceSafe() [2/2]

Real STK::Stat::covarianceSafe ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
ExprBase< Weights > const W,
bool  unbiased = false 
)

The covariance between X and Y will be computed using the formula.

See also
Mark Galassi, Jim Davies, James Theiler, Brian Gough, Gerard Jungman, Michael Booth, and Fabrice Rossi. http://www.gnu.org/software/gsl/manual GNU Scientific Library - Reference manual, Version 1.15, 2011. http://www.gnu.org/software/gsl/manual/html_node/Weighted-Samples.html Sec. 21.7 Weighted Samples

$ cov(X,Y)=\frac{\sum_{i=1}^{N}w_i}{\left(\sum_{i=1}^{N}w_i\right)^2-\sum_{i=1}^{N}w_i^2}
\sum_{i=1}^N w_i \left(  x_{i}-\bar{x} \right)  \left( y_{i}-\bar{y} \right). $

Parameters
X,Y,Wvariables and weights
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 208 of file STK_Stat_Covariance.h.

213{
217 typedef typename hidden::Traits<XArray>::Type XType;
218 typedef typename hidden::Traits<YArray>::Type YType;
219 typedef typename hidden::Promote<XType, YType>::result_type Type;
220
221#ifdef STK_BOUNDS_CHECK
222 if (X.range() != Y.range())
223 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceWithFixedMeanSafe(X,Y,W,xMu,yMu,unbiased),ranges are not the sames);
224#endif
225
226 Type xMean = mean(X.asDerived(),W.asDerived()), yMean = mean(Y.asDerived(),W.asDerived());
227 // compute covariance
228 Type xsum = 0.0, ysum = 0.0, xdev, ydev, sumWeights = 0.0, sum2Weights = 0.0, cov = 0.0;
229 for (int i=X.begin(); i<X.end(); i++)
230 {
231 if (Arithmetic<Type>::isFinite(X[i])
232 && Arithmetic<Type>::isFinite(Y[i])
233 && Arithmetic<Type>::isFinite(W[i])
234 )
235 {
236 xsum += (xdev = X[i] - xMean); // deviation from the mean
237 ysum += (ydev = Y[i] - yMean); // deviation from the mean
238 Type Wi = std::abs((Real)W[i]);
239 cov += Wi * (xdev*ydev); // cross product
240 sumWeights += Wi; // sum absolute weights
241 sum2Weights += Wi * Wi; // sum squared weights
242 }
243 }
244 // compute the variance
245 return (unbiased) ?
246 (sumWeights*sumWeights > sum2Weights) ?
247 (cov - xsum*ysum)/(sumWeights - sum2Weights/sumWeights) : 0.
248// (cov - xsum*ysum/sumWeights)/(sumWeights - sum2Weights/sumWeights) : 0.
249 : (sumWeights) ? (cov - xsum*ysum)/(sumWeights) : 0.;
250}

References STK::Stat::covarianceWithFixedMeanSafe(), STK::Arithmetic< Type >::isFinite(), STK::Stat::mean(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

◆ covarianceWithFixedMean() [1/4]

template<class Array , class Mean >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covarianceWithFixedMean ( ExprBase< Array > const V,
ExprBase< Mean > const mean,
bool  unbiased = false 
)

Compute the covariance matrix using the columns of the data set V.

\[ \hat{\Sigma} = \frac{1}{n}
                     \sum_{i=1}^n (V_i-mu) (V_i-mu)^T.
\]

Parameters
Vvariable
meanmean of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 589 of file STK_Stat_Covariance.h.

590{
591#ifdef STK_BOUNDS_CHECK
592 if (V.cols()!=mean.range()) STKRUNTIME_ERROR_NO_ARG(covarianceWithFixedMean,V.cols()!=mean.range());
593#endif
594 CArraySquare<typename Array::Type, Array::sizeCols_, Array::orient_> cov_(V.cols());
595 for (int j= cov_.begin(); j< cov_.end(); j++)
596 {
597 cov_(j, j) = varianceWithFixedMean(V.col(j), mean[j], unbiased);
598 for (int i= cov_.begin(); i<j; i++)
599 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.col(i), V.col(j), mean[i], mean[j], unbiased));}
600 }
601 return cov_;
602}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), STKRUNTIME_ERROR_NO_ARG, and STK::Stat::varianceWithFixedMean().

◆ covarianceWithFixedMean() [2/4]

template<class Array , class Weights , class Mean >
CArraySquare< typename Array::Type, Array::sizeCols_, Array::orient_ > STK::Stat::covarianceWithFixedMean ( ExprBase< Array > const V,
Weights const W,
ExprBase< Mean > const mean,
bool  unbiased = false 
)

Compute the weighted covariance matrix using the column of the data set V.

\[ \hat{\sigma}^2 = \frac{1}{n}
                     \sum_{i=1}^n  w_{i} w_j (V_i-\hat{\mu}) (V_i-\hat{\mu})^T.
\]

Parameters
V,Wthe variables and the weights
meanthe (weighted) mean of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 643 of file STK_Stat_Covariance.h.

644{
645#ifdef STK_BOUNDS_CHECK
646 if (V.cols()!=mean.range()) STKRUNTIME_ERROR_NO_ARG(covarianceWithFixedMean,V.cols()!=mean.range());
647 if (W.range()!=mean.range()) STKRUNTIME_ERROR_NO_ARG(covarianceWithFixedMean,W.range()!=mean.range());
648#endif
649
650 CArraySquare<typename Array::Type, Array::sizeCols_, Array::orient_> cov_(V.cols(), Array::orient_);
651 for (int j= cov_.begin(); j< cov_.end(); j++)
652 {
653 cov_(j, j) = varianceWithFixedMean(V.col(j), W, mean[j], unbiased);
654 for (int i= cov_.begin(); i<j; i++)
655 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.col(i), V.col(j), W, mean[i], mean[j], unbiased));}
656 }
657 return cov_;
658}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), STKRUNTIME_ERROR_NO_ARG, and STK::Stat::varianceWithFixedMean().

◆ covarianceWithFixedMean() [3/4]

Real STK::Stat::covarianceWithFixedMean ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
ExprBase< Weights > const W,
typename hidden::Traits< XArray >::Type const xMean,
typename hidden::Traits< YArray >::Type const yMean,
bool  unbiased = false 
)

The covariance between X and Y will be computed using the formula.

See also
Mark Galassi, Jim Davies, James Theiler, Brian Gough, Gerard Jungman, Michael Booth, and Fabrice Rossi. http://www.gnu.org/software/gsl/manual GNU Scientific Library - Reference manual, Version 1.15, 2011. http://www.gnu.org/software/gsl/manual/html_node/Weighted-Samples.html Sec. 21.7 Weighted Samples

$ cov(X,Y)=\frac{\sum_{i=1}^{N}w_i}{\left(\sum_{i=1}^{N}w_i\right)^2-\sum_{i=1}^{N}w_i^2}
\sum_{i=1}^N w_i \left(  x_{i}-\bar{x} \right)  \left( y_{i}-\bar{y} \right). $

Parameters
X,Y,Wvariables and weights
xMean,yMeanmeans of the varaibles
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 371 of file STK_Stat_Covariance.h.

378{
382 typedef typename hidden::Traits<XArray>::Type XType;
383 typedef typename hidden::Traits<YArray>::Type YType;
384 typedef typename hidden::Promote<XType, YType>::result_type Type;
385
386#ifdef STK_BOUNDS_CHECK
387 if (X.range() != Y.range())
388 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceWithFixedMeanSafe(X,Y,W,xMu,yMu,unbiased),ranges are not the sames);
389#endif
390
391 Type xsum = 0.0, ysum = 0.0, xdev, ydev, sumWeights = 0.0, sum2Weights = 0.0, cov = 0.0;
392 for (int i=X.begin(); i<X.end(); i++)
393 {
394 xsum += (xdev = X[i] - xMean);
395 ysum += (ydev = Y[i] - yMean);
396 Type Wi = std::abs((Real)W[i]);
397 cov += Wi * (xdev*ydev);
398 sumWeights += Wi;
399 sum2Weights += Wi * Wi;
400 }
401 // compute the variance
402 return (unbiased) ?
403 (sumWeights*sumWeights > sum2Weights) ?
404 (cov - xsum*ysum)/(sumWeights - sum2Weights/sumWeights) : 0.
405// (cov - xsum*ysum/sumWeights)/(sumWeights - sum2Weights/sumWeights) : 0.
406 : (sumWeights) ? (cov - xsum*ysum)/(sumWeights) : 0.;
407}

References STK::Stat::covarianceWithFixedMeanSafe(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

◆ covarianceWithFixedMean() [4/4]

template<class XArray , class YArray >
Real STK::Stat::covarianceWithFixedMean ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
typename hidden::Traits< XArray >::Type const xMu,
typename hidden::Traits< YArray >::Type const yMu,
bool  unbiased = false 
)

Compute the covariance of the variables X and Y with fixed means.

\[ \hat{cov}(X,Y) = \frac{1}{n}
                \sum_{i=1}^n (X_i - \mu_X)(Y_i - \mu_Y).
\]

Parameters
X,Yvariables
xMu,yMumeans of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 264 of file STK_Stat_Covariance.h.

270{
273 typedef typename hidden::Traits<XArray>::Type XType;
274 typedef typename hidden::Traits<YArray>::Type YType;
275 typedef typename hidden::Promote<XType, YType>::result_type Type;
276
277#ifdef STK_BOUNDS_CHECK
278 if (X.range() != Y.range())
279 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceWithFixedMean(X,Y,xMu,yMu,unbiased),ranges are not the sames);
280#endif
281
282 int nobs = X.size();
283 Type xsum = 0.0, ysum = 0.0, cov = 0.0, xdev, ydev;
284 for (int i=X.begin(); i<X.end(); i++)
285 {
286 xsum += (xdev = X[i] - xMu);
287 ysum += (ydev = Y[i] - yMu);
288 cov += (xdev*ydev);
289 }
290 // compute the covariance with corrected bias
291 return (unbiased) ?
292 (nobs > 1) ? ((cov - (xsum*ysum)/Type(nobs))/Type(nobs -1))
293 : Arithmetic<Real>::infinity()
294
295 :
296 (nobs > 0) ? (cov - (xsum*ysum)/Type(nobs))/Type(nobs)
297 : Arithmetic<Real>::infinity();
298}

References STK::Stat::covarianceWithFixedMean(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

Referenced by STK::Stat::covariance(), STK::Stat::covariance(), STK::Stat::covarianceByRow(), STK::Stat::covarianceByRow(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMeanByRow(), STK::Stat::covarianceWithFixedMeanByRow(), STK::Stat::Multivariate< Array, Real >::run(), and STK::Stat::Multivariate< Array, Real >::run().

◆ covarianceWithFixedMeanByRow() [1/2]

template<class Array , class Mean >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceWithFixedMeanByRow ( ExprBase< Array > const V,
ExprBase< Mean > const mean,
bool  unbiased = false 
)

Compute the covariance matrix using the rows of the data set V.

\[ \hat{\Sigma} = \frac{1}{n}
                     \sum_{i=1}^n (V_i-mu)^T (V_i-mu).
\]

Parameters
Vvariable
meanmean of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 616 of file STK_Stat_Covariance.h.

617{
618#ifdef STK_BOUNDS_CHECK
619 if (V.rows()!=mean.range()) STKRUNTIME_ERROR_NO_ARG(covarianceWithFixedMean,V.rows()!=mean.range());
620#endif
621 CArraySquare<typename Array::Type, Array::sizeRows_, Array::orient_> cov_(V.rows());
622 for (int j= cov_.begin(); j< cov_.end(); j++)
623 {
624 cov_(j, j) = varianceWithFixedMean(V.row(j), mean[j], unbiased);
625 for (int i= cov_.begin(); i<j; i++)
626 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.row(i), V.row(j), mean[i], mean[j], unbiased));}
627 }
628 return cov_;
629}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), STKRUNTIME_ERROR_NO_ARG, and STK::Stat::varianceWithFixedMean().

◆ covarianceWithFixedMeanByRow() [2/2]

template<class Array , class Weights , class Mean >
CArraySquare< typename Array::Type, Array::sizeRows_, Array::orient_ > STK::Stat::covarianceWithFixedMeanByRow ( ExprBase< Array > const V,
Weights const W,
ExprBase< Mean > const mean,
bool  unbiased = false 
)

Compute the weighted covariance matrix using the rows of the data set V.

\[ \hat{\sigma}^2 = \frac{1}{n}
                     \sum_{i=1}^n  w_{i} w_j (V_i-\hat{\mu})^T (V_i-\hat{\mu}).
\]

Parameters
V,Wthe variables and the weights
meanthe (weighted) mean of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 672 of file STK_Stat_Covariance.h.

673{
674#ifdef STK_BOUNDS_CHECK
675 if (V.rows()!=mean.range()) STKRUNTIME_ERROR_NO_ARG(covarianceWithFixedMean,V.cols()!=mean.range());
676 if (W.range()!=mean.range()) STKRUNTIME_ERROR_NO_ARG(covarianceWithFixedMean,W.range()!=mean.range());
677#endif
678
679 CArraySquare<typename Array::Type, Array::sizeRows_, Array::orient_> cov_(V.rows());
680 for (int j= cov_.begin(); j< cov_.end(); j++)
681 {
682 cov_(j, j) = varianceWithFixedMean(V.row(j), W, mean[j], unbiased);
683 for (int i= cov_.begin(); i<j; i++)
684 { cov_(j,i) = ( cov_(i, j) = covarianceWithFixedMean(V.row(i), V.row(j), W, mean[i], mean[j], unbiased));}
685 }
686 return cov_;
687}

References STK::Stat::covarianceWithFixedMean(), STK::Stat::mean(), STKRUNTIME_ERROR_NO_ARG, and STK::Stat::varianceWithFixedMean().

◆ covarianceWithFixedMeanSafe() [1/2]

Real STK::Stat::covarianceWithFixedMeanSafe ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
ExprBase< Weights > const W,
typename hidden::Traits< XArray >::Type const xMu,
typename hidden::Traits< YArray >::Type const yMu,
bool  unbiased = false 
)

The covariance between X and Y will be computed using the formula.

See also
Mark Galassi, Jim Davies, James Theiler, Brian Gough, Gerard Jungman, Michael Booth, and Fabrice Rossi. http://www.gnu.org/software/gsl/manual GNU Scientific Library - Reference manual, Version 1.15, 2011. http://www.gnu.org/software/gsl/manual/html_node/Weighted-Samples.html Sec. 21.7 Weighted Samples

$ cov(X,Y)=\frac{\sum_{i=1}^{N}w_i}{\left(\sum_{i=1}^{N}w_i\right)^2-\sum_{i=1}^{N}w_i^2}
\sum_{i=1}^N w_i \left(  x_{i}-\mu_{x} \right)  \left( y_{i}-\mu_{y} \right). $

Parameters
X,Y,Wvariables and weights
xMu,yMumeans of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 428 of file STK_Stat_Covariance.h.

435{
439 typedef typename hidden::Traits<XArray>::Type XType;
440 typedef typename hidden::Traits<YArray>::Type YType;
441 typedef typename hidden::Promote<XType, YType>::result_type Type;
442
443#ifdef STK_BOUNDS_CHECK
444 if (X.range() != Y.range())
445 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceWithFixedMeanSafe(X,Y,W,xMu,yMu,unbiased),ranges are not the sames);
446#endif
447
448 // compute covariance
449 Type xsum = 0.0, ysum = 0.0, xdev, ydev, sumWeights = 0.0, sum2Weights = 0.0, cov = 0.0;
450 for (int i=X.begin(); i<X.end(); i++)
451 {
452 if (Arithmetic<Type>::isFinite(X[i])
453 && Arithmetic<Type>::isFinite(Y[i])
454 && Arithmetic<Type>::isFinite(W[i])
455 )
456 {
457 xsum += (xdev = X[i] - xMu); // deviation from the mean
458 ysum += (ydev = Y[i] - yMu); // deviation from the mean
459 Type Wi = std::abs((Real)W[i]);
460 cov += Wi * (xdev*ydev); // cross product
461 sumWeights += Wi; // sum absolute weights
462 sum2Weights += Wi * Wi; // sum squared weights
463 }
464 }
465 // compute the variance
466 return (unbiased) ?
467 (sumWeights*sumWeights > sum2Weights) ?
468 (cov - xsum*ysum)/(sumWeights - sum2Weights/sumWeights) : 0.
469// (cov - xsum*ysum/sumWeights)/(sumWeights - sum2Weights/sumWeights) : 0.
470 : (sumWeights) ? (cov - xsum*ysum)/(sumWeights) : 0.;
471}

References STK::Stat::covarianceWithFixedMeanSafe(), STK::Arithmetic< Type >::isFinite(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

◆ covarianceWithFixedMeanSafe() [2/2]

template<class XArray , class YArray >
Real STK::Stat::covarianceWithFixedMeanSafe ( ExprBase< XArray > const X,
ExprBase< YArray > const Y,
typename hidden::Traits< XArray >::Type const xMu,
typename hidden::Traits< YArray >::Type const yMu,
bool  unbiased = false 
)

Compute the covariance of the variables X and Y with fixed means.

\[ \hat{cov}(X,Y) = \frac{1}{n}
                \sum_{i=1}^n (X_i - \mu_X)(Y_i - \mu_Y).
\]

Parameters
X,Yvariables
xMu,yMumeans of the variables
unbiasedtrue if we want an unbiased estimate of the variance, false otherwise (default is false)

Definition at line 312 of file STK_Stat_Covariance.h.

318{
321 typedef typename hidden::Traits<XArray>::Type XType;
322 typedef typename hidden::Traits<YArray>::Type YType;
323 typedef typename hidden::Promote<XType, YType>::result_type Type;
324
325 #ifdef STK_BOUNDS_CHECK
326 if (X.range() != Y.range())
327 STKRUNTIME_ERROR_NO_ARG(Error in Stat::covarianceWithFixedMeanSafe(X,Y,xMu,yMu,unbiased),ranges are not the sames);
328 #endif
329
330 int nobs = X.size();
331 Type xsum = 0.0, ysum = 0.0, cov = 0.0, xdev, ydev;
332 for (int i=X.begin(); i<X.end(); i++)
333 {
334 if (Arithmetic<Type>::isFinite(X[i]) && Arithmetic<Type>::isFinite(Y[i]))
335 {
336 xsum += (xdev = X[i] - xMu); // deviation from the mean
337 ysum += (ydev = Y[i] - yMu); // deviation from the mean
338 cov += (xdev*ydev); // squared value
339 }
340 else nobs--;
341 }
342 // compute the covariance with corrected bias
343 return (unbiased) ?
344 (nobs > 1) ? ((cov - (xsum*ysum)/Type(nobs))/Type(nobs -1))
345 : Arithmetic<Real>::infinity()
346
347 :
348 (nobs > 0) ? (cov - (xsum*ysum)/Type(nobs))/Type(nobs)
349 : Arithmetic<Real>::infinity();
350}

References STK::Stat::covarianceWithFixedMeanSafe(), STK::Arithmetic< Type >::isFinite(), STK_STATIC_ASSERT_ONE_DIMENSION_ONLY, and STKRUNTIME_ERROR_NO_ARG.

Referenced by STK::Stat::covariance(), STK::Stat::covarianceSafe(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMeanSafe(), and STK::Stat::covarianceWithFixedMeanSafe().

◆ max()

template<class Derived >
hidden::FunctorTraits< Derived, MaxOp >::Row STK::Stat::max ( Derived const A)

Compute the maximal(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual max of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the maximal values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe array
Returns
the maximal value(s) of A or NA if there is no available value.

Definition at line 753 of file STK_Stat_Functors.h.

Referenced by main(), STK::Stat::ConfusionMatrix< TrueArray, PredArray >::operator()(), STK::Stat::MaxOp< Derived >::operator()(), STK::Stat::MaxSafeOp< Derived >::operator()(), STK::Stat::MaxOp< Derived >::operator()(), STK::Stat::MaxSafeOp< Derived >::operator()(), STK::Stat::Multivariate< Array, Real >::run(), and STK::Stat::Multivariate< Array, Real >::run().

◆ maxSafe()

template<class Derived >
hidden::FunctorTraits< Derived, MaxSafeOp >::Row STK::Stat::maxSafe ( Derived const A)

Compute safely the maximal(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual max of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the maximal values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe array
Returns
the maximal value(s) of A or NA if there is no available value.

Definition at line 765 of file STK_Stat_Functors.h.

◆ mean()

template<class Derived >
hidden::FunctorTraits< Derived, MeanOp >::Row STK::Stat::mean ( Derived const A)

Compute the mean(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual mean of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the mean values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
Returns
the mean(s) or NA if there is no available value

Definition at line 799 of file STK_Stat_Functors.h.

Referenced by STK::Stat::center(), STK::Stat::center(), STK::Stat::centerByCol(), STK::Stat::centerByCol(), STK::Stat::centerByRow(), STK::Stat::centerByRow(), STK::IGaussianModel< Array >::compMean(), STK::IGaussianModel< Array >::compWeightedMean(), STK::Stat::covariance(), STK::Stat::covariance(), STK::Stat::covariance(), STK::Stat::covariance(), STK::Stat::covarianceByRow(), STK::Stat::covarianceByRow(), STK::Stat::covarianceSafe(), STK::Stat::covarianceSafe(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMeanByRow(), STK::Stat::covarianceWithFixedMeanByRow(), main(), main(), STK::Stat::Multivariate< Array, Real >::run(), STK::Stat::Multivariate< Array, Real >::run(), STK::Stat::standardize(), STK::Stat::standardize(), STK::Stat::standardizeByCol(), STK::Stat::standardizeByCol(), STK::Stat::standardizeByRow(), STK::Stat::standardizeByRow(), STK::Stat::uncenter(), STK::Stat::uncenterByCol(), STK::Stat::uncenterByRow(), STK::Stat::unstandardize(), STK::Stat::unstandardizeByCol(), STK::Stat::unstandardizeByRow(), STK::Stat::varianceWithFixedMean(), STK::Stat::varianceWithFixedMean(), STK::Stat::varianceWithFixedMeanByCol(), STK::Stat::varianceWithFixedMeanByCol(), STK::Stat::varianceWithFixedMeanByRow(), STK::Stat::varianceWithFixedMeanByRow(), STK::Stat::varianceWithFixedMeanSafe(), STK::Stat::varianceWithFixedMeanSafe(), STK::Stat::varianceWithFixedMeanSafeByCol(), STK::Stat::varianceWithFixedMeanSafeByCol(), STK::Stat::varianceWithFixedMeanSafeByRow(), and STK::Stat::varianceWithFixedMeanSafeByRow().

◆ meanSafe()

template<class Derived >
hidden::FunctorTraits< Derived, MeanSafeOp >::Row STK::Stat::meanSafe ( Derived const A)

Compute safely the mean(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual mean of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the mean values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
Returns
the mean(s) or NA if there is no available value

Definition at line 811 of file STK_Stat_Functors.h.

Referenced by main().

◆ min()

template<class Derived >
hidden::FunctorTraits< Derived, MinOp >::Row STK::Stat::min ( Derived const A)

Compute the minimal(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual min of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the minimal values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe array
Returns
the minimal value(s) of A or NA if there is no available value.

Definition at line 726 of file STK_Stat_Functors.h.

Referenced by main(), STK::Stat::ConfusionMatrix< TrueArray, PredArray >::operator()(), STK::Stat::MinOp< Derived >::operator()(), STK::Stat::MinSafeOp< Derived >::operator()(), STK::Stat::MinOp< Derived >::operator()(), STK::Stat::MinSafeOp< Derived >::operator()(), STK::Stat::Multivariate< Array, Real >::run(), and STK::Stat::Multivariate< Array, Real >::run().

◆ minSafe()

template<class Derived >
hidden::FunctorTraits< Derived, MinSafeOp >::Row STK::Stat::minSafe ( Derived const A)

Compute safely the minimal(s) [weighted] value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual minimal value of the vector. If A is a two-dimensional array, the function will return (by value) an STK::Array2DPoint with the minimal values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe array
Returns
the minimal value(s) of A or NA if there is no available value.

Definition at line 740 of file STK_Stat_Functors.h.

◆ standardize() [1/2]

template<class Array , class RowVector >
void STK::Stat::standardize ( Array &  m,
RowVector &  mean,
RowVector &  std,
bool  unbiased = false 
)

Compute the mean and the standard deviation by columns of the variable m and standardize it.

Parameters
m,mean,stdarrays with the data, the means and the standard deviations
unbiasedfalse if the standard deviation as to be corrected, true otherwise

Definition at line 129 of file STK_Stat_Transform.h.

130{
131 typedef typename Array::Type Type;
132 enum
133 {
134 sizeRows_ = Array::sizeRows_,
135 sizeCols_ = Array::sizeCols_
136 };
137 // center V
138 centerByCol(m, mean);
139 // compute variance with mean 0 as V is centered
140 std = Stat::varianceWithFixedMeanByCol(m, Const::Vector<Type, sizeCols_>(m.cols())*Type(0), unbiased);
141 for (int j= m.beginCols(); j< m.endCols(); j++)
142 {
143 Real dev= std[j];
144 std[j] = (dev = Arithmetic<Real>::isFinite(dev) ? std::sqrt((double)dev) : 0.);
145 if (dev) { m.col(j) /= dev;}
146 }
147}
void centerByCol(Array &m, RowVector &mean)

References STK::Stat::centerByCol(), STK::Arithmetic< Type >::isFinite(), STK::Stat::mean(), and STK::Stat::varianceWithFixedMeanByCol().

Referenced by main(), STK::IAAModel< Array >::standardize(), STK::IAAModel< Array >::standardize(), STK::Stat::standardizeByCol(), and STK::Stat::standardizeByCol().

◆ standardize() [2/2]

template<class Array , class Weights , class RowVector >
void STK::Stat::standardize ( Array &  m,
Weights const W,
RowVector &  mean,
RowVector &  std,
bool  unbiased = false 
)

Compute the weighted means and standard deviations by columns of the variable m and standardize it.

Parameters
m,W,mean,stdcontainers with the Data, the weights, the means and standard deviations
unbiasedfalse if the standard deviation as to be corrected, true otherwise

Definition at line 188 of file STK_Stat_Transform.h.

189{
190 typedef typename Array::Type Type;
191 enum
192 {
193 sizeRows_ = Array::sizeRows_,
194 sizeCols_ = Array::sizeCols_
195 };
196 // center m
197 centerByCol(m, W, mean);
198 // compute variance with mean 0 as m is centered
199 std = Stat::varianceWithFixedMeanByCol(m, W, Const::Vector<Type, sizeCols_>(m.cols())*Type(0), unbiased);
200 // center
201 for (int j= m.beginCols(); j< m.endCols(); j++)
202 {
203 // compute standard deviation
204 Real dev = std[j];
205 // take square root and save result
206 std[j] = (dev = Arithmetic<Real>::isFinite(dev) ? std::sqrt((double)dev) : 0.);
207 // standardize data if necessary
208 if (dev) { m.col(j) /= std[j];}
209 }
210}

References STK::Stat::centerByCol(), STK::Arithmetic< Type >::isFinite(), STK::Stat::mean(), and STK::Stat::varianceWithFixedMeanByCol().

◆ standardizeByRow() [1/2]

template<class Array , class ColVector >
void STK::Stat::standardizeByRow ( Array &  m,
ColVector &  mean,
ColVector &  std,
bool  unbiased = false 
)

Compute the mean and the standard deviation by rows of the variable m and standardize it.

Parameters
m,mean,stdarrays with the data, the means and the standard deviations
unbiasedfalse if the standard deviation as to be corrected, true otherwise

Definition at line 160 of file STK_Stat_Transform.h.

161{
162 typedef typename Array::Type Type;
163 enum
164 {
165 sizeRows_ = Array::sizeRows_,
166 sizeCols_ = Array::sizeCols_
167 };
168 // center V
169 centerByRow(m, mean);
170 // compute variance with mean 0 as V is centered
171 std = Stat::varianceWithFixedMeanByRow(m, Const::Point<Type, sizeRows_>(m.rows())*Type(0), unbiased);
172 for (int i= m.beginRows(); i< m.endRows(); i++)
173 {
174 Real dev = std[i];
175 std[i] = (dev = Arithmetic<Real>::isFinite(dev) ? std::sqrt((double)dev) : 0.);
176 if (dev) { m.row(i) /= dev;}
177 }
178}
void centerByRow(Array &m, ColVector &mean)
Compute the mean by row of the variables in the container V and center V.

References STK::Stat::centerByRow(), STK::Arithmetic< Type >::isFinite(), STK::Stat::mean(), and STK::Stat::varianceWithFixedMeanByRow().

◆ standardizeByRow() [2/2]

template<class Array , class Weights , class ColVector >
void STK::Stat::standardizeByRow ( Array &  m,
Weights const W,
ColVector &  mean,
ColVector &  std,
bool  unbiased = false 
)

Compute the weighted means and standard deviations by rows of the variable m and standardize it.

Parameters
m,W,mean,stdcontainers with the Data, the weights, the means and standard deviations
unbiasedfalse if the standard deviation as to be corrected, true otherwise

Definition at line 223 of file STK_Stat_Transform.h.

224{
225 typedef typename Array::Type Type;
226 enum
227 {
228 sizeRows_ = Array::sizeRows_,
229 sizeCols_ = Array::sizeCols_
230 };
231 // center m
232 centerByRow(m, W, mean);
233 // compute variance with mean 0 as m is centered
234 std = Stat::varianceWithFixedMeanByRow(m, W, Const::Point<Type, sizeRows_>(m.rows())*Type(0), unbiased);
235 // center
236 for (int i= m.beginRows(); i< m.endRows(); i++)
237 {
238 // compute standard deviation
239 Real dev = std[i];
240 // take square root and save result
241 std[i] = (dev = Arithmetic<Real>::isFinite(dev) ? std::sqrt((double)dev) : 0.);
242 // standardize data if necessary
243 if (dev) { m.row(i) /= std[i];}
244 }
245}

References STK::Stat::centerByRow(), STK::Arithmetic< Type >::isFinite(), STK::Stat::mean(), and STK::Stat::varianceWithFixedMeanByRow().

◆ sum()

template<class Derived >
hidden::FunctorTraits< Derived, SumOp >::Row STK::Stat::sum ( Derived const A)

Compute the sum of A.

If A is a row-vector or a column-vector then the function will return the usual sum of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the sum of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min
Parameters
Athe data
Returns
the mean(s) or NA if there is no available value

Definition at line 776 of file STK_Stat_Functors.h.

Referenced by STK::Stat::Univariate< TContainer1D, Real >::compStatistics(), STK::Stat::Univariate< TContainer1D, Real >::compWeightedStatistics(), STK::IMixtureComposer::cStep(), STK::MixtureSemiLearner::cStep(), STK::Stat::SumOp< Derived >::operator()(), STK::Stat::SumSafeOp< Derived >::operator()(), STK::Stat::MeanOp< Derived >::operator()(), STK::Stat::MeanSafeOp< Derived >::operator()(), STK::Stat::VarianceOp< Derived >::operator()(), STK::Stat::VarianceSafeOp< Derived >::operator()(), STK::Stat::SumOp< Derived >::operator()(), STK::Stat::SumSafeOp< Derived >::operator()(), STK::Stat::MeanOp< Derived >::operator()(), STK::Stat::MeanSafeOp< Derived >::operator()(), STK::Stat::VarianceOp< Derived >::operator()(), STK::Stat::VarianceSafeOp< Derived >::operator()(), STK::Stat::VarianceWithFixedMeanOp< Derived >::operator()(), STK::Stat::VarianceWithFixedMeanSafeOp< Derived >::operator()(), STK::Stat::VarianceWithFixedMeanOp< Derived >::operator()(), and STK::Stat::VarianceWithFixedMeanSafeOp< Derived >::operator()().

◆ sumSafe()

template<class Derived >
hidden::FunctorTraits< Derived, SumSafeOp >::Row STK::Stat::sumSafe ( Derived const A)

Compute safely the mean(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual mean of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the mean values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
Returns
the mean(s) or NA if there is no available value

Definition at line 787 of file STK_Stat_Functors.h.

◆ unstandardize() [1/2]

template<class Array , class RowVector >
void STK::Stat::unstandardize ( Array &  m,
RowVector const mean,
RowVector const std 
)

undo the standardization by columns of the standardized variable m

Parameters
m,mean,stdthe containers with the Data, the means and the standard deviations

Definition at line 331 of file STK_Stat_Transform.h.

332{
333 unstandardizeByCol(m, std); // unstandardize
334 uncenterByCol(m, mean); // uncenter
335}
void unstandardizeByCol(Array &m, RowVector const &std)
void uncenterByCol(Array &m, RowVector const &mean)

References STK::Stat::mean(), STK::Stat::uncenterByCol(), and STK::Stat::unstandardizeByCol().

◆ unstandardize() [2/2]

template<class Array , class RowVector >
void STK::Stat::unstandardize ( Array &  m,
RowVector const std 
)

undo the standardization by columns of the standardized variable m.

Parameters
m,stdthe container with the data and the standard deviations

Definition at line 292 of file STK_Stat_Transform.h.

293{
294 typedef typename Array::Type Type;
295 enum
296 {
297 sizeRows_ = Array::sizeRows_,
298 sizeCols_ = Array::sizeCols_
299 };
300 if (m.cols() != std.range())
301 STKRUNTIME_ERROR_NO_ARG(Error in Stat::unstandardize(m, std),ranges are not the sames);
302 m = m.prod(Const::Vector<Type, sizeRows_>(m.rows())*std);
303}

References STKRUNTIME_ERROR_NO_ARG, and STK::Stat::unstandardize().

Referenced by main(), STK::IAAModel< Array >::standardize(), STK::IAAModel< Array >::standardize(), STK::Stat::unstandardize(), STK::Stat::unstandardizeByCol(), STK::Stat::unstandardizeByCol(), and STK::IAAModel< Array >::unstandardizeResults().

◆ unstandardizeByRow() [1/2]

template<class Array , class ColVector >
void STK::Stat::unstandardizeByRow ( Array &  m,
ColVector const mean,
ColVector const std 
)

undo the standardization by rows of the standardized variable m

Parameters
m,mean,stdthe containers with the Data, the means and the standard deviations

Definition at line 345 of file STK_Stat_Transform.h.

346{
347 unstandardizeByRow(m, std); // unstandardize
348 uncenterByRow(m, mean); // uncenter
349}
void unstandardizeByRow(Array &m, ColVector const &std)
undo the standardization by rows of the standardized variable m.
void uncenterByRow(Array &m, ColVector const &mean)
Add the means to the rows of the container m.

References STK::Stat::mean(), STK::Stat::uncenterByRow(), and STK::Stat::unstandardizeByRow().

◆ unstandardizeByRow() [2/2]

template<class Array , class ColVector >
void STK::Stat::unstandardizeByRow ( Array &  m,
ColVector const std 
)

undo the standardization by rows of the standardized variable m.

Parameters
m,stdthe container with the data and the standard deviations

Definition at line 313 of file STK_Stat_Transform.h.

314{
315 typedef typename Array::Type Type;
316 enum
317 {
318 sizeRows_ = Array::sizeRows_,
319 sizeCols_ = Array::sizeCols_
320 };
321 if (m.rows() != std.range())
322 STKRUNTIME_ERROR_NO_ARG(Error in Stat::unstandardizeByRow(m, std),ranges are not the sames);
323 m = m.prod(std * Const::Point<Type, sizeCols_>(m.cols()));
324}

References STKRUNTIME_ERROR_NO_ARG, and STK::Stat::unstandardizeByRow().

Referenced by STK::Stat::unstandardizeByRow(), and STK::Stat::unstandardizeByRow().

◆ variance()

template<class Derived >
hidden::FunctorTraits< Derived, VarianceOp >::Row STK::Stat::variance ( Derived const A,
bool  unbiased = false 
)

Compute the variance(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual variance of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the variance values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
unbiasedthe unbiased variance(s) if true or the Maximum-likelihood variance otherwise (the default)
Returns
the variance(s) or NA if there is no available value

Definition at line 826 of file STK_Stat_Functors.h.

827{ return typename hidden::FunctorTraits<Derived, VarianceOp>::ColOp(A)(unbiased);}

Referenced by main().

◆ varianceSafe()

template<class Derived >
hidden::FunctorTraits< Derived, VarianceSafeOp >::Row STK::Stat::varianceSafe ( Derived const A,
bool  unbiased = false 
)

Compute safely the variance(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual variance of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the variance values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
unbiasedthe unbiased variance(s) if true or the Maximum-likelihood variance otherwise (the default)
Returns
the variance(s) or NA if there is no available value

Definition at line 869 of file STK_Stat_Functors.h.

870{ return typename hidden::FunctorTraits<Derived, VarianceSafeOp>::ColOp(A)(unbiased);}

◆ varianceWithFixedMean()

template<class Derived , class MeanType >
hidden::FunctorTraits< Derived, VarianceWithFixedMeanOp >::Row STK::Stat::varianceWithFixedMean ( Derived const A,
MeanType const mean,
bool  unbiased 
)

Compute the VarianceWithFixedMean(s) value(s) of A.

If A is a row-vector or a column-vector then the function will return the usual variance of the vector. If A is a two-dimensional array, the function will return (by value) an Array2DPoint with the variance values of each columns.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
meanThe mean (s) to use
unbiasedthe unbiased variance(s) if true or the Maximum-likelihood variance otherwise (the default)
Returns
the variance(s) or NA if there is no available value

Definition at line 912 of file STK_Stat_Functors.h.

913{ return typename hidden::FunctorTraits<Derived, VarianceWithFixedMeanOp>::ColOp(A)(mean, unbiased);}

References STK::Stat::mean().

Referenced by STK::JointGaussianModel< Array, WColVector >::computeParameters(), STK::JointGaussianModel< Array, WColVector >::computeParameters(), STK::Stat::covariance(), STK::Stat::covariance(), STK::Stat::covarianceByRow(), STK::Stat::covarianceByRow(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMean(), STK::Stat::covarianceWithFixedMeanByRow(), STK::Stat::covarianceWithFixedMeanByRow(), main(), STK::DiagGaussian_sjk< Array >::randomInit(), STK::DiagGaussian_sjsk< Array >::randomInit(), STK::HDGaussian_AjkBkQkD< Array >::randomInit(), STK::Stat::Multivariate< Array, Real >::run(), STK::DiagGaussian_sjk< Array >::run(), STK::DiagGaussian_sjsk< Array >::run(), STK::HDGaussian_AjkBkQkD< Array >::run(), and STK::Stat::Multivariate< Array, Real >::run().

◆ varianceWithFixedMeanSafe()

template<class Derived , class MeanType >
hidden::FunctorTraits< Derived, VarianceWithFixedMeanSafeOp >::Row STK::Stat::varianceWithFixedMeanSafe ( Derived const A,
MeanType const mean,
bool  unbiased 
)

Compute safely the VarianceWithFixedMean(s) value(s) of A.

See also
STK::Stat::mean, STK::Stat::max, STK::Stat::variance, STK::Stat::varianceWithFixedMean, STK::Stat::min, STK::Stat::sum.
Parameters
Athe data
meanThe mean (s) to use
unbiasedthe unbiased variance(s) if true or the Maximum-likelihood variance otherwise (the default)
Returns
the variance(s) or NA if there is no available value

Definition at line 951 of file STK_Stat_Functors.h.

952{ return typename hidden::FunctorTraits<Derived, VarianceWithFixedMeanSafeOp>::ColOp(A)(mean, unbiased);}

References STK::Stat::mean().