![]() |
Hex
1.0
Hydrogen-electron collision solver
|
Symmetric diagonal matrix. More...
#include <matrix.h>
Public Member Functions | |
SymDiaMatrix () | |
Empty constructor. More... | |
SymDiaMatrix (int n) | |
Size constructor. More... | |
SymDiaMatrix (int n, const iArrayView id, const cArrayView v) | |
Data constructor. More... | |
SymDiaMatrix (SymDiaMatrix const &A) | |
Copy constructor. More... | |
SymDiaMatrix (SymDiaMatrix &&A) | |
Move constructor. More... | |
template<class Functor > | |
SymDiaMatrix & | populate (unsigned d, Functor f) |
Plain symmetrical populator. More... | |
~SymDiaMatrix () | |
void | drop () |
Free all fields, set dimensions to zero. More... | |
int | diag (int i) const |
Diagonal index. More... | |
size_t | size () const |
Matrix dimension. More... | |
int | bandwidth () const |
Bandwidth. More... | |
bool | is_compatible (SymDiaMatrix const &B) const |
Check compatibility of matrices. More... | |
SymDiaMatrix const & | operator= (SymDiaMatrix &&A) |
SymDiaMatrix const & | operator= (SymDiaMatrix const &A) |
SymDiaMatrix const & | operator+= (SymDiaMatrix const &B) |
SymDiaMatrix const & | operator-= (SymDiaMatrix const &B) |
cArray | dot (const cArrayView B, MatrixTriangle triangle=both, bool parallelize=false) const |
Dot product. More... | |
cArray | lowerSolve (const cArrayView b) const |
Back-substitution (lower). More... | |
cArray | upperSolve (const cArrayView b) const |
Back-substitution (upper). More... | |
SymDiaMatrix | kron (SymDiaMatrix const &B) const |
Kronecker product. More... | |
void | link (std::string name) |
Link matrix to a disk file. More... | |
std::string | linkedto () const |
Return the name of the linked disk file. More... | |
cArray | toPaddedRows () const |
Zero-pad rows. More... | |
cArray | toPaddedCols () const |
Zero-pad columns. More... | |
CooMatrix | tocoo (MatrixTriangle triangle=both) const |
Convert matrix part to CooMatrix. More... | |
RowMatrix< Complex > | torow (MatrixTriangle triangle=both) const |
Convert matrix part to RowMatrix. More... | |
iArray const & | diag () const |
Diagonal indices. More... | |
iArray & | diag () |
cArrayView | main_diagonal () const |
Main diagonal. More... | |
cArrayView | main_diagonal () |
cArray const & | data () const |
Data pointer. More... | |
cArray & | data () |
Complex const * | dptr (int i) const |
Pointer to diagonal data. More... | |
Complex * | dptr (int i) |
bool | hdfload () |
Load from file. More... | |
bool | hdfload (std::string name) |
bool | hdfsave () const |
Save data to file. More... | |
bool | hdfsave (std::string name, bool docompress=false, int consec=10) const |
Friends | |
std::ostream & | operator<< (std::ostream &out, SymDiaMatrix const &A) |
Output to a text stream. More... | |
A multi-diagonal matrix is a sparse matrix of the structure
\[ A = \pmatrix { \ast & \dots & \ast & & & \cr \vdots & \ddots & & \ddots & & \cr \ast & & \ddots & & \ddots & \cr & \ddots & & \ddots & & \ast \cr & & \ddots & & \ddots & \vdots \cr & & & \ast & \dots & \ast \cr } \ , \]
i.e. it is banded, with nonzero elements only near to the diagonal. Some of the diagonals may be identically zero. This class holds all nonzero main and upper diagonals (lower diagonals are not necessary in the symmetric case). The diagonal storage has several advantages:
SymDiaMatrix::SymDiaMatrix | ( | ) |
SymDiaMatrix::SymDiaMatrix | ( | int | n | ) |
SymDiaMatrix::SymDiaMatrix | ( | int | n, |
const iArrayView | id, | ||
const cArrayView | v | ||
) |
n | Size of the matrix. |
id | Identifyiers of the diagonals (positive integers expected). |
v | Stacked (and padded if ncessary) diagonals. |
SymDiaMatrix::SymDiaMatrix | ( | SymDiaMatrix const & | A | ) |
SymDiaMatrix::SymDiaMatrix | ( | SymDiaMatrix && | A | ) |
|
inline |
|
inline |
Return the bandwidth of the matrix, i.e. number of all (upper, main an lower) diagonals that would have to be stored in a full banded-matrix format.
|
inline |
Return direct-access data pointer.
|
inline |
|
inline |
Return array of indices of the stored diagonals. These are always only main and upper diagonals, so all numbers are non-negative. The array is sorted and begins with zero. Its length is always larger than zero.
|
inline |
|
inline |
Return diagonal index for of i-th stored diagonal. The zero-th stored diagonal is always the main diagonal (= 0), but it doesn't have to hold for next diagonals.
cArray SymDiaMatrix::dot | ( | const cArrayView | B, |
MatrixTriangle | triangle = both , |
||
bool | parallelize = false |
||
) | const |
This is a key member of the structure, defining e.g. the speed of conjugate gradients and evaluation of the scattering amplitudes.
B | Dense matrix. It is supposed to be stored by columns and to have dimensions n times k, where n is the column count of (*this) matrix. Also, though only a view of the array is required, it is assumed that B is actually rArray, i.e. that it is aligned with the alignment 2*sizeof(T). |
triangle | Whether to use only the upper or only the lower or both triangles of the othwerwise symmetric matrix. |
parallelize | Whether to use OpenMP to parallelize the SpMV operation. |
|
inline |
i | Index of the diagonal in the "idiag_" array. The maximal value is thus less then the number stored diagonals. |
|
inline |
|
inline |
|
inline |
Load the matrix from a HDF5 file created by the routine hdfsave.
bool SymDiaMatrix::hdfload | ( | std::string | name | ) |
|
inline |
Save the matrix to the HDF5 file as a set of four datasets:
bool SymDiaMatrix::hdfsave | ( | std::string | name, |
bool | docompress = false , |
||
int | consec = 10 |
||
) | const |
bool SymDiaMatrix::is_compatible | ( | SymDiaMatrix const & | B | ) | const |
Check that the matrix B has the same dimensions as *this matrix and also that they keep the same diagonals. Such matrices can be very effectively summed and subtracted – just by doing the operation on the stored element arrays.
SymDiaMatrix SymDiaMatrix::kron | ( | SymDiaMatrix const & | B | ) | const |
Compute Kronecker product with other matrix.
void SymDiaMatrix::link | ( | std::string | name | ) |
|
inline |
cArray SymDiaMatrix::lowerSolve | ( | const cArrayView | b | ) | const |
Assume the matrix is normalized lower-triangular (i.e. has unit main diagonal and zero upper triangle) and do the triangular solve.
b | Right hand side of the triangular system. |
|
inline |
Return direct-access view of the main diagonal.
|
inline |
SymDiaMatrix const & SymDiaMatrix::operator+= | ( | SymDiaMatrix const & | B | ) |
SymDiaMatrix const & SymDiaMatrix::operator-= | ( | SymDiaMatrix const & | B | ) |
SymDiaMatrix const & SymDiaMatrix::operator= | ( | SymDiaMatrix && | A | ) |
SymDiaMatrix const & SymDiaMatrix::operator= | ( | SymDiaMatrix const & | A | ) |
|
inline |
Given a functor of the signature
the function will call the functor with row and column number of every element that is to be set.
d | How many upper diagonals to populate. The main diagonal will be populated always. |
f | The functor that will compute the matrix elements. |
|
inline |
Return row/column count. The matrix is symmetric and so both counts are equal.
CooMatrix SymDiaMatrix::tocoo | ( | MatrixTriangle | triangle = both | ) | const |
cArray SymDiaMatrix::toPaddedCols | ( | ) | const |
Pad rows with zeros as in toPaddedRows, then concatenate columns and return as a single array.
cArray SymDiaMatrix::toPaddedRows | ( | ) | const |
Pad rows with zeros as below:
\[ \left( \matrix { \ast & \ast & & & \cr \ast & \ast & \ast & & \cr \ast & \ast & \ast & \ast & \cr \ast & \ast & \ast & \ast & \ast \cr & \ast & \ast & \ast & \ast \cr & & \ast & \ast & \ast \cr & & & \ast & \ast \cr } \right) \longrightarrow \matrix { 0 & 0 & 0 \cr & 0 & 0 \cr & & 0 \cr & & . \cr & & . \cr & & . \cr & & . \cr } \left( \matrix { \ast & \ast & & & \cr \ast & \ast & \ast & & \cr \ast & \ast & \ast & \ast & \cr \ast & \ast & \ast & \ast & \ast \cr & \ast & \ast & \ast & \ast \cr & & \ast & \ast & \ast \cr & & & \ast & \ast \cr } \right) \matrix { . & & \cr . & & \cr . & & \cr . & & \cr 0 & & \cr 0 & 0 & \cr 0 & 0 & 0 \cr } \longrightarrow \matrix { 0 & 0 & 0 & \ast & \ast \cr 0 & 0 & \ast & \ast & \ast \cr 0 & \ast & \ast & \ast & \ast \cr \ast & \ast & \ast & \ast & \ast \cr \ast & \ast & \ast & \ast & 0 \cr \ast & \ast & \ast & 0 & 0 \cr \ast & \ast & 0 & 0 & 0 \cr } \]
Then concatenate rows and return as a single array.
RowMatrix< Complex > SymDiaMatrix::torow | ( | MatrixTriangle | triangle = both | ) | const |
cArray SymDiaMatrix::upperSolve | ( | const cArrayView | b | ) | const |
Assume the matrix is normalized upper-triangular (i.e. has unit main diagonal and zero lower triangle) and do the triangular solve.
b | Right hand side of the triangular system. |
|
friend |