INMOST
A toolkit for distributed mathematical modeling
INMOST::AbstractMatrixReadOnly< Var > Class Template Referenceabstract
Inheritance diagram for INMOST::AbstractMatrixReadOnly< Var >:
Collaboration diagram for INMOST::AbstractMatrixReadOnly< Var >:

Public Types

typedef unsigned enumerator
 

Public Member Functions

virtual enumerator Rows () const =0
 Obtain number of rows. More...
 
virtual enumerator Cols () const =0
 Obtain number of columns. More...
 
virtual Var compute (enumerator i, enumerator j) const =0
 Compute element of the matrix by row and column indices without right to change the element. More...
 
bool CheckNans () const
 Check all matrix entries for not a number. More...
 
bool CheckInfs () const
 Check all matrix entries for infinity. More...
 
bool CheckNansInfs () const
 Check all matrix entries for not a number and infinity. More...
 
Var Det () const
 Matrix determinant.
 
bool SVD (AbstractMatrix< Var > &U, AbstractMatrix< Var > &Sigma, AbstractMatrix< Var > &V, bool order_singular_values=true, bool nonnegative=true) const
 Singular value decomposition. More...
 
bool cSVD (AbstractMatrix< Var > &U, AbstractMatrix< Var > &Sigma, AbstractMatrix< Var > &V) const
 Singular value decomposition. More...
 
ConstMatrixTranspose< Var > Transpose () const
 Transpose current matrix. More...
 
ConstMatrixConjugateTranspose< Var > ConjugateTranspose () const
 Transpose and conjugate current matrix. More...
 
ConstMatrixConjugate< Var > Conjugate () const
 Conjugate current matrix. More...
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > CrossProduct (const AbstractMatrixReadOnly< typeB > &other) const
 Cross-product operation for a vector. More...
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > Transform (const AbstractMatrixReadOnly< typeB > &other) const
 Transformation matrix from current vector to provided vector using shortest arc rotation. More...
 
template<typename typeB >
MatrixDifference< Var, typeB > operator- (const AbstractMatrixReadOnly< typeB > &other) const
 Subtract a matrix. More...
 
template<typename typeB >
MatrixSum< Var, typeB > operator+ (const AbstractMatrixReadOnly< typeB > &other) const
 Add a matrix. More...
 
template<typename typeB >
MatrixMul< Var, typeB, typename Promote< Var, typeB >::type > operator* (const AbstractMatrixReadOnly< typeB > &other) const
 Multiply the matrix by another matrix. More...
 
MatrixUnaryMinus< Var > operator- () const
 Unary minus. Change sign of each element of the matrix.
 
template<typename typeB >
KroneckerProduct< Var, typeB > Kronecker (const AbstractMatrixReadOnly< typeB > &other) const
 Kronecker product, latex symbol \otimes. More...
 
Matrix< Var > Invert (int *ierr=NULL) const
 Inverts matrix using Crout-LU decomposition with full pivoting for maximum element. More...
 
Matrix< Var > CholeskyInvert (int *ierr=NULL) const
 Inverts symmetric positive-definite matrix using Cholesky decomposition.
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > Solve (const AbstractMatrixReadOnly< typeB > &B, int *ierr=NULL) const
 Finds X in A*X=B, where A and B are general matrices. More...
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > CholeskySolve (const AbstractMatrixReadOnly< typeB > &B, int *ierr=NULL) const
 Finds X in A*X=B, where A is a square symmetric positive definite matrix. More...
 
Var Trace () const
 Calculate sum of the diagonal elements of the matrix. More...
 
void Print (INMOST_DATA_REAL_TYPE threshold=1.0e-10, std::ostream &sout=std::cout) const
 Output matrix to screen. More...
 
bool isSymmetric (double eps=1.0e-7) const
 Check if the matrix is symmetric. More...
 
template<typename typeB >
Promote< Var, typeB >::type DotProduct (const AbstractMatrixReadOnly< typeB > &other) const
 Computes dot product by summing up multiplication of entries with the same indices in the current and the provided matrix. More...
 
template<typename typeB >
Promote< Var, typeB >::type operator^ (const AbstractMatrixReadOnly< typeB > &other) const
 Computes dot product by summing up multiplication of entries with the same indices in the current and the provided matrix. More...
 
SelfPromote< Var >::type FrobeniusNorm () const
 Computes frobenius norm of the matrix. More...
 
Var MaxNorm () const
 Computes maximum absolute value of the matrix. More...
 
Matrix< Var > PseudoInvert (INMOST_DATA_REAL_TYPE tol=0, int *ierr=NULL) const
 Calculates Moore-Penrose pseudo-inverse of the matrix. More...
 
Matrix< Var > Power (INMOST_DATA_REAL_TYPE n, int *ierr=NULL) const
 Calcuate A^n, where n is some real value. More...
 
Matrix< Var > Root (INMOST_DATA_ENUM_TYPE iter=25, INMOST_DATA_REAL_TYPE tol=1.0e-7, int *ierr=NULL) const
 Calculate square root of A matrix by Babylonian method. More...
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > PseudoSolve (const AbstractMatrixReadOnly< typeB > &B, INMOST_DATA_REAL_TYPE tol=0, int *ierr=NULL) const
 Solves the system of equations of the form A*X=B, with A and B matrices. More...
 
Matrix< Var > ExtractSubMatrix (enumerator ibeg, enumerator iend, enumerator jbeg, enumerator jend) const
 Extract submatrix of a matrix. More...
 
ConstMatrixRepack< Var > Repack (enumerator rows, enumerator cols) const
 Change representation of the matrix into matrix of another size. More...
 
ConstSubMatrix< Var > operator() (enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col) const
 Extract submatrix of a matrix for in-place manipulation of elements. More...
 
ConstBlockOfMatrix< Var > BlockOf (enumerator nrows, enumerator ncols, enumerator offset_row, enumerator offset_col) const
 Define matrix as a part of a matrix of larger size with in-place manipulation of elements. More...
 
template<typename typeB >
MatrixMulCoef< Var, typeB, typename Promote< Var, typeB >::type > operator* (const typeB &coef) const
 Multiply the matrix by a coefficient. More...
 
template<class A >
MatrixMulShellCoef< Var, shell_expression< A >, typename Promote< Var, variable >::type > operator* (shell_expression< A > const &coef) const
 Multiply the matrix by a coefficient. More...
 
template<typename typeB >
MatrixDivCoef< Var, typeB, typename Promote< Var, typeB >::type > operator/ (const typeB &coef) const
 Divide the matrix by a coefficient of a different type. More...
 
template<class A >
MatrixDivShellCoef< Var, shell_expression< A >, typename Promote< Var, variable >::type > operator/ (shell_expression< A > const &coef) const
 Divide the matrix by a coefficient of a different type. More...
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > operator/ (const AbstractMatrixReadOnly< typeB > &other) const
 Performs B^{-1}*A, multiplication by inverse matrix from left. More...
 
ConstMatrixConcatCols< Var > ConcatCols (const AbstractMatrixReadOnly< Var > &B) const
 Concatenate B matrix as columns of current matrix. More...
 
template<typename VarB >
ConstMatrixConcatCols2< Var, VarB, typename Promote< Var, VarB >::type > ConcatCols (const AbstractMatrixReadOnly< VarB > &B) const
 Concatenate B matrix as columns of current matrix. More...
 
ConstMatrixConcatRows< Var > ConcatRows (const AbstractMatrixReadOnly< Var > &B) const
 Concatenate B matrix as rows of current matrix. More...
 
template<typename VarB >
ConstMatrixConcatRows2< Var, VarB, typename Promote< Var, VarB >::type > ConcatRows (const AbstractMatrixReadOnly< VarB > &B) const
 Concatenate B matrix as rows of current matrix. More...
 
virtual ~AbstractMatrixReadOnly ()
 Destructor.
 
virtual __INLINE INMOST_DATA_ENUM_TYPE GetMatrixCount () const
 Retrieve number of indices of derivatives.
 
__INLINE Matrix< Promote< variable, variable >::type > CholeskySolve (const AbstractMatrixReadOnly< variable > &B, int *ierr) const
 
__INLINE Matrix< Promote< INMOST_DATA_REAL_TYPE, variable >::type > CholeskySolve (const AbstractMatrixReadOnly< variable > &B, int *ierr) const
 
__INLINE Matrix< Promote< variable, INMOST_DATA_REAL_TYPE >::type > CholeskySolve (const AbstractMatrixReadOnly< INMOST_DATA_REAL_TYPE > &B, int *ierr) const
 

Additional Inherited Members

- Static Protected Attributes inherited from INMOST::AbstractMatrixBase
static thread_private< Sparse::RowMergermerger
 

Detailed Description

template<typename Var>
class INMOST::AbstractMatrixReadOnly< Var >

Definition at line 190 of file inmost_dense.h.

Member Function Documentation

◆ BlockOf()

template<typename Var >
ConstBlockOfMatrix< Var > INMOST::AbstractMatrixReadOnly< Var >::BlockOf ( enumerator  nrows,
enumerator  ncols,
enumerator  offset_row,
enumerator  offset_col 
) const

Define matrix as a part of a matrix of larger size with in-place manipulation of elements.

Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix. Then the method returns B = {a_ij}, if i in [offset_row,offset_row+n), and j in [offset_col,offset_col+m) and B = {0} otherwise.

Parameters
nrowsNumber of rows in larger matrix.
ncolsNumber of columns in larger matrix.
offset_rowOffset for row number.
offset_colOffset for column number.
Returns
Submatrix of the original matrix.

Definition at line 4888 of file inmost_dense.h.

◆ CheckInfs()

template<typename Var >
bool INMOST::AbstractMatrixReadOnly< Var >::CheckInfs ( ) const
inline

Check all matrix entries for infinity.

Also checks derivatives for matrices of variables.

Definition at line 225 of file inmost_dense.h.

◆ CheckNans()

template<typename Var >
bool INMOST::AbstractMatrixReadOnly< Var >::CheckNans ( ) const
inline

Check all matrix entries for not a number.

Also checks derivatives for matrices of variables.

Definition at line 222 of file inmost_dense.h.

◆ CheckNansInfs()

template<typename Var >
bool INMOST::AbstractMatrixReadOnly< Var >::CheckNansInfs ( ) const
inline

Check all matrix entries for not a number and infinity.

Also checks derivatives for matrices of variables.

Definition at line 228 of file inmost_dense.h.

◆ CholeskySolve()

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::CholeskySolve ( const AbstractMatrixReadOnly< typeB > &  B,
int *  ierr = NULL 
) const

Finds X in A*X=B, where A is a square symmetric positive definite matrix.

Uses Cholesky decomposition algorithm.

Parameters
BRight hand side matrix.
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr equal to a positive integer that represents the row on which the failure occurred (starting from 1), in case of no failure *ierr = 0.
See also
Matrix::PseudoInvert.
Returns
Inverse matrix,

Definition at line 4120 of file inmost_dense.h.

◆ Cols()

template<typename Var >
virtual enumerator INMOST::AbstractMatrixReadOnly< Var >::Cols ( ) const
pure virtual

Obtain number of columns.

Returns
Number of columns.

Implemented in INMOST::KroneckerProduct< VarA, VarB >, INMOST::MatrixDivShellCoef< VarA, VarB, VarR >, INMOST::MatrixMulShellCoef< VarA, VarB, VarR >, INMOST::MatrixDivCoef< VarA, VarB, VarR >, INMOST::MatrixMulCoef< VarA, VarB, VarR >, INMOST::MatrixMul< variable, variable, Promote< variable, variable >::type >, INMOST::MatrixMul< variable, INMOST_DATA_REAL_TYPE, Promote< variable, INMOST_DATA_REAL_TYPE >::type >, INMOST::MatrixMul< INMOST_DATA_REAL_TYPE, variable, Promote< INMOST_DATA_REAL_TYPE, variable >::type >, INMOST::MatrixMul< VarA, VarB, VarR >, INMOST::MatrixRepack< Var >, INMOST::ConstMatrixRepack< Var >, INMOST::MatrixConcatCols< Var >, INMOST::ConstMatrixConcatCols2< VarA, VarB, VarR >, INMOST::ConstMatrixConcatCols< Var >, INMOST::MatrixConcatRows< Var >, INMOST::ConstMatrixConcatRows2< VarA, VarB, VarR >, INMOST::ConstMatrixConcatRows< Var >, INMOST::MatrixTranspose< Var >, INMOST::MatrixUnaryMinus< Var >, INMOST::ConstMatrixConjugate< Var >, INMOST::ConstMatrixConjugateTranspose< Var >, INMOST::ConstMatrixTranspose< Var >, INMOST::MatrixDifference< VarA, VarB >, INMOST::MatrixSum< VarA, VarB >, INMOST::MatrixDiag< Var >, INMOST::MatrixCol< Var >, INMOST::MatrixRow< Var >, INMOST::MatrixUnit< Var >, INMOST::ConstBlockOfMatrix< Var >, INMOST::BlockOfMatrix< Var >, INMOST::ConstSubMatrix< Var >, INMOST::SubMatrix< Var >, INMOST::Matrix< Var, storage_type >, INMOST::Matrix< INMOST::Promote< INMOST::multivar_expression, INMOST_DATA_REAL_TYPE >::type >, INMOST::Matrix< VarR >, INMOST::Matrix< INMOST::Promote< INMOST_DATA_REAL_TYPE, INMOST::multivar_expression >::type >, INMOST::Matrix< INMOST_DATA_REAL_TYPE >, INMOST::Matrix< INMOST::Promote< INMOST::multivar_expression, INMOST::multivar_expression >::type >, and INMOST::SymmetricMatrix< Var, storage_type >.

◆ compute()

template<typename Var >
virtual Var INMOST::AbstractMatrixReadOnly< Var >::compute ( enumerator  i,
enumerator  j 
) const
pure virtual

Compute element of the matrix by row and column indices without right to change the element.

Parameters
iColumn index.
jRow index.
Returns
Reference to constant element.

Implemented in INMOST::KroneckerProduct< VarA, VarB >, INMOST::MatrixMul< variable, variable, Promote< variable, variable >::type >, INMOST::MatrixMul< variable, INMOST_DATA_REAL_TYPE, Promote< variable, INMOST_DATA_REAL_TYPE >::type >, INMOST::MatrixMul< INMOST_DATA_REAL_TYPE, variable, Promote< INMOST_DATA_REAL_TYPE, variable >::type >, INMOST::MatrixRepack< Var >, INMOST::ConstMatrixRepack< Var >, INMOST::MatrixConcatCols< Var >, INMOST::ConstMatrixConcatCols< Var >, INMOST::MatrixConcatRows< Var >, INMOST::ConstMatrixConcatRows< Var >, INMOST::MatrixTranspose< Var >, INMOST::MatrixUnaryMinus< Var >, INMOST::ConstMatrixConjugate< Var >, INMOST::ConstMatrixConjugateTranspose< Var >, INMOST::ConstMatrixTranspose< Var >, INMOST::MatrixDifference< VarA, VarB >, INMOST::MatrixSum< VarA, VarB >, INMOST::MatrixDiag< Var >, INMOST::MatrixCol< Var >, INMOST::MatrixRow< Var >, INMOST::MatrixUnit< Var >, INMOST::ConstBlockOfMatrix< Var >, INMOST::BlockOfMatrix< Var >, INMOST::ConstSubMatrix< Var >, INMOST::SubMatrix< Var >, INMOST::Matrix< Var, storage_type >, INMOST::Matrix< INMOST::Promote< INMOST::multivar_expression, INMOST_DATA_REAL_TYPE >::type >, INMOST::Matrix< VarR >, INMOST::Matrix< INMOST::Promote< INMOST_DATA_REAL_TYPE, INMOST::multivar_expression >::type >, INMOST::Matrix< INMOST_DATA_REAL_TYPE >, INMOST::Matrix< INMOST::Promote< INMOST::multivar_expression, INMOST::multivar_expression >::type >, and INMOST::SymmetricMatrix< Var, storage_type >.

◆ ConcatCols() [1/2]

template<typename Var >
ConstMatrixConcatCols<Var> INMOST::AbstractMatrixReadOnly< Var >::ConcatCols ( const AbstractMatrixReadOnly< Var > &  B) const
inline

Concatenate B matrix as columns of current matrix.

Assumes that number of rows of current matrix is equal to number of rows of B matrix.

Parameters
BMatrix to be concatenated to current matrix.
Returns
Result of concatenation of current matrix and parameter.
See also
Matrix::ConcatRows

Definition at line 592 of file inmost_dense.h.

◆ ConcatCols() [2/2]

template<typename Var >
template<typename VarB >
ConstMatrixConcatCols2<Var, VarB, typename Promote<Var, VarB>::type> INMOST::AbstractMatrixReadOnly< Var >::ConcatCols ( const AbstractMatrixReadOnly< VarB > &  B) const
inline

Concatenate B matrix as columns of current matrix.

Assumes that number of rows of current matrix is equal to number of rows of B matrix.

Parameters
BMatrix to be concatenated to current matrix.
Returns
Result of concatenation of current matrix and parameter.
See also
Matrix::ConcatRows

Definition at line 604 of file inmost_dense.h.

◆ ConcatRows() [1/2]

template<typename Var >
ConstMatrixConcatRows<Var> INMOST::AbstractMatrixReadOnly< Var >::ConcatRows ( const AbstractMatrixReadOnly< Var > &  B) const
inline

Concatenate B matrix as rows of current matrix.

Assumes that number of colums of current matrix is equal to number of columns of B matrix.

Parameters
BMatrix to be concatenated to current matrix.
Returns
Result of concatenation of current matrix and parameter.
See also
Matrix::ConcatCols

Definition at line 614 of file inmost_dense.h.

◆ ConcatRows() [2/2]

template<typename Var >
template<typename VarB >
ConstMatrixConcatRows2<Var,VarB,typename Promote<Var,VarB>::type> INMOST::AbstractMatrixReadOnly< Var >::ConcatRows ( const AbstractMatrixReadOnly< VarB > &  B) const
inline

Concatenate B matrix as rows of current matrix.

Assumes that number of colums of current matrix is equal to number of columns of B matrix.

Parameters
BMatrix to be concatenated to current matrix.
Returns
Result of concatenation of current matrix and parameter.
See also
Matrix::ConcatCols

Definition at line 626 of file inmost_dense.h.

◆ Conjugate()

template<typename Var >
ConstMatrixConjugate<Var> INMOST::AbstractMatrixReadOnly< Var >::Conjugate ( ) const
inline

Conjugate current matrix.

Returns
Conjugated matrix.

Definition at line 286 of file inmost_dense.h.

◆ ConjugateTranspose()

template<typename Var >
ConstMatrixConjugateTranspose<Var> INMOST::AbstractMatrixReadOnly< Var >::ConjugateTranspose ( ) const
inline

Transpose and conjugate current matrix.

Returns
Transposed and conjugated matrix.

Definition at line 283 of file inmost_dense.h.

◆ CrossProduct()

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::CrossProduct ( const AbstractMatrixReadOnly< typeB > &  other) const

Cross-product operation for a vector.

Both right hand side and left hand side should be a vector

Parameters
otherThe right hand side of cross product.
Returns
The cross product of current and right hand side vector.
Warning
Works for 3x1 vector and 3xm m-vectors as right hand side.

Definition at line 3676 of file inmost_dense.h.

◆ cSVD()

template<typename Var >
bool INMOST::AbstractMatrixReadOnly< Var >::cSVD ( AbstractMatrix< Var > &  U,
AbstractMatrix< Var > &  Sigma,
AbstractMatrix< Var > &  V 
) const

Singular value decomposition.

Specific for complex matrices.

Definition at line 5209 of file inmost_dense.h.

◆ DotProduct()

template<typename Var >
template<typename typeB >
Promote<Var, typeB>::type INMOST::AbstractMatrixReadOnly< Var >::DotProduct ( const AbstractMatrixReadOnly< typeB > &  other) const
inline

Computes dot product by summing up multiplication of entries with the same indices in the current and the provided matrix.

Parameters
otherProvided matrix.
Returns
Dot product of two matrices.

Definition at line 428 of file inmost_dense.h.

◆ ExtractSubMatrix()

template<typename Var >
Matrix< Var > INMOST::AbstractMatrixReadOnly< Var >::ExtractSubMatrix ( enumerator  ibeg,
enumerator  iend,
enumerator  jbeg,
enumerator  jend 
) const

Extract submatrix of a matrix.

Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix. Then the method returns B = {a_ij}, i in [ibeg,iend), and j in [jbeg,jend).

Parameters
ibegStarting row in the original matrix.
iendLast row (excluded) in the original matrix.
jbegStarting column in the original matrix.
jendLast column (excluded) in the original matrix.
Returns
Submatrix of the original matrix.

Definition at line 4745 of file inmost_dense.h.

◆ FrobeniusNorm()

template<typename Var >
SelfPromote<Var>::type INMOST::AbstractMatrixReadOnly< Var >::FrobeniusNorm ( ) const
inline

Computes frobenius norm of the matrix.

Returns
Frobenius norm of the matrix.

Definition at line 449 of file inmost_dense.h.

◆ Invert()

template<typename Var >
Matrix< Var > INMOST::AbstractMatrixReadOnly< Var >::Invert ( int *  ierr = NULL) const

Inverts matrix using Crout-LU decomposition with full pivoting for maximum element.

Works well for non-singular matrices, for singular matrices look see Matrix::PseudoInvert.

Parameters
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr equal to a positive integer that represents the row on which the failure occurred (starting from 1), in case of no failure *ierr = 0.
Returns
Inverse matrix.
See also
Matrix::PseudoInvert.
Todo:
(test) Activate and test implementation with Solve.

Definition at line 4100 of file inmost_dense.h.

◆ isSymmetric()

template<typename Var >
bool INMOST::AbstractMatrixReadOnly< Var >::isSymmetric ( double  eps = 1.0e-7) const
inline

Check if the matrix is symmetric.

Returns
Returns true if the matrix is symmetric, otherwise false.

Definition at line 411 of file inmost_dense.h.

◆ Kronecker()

template<typename Var >
template<typename typeB >
KroneckerProduct< Var, typeB > INMOST::AbstractMatrixReadOnly< Var >::Kronecker ( const AbstractMatrixReadOnly< typeB > &  other) const

Kronecker product, latex symbol \otimes.

Parameters
otherMatrix on the right of the Kronecker product.
Returns
Result of the Kronecker product of current and another matrix.

Definition at line 4093 of file inmost_dense.h.

◆ MaxNorm()

template<typename Var >
Var INMOST::AbstractMatrixReadOnly< Var >::MaxNorm ( ) const
inline

Computes maximum absolute value of the matrix.

Returns
Maximum norm of the matrix.

Definition at line 459 of file inmost_dense.h.

◆ operator()()

template<typename Var >
ConstSubMatrix< Var > INMOST::AbstractMatrixReadOnly< Var >::operator() ( enumerator  first_row,
enumerator  last_row,
enumerator  first_col,
enumerator  last_col 
) const

Extract submatrix of a matrix for in-place manipulation of elements.

Let A = {a_ij}, i in [0,n), j in [0,m) be original matrix. Then the method returns B = {a_ij}, i in [ibeg,iend), and j in [jbeg,jend).

Parameters
first_rowStarting row in the original matrix.
last_rowLast row (excluded) in the original matrix.
first_colStarting column in the original matrix.
last_colLast column (excluded) in the original matrix.
Returns
Submatrix of the original matrix.

Definition at line 4877 of file inmost_dense.h.

◆ operator*() [1/3]

template<typename Var >
template<typename typeB >
MatrixMul< Var, typeB, typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::operator* ( const AbstractMatrixReadOnly< typeB > &  other) const

Multiply the matrix by another matrix.

Parameters
otherMatrix to be multiplied from right.
Returns
Matrix multipled by another matrix.

Definition at line 3917 of file inmost_dense.h.

◆ operator*() [2/3]

template<typename Var >
template<typename typeB >
MatrixMulCoef< Var, typeB, typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::operator* ( const typeB &  coef) const

Multiply the matrix by a coefficient.

Parameters
coefCoefficient.
Returns
Matrix multiplied by the coefficient.

Definition at line 4036 of file inmost_dense.h.

◆ operator*() [3/3]

template<typename Var >
template<typename A >
MatrixMulShellCoef< Var, shell_expression< A >, typename Promote< Var, variable >::type > INMOST::AbstractMatrixReadOnly< Var >::operator* ( shell_expression< A > const &  coef) const

Multiply the matrix by a coefficient.

Parameters
coefCoefficient.
Returns
Matrix multiplied by the coefficient.

Definition at line 4045 of file inmost_dense.h.

◆ operator+()

template<typename Var >
template<typename typeB >
MatrixSum< Var, typeB > INMOST::AbstractMatrixReadOnly< Var >::operator+ ( const AbstractMatrixReadOnly< typeB > &  other) const

Add a matrix.

Parameters
otherThe matrix to be added.
Returns
Sum of current matrix with another matrix.

Definition at line 3883 of file inmost_dense.h.

◆ operator-()

template<typename Var >
template<typename typeB >
MatrixDifference< Var, typeB > INMOST::AbstractMatrixReadOnly< Var >::operator- ( const AbstractMatrixReadOnly< typeB > &  other) const

Subtract a matrix.

Parameters
otherThe matrix to be subtracted.
Returns
Difference of current matrix with another matrix.

Definition at line 3849 of file inmost_dense.h.

◆ operator/() [1/3]

template<typename Var >
template<typename typeB >
Matrix<typename Promote<Var, typeB>::type> INMOST::AbstractMatrixReadOnly< Var >::operator/ ( const AbstractMatrixReadOnly< typeB > &  other) const
inline

Performs B^{-1}*A, multiplication by inverse matrix from left.

Throws exception if matrix is not invertible. See Mesh::PseudoSolve for singular matrices.

Parameters
otherMatrix to be inverted and multiplied from left.
Returns
Multiplication of current matrix by inverse of another
See also
Matrix::PseudoInvert.

Definition at line 580 of file inmost_dense.h.

◆ operator/() [2/3]

template<typename Var >
template<typename typeB >
MatrixDivCoef< Var, typeB, typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::operator/ ( const typeB &  coef) const

Divide the matrix by a coefficient of a different type.

Parameters
coefCoefficient.
Returns
Matrix divided by the coefficient.

Definition at line 4074 of file inmost_dense.h.

◆ operator/() [3/3]

template<typename Var >
template<typename A >
MatrixDivShellCoef< Var, shell_expression< A >, typename Promote< Var, variable >::type > INMOST::AbstractMatrixReadOnly< Var >::operator/ ( shell_expression< A > const &  coef) const

Divide the matrix by a coefficient of a different type.

Parameters
coefCoefficient.
Returns
Matrix divided by the coefficient.

Definition at line 4065 of file inmost_dense.h.

◆ operator^()

template<typename Var >
template<typename typeB >
Promote<Var, typeB>::type INMOST::AbstractMatrixReadOnly< Var >::operator^ ( const AbstractMatrixReadOnly< typeB > &  other) const
inline

Computes dot product by summing up multiplication of entries with the same indices in the current and the provided matrix.

Parameters
otherProvided matrix.
Returns
Dot product of two matrices.

Definition at line 443 of file inmost_dense.h.

◆ Power()

template<typename Var >
Matrix< Var > INMOST::AbstractMatrixReadOnly< Var >::Power ( INMOST_DATA_REAL_TYPE  n,
int *  ierr = NULL 
) const

Calcuate A^n, where n is some real value.

Warning
Currently works for symmetric matrices only!
Parameters
nReal value.
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr = 1, in case of no failure *ierr = 0. The error may be caused by error in SVD algorithm.
Returns
Matrix in power of n.

Definition at line 4817 of file inmost_dense.h.

◆ Print()

template<typename Var >
void INMOST::AbstractMatrixReadOnly< Var >::Print ( INMOST_DATA_REAL_TYPE  threshold = 1.0e-10,
std::ostream &  sout = std::cout 
) const
inline

Output matrix to screen.

Does not output derivatices.

Parameters
thresholdElements smaller then the number are considered to be zero.

Definition at line 389 of file inmost_dense.h.

◆ PseudoInvert()

template<typename Var >
Matrix< Var > INMOST::AbstractMatrixReadOnly< Var >::PseudoInvert ( INMOST_DATA_REAL_TYPE  tol = 0,
int *  ierr = NULL 
) const

Calculates Moore-Penrose pseudo-inverse of the matrix.

Parameters
tolThershold for singular values. Singular values smaller then threshold are considered to be zero.
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr = 1, in case of no failure *ierr = 0.
Returns
A pseudo-inverse of the matrix.

Definition at line 4762 of file inmost_dense.h.

◆ PseudoSolve()

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::PseudoSolve ( const AbstractMatrixReadOnly< typeB > &  B,
INMOST_DATA_REAL_TYPE  tol = 0,
int *  ierr = NULL 
) const

Solves the system of equations of the form A*X=B, with A and B matrices.

Uses Moore-Penrose pseudo-inverse of the matrix A and calculates X = A^+*B.

Parameters
BMatrix at the right hand side.
tolThershold for singular values. Singular values smaller then threshold are considered to be zero.
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr = 1, in case of no failure *ierr = 0.
Returns
A pair of the solution matrix X and boolean. If boolean is true, then the matrix was inverted successfully.

Definition at line 4845 of file inmost_dense.h.

◆ Repack()

template<typename Var >
ConstMatrixRepack<Var> INMOST::AbstractMatrixReadOnly< Var >::Repack ( enumerator  rows,
enumerator  cols 
) const
inline

Change representation of the matrix into matrix of another size.

Useful to change representation from matrix into vector and back. Replaces original number of columns and rows with a new one.

Returns
Matrix with same entries and provided number of rows and columns.

Definition at line 519 of file inmost_dense.h.

◆ Root()

template<typename Var >
Matrix< Var > INMOST::AbstractMatrixReadOnly< Var >::Root ( INMOST_DATA_ENUM_TYPE  iter = 25,
INMOST_DATA_REAL_TYPE  tol = 1.0e-7,
int *  ierr = NULL 
) const

Calculate square root of A matrix by Babylonian method.

Parameters
iterNumber of iterations.
tolConvergence tolerance.
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr = 1, in case of no failure *ierr = 0.
Returns
Square root of a matrix.

Definition at line 4792 of file inmost_dense.h.

◆ Rows()

template<typename Var >
virtual enumerator INMOST::AbstractMatrixReadOnly< Var >::Rows ( ) const
pure virtual

Obtain number of rows.

Returns
Number of rows.

Implemented in INMOST::KroneckerProduct< VarA, VarB >, INMOST::MatrixDivShellCoef< VarA, VarB, VarR >, INMOST::MatrixMulShellCoef< VarA, VarB, VarR >, INMOST::MatrixDivCoef< VarA, VarB, VarR >, INMOST::MatrixMulCoef< VarA, VarB, VarR >, INMOST::MatrixMul< variable, variable, Promote< variable, variable >::type >, INMOST::MatrixMul< variable, INMOST_DATA_REAL_TYPE, Promote< variable, INMOST_DATA_REAL_TYPE >::type >, INMOST::MatrixMul< INMOST_DATA_REAL_TYPE, variable, Promote< INMOST_DATA_REAL_TYPE, variable >::type >, INMOST::MatrixMul< VarA, VarB, VarR >, INMOST::MatrixRepack< Var >, INMOST::ConstMatrixRepack< Var >, INMOST::MatrixConcatCols< Var >, INMOST::ConstMatrixConcatCols2< VarA, VarB, VarR >, INMOST::ConstMatrixConcatCols< Var >, INMOST::MatrixConcatRows< Var >, INMOST::ConstMatrixConcatRows2< VarA, VarB, VarR >, INMOST::ConstMatrixConcatRows< Var >, INMOST::MatrixTranspose< Var >, INMOST::MatrixUnaryMinus< Var >, INMOST::ConstMatrixConjugate< Var >, INMOST::ConstMatrixConjugateTranspose< Var >, INMOST::ConstMatrixTranspose< Var >, INMOST::MatrixDifference< VarA, VarB >, INMOST::MatrixSum< VarA, VarB >, INMOST::MatrixDiag< Var >, INMOST::MatrixCol< Var >, INMOST::MatrixRow< Var >, INMOST::MatrixUnit< Var >, INMOST::ConstBlockOfMatrix< Var >, INMOST::BlockOfMatrix< Var >, INMOST::ConstSubMatrix< Var >, INMOST::SubMatrix< Var >, INMOST::Matrix< Var, storage_type >, INMOST::Matrix< INMOST::Promote< INMOST::multivar_expression, INMOST_DATA_REAL_TYPE >::type >, INMOST::Matrix< VarR >, INMOST::Matrix< INMOST::Promote< INMOST_DATA_REAL_TYPE, INMOST::multivar_expression >::type >, INMOST::Matrix< INMOST_DATA_REAL_TYPE >, INMOST::Matrix< INMOST::Promote< INMOST::multivar_expression, INMOST::multivar_expression >::type >, and INMOST::SymmetricMatrix< Var, storage_type >.

◆ Solve()

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::Solve ( const AbstractMatrixReadOnly< typeB > &  B,
int *  ierr = NULL 
) const

Finds X in A*X=B, where A and B are general matrices.

Converts system into A^T*A*X=A^T*B. Inverts matrix A^T*A using Crout-LU decomposition with full pivoting for maximum element. Works well for non-singular matrices, for singular matrices see Matrix::PseudoInvert.

Parameters
BRight hand side matrix.
ierrReturns error on fail. If ierr is NULL, then throws an exception. If *ierr == -1 on input, then prints out information in case of failure. In case of failure *ierr equal to a positive integer that represents the row on which the failure occurred (starting from 1), in case of no failure *ierr = 0.
Returns
Inverse matrix,
See also
Matrix::PseudoInvert.
Warning
Number of rows in matrices A and B should match.
Todo:
  1. Test implementation.
  2. Maximum product transversal + block pivoting instead of pivoting by maximum element.

Definition at line 4588 of file inmost_dense.h.

◆ SVD()

template<typename Var >
bool INMOST::AbstractMatrixReadOnly< Var >::SVD ( AbstractMatrix< Var > &  U,
AbstractMatrix< Var > &  Sigma,
AbstractMatrix< Var > &  V,
bool  order_singular_values = true,
bool  nonnegative = true 
) const

Singular value decomposition.

Reconstruct matrix: A = U*Sigma*V.Transpose(). Source http://www.public.iastate.edu/~dicook/JSS/paper/code/svd.c With few corrections for robustness.

Parameters
ULeft unitary matrix, U^T U = I.
SigmaDiagonal matrix with singular values.
VRight unitary matrix, not transposed.
order_singular_valuesReturn singular values in descending order.
nonnegativeChange sign of singular values.
Warning
Somehow result differ in auto-differential mode.
Todo:
Test implementation for auto-differentiation.

Definition at line 5514 of file inmost_dense.h.

◆ Trace()

template<typename Var >
Var INMOST::AbstractMatrixReadOnly< Var >::Trace ( ) const
inline

Calculate sum of the diagonal elements of the matrix.

Returns
Trace of the matrix.

Definition at line 379 of file inmost_dense.h.

◆ Transform()

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrixReadOnly< Var >::Transform ( const AbstractMatrixReadOnly< typeB > &  other) const

Transformation matrix from current vector to provided vector using shortest arc rotation.

Parameters
otherVector to transform to.
Returns
A sqaure (rotation) matrix that transforms current vector into right hand side vector.
Warning
Works only for 3x1 vectors.

Definition at line 3735 of file inmost_dense.h.

◆ Transpose()

template<typename Var >
ConstMatrixTranspose<Var> INMOST::AbstractMatrixReadOnly< Var >::Transpose ( ) const
inline

Transpose current matrix.

Returns
Transposed matrix.

Definition at line 280 of file inmost_dense.h.


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