fml  0.1-0
Fused Matrix Library
fml::stats Namespace Reference

Statistics kernels. More...

Functions

template<typename REAL >
void cor (cpumat< REAL > &x, cpumat< REAL > &cov)
 Covariance. More...
 
template<typename REAL >
void cor (cpumat< REAL > &x, cpumat< REAL > &y, cpumat< REAL > &cov)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
REAL cor (const cpuvec< REAL > &x, const cpuvec< REAL > &y)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void cov (cpumat< REAL > &x, cpumat< REAL > &cov)
 Covariance. More...
 
template<typename REAL >
void cov (cpumat< REAL > &x, cpumat< REAL > &y, cpumat< REAL > &cov)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
REAL cov (const cpuvec< REAL > &x, const cpuvec< REAL > &y)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void pca (const bool rm_mean, const bool rm_sd, cpumat< REAL > &x, cpuvec< REAL > &sdev, cpumat< REAL > &rot)
 Principal components analysis. More...
 
template<typename REAL >
void pca (const bool rm_mean, const bool rm_sd, cpumat< REAL > &x, cpuvec< REAL > &sdev)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void cor (gpumat< REAL > &x, gpumat< REAL > &cov)
 Covariance. More...
 
template<typename REAL >
void cor (gpumat< REAL > &x, gpumat< REAL > &y, gpumat< REAL > &cov)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
REAL cor (const gpuvec< REAL > &x, const gpuvec< REAL > &y)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void cov (gpumat< REAL > &x, gpumat< REAL > &cov)
 Covariance. More...
 
template<typename REAL >
void cov (gpumat< REAL > &x, gpumat< REAL > &y, gpumat< REAL > &cov)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
REAL cov (const gpuvec< REAL > &x, const gpuvec< REAL > &y)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void pca (const bool rm_mean, const bool rm_sd, gpumat< REAL > &x, gpuvec< REAL > &sdev, gpumat< REAL > &rot)
 Principal components analysis. More...
 
template<typename REAL >
void pca (const bool rm_mean, const bool rm_sd, gpumat< REAL > &x, gpuvec< REAL > &sdev)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void cor (mpimat< REAL > &x, mpimat< REAL > &cov)
 Covariance. More...
 
template<typename REAL >
void cor (mpimat< REAL > &x, mpimat< REAL > &y, mpimat< REAL > &cov)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void cov (mpimat< REAL > &x, mpimat< REAL > &cov)
 Covariance. More...
 
template<typename REAL >
void cov (mpimat< REAL > &x, mpimat< REAL > &y, mpimat< REAL > &cov)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 
template<typename REAL >
void pca (const bool rm_mean, const bool rm_sd, mpimat< REAL > &x, cpuvec< REAL > &sdev, mpimat< REAL > &rot)
 Principal components analysis. More...
 
template<typename REAL >
void pca (const bool rm_mean, const bool rm_sd, mpimat< REAL > &x, cpuvec< REAL > &sdev)
 This is an overloaded member function, provided for convenience. It differs from the above function only in what argument(s) it accepts.
 

Detailed Description

Statistics kernels.

Function Documentation

◆ cor() [1/3]

template<typename REAL >
void fml::stats::cor ( cpumat< REAL > &  x,
cpumat< REAL > &  cov 
)

Covariance.

Computes the lower triangle of the Pearson correlation matrix. Centering is done in-place.

Parameters
[in,out]x,yInput data. For the matrix variants, data is mean-centered on return.
[out]covThe correlation matrix.

Exceptions
In the matrix-matrix and vector-vector variants, if the object dimensions/sizes are non-conformable, a runtime_error exception is thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ cor() [2/3]

template<typename REAL >
void fml::stats::cor ( gpumat< REAL > &  x,
gpumat< REAL > &  cov 
)

Covariance.

Computes the lower triangle of the Pearson correlation matrix. Centering is done in-place.

Parameters
[in,out]x,yInput data. For the matrix variants, data is mean-centered on return.
[out]covThe correlation matrix.

Exceptions
In the matrix-matrix and vector-vector variants, if the object dimensions/sizes are non-conformable, a runtime_error exception is thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ cor() [3/3]

template<typename REAL >
void fml::stats::cor ( mpimat< REAL > &  x,
mpimat< REAL > &  cov 
)

Covariance.

Computes the lower triangle of the Pearson correlation matrix. Centering is done in-place.

Parameters
[in,out]x,yInput data. For the matrix variants, data is mean-centered on return.
[out]covThe correlation matrix.

Exceptions
In the matrix-matrix and vector-vector variants, if the object dimensions/sizes are non-conformable, a runtime_error exception is thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ cov() [1/3]

template<typename REAL >
void fml::stats::cov ( cpumat< REAL > &  x,
cpumat< REAL > &  cov 
)

Covariance.

Computes the lower triangle of the variance-covariance matrix. Centering is done in-place.

Parameters
[in,out]x,yInput data. For the matrix variants, data is mean-centered on return.
[out]covThe covariance matrix.

Exceptions
In the matrix-matrix and vector-vector variants, if the object dimensions/sizes are non-conformable, a runtime_error exception is thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ cov() [2/3]

template<typename REAL >
void fml::stats::cov ( gpumat< REAL > &  x,
gpumat< REAL > &  cov 
)

Covariance.

Computes the lower triangle of the variance-covariance matrix. Centering is done in-place.

Parameters
[in,out]x,yInput data. For the matrix variants, data is mean-centered on return.
[out]covThe covariance matrix.

Exceptions
In the matrix-matrix and vector-vector variants, if the object dimensions/sizes are non-conformable, a runtime_error exception is thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ cov() [3/3]

template<typename REAL >
void fml::stats::cov ( mpimat< REAL > &  x,
mpimat< REAL > &  cov 
)

Covariance.

Computes the lower triangle of the variance-covariance matrix. Centering is done in-place.

Parameters
[in,out]x,yInput data. For the matrix variants, data is mean-centered on return.
[out]covThe covariance matrix.

Exceptions
In the matrix-matrix and vector-vector variants, if the object dimensions/sizes are non-conformable, a runtime_error exception is thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ pca() [1/3]

template<typename REAL >
void fml::stats::pca ( const bool  rm_mean,
const bool  rm_sd,
cpumat< REAL > &  x,
cpuvec< REAL > &  sdev,
cpumat< REAL > &  rot 
)

Principal components analysis.

Parameters
[in]rm_mean,rm_sdShould the column means/sds be removed?
[in,out]xInput data. Values are overwritten.
[out]sdevStandard deviations of the principal components.
[out]rotThe variable loadings.

Implementation Details
Uses linalg::svd().

Memory Allocations
If the dimensions of the outputs are inappropriately sized, they will automatically be re-allocated.

Exceptions
If a reallocation is triggered and fails, a bad_alloc exception will be thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ pca() [2/3]

template<typename REAL >
void fml::stats::pca ( const bool  rm_mean,
const bool  rm_sd,
gpumat< REAL > &  x,
gpuvec< REAL > &  sdev,
gpumat< REAL > &  rot 
)

Principal components analysis.

Parameters
[in]rm_mean,rm_sdShould the column means/sds be removed?
[in,out]xInput data. Values are overwritten.
[out]sdevStandard deviations of the principal components.
[out]rotThe variable loadings.

Implementation Details
Uses linalg::svd().

Memory Allocations
If the dimensions of the outputs are inappropriately sized, they will automatically be re-allocated.

Exceptions
If a reallocation is triggered and fails, a bad_alloc exception will be thrown.

Template Parameters
REALshould be 'float' or 'double'.

◆ pca() [3/3]

template<typename REAL >
void fml::stats::pca ( const bool  rm_mean,
const bool  rm_sd,
mpimat< REAL > &  x,
cpuvec< REAL > &  sdev,
mpimat< REAL > &  rot 
)

Principal components analysis.

Parameters
[in]rm_mean,rm_sdShould the column means/sds be removed?
[in,out]xInput data. Values are overwritten.
[out]sdevStandard deviations of the principal components.
[out]rotThe variable loadings.

Implementation Details
Uses linalg::svd().

Memory Allocations
If the dimensions of the outputs are inappropriately sized, they will automatically be re-allocated.

Exceptions
If a reallocation is triggered and fails, a bad_alloc exception will be thrown.

Template Parameters
REALshould be 'float' or 'double'.