INMOST
A toolkit for distributed mathematical modeling
|
This class allows to address a matrix as a block of an empty matrix of larger size. More...
#include <inmost_dense.h>
Public Types | |
typedef AbstractMatrixReadOnly< Var >::enumerator | enumerator |
Public Types inherited from INMOST::AbstractMatrixReadOnly< Var > | |
typedef unsigned | enumerator |
Public Member Functions | |
__INLINE enumerator | Rows () const |
Number of rows in submatrix. More... | |
__INLINE enumerator | Cols () const |
Number of columns in submatrix. More... | |
ConstBlockOfMatrix (const AbstractMatrixReadOnly< Var > &rM, enumerator num_rows, enumerator num_cols, enumerator offset_row, enumerator offset_col) | |
Create submatrix for a matrix. More... | |
ConstBlockOfMatrix (const ConstBlockOfMatrix &b) | |
__INLINE Var | compute (enumerator i, enumerator j) const |
Access element of the matrix by row and column indices without right to change the element. More... | |
::INMOST::Matrix< Var > | MakeMatrix () |
Convert block of matrix into matrix. More... | |
Public Member Functions inherited from INMOST::AbstractMatrixReadOnly< Var > | |
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::RowMerger > | merger |
This class allows to address a matrix as a block of an empty matrix of larger size.
Definition at line 2394 of file inmost_dense.h.
|
inline |
Create submatrix for a matrix.
rM | Reference to the matrix that stores elements. |
num_rows | Number of rows in the larger matrix. |
num_cols | Number of columns in the larger matrix. |
first_row | Offset for row index in the larger matrix. |
first_column | Offset for column index in the larger matrix. |
Definition at line 2418 of file inmost_dense.h.
|
inlinevirtual |
Number of columns in submatrix.
Implements INMOST::AbstractMatrixReadOnly< Var >.
Definition at line 2411 of file inmost_dense.h.
|
inlinevirtual |
Access element of the matrix by row and column indices without right to change the element.
i | Column index. |
j | Row index. |
Implements INMOST::AbstractMatrixReadOnly< Var >.
Definition at line 2426 of file inmost_dense.h.
|
inline |
Convert block of matrix into matrix.
Note, that modifying returned matrix does not affect elements of the matrix or original matrix used to create block of matrix.
Definition at line 2440 of file inmost_dense.h.
|
inlinevirtual |
Number of rows in submatrix.
Implements INMOST::AbstractMatrixReadOnly< Var >.
Definition at line 2408 of file inmost_dense.h.