bpp-core  2.2.0
bpp::MatrixTools Class Reference

Functions dealing with matrices. More...

#include <Bpp/Numeric/Matrix/MatrixTools.h>

Public Member Functions

 MatrixTools ()
 
 ~MatrixTools ()
 

Static Public Member Functions

template<class MatrixA , class MatrixO >
static void copy (const MatrixA &A, MatrixO &O)
 Copy operation. This function supplies the lack of inheritence of the assigment operator :D . More...
 
template<class Matrix >
static void getId (size_t n, Matrix &O)
 Get a identity matrix of a given size. More...
 
template<class Scalar >
static void diag (const std::vector< Scalar > &D, Matrix< Scalar > &O)
 
template<class Scalar >
static void diag (const Scalar x, size_t n, Matrix< Scalar > &O)
 
template<class Scalar >
static void diag (const Matrix< Scalar > &M, std::vector< Scalar > &O) throw (DimensionException)
 
template<class Matrix , class Scalar >
static void fill (Matrix &M, Scalar x)
 Set all elements in M to value x. More...
 
template<class Matrix , class Scalar >
static void scale (Matrix &A, Scalar a, Scalar b=0)
 Multiply all elements of a matrix by a given value, and add a constant. More...
 
template<class Scalar >
static void mult (const Matrix< Scalar > &A, const Matrix< Scalar > &B, Matrix< Scalar > &O) throw (DimensionException)
 
template<class Scalar >
static void mult (const Matrix< Scalar > &A, const std::vector< Scalar > &D, const Matrix< Scalar > &B, Matrix< Scalar > &O) throw (DimensionException)
 Compute A . D . B where D is a diagonal matrix in O(n^3). More...
 
template<class Scalar >
static void mult (const Matrix< Scalar > &A, const std::vector< Scalar > &D, const std::vector< Scalar > &U, const std::vector< Scalar > &L, const Matrix< Scalar > &B, Matrix< Scalar > &O) throw (DimensionException)
 Compute A . (U+D+L) . B where D is a diagonal matrix, U (resp. L) is a matrix in which the only non-zero terms are on the diagonal that is over (resp. under) the main diagonal, in O(n^3). More...
 
template<class MatrixA , class MatrixB >
static void add (MatrixA &A, const MatrixB &B) throw (DimensionException)
 Add matrix B to matrix A. More...
 
template<class MatrixA , class MatrixB , class Scalar >
static void add (MatrixA &A, Scalar &x, const MatrixB &B) throw (DimensionException)
 Add matrix x.B to matrix A. More...
 
template<class Matrix >
static void pow (const Matrix &A, size_t p, Matrix &O) throw (DimensionException)
 Compute the power of a given matrix. More...
 
template<class Scalar >
static void pow (const Matrix< Scalar > &A, double p, Matrix< Scalar > &O) throw (DimensionException)
 Compute the power of a given matrix, using eigen value decomposition. More...
 
template<class Scalar >
static void exp (const Matrix< Scalar > &A, Matrix< Scalar > &O) throw (DimensionException)
 Perform matrix exponentiation using diagonalization. More...
 
template<class Matrix , class Scalar >
static void Taylor (const Matrix &A, size_t p, std::vector< RowMatrix< Scalar > > &vO) throw (DimensionException)
 Compute a vector of the first powers of a given matrix. More...
 
template<class Matrix >
static std::vector< size_t > whichMax (const Matrix &m)
 
template<class Matrix >
static std::vector< size_t > whichMin (const Matrix &m)
 
template<class Real >
static Real max (const Matrix< Real > &m)
 
template<class Real >
static Real min (const Matrix< Real > &m)
 
template<class Matrix >
static void print (const Matrix &m, std::ostream &out=std::cout)
 Print a matrix to a stream. More...
 
template<class Matrix >
static void print (const Matrix &m, bpp::OutputStream &out, char pIn='(', char pOut=')')
 Print a matrix to a stream. More...
 
template<class Matrix >
static void printForR (const Matrix &m, const std::string &variableName="x", std::ostream &out=std::cout)
 Print a matrix to a stream, so that it is read by R. More...
 
template<class Real >
static void print (const std::vector< Real > &v, std::ostream &out=std::cout)
 Print a vector to a stream. More...
 
template<class Matrix >
static bool isSquare (const Matrix &A)
 
template<class Scalar >
static Scalar inv (const Matrix< Scalar > &A, Matrix< Scalar > &O) throw (DimensionException, ZeroDivisionException)
 
template<class Scalar >
static double det (const Matrix< Scalar > &A) throw (DimensionException)
 Get determinant of a square matrix. More...
 
template<class MatrixA , class MatrixO >
static void transpose (const MatrixA &A, MatrixO &O)
 
template<class Scalar >
static void covar (const Matrix< Scalar > &A, Matrix< Scalar > &O)
 Compute the variance-covariance matrix of an input matrix. More...
 
template<class Scalar >
static void kroneckerMult (const Matrix< Scalar > &A, const Matrix< Scalar > &B, Matrix< Scalar > &O)
 Compute the Kronecker product of two row matrices. More...
 
template<class Scalar >
static void hadamardMult (const Matrix< Scalar > &A, const Matrix< Scalar > &B, Matrix< Scalar > &O)
 Compute the Hadamard product of two row matrices with same dimensions. More...
 
template<class Scalar >
static void hadamardMult (const Matrix< Scalar > &A, const std::vector< Scalar > &B, Matrix< Scalar > &O, bool row=true)
 Compute the "Hadamard" product of a row matrix and a vector containing weights, according to rows or columns. More...
 
template<class Scalar >
static void kroneckerSum (const Matrix< Scalar > &A, const Matrix< Scalar > &B, Matrix< Scalar > &O)
 Compute the Kronecker sum of two row matrices. More...
 
template<class Scalar >
static void kroneckerSum (const std::vector< Matrix< Scalar > *> &vA, Matrix< Scalar > &O)
 Compute the Kronecker sum of n row matrices. More...
 
template<class Scalar >
static void toVVdouble (const Matrix< Scalar > &M, std::vector< std::vector< Scalar > > &vO)
 Convert to a vector of vector. More...
 
template<class Scalar >
static Scalar sumElements (const Matrix< Scalar > &M)
 Sum all elements in M. More...
 
template<class Scalar >
static Scalar lap (Matrix< Scalar > &assignCost, std::vector< int > &rowSol, std::vector< int > &colSol, std::vector< Scalar > &u, std::vector< Scalar > &v) throw (Exception)
 Linear Assignment Problem. More...
 

Detailed Description

Functions dealing with matrices.

Definition at line 57 of file MatrixTools.h.

Constructor & Destructor Documentation

◆ MatrixTools()

bpp::MatrixTools::MatrixTools ( )
inline

Definition at line 60 of file MatrixTools.h.

◆ ~MatrixTools()

bpp::MatrixTools::~MatrixTools ( )
inline

Definition at line 61 of file MatrixTools.h.

Member Function Documentation

◆ add() [1/2]

template<class MatrixA , class MatrixB >
static void bpp::MatrixTools::add ( MatrixA &  A,
const MatrixB &  B 
)
throw (DimensionException
)
inlinestatic

Add matrix B to matrix A.

Parameters
A[in] Matrix A
B[in] Matrix B
Exceptions
DimensionExceptionIf A and B have note the same size.

Definition at line 300 of file MatrixTools.h.

Referenced by covar().

◆ add() [2/2]

template<class MatrixA , class MatrixB , class Scalar >
static void bpp::MatrixTools::add ( MatrixA &  A,
Scalar &  x,
const MatrixB &  B 
)
throw (DimensionException
)
inlinestatic

Add matrix x.B to matrix A.

Parameters
A[in,out] Matrix A
x[in] Scalar x
B[in] Matrix B
Exceptions
DimensionExceptionIf A and B have note the same size.

Definition at line 326 of file MatrixTools.h.

◆ copy()

template<class MatrixA , class MatrixO >
static void bpp::MatrixTools::copy ( const MatrixA &  A,
MatrixO &  O 
)
inlinestatic

Copy operation. This function supplies the lack of inheritence of the assigment operator :D .

Parameters
A[in] Original matrix.
O[out] A copy of the given matrix.

Definition at line 71 of file MatrixTools.h.

Referenced by pow(), and Taylor().

◆ covar()

template<class Scalar >
static void bpp::MatrixTools::covar ( const Matrix< Scalar > &  A,
Matrix< Scalar > &  O 
)
inlinestatic

Compute the variance-covariance matrix of an input matrix.

The input matrix represent a n-sample of a random vector of dimension r. It is assumed to have r rows and n columns. The variance matrix is then computed as

\[ V = A\cdot A^T - \mu\cdot\mu^T\]

, where $\mu$ is the mean vector of the sample. the output matrix is a square matrix of size r.

Parameters
A[in] The intput matrix.
O[out] The resulting variance covariance matrix.

Definition at line 731 of file MatrixTools.h.

References add(), bpp::Matrix< Scalar >::getNumberOfColumns(), bpp::Matrix< Scalar >::getNumberOfRows(), mult(), bpp::Matrix< Scalar >::resize(), scale(), and transpose().

Referenced by bpp::AdaptiveKernelDensityEstimation::init_().

◆ det()

template<class Scalar >
static double bpp::MatrixTools::det ( const Matrix< Scalar > &  A)
throw (DimensionException
)
inlinestatic

Get determinant of a square matrix.

This implementation is in $o(n^3)$ and uses the LU decomposition method.

Parameters
A[in] The input matrix.
Returns
The determinant of A.
Exceptions
DimensionExceptionIf A is not a square matrix.

Definition at line 694 of file MatrixTools.h.

References isSquare().

Referenced by bpp::AdaptiveKernelDensityEstimation::init_().

◆ diag() [1/3]

template<class Scalar >
static void bpp::MatrixTools::diag ( const std::vector< Scalar > &  D,
Matrix< Scalar > &  O 
)
inlinestatic
Parameters
D[in] A vector of diagonal elements.
O[out] A diagonal matrix with diagonal elements taken from a vector.

Definition at line 106 of file MatrixTools.h.

References bpp::Matrix< Scalar >::resize().

◆ diag() [2/3]

template<class Scalar >
static void bpp::MatrixTools::diag ( const Scalar  x,
size_t  n,
Matrix< Scalar > &  O 
)
inlinestatic
Parameters
x[in] A scalar
n[in] the dimension of the output matrix
O[out] A diagonal matrix with diagonal elements equal to x

Definition at line 122 of file MatrixTools.h.

References bpp::Matrix< Scalar >::resize().

◆ diag() [3/3]

template<class Scalar >
static void bpp::MatrixTools::diag ( const Matrix< Scalar > &  M,
std::vector< Scalar > &  O 
)
throw (DimensionException
)
inlinestatic
Parameters
M[in] The matrix.
O[out] The diagonal elements of a square matrix as a vector.
Exceptions
DimensionExceptionIf M is not a square matrix.

Definition at line 137 of file MatrixTools.h.

◆ exp()

template<class Scalar >
static void bpp::MatrixTools::exp ( const Matrix< Scalar > &  A,
Matrix< Scalar > &  O 
)
throw (DimensionException
)
inlinestatic

Perform matrix exponentiation using diagonalization.

Warning
This method currently relies only on diagonalization, so it won't work if your matrix is not diagonalizable. The function may be extended later to deal with other cases.
Parameters
A[in] The matrix.
O[out] $\prod_{i=1}^p m$.
Exceptions
DimensionExceptionIf m is not a square matrix.

Definition at line 416 of file MatrixTools.h.

References bpp::VectorTools::exp(), bpp::EigenValue< Real >::getRealEigenValues(), bpp::EigenValue< Real >::getV(), inv(), and mult().

◆ fill()

template<class Matrix , class Scalar >
static void bpp::MatrixTools::fill ( Matrix M,
Scalar  x 
)
inlinestatic

Set all elements in M to value x.

Parameters
MA matrix.
xThe value to use.

Definition at line 152 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ getId()

template<class Matrix >
static void bpp::MatrixTools::getId ( size_t  n,
Matrix O 
)
inlinestatic

Get a identity matrix of a given size.

Parameters
nthe size of the matrix.
O[out] A identity matrix of size n.

Definition at line 90 of file MatrixTools.h.

References bpp::Matrix< Scalar >::resize().

Referenced by inv().

◆ hadamardMult() [1/2]

template<class Scalar >
static void bpp::MatrixTools::hadamardMult ( const Matrix< Scalar > &  A,
const Matrix< Scalar > &  B,
Matrix< Scalar > &  O 
)
inlinestatic

Compute the Hadamard product of two row matrices with same dimensions.

Parameters
A[in] The first row matrix.
B[in] The second row matrix.
O[out] The Hadamard product.

Definition at line 796 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), bpp::Matrix< Scalar >::getNumberOfRows(), and bpp::Matrix< Scalar >::resize().

Referenced by bpp::DualityDiagram::compute_(), and bpp::CorrespondenceAnalysis::CorrespondenceAnalysis().

◆ hadamardMult() [2/2]

template<class Scalar >
static void bpp::MatrixTools::hadamardMult ( const Matrix< Scalar > &  A,
const std::vector< Scalar > &  B,
Matrix< Scalar > &  O,
bool  row = true 
)
inlinestatic

Compute the "Hadamard" product of a row matrix and a vector containing weights, according to rows or columns.

Parameters
A[in] The row matrix.
B[in] The vector of row or column weights.
O[out] The 'Hadamard' product.
rowBoolean. If row is set to 'true', the vector contains weights for rows. Otherwise the vector contains weights for columns.

Definition at line 823 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), bpp::Matrix< Scalar >::getNumberOfRows(), and bpp::Matrix< Scalar >::resize().

◆ inv()

template<class Scalar >
static Scalar bpp::MatrixTools::inv ( const Matrix< Scalar > &  A,
Matrix< Scalar > &  O 
)
throw (DimensionException,
ZeroDivisionException
)
inlinestatic
Parameters
A[in] The matrix to inverse.
O[out] The inverse matrix of A.
Returns
x the minimum absolute value of the diagonal of the LU decomposition
Exceptions
DimensionExceptionIf A is not a square matrix.

Definition at line 675 of file MatrixTools.h.

References getId(), and isSquare().

Referenced by exp(), and pow().

◆ isSquare()

template<class Matrix >
static bool bpp::MatrixTools::isSquare ( const Matrix A)
inlinestatic
Returns
True if the matrix is a square matrix.
Parameters
AA matrix.

Definition at line 666 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

Referenced by det(), and inv().

◆ kroneckerMult()

template<class Scalar >
static void bpp::MatrixTools::kroneckerMult ( const Matrix< Scalar > &  A,
const Matrix< Scalar > &  B,
Matrix< Scalar > &  O 
)
inlinestatic

Compute the Kronecker product of two row matrices.

Parameters
A[in] The first row matrix.
B[in] The second row matrix.
O[out] The product $A \otimes B$.

Definition at line 765 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), bpp::Matrix< Scalar >::getNumberOfRows(), and bpp::Matrix< Scalar >::resize().

◆ kroneckerSum() [1/2]

template<class Scalar >
static void bpp::MatrixTools::kroneckerSum ( const Matrix< Scalar > &  A,
const Matrix< Scalar > &  B,
Matrix< Scalar > &  O 
)
inlinestatic

Compute the Kronecker sum of two row matrices.

Parameters
A[in] The first row matrix.
B[in] The second row matrix.
O[out] The product $A \oplus B$.

Definition at line 861 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), bpp::Matrix< Scalar >::getNumberOfRows(), and bpp::Matrix< Scalar >::resize().

◆ kroneckerSum() [2/2]

template<class Scalar >
static void bpp::MatrixTools::kroneckerSum ( const std::vector< Matrix< Scalar > *> &  vA,
Matrix< Scalar > &  O 
)
inlinestatic

Compute the Kronecker sum of n row matrices.

Parameters
vA[in] A vector of row matrices of any size.
O[out] The product $\bigoplus_i A_i$.

Definition at line 891 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), bpp::Matrix< Scalar >::getNumberOfRows(), and bpp::Matrix< Scalar >::resize().

◆ lap()

template<class Scalar >
static Scalar bpp::MatrixTools::lap ( Matrix< Scalar > &  assignCost,
std::vector< int > &  rowSol,
std::vector< int > &  colSol,
std::vector< Scalar > &  u,
std::vector< Scalar > &  v 
)
throw (Exception
)
inlinestatic

Linear Assignment Problem.

The algorithm coded here is described in

  • A Shortest Augmenting Path Algorithm for Dense and Sparse Linear Assignment Problems, Computing 38, 325-340, 1987 by R. Jonker and A. Volgenant, University of Amsterdam.
Parameters
assignCost[input/output] Cost matrix
rowSol[output] Column assigned to row in solution
colSol[output] Row assigned to column in solution
u[output] Dual variables, row reduction numbers
v[output] Dual variables, column reduction numbers
Returns
The optimal cost.

Definition at line 976 of file MatrixTools.h.

References min().

◆ max()

template<class Real >
static Real bpp::MatrixTools::max ( const Matrix< Real > &  m)
inlinestatic
Returns
The maximum value in the matrix.
Parameters
mThe matrix.

Definition at line 521 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ min()

template<class Real >
static Real bpp::MatrixTools::min ( const Matrix< Real > &  m)
inlinestatic
Returns
The minimum value in the matrix.
Parameters
mThe matrix.

Definition at line 547 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

Referenced by lap().

◆ mult() [1/3]

template<class Scalar >
static void bpp::MatrixTools::mult ( const Matrix< Scalar > &  A,
const Matrix< Scalar > &  B,
Matrix< Scalar > &  O 
)
throw (DimensionException
)
inlinestatic
Parameters
A[in] First matrix.
B[in] Second matrix.
O[out] The dot product of two matrices.

Definition at line 190 of file MatrixTools.h.

Referenced by bpp::DualityDiagram::compute_(), covar(), exp(), bpp::AdaptiveKernelDensityEstimation::init_(), bpp::AdaptiveKernelDensityEstimation::kDensity(), pow(), and Taylor().

◆ mult() [2/3]

template<class Scalar >
static void bpp::MatrixTools::mult ( const Matrix< Scalar > &  A,
const std::vector< Scalar > &  D,
const Matrix< Scalar > &  B,
Matrix< Scalar > &  O 
)
throw (DimensionException
)
inlinestatic

Compute A . D . B where D is a diagonal matrix in O(n^3).

Since D is a diagonal matrix, this function is more efficient than doing mult(mult(A, diag(D)), B), which involves two 0(n^3) operations.

Parameters
A[in] The first matrix.
D[in] The diagonal matrix (only diagonal elements in a vector)
B[in] The second matrix.
O[out] The result matrix.
Exceptions
DimensionExceptionIf matrices have not the appropriate size.

Definition at line 224 of file MatrixTools.h.

◆ mult() [3/3]

template<class Scalar >
static void bpp::MatrixTools::mult ( const Matrix< Scalar > &  A,
const std::vector< Scalar > &  D,
const std::vector< Scalar > &  U,
const std::vector< Scalar > &  L,
const Matrix< Scalar > &  B,
Matrix< Scalar > &  O 
)
throw (DimensionException
)
inlinestatic

Compute A . (U+D+L) . B where D is a diagonal matrix, U (resp. L) is a matrix in which the only non-zero terms are on the diagonal that is over (resp. under) the main diagonal, in O(n^3).

Since D is a diagonal matrix, this function is more efficient than doing mult(mult(A, diag(D)), B), which involves two 0(n^3) operations.

Parameters
A[in] The first matrix.
D[in] The diagonal matrix (only diagonal elements in a vector)
U[in] The upper diagonal matrix (only upper diagonal elements in a vector)
L[in] The lower diagonal matrix (only lower diagonal elements in a vector)
B[in] The second matrix.
O[out] The result matrix.
Exceptions
DimensionExceptionIf matrices have not the appropriate size.

Definition at line 263 of file MatrixTools.h.

◆ pow() [1/2]

template<class Matrix >
static void bpp::MatrixTools::pow ( const Matrix A,
size_t  p,
Matrix O 
)
throw (DimensionException
)
inlinestatic

Compute the power of a given matrix.

Parameters
A[in] The matrix.
pThe number of multiplications.
O[out] $ A^p $ computed recursively: $ A^{2n} = (A^n)^2 $ $ A^{2n+1} = A*(A^n)^2 $
If p = 0, sends the identity matrix.
Exceptions
DimensionExceptionIf m is not a square matrix.

Definition at line 355 of file MatrixTools.h.

References copy(), and mult().

Referenced by bpp::FullHmmTransitionMatrix::getEquilibriumFrequencies().

◆ pow() [2/2]

template<class Scalar >
static void bpp::MatrixTools::pow ( const Matrix< Scalar > &  A,
double  p,
Matrix< Scalar > &  O 
)
throw (DimensionException
)
inlinestatic

Compute the power of a given matrix, using eigen value decomposition.

Parameters
A[in] The matrix.
pThe power of the matrix.
O[out] $\prod_{i=1}^p m$. If p = 0, sends the identity matrix.
Exceptions
DimensionExceptionIf m is not a square matrix.

Definition at line 394 of file MatrixTools.h.

References bpp::EigenValue< Real >::getRealEigenValues(), bpp::EigenValue< Real >::getV(), inv(), mult(), and bpp::VectorTools::pow().

◆ print() [1/3]

template<class Matrix >
static void bpp::MatrixTools::print ( const Matrix m,
std::ostream &  out = std::cout 
)
inlinestatic

Print a matrix to a stream.

Parameters
mThe matrix to print.
outThe stream to use.

Definition at line 573 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ print() [2/3]

template<class Matrix >
static void bpp::MatrixTools::print ( const Matrix m,
bpp::OutputStream out,
char  pIn = '(',
char  pOut = ')' 
)
inlinestatic

Print a matrix to a stream.

Parameters
mThe matrix to print.
outThe stream to use.
pInleft delimiter (default: "(")
pOutright delimiter (default: ")")

Definition at line 598 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ print() [3/3]

template<class Real >
static void bpp::MatrixTools::print ( const std::vector< Real > &  v,
std::ostream &  out = std::cout 
)
inlinestatic

Print a vector to a stream.

Parameters
vThe vector to print.
outThe stream to use.

Definition at line 649 of file MatrixTools.h.

◆ printForR()

template<class Matrix >
static void bpp::MatrixTools::printForR ( const Matrix m,
const std::string &  variableName = "x",
std::ostream &  out = std::cout 
)
inlinestatic

Print a matrix to a stream, so that it is read by R.

Parameters
mThe matrix to print.
variableNameThe name of the R variable handeling the matrix
outThe stream to use.

Definition at line 625 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ scale()

template<class Matrix , class Scalar >
static void bpp::MatrixTools::scale ( Matrix A,
Scalar  a,
Scalar  b = 0 
)
inlinestatic

Multiply all elements of a matrix by a given value, and add a constant.

Performs $\forall i \forall j m_{i,j} = a.m_{i,j}+b$.

Parameters
AA matrix.
aMultiplicator.
bConstant.

Definition at line 173 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

Referenced by bpp::CorrespondenceAnalysis::CorrespondenceAnalysis(), covar(), bpp::AdaptiveKernelDensityEstimation::init_(), and bpp::AdaptiveKernelDensityEstimation::kDensity().

◆ sumElements()

template<class Scalar >
static Scalar bpp::MatrixTools::sumElements ( const Matrix< Scalar > &  M)
inlinestatic

Sum all elements in M.

Parameters
MA matrix.
Returns
The sum of all elements.

Definition at line 946 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ Taylor()

template<class Matrix , class Scalar >
static void bpp::MatrixTools::Taylor ( const Matrix A,
size_t  p,
std::vector< RowMatrix< Scalar > > &  vO 
)
throw (DimensionException
)
inlinestatic

Compute a vector of the first powers of a given matrix.

Parameters
A[in] The matrix.
pThe number of powers.
vO[out] the vector of the powers (from 0 to p)
Exceptions
DimensionExceptionIf m is not a square matrix.

Definition at line 438 of file MatrixTools.h.

References copy(), and mult().

◆ toVVdouble()

template<class Scalar >
static void bpp::MatrixTools::toVVdouble ( const Matrix< Scalar > &  M,
std::vector< std::vector< Scalar > > &  vO 
)
inlinestatic

Convert to a vector of vector.

Parameters
M[in] A matrix object.
vO[out] The output vector of vector (will be resized accordingly).

Definition at line 925 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ transpose()

template<class MatrixA , class MatrixO >
static void bpp::MatrixTools::transpose ( const MatrixA &  A,
MatrixO &  O 
)
inlinestatic
Parameters
A[in] The matrix to transpose.
O[out] The transposition of A.

Definition at line 706 of file MatrixTools.h.

Referenced by bpp::DualityDiagram::compute_(), and covar().

◆ whichMax()

template<class Matrix >
static std::vector<size_t> bpp::MatrixTools::whichMax ( const Matrix m)
inlinestatic
Returns
The position of the maximum value in the matrix.
Parameters
mThe matrix.

Definition at line 458 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().

◆ whichMin()

template<class Matrix >
static std::vector<size_t> bpp::MatrixTools::whichMin ( const Matrix m)
inlinestatic
Returns
The position of the minimum value in the matrix.
Parameters
mThe matrix.

Definition at line 490 of file MatrixTools.h.

References bpp::Matrix< Scalar >::getNumberOfColumns(), and bpp::Matrix< Scalar >::getNumberOfRows().


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