INMOST
A toolkit for distributed mathematical modeling
INMOST::AbstractMatrix< Var > Class Template Referenceabstract

Abstract class for a matrix, used to abstract away all the data storage and access and provide common implementation of the algorithms. More...

#include <inmost_dense.h>

Inheritance diagram for INMOST::AbstractMatrix< Var >:
Collaboration diagram for INMOST::AbstractMatrix< Var >:

Public Types

typedef AbstractMatrixReadOnly< Var >::enumerator enumerator
 
- Public Types inherited from INMOST::AbstractMatrixReadOnly< Var >
typedef unsigned enumerator
 

Public Member Functions

 AbstractMatrix ()
 Construct empty matrix.
 
AbstractMatrixoperator= (AbstractMatrix const &other)
 Assign matrix of the same type. More...
 
template<typename typeB >
AbstractMatrixoperator= (AbstractMatrixReadOnly< typeB > const &other)
 Assign matrix of another type. More...
 
AbstractMatrixoperator= (Var const &b)
 Assign value to all entries of the matrix. More...
 
virtual Var & operator() (enumerator i, enumerator j)=0
 Access element of the matrix by row and column indices. More...
 
virtual const Var & get (enumerator i, enumerator j) const =0
 Access element of the matrix by row and column indices. More...
 
virtual void Resize (enumerator rows, enumerator cols)=0
 Resize the matrix into different size. More...
 
void Zero ()
 Set all the elements of the matrix to zero.
 
virtual void Swap (AbstractMatrix< Var > &b)
 Exchange contents of two matrices.
 
template<typename typeB >
AbstractMatrixoperator-= (const AbstractMatrixReadOnly< typeB > &other)
 Subtract a matrix and store result in the current. More...
 
template<typename typeB >
AbstractMatrixoperator-= (const AbstractMatrix< typeB > &other)
 Subtract a matrix and store result in the current. More...
 
template<typename typeB >
AbstractMatrixoperator+= (const AbstractMatrixReadOnly< typeB > &other)
 Add a matrix and store result in the current. More...
 
template<typename typeB >
AbstractMatrixoperator+= (const AbstractMatrix< typeB > &other)
 Add a matrix and store result in the current. More...
 
template<typename typeB >
AbstractMatrixoperator*= (const AbstractMatrixReadOnly< typeB > &B)
 Multiply matrix with another matrix in-place. More...
 
template<typename typeB >
AbstractMatrixoperator*= (typeB coef)
 Multiply the matrix by the coefficient of the same type and store the result. More...
 
template<typename typeB >
AbstractMatrixoperator/= (typeB coef)
 Divide the matrix by the coefficient of the same type and store the result. More...
 
SubMatrix< Var > operator() (enumerator first_row, enumerator last_row, enumerator first_col, enumerator last_col)
 Extract submatrix of a matrix for in-place manipulation of elements. More...
 
BlockOfMatrix< Var > BlockOf (enumerator nrows, enumerator ncols, enumerator offset_row, enumerator offset_col)
 Define matrix as a part of a matrix of larger size with in-place manipulation of elements. More...
 
MatrixTranspose< Var > Transpose ()
 Transpose the current matrix with access to elements. More...
 
MatrixRepack< Var > Repack (enumerator rows, enumerator cols)
 Change representation of the matrix into matrix of another size. More...
 
MatrixConcatCols< Var > ConcatCols (AbstractMatrix< Var > &B)
 Concatenate B matrix as columns of current matrix. More...
 
MatrixConcatRows< Var > ConcatRows (AbstractMatrix< Var > &B)
 Concatenate B matrix as rows of current matrix. 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...
 
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > Transform (const AbstractMatrix< typeB > &other) const
 Transformation matrix from current vector to provided vector using shortest arc rotation. 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 >
Matrix< typename Promote< Var, typeB >::type > CrossProduct (const AbstractMatrix< typeB > &other) const
 Cross-product operation for a vector. 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 >
Promote< Var, typeB >::type DotProduct (const AbstractMatrix< 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 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 AbstractMatrix< 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 frobenious norm of the matrix. More...
 
Var MaxNorm () const
 Computes maximum absolute value of the 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...
 
void MPT (INMOST_DATA_ENUM_TYPE *Perm, INMOST_DATA_REAL_TYPE *SL=NULL, INMOST_DATA_REAL_TYPE *SR=NULL) const
 Maximum product transversal. More...
 
__INLINE INMOST_DATA_ENUM_TYPE GetMatrixCount () const
 Retrieve number of indices of derivatives.
 
virtual ~AbstractMatrix ()
 Destructor.
 
template<typename typeB >
AbstractMatrix< Var > & operator-= (const AbstractMatrix< typeB > &other)
 
template<typename typeB >
AbstractMatrix< Var > & operator-= (const AbstractMatrixReadOnly< typeB > &other)
 
template<typename typeB >
AbstractMatrix< Var > & operator+= (const AbstractMatrix< typeB > &other)
 
template<typename typeB >
AbstractMatrix< Var > & operator+= (const AbstractMatrixReadOnly< typeB > &other)
 
template<typename typeB >
AbstractMatrix< Var > & operator*= (const AbstractMatrixReadOnly< typeB > &B)
 
template<typename typeB >
AbstractMatrix< Var > & operator*= (typeB coef)
 
template<typename typeB >
AbstractMatrix< Var > & operator/= (typeB coef)
 
- Public Member Functions inherited from INMOST::AbstractMatrixReadOnly< Var >
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.
 
__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::AbstractMatrix< Var >

Abstract class for a matrix, used to abstract away all the data storage and access and provide common implementation of the algorithms.

Definition at line 646 of file inmost_dense.h.

Member Function Documentation

◆ BlockOf()

template<typename Var >
BlockOfMatrix< Var > INMOST::AbstractMatrix< Var >::BlockOf ( enumerator  nrows,
enumerator  ncols,
enumerator  offset_row,
enumerator  offset_col 
)

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 4883 of file inmost_dense.h.

◆ CheckInfs()

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

Check all matrix entries for infinity.

Also checks derivatives for matrices of variables.

Definition at line 827 of file inmost_dense.h.

◆ CheckNans()

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

Check all matrix entries for not a number.

Also checks derivatives for matrices of variables.

Definition at line 818 of file inmost_dense.h.

◆ CheckNansInfs()

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

Check all matrix entries for not a number and infinity.

Also checks derivatives for matrices of variables.

Definition at line 836 of file inmost_dense.h.

◆ ConcatCols()

template<typename Var >
MatrixConcatCols<Var> INMOST::AbstractMatrix< Var >::ConcatCols ( AbstractMatrix< Var > &  B)
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 808 of file inmost_dense.h.

◆ ConcatRows()

template<typename Var >
MatrixConcatRows<Var> INMOST::AbstractMatrix< Var >::ConcatRows ( AbstractMatrix< Var > &  B)
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 815 of file inmost_dense.h.

◆ CrossProduct() [1/2]

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrix< Var >::CrossProduct ( const AbstractMatrix< 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 3696 of file inmost_dense.h.

◆ CrossProduct() [2/2]

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrix< 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 3715 of file inmost_dense.h.

◆ DotProduct() [1/2]

template<typename Var >
template<typename typeB >
Promote<Var, typeB>::type INMOST::AbstractMatrix< Var >::DotProduct ( const AbstractMatrix< 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 879 of file inmost_dense.h.

◆ DotProduct() [2/2]

template<typename Var >
template<typename typeB >
Promote<Var, typeB>::type INMOST::AbstractMatrix< 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 895 of file inmost_dense.h.

◆ FrobeniusNorm()

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

Computes frobenious norm of the matrix.

Returns
Frobenius norm of the matrix.

Definition at line 925 of file inmost_dense.h.

◆ get()

◆ isSymmetric()

template<typename Var >
bool INMOST::AbstractMatrix< 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 977 of file inmost_dense.h.

◆ MaxNorm()

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

Computes maximum absolute value of the matrix.

Returns
Maximum norm of the matrix.

Definition at line 935 of file inmost_dense.h.

◆ MPT()

template<typename Var >
void INMOST::AbstractMatrix< Var >::MPT ( INMOST_DATA_ENUM_TYPE *  Perm,
INMOST_DATA_REAL_TYPE *  SL = NULL,
INMOST_DATA_REAL_TYPE *  SR = NULL 
) const

Maximum product transversal.

Computes unsymmetric reordering that maximizes product on diagonal. Returns reordering matrix P and scaling matrix S that transforms matrix into I-dominant matrix.

Parameters
PermArray for reordering, size of columns of the matrix.
SLDiagonal for rescaling matrix from left, size of columns of the matrix.
SRDiagonal for rescaling matrix from right, size of rows of the matrix.
Todo:
  1. Test rescaling.
  2. Test on non-square matrices.

Definition at line 4986 of file inmost_dense.h.

◆ operator()() [1/2]

template<typename Var >
SubMatrix< Var > INMOST::AbstractMatrix< Var >::operator() ( enumerator  first_row,
enumerator  last_row,
enumerator  first_col,
enumerator  last_col 
)

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 4872 of file inmost_dense.h.

◆ operator()() [2/2]

◆ operator*=() [1/2]

template<typename Var >
template<typename typeB >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator*= ( const AbstractMatrixReadOnly< typeB > &  B)

Multiply matrix with another matrix in-place.

Parameters
BAnother matrix to the right in multiplication.
Returns
Reference to current matrix.

◆ operator*=() [2/2]

template<typename Var >
template<typename typeB >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator*= ( typeB  coef)

Multiply the matrix by the coefficient of the same type and store the result.

Parameters
coefCoefficient.
Returns
Reference to the current matrix.

◆ operator+=() [1/2]

template<typename Var >
template<typename typeB >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator+= ( const AbstractMatrix< typeB > &  other)

Add a matrix and store result in the current.

Parameters
otherThe matrix to be added.
Returns
Reference to the current matrix.

◆ operator+=() [2/2]

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

Add a matrix and store result in the current.

Parameters
otherThe matrix to be added.
Returns
Reference to the current matrix.

◆ operator-=() [1/2]

template<typename Var >
template<typename typeB >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator-= ( const AbstractMatrix< typeB > &  other)

Subtract a matrix and store result in the current.

Parameters
otherThe matrix to be subtracted.
Returns
Reference to the current matrix.

◆ operator-=() [2/2]

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

Subtract a matrix and store result in the current.

Parameters
otherThe matrix to be subtracted.
Returns
Reference to the current matrix.

◆ operator/=()

template<typename Var >
template<typename typeB >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator/= ( typeB  coef)

Divide the matrix by the coefficient of the same type and store the result.

Parameters
coefCoefficient.
Returns
Reference to the current matrix.

◆ operator=() [1/3]

template<typename Var >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator= ( AbstractMatrix< Var > const &  other)
inline

Assign matrix of the same type.

Parameters
otherAnother matrix of the same type.
Returns
Reference to matrix.

Definition at line 685 of file inmost_dense.h.

◆ operator=() [2/3]

template<typename Var >
template<typename typeB >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator= ( AbstractMatrixReadOnly< typeB > const &  other)
inline

Assign matrix of another type.

Parameters
otherAnother matrix of different type.
Returns
Reference to matrix.

Definition at line 697 of file inmost_dense.h.

◆ operator=() [3/3]

template<typename Var >
AbstractMatrix& INMOST::AbstractMatrix< Var >::operator= ( Var const &  b)
inline

Assign value to all entries of the matrix.

Parameters
bAssigned value.
Returns
Reference to matrix.

Definition at line 708 of file inmost_dense.h.

◆ operator^() [1/2]

template<typename Var >
template<typename typeB >
Promote<Var, typeB>::type INMOST::AbstractMatrix< Var >::operator^ ( const AbstractMatrix< 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 910 of file inmost_dense.h.

◆ operator^() [2/2]

template<typename Var >
template<typename typeB >
Promote<Var, typeB>::type INMOST::AbstractMatrix< 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 919 of file inmost_dense.h.

◆ Print()

template<typename Var >
void INMOST::AbstractMatrix< 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 955 of file inmost_dense.h.

◆ Repack()

template<typename Var >
MatrixRepack<Var> INMOST::AbstractMatrix< Var >::Repack ( enumerator  rows,
enumerator  cols 
)
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 801 of file inmost_dense.h.

◆ Resize()

◆ Trace()

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

Calculate sum of the diagonal elements of the matrix.

Returns
Trace of the matrix.

Definition at line 945 of file inmost_dense.h.

◆ Transform() [1/2]

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrix< Var >::Transform ( const AbstractMatrix< 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 3774 of file inmost_dense.h.

◆ Transform() [2/2]

template<typename Var >
template<typename typeB >
Matrix< typename Promote< Var, typeB >::type > INMOST::AbstractMatrix< 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 3812 of file inmost_dense.h.

◆ Transpose()

template<typename Var >
MatrixTranspose<Var> INMOST::AbstractMatrix< Var >::Transpose ( )
inline

Transpose the current matrix with access to elements.

Returns
Transposed matrix.

Definition at line 796 of file inmost_dense.h.


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