HDK
|
#include <UT_SparseMatrix.h>
Classes | |
struct | ColumnValue |
struct | Triplet |
Public Types | |
typedef UT_Array< UT_ExintArray > | LevelSet |
Public Member Functions | |
UT_SparseMatrixCSRT () | |
Constructs an uninitialized sparse matrix. More... | |
UT_SparseMatrixCSRT (exint rows, exint columns) | |
Constructs a sparse matrix initialized with the given rows and columns. More... | |
UT_SparseMatrixCSRT (const UT_SparseMatrixCSRT &matrix) | |
void | init (exint rows, exint columns) |
void | zero () |
Clears all data in the matrix. More... | |
void | swap (UT_SparseMatrixCSRT &other) |
Swaps size and data elements with the given sparse matrix. More... | |
exint | getNumRows () const |
The number of rows in the matrix. More... | |
exint | getNumCols () const |
The number of columns in the matrix. More... | |
exint | getNumNonZeros () const |
The total number of non-zero elements in the matrix. More... | |
exint | getNumCellsInRow (exint row) const |
The number of non-zero elements in the given row. More... | |
void | setValues (UT_Array< Triplet > &triplets, T tolerance=0.0) |
void | setValuesFromRows (UT_Array< UT_Array< ColumnValue >> &rowValues, T tolerance=0.0) |
void | addValues (UT_Array< Triplet > &triplets, bool alreadySorted=false, T tolerance=1e-5) |
void | extractDiagonal (UT_VectorT< T > &result, bool invert=false, T tolerance=0, T defaultValue=0) const |
void | multVec (const UT_VectorT< T > &v, UT_VectorT< T > &result) const |
Matrix-vector multiplication. Computes result = M * v. More... | |
UT_PreciseT< T > | multVecAndDot (const UT_VectorT< T > &v, UT_VectorT< T > &result) const |
void | subtractMultVec (const UT_VectorT< T > &x, const UT_VectorT< T > &y, UT_VectorT< T > &result) const |
void | multMatrix (const UT_SparseMatrixCSRT< T > &B, UT_SparseMatrixCSRT &result, T tolerance=0.0) const |
void | transpose (UT_SparseMatrixCSRT< T > &result) const |
Computes the transpose of this matrix. result can refer to this matrix. More... | |
void | scale (const T &scalar) |
void | scaleRows (const UT_VectorT< T > &scalars) |
void | normalizeRows (T tolerance=1e-5) |
void | rowNorms (UT_VectorT< T > &result, int p=2, bool invert=false, T tolerance=0, T defaultValue=0) const |
void | buildDependencyLevelSet (LevelSet &levelSet, bool buildForUpper=false) const |
void | solveLowerTriangular (UT_VectorT< T > &x, const UT_VectorT< T > &b, bool unitDiagonal=false, T tolerance=1e-5, bool negate=false) const |
void | solveUpperTriangular (UT_VectorT< T > &x, const UT_VectorT< T > &b, bool unitDiagonal=false, T tolerance=1e-5, bool negate=false) const |
void | solveLowerTriangularTranspose (UT_VectorT< T > &x, const UT_VectorT< T > &b, bool unitDiagonal=false, T tolerance=1e-5, bool negate=false) const |
void | solveUpperTriangularTranspose (UT_VectorT< T > &x, const UT_VectorT< T > &b, bool unitDiagonal=false, T tolerance=1e-5, bool negate=false) const |
int | incompleteCholeskyFactorization (T tol=1e-5, bool useModified=false, T tau=0.97, T mindiagratio=0.25) |
void | incompleteLUFactorization (UT_SparseMatrixCSRT< T > &L, UT_SparseMatrixCSRT< T > &U, T tolerance, const LevelSet &workUnits, UT_Interrupt *boss=NULL) const |
void | extractSubMatrix (UT_SparseMatrixCSRT< T > &out, exint rowstart, exint rowend, exint colstart, exint colend) const |
template<typename CellFunction > | |
void | visit (const CellFunction &inspect, bool serial=false) |
template<typename CellFunction > | |
void | removeIf (const CellFunction &condition) |
bool | hasNan () const |
Returns true if any non-zero element in this matrix has value NaN. More... | |
void | testForNan () const |
Asserts that the matrix contains no element with value NaN. More... | |
void | printFull (std::ostream &os) const |
Writes the full matrix, including zeros, to the given output stream. More... | |
void | save (std::ostream &os) const |
Serializes the matrix to the given output stream. More... | |
void | load (UT_IStream &is) |
Deserializes the matrix from the given input stream. More... | |
UT_SparseMatrixCSRT< T > & | operator= (const UT_SparseMatrixCSRT< T > &m) |
Deep copies all data from the given matrix. More... | |
UT_SparseMatrixCSRT< T > & | operator*= (const T &scalar) |
Equivalent to A.scale(scalar) More... | |
Friends | |
class | UT_SparseMatrixBuilderT< T > |
Sparse matrix class that uses the compressed sparse row storage scheme. This is a general-purpose sparse matrix class that can contain an arbitrary number of non-zero elements per row.
Definition at line 662 of file UT_SparseMatrix.h.
typedef UT_Array<UT_ExintArray> UT_SparseMatrixCSRT< T >::LevelSet |
Analyzes the matrix and constructs a directed acyclic graph as a level set, where each level contains rows that are independent of each other. More specifically, level i of the levelSet contains rows that are only dependent on rows present in levels [0, i-1]. The size of the levelSet is the number of levels. If buildForUpper is true, this constructs a level set for the upper triangular part of the matrix. Otherwise, this constructs a level set for the lower triangular part of the matrix.
Definition at line 814 of file UT_SparseMatrix.h.
UT_SparseMatrixCSRT< T >::UT_SparseMatrixCSRT | ( | ) |
Constructs an uninitialized sparse matrix.
UT_SparseMatrixCSRT< T >::UT_SparseMatrixCSRT | ( | exint | rows, |
exint | columns | ||
) |
Constructs a sparse matrix initialized with the given rows and columns.
|
explicit |
Constructs a sparse matrix that is a deep copy of the given sparse matrix.
void UT_SparseMatrixCSRT< T >::addValues | ( | UT_Array< Triplet > & | triplets, |
bool | alreadySorted = false , |
||
T | tolerance = 1e-5 |
||
) |
Adds the given cells to the matrix. Triplets with the same coordinates as cells already present in the matrix will be added together. The array of triplets will be sorted unless alreadySorted is true. Cells whose magnitude is below the tolerance value after the addition will be discarded
void UT_SparseMatrixCSRT< T >::buildDependencyLevelSet | ( | LevelSet & | levelSet, |
bool | buildForUpper = false |
||
) | const |
void UT_SparseMatrixCSRT< T >::extractDiagonal | ( | UT_VectorT< T > & | result, |
bool | invert = false , |
||
T | tolerance = 0 , |
||
T | defaultValue = 0 |
||
) | const |
Writes the diagonal values of the matrix into the given vector If invert is true, the values will be inverted (i.e. 1/x). If a row diagonal's value <= tolerance, then defaultValue is placed in its result row.
void UT_SparseMatrixCSRT< T >::extractSubMatrix | ( | UT_SparseMatrixCSRT< T > & | out, |
exint | rowstart, | ||
exint | rowend, | ||
exint | colstart, | ||
exint | colend | ||
) | const |
Initializes out to a matrix with size (rowend-rowstart, colend-colstart) and fills it with the elements defined by the submatrix in rows [rowstart, rowend), [colstart, colend).
|
inline |
The number of non-zero elements in the given row.
Definition at line 734 of file UT_SparseMatrix.h.
|
inline |
The number of columns in the matrix.
Definition at line 728 of file UT_SparseMatrix.h.
|
inline |
The total number of non-zero elements in the matrix.
Definition at line 731 of file UT_SparseMatrix.h.
|
inline |
The number of rows in the matrix.
Definition at line 725 of file UT_SparseMatrix.h.
bool UT_SparseMatrixCSRT< T >::hasNan | ( | ) | const |
Returns true if any non-zero element in this matrix has value NaN.
int UT_SparseMatrixCSRT< T >::incompleteCholeskyFactorization | ( | T | tol = 1e-5 , |
bool | useModified = false , |
||
T | tau = 0.97 , |
||
T | mindiagratio = 0.25 |
||
) |
Factorizes this matrix into a triangular matrix M such that M*M^T ~= A This keeps the non-zero structure of the matrix. The matrix must be symmetric. This is unstable depending on the diagonal; Returns -1 if the diagonal contains zeros within the given tolerance, and -2 if not positive definite (diagonal has negative numbers). Returns 0 on success. If useModified = true, then this uses the Modified Incomplete Cholesky Factorization, where discarded entries are accounted for by adjusting the diagonal. tau is a tuning constant that scales the adjustment from the missing entries. mindiagratio is the relative tolerance at which the original diagonal will be used if the matrix becomes negative definite as a result of the computation.
void UT_SparseMatrixCSRT< T >::incompleteLUFactorization | ( | UT_SparseMatrixCSRT< T > & | L, |
UT_SparseMatrixCSRT< T > & | U, | ||
T | tolerance, | ||
const LevelSet & | workUnits, | ||
UT_Interrupt * | boss = NULL |
||
) | const |
Factorizes this matrix into L and U where LU approximates this matrix. The non-zero structure of LU matches this matrix; that is, this produces ILU(0) of this matrix. A lower triangular level set is needed for multithreading. e.g. LevelSet levelSet; A.buildDependencyLevelSet(levelSet); A.incompleteLUFactorization(L, U, 1e-5, levelSet);
void UT_SparseMatrixCSRT< T >::init | ( | exint | rows, |
exint | columns | ||
) |
Sets the size of the matrix with the given rows and columns and clears the matrix.
void UT_SparseMatrixCSRT< T >::load | ( | UT_IStream & | is | ) |
Deserializes the matrix from the given input stream.
void UT_SparseMatrixCSRT< T >::multMatrix | ( | const UT_SparseMatrixCSRT< T > & | B, |
UT_SparseMatrixCSRT< T > & | result, | ||
T | tolerance = 0.0 |
||
) | const |
Matrix-matrix multiplication. Computes result = M * B. Resultant elements whose values are smaller than the given tolerance are deleted.
void UT_SparseMatrixCSRT< T >::multVec | ( | const UT_VectorT< T > & | v, |
UT_VectorT< T > & | result | ||
) | const |
Matrix-vector multiplication. Computes result = M * v.
UT_PreciseT<T> UT_SparseMatrixCSRT< T >::multVecAndDot | ( | const UT_VectorT< T > & | v, |
UT_VectorT< T > & | result | ||
) | const |
Matrix-vector multiplication.Computes result = M * v and returns the dot product of v and the result.
void UT_SparseMatrixCSRT< T >::normalizeRows | ( | T | tolerance = 1e-5 | ) |
Multiples each row such that the row diagonal becomes 1. If the row diagonal is below the tolerance, the row is unchanged.
|
inline |
Equivalent to A.scale(scalar)
Definition at line 957 of file UT_SparseMatrix.h.
UT_SparseMatrixCSRT<T>& UT_SparseMatrixCSRT< T >::operator= | ( | const UT_SparseMatrixCSRT< T > & | m | ) |
Deep copies all data from the given matrix.
void UT_SparseMatrixCSRT< T >::printFull | ( | std::ostream & | os | ) | const |
Writes the full matrix, including zeros, to the given output stream.
|
inline |
Condition is called for each cell exactly once. If the condition returns true, then the cell is removed from the matrix. The function signature of the condition should be (exint row, exint col, T value) -> bool. condition can modify the cell vakyes if it has the signature (exint, exint, &T)
Definition at line 925 of file UT_SparseMatrix.h.
void UT_SparseMatrixCSRT< T >::rowNorms | ( | UT_VectorT< T > & | result, |
int | p = 2 , |
||
bool | invert = false , |
||
T | tolerance = 0 , |
||
T | defaultValue = 0 |
||
) | const |
Calculates the norms of each row. result must have the same size as the number of rows in this matrix. p determines the exponent of each column. i.e. result[i] = SUM( |Aij|^p ) for j in row(i) If invert is true, then inverts each sum i.e. result[i] = 1.0 / SUM( |Aij|^p ) for j in row(i) If the row norm is smaller than the tolerance, defaultValue is placed in its result entry instead.
void UT_SparseMatrixCSRT< T >::save | ( | std::ostream & | os | ) | const |
Serializes the matrix to the given output stream.
void UT_SparseMatrixCSRT< T >::scale | ( | const T & | scalar | ) |
Element-wise scalar multiplication. Every element in this matrix is multiplied by the given scalar.
void UT_SparseMatrixCSRT< T >::scaleRows | ( | const UT_VectorT< T > & | scalars | ) |
Element-wise scalar multiplication. Every element in row i of the matrix is multiplied by scalars[i]
void UT_SparseMatrixCSRT< T >::setValues | ( | UT_Array< Triplet > & | triplets, |
T | tolerance = 0.0 |
||
) |
Clears then fills the matrix with the given array of triplets. This sorts the triplets array. Triplets with the same indices will be summed together. A cell whose value is smaller than the tolerance will be deleted.
void UT_SparseMatrixCSRT< T >::setValuesFromRows | ( | UT_Array< UT_Array< ColumnValue >> & | rowValues, |
T | tolerance = 0.0 |
||
) |
Clears then fills the matrix with the given array of ColumnValue arrays. rowValues must have the same size as the number of rows in this matrix. This data layout makes it easier to multithread the construction of the sparse matrix entries.
void UT_SparseMatrixCSRT< T >::solveLowerTriangular | ( | UT_VectorT< T > & | x, |
const UT_VectorT< T > & | b, | ||
bool | unitDiagonal = false , |
||
T | tolerance = 1e-5 , |
||
bool | negate = false |
||
) | const |
Assumes this is a lower triangular matrix. Solves the equation Ax = b. Values may be present in the upper triangular part of the matrix but they will be ignored. If unitDiagonal is true, then 1 will be used instead of the current diagonal values. If a diagonal of A is zero within the tolerance, 0 will be placed in the corresponding row for x. If negate is true, then x is negated x can refer to the same vector as b to solve in-place.
void UT_SparseMatrixCSRT< T >::solveLowerTriangularTranspose | ( | UT_VectorT< T > & | x, |
const UT_VectorT< T > & | b, | ||
bool | unitDiagonal = false , |
||
T | tolerance = 1e-5 , |
||
bool | negate = false |
||
) | const |
Solves the lower triangular system A^T x = b. A is assumed to be an upper triangular matrix. Values may be present in the lower triangular part of the matrix but they will be ignored. If unitDiagonal is true, then 1 will be used instead of the current diagonal values. If a diagonal of A is zero within the tolerance, 0 will be placed in the corresponding row for x. If negate is true, then x is negated.
void UT_SparseMatrixCSRT< T >::solveUpperTriangular | ( | UT_VectorT< T > & | x, |
const UT_VectorT< T > & | b, | ||
bool | unitDiagonal = false , |
||
T | tolerance = 1e-5 , |
||
bool | negate = false |
||
) | const |
Assumes this is an upper triangular matrix. Solves the equation Ax = b. Values may be present in the lower triangular part of the matrix but they will be ignored. If unitDiagonal is true, then 1 will be used instead of the current diagonal values. If a diagonal of A is zero within the tolerance, 0 will be placed in the corresponding row for x. If negate is true, then x is negated x can refer to the same vector as b to solve in-place.
void UT_SparseMatrixCSRT< T >::solveUpperTriangularTranspose | ( | UT_VectorT< T > & | x, |
const UT_VectorT< T > & | b, | ||
bool | unitDiagonal = false , |
||
T | tolerance = 1e-5 , |
||
bool | negate = false |
||
) | const |
Solves the upper triangular system A^T x = b. A is assumed to be a lower triangular matrix. Values may be present in the upper triangular part of the matrix but they will be ignored. If unitDiagonal is true, then 1 will be used instead of the current diagonal values. If a diagonal of A is zero within the tolerance, 0 will be placed in the corresponding row for x. If negate is true, then x is negated.
void UT_SparseMatrixCSRT< T >::subtractMultVec | ( | const UT_VectorT< T > & | x, |
const UT_VectorT< T > & | y, | ||
UT_VectorT< T > & | result | ||
) | const |
Matrix-vector multiplication. Computes result = y - Ax. result can refer to the same vector as y.
void UT_SparseMatrixCSRT< T >::swap | ( | UT_SparseMatrixCSRT< T > & | other | ) |
Swaps size and data elements with the given sparse matrix.
void UT_SparseMatrixCSRT< T >::testForNan | ( | ) | const |
Asserts that the matrix contains no element with value NaN.
void UT_SparseMatrixCSRT< T >::transpose | ( | UT_SparseMatrixCSRT< T > & | result | ) | const |
Computes the transpose of this matrix. result can refer to this matrix.
|
inline |
Calls inspect(row, col, value) on each cell. Each (row, col) is visited exactly once. This is done in parallel over the rows by default, but can be done serially by setting serial = true. inspect can modify the cell values if it has the signature (exint, exint, &T).
Definition at line 908 of file UT_SparseMatrix.h.
void UT_SparseMatrixCSRT< T >::zero | ( | ) |
Clears all data in the matrix.
|
friend |
Definition at line 1015 of file UT_SparseMatrix.h.