#include <getmatrix.h>
Matrix in CRS format (compressed row storage; also named CSR), see an introduction.
Definition at line 131 of file getmatrix.h.
◆ CRS_Matrix() [1/2]
CRS_Matrix::CRS_Matrix |
( |
| ) |
|
◆ CRS_Matrix() [2/2]
◆ ~CRS_Matrix()
CRS_Matrix::~CRS_Matrix |
( |
| ) |
|
|
overridevirtual |
◆ Compare2Old()
bool CRS_Matrix::Compare2Old |
( |
int |
nnode, |
|
|
int const |
id[], |
|
|
int const |
ik[], |
|
|
double const |
sk[] |
|
) |
| const |
Compare this
CRS matrix with an external CRS matrix stored in C-Style.
The method prints statements on differences found.
- Parameters
-
[in] | nnode | row number of external matrix |
[in] | id | start indices of matrix rows of external matrix |
[in] | ik | column indices of external matrix |
[in] | sk | non-zero values of external matrix |
- Returns
- true iff all data are identical.
Definition at line 146 of file getmatrix.cpp.
◆ Debug()
void CRS_Matrix::Debug |
( |
| ) |
const |
|
overridevirtual |
◆ Defect()
void CRS_Matrix::Defect |
( |
std::vector< double > & |
w, |
|
|
std::vector< double > const & |
f, |
|
|
std::vector< double > const & |
u |
|
) |
| const |
|
overridevirtual |
Calculates the defect/residuum w := f - K*u.
- Parameters
-
[in,out] | w | resulting vector (preallocated) |
[in] | f | load vector |
[in] | u | vector |
Implements Matrix.
Definition at line 66 of file getmatrix.cpp.
◆ fetch()
int CRS_Matrix::fetch |
( |
int |
row, |
|
|
int |
col |
|
) |
| const |
|
inlineoverridevirtual |
Finds in a CRS matrix the access index for an entry at row row
and column col
.
- Parameters
-
[in] | row | row index |
[in] | col | column index |
- Returns
- index for element (
row
, col
). If no appropriate entry exists then -1 will be returned.
- Warning
- assert() stops the function in case that matrix element (
row
, col
) doesn't exist.
Implements Matrix.
Definition at line 103 of file getmatrix.cpp.
◆ GetDiag()
void CRS_Matrix::GetDiag |
( |
std::vector< double > & |
d | ) |
const |
|
overridevirtual |
Extracts the diagonal elements of the sparse matrix.
- Parameters
-
[in,out] | d | (prellocated) vector of diagonal elements |
Implements Matrix.
Definition at line 86 of file getmatrix.cpp.
◆ Mult()
void CRS_Matrix::Mult |
( |
std::vector< double > & |
w, |
|
|
std::vector< double > const & |
u |
|
) |
| const |
|
overridevirtual |
**
Performs the matrix-vector product w := K*u.
- Parameters
-
[in,out] | w | resulting vector (preallocated) |
[in] | u | vector |
Implements Matrix.
Definition at line 49 of file getmatrix.cpp.
◆ DefectRestrict
void DefectRestrict |
( |
CRS_Matrix const & |
SK, |
|
|
BisectInterpolation const & |
P, |
|
|
std::vector< double > & |
fc, |
|
|
std::vector< double > & |
ff, |
|
|
std::vector< double > & |
uf |
|
) |
| |
|
friend |
Calculates the defect and projects it to the next coarser level \( f_C := P^T \cdot (f_F - SK\cdot u_F) \).
- Parameters
-
[in] | SK | matrix on fine mesh |
[in] | P | prolongation operator |
[in,out] | fc | resulting coarse mesh vector (preallocated) |
[in] | ff | r.h.s. on fine mesh |
[in] | uf | status vector on fine mesh |
◆ _id
std::vector<int> CRS_Matrix::_id |
|
protected |
start indices of matrix rows
Definition at line 225 of file getmatrix.h.
◆ _ik
std::vector<int> CRS_Matrix::_ik |
|
protected |
◆ _nnz
number of non-zero entries
Definition at line 224 of file getmatrix.h.
◆ _sk
std::vector<double> CRS_Matrix::_sk |
|
protected |
The documentation for this class was generated from the following files: