jacobi_oo_STL
|
#include <getmatrix.h>
Public Member Functions | |
CRS_Matrix (Mesh const &mesh) | |
~CRS_Matrix () | |
void | Derive_Matrix_Pattern () |
void | CalculateLaplace (std::vector< double > &f) |
void | ApplyDirichletBC (std::vector< double > const &u, std::vector< double > &f) |
void | GetDiag (std::vector< double > &d) const |
void | Mult (std::vector< double > &w, std::vector< double > const &u) const |
void | Defect (std::vector< double > &w, std::vector< double > const &f, std::vector< double > const &u) const |
int | Nrows () const |
void | Debug () const |
int | fetch (int row, int col) const |
void | AddElem_3 (int const ial[3], double const ske[3][3], double const fe[3], std::vector< double > &f) |
bool | Compare2Old (int nnode, int const id[], int const ik[], double const sk[]) const |
Square matrix in CRS format (compressed row storage; also named CSR), see an introduction.
Definition at line 42 of file getmatrix.h.
|
explicit |
Intializes the CRS matrix structure from the given discetization in mesh
.
The sparse matrix pattern is generated but the values are 0.
[in] | mesh | given discretization |
mesh
is stored inside this class. Therefore, changing mesh
outside requires also to call method Derive_Matrix_Pattern
explicitely.Definition at line 88 of file getmatrix.cpp.
|
inline |
Destructor.
Definition at line 63 of file getmatrix.h.
void CRS_Matrix::AddElem_3 | ( | int const | ial[3], |
double const | ske[3][3], | ||
double const | fe[3], | ||
std::vector< double > & | f | ||
) |
Adds the element stiffness matrix ske
and the element load vector fe
of one triangular element with linear shape functions to the appropriate positions in the stiffness matrix, stored as CSR matrix K(sk
,id
, ik
).
[in] | ial | node indices of the three element vertices |
[in] | ske | element stiffness matrix |
[in] | fe | element load vector |
[in,out] | f | distributed local vector storing the right hand side |
Definition at line 325 of file getmatrix.cpp.
void CRS_Matrix::ApplyDirichletBC | ( | std::vector< double > const & | u, |
std::vector< double > & | f | ||
) |
Applies Dirichlet boundary conditions to stiffness matrix and to load vector f
. The penalty method is used for incorporating the given values u
.
[in] | u | (global) vector with Dirichlet data |
[in,out] | f | load vector |
Definition at line 202 of file getmatrix.cpp.
void CRS_Matrix::CalculateLaplace | ( | std::vector< double > & | f | ) |
Calculates the entries of f.e. stiffness matrix and load/rhs vector f
for the Laplace operator in 2D. No memory is allocated.
[in,out] | f | (preallocated) rhs/load vector |
Definition at line 169 of file getmatrix.cpp.
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.
[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 |
Definition at line 233 of file getmatrix.cpp.
void CRS_Matrix::Debug | ( | ) | const |
Show the matrix entries.
Definition at line 148 of file getmatrix.cpp.
void CRS_Matrix::Defect | ( | std::vector< double > & | w, |
std::vector< double > const & | f, | ||
std::vector< double > const & | u | ||
) | const |
Calculates the defect/residuum w := f - K*u.
[in,out] | w | resulting vector (preallocated) |
[in] | f | load vector |
[in] | u | vector |
Definition at line 286 of file getmatrix.cpp.
void CRS_Matrix::Derive_Matrix_Pattern | ( | ) |
Generates the sparse matrix pattern and overwrites the existing pattern.
The sparse matrix pattern is generated but the values are 0.
Definition at line 95 of file getmatrix.cpp.
int CRS_Matrix::fetch | ( | int | row, |
int | col | ||
) | const |
Finds in a CRS matrix the access index for an entry at row row
and column col
.
[in] | row | row index |
[in] | col | column index |
row
, col
). If no appropriate entry exists then -1 will be returned.row
, col
) doesn't exist. Definition at line 304 of file getmatrix.cpp.
void CRS_Matrix::GetDiag | ( | std::vector< double > & | d | ) | const |
Extracts the diagonal elemenst of the sparse matrix.
[in,out] | d | (prellocated) vector of diagonal elements |
Definition at line 220 of file getmatrix.cpp.
void CRS_Matrix::Mult | ( | std::vector< double > & | w, |
std::vector< double > const & | u | ||
) | const |
Performs the matrix-vector product w := K*u.
[in,out] | w | resulting vector (preallocated) |
[in] | u | vector |
Definition at line 269 of file getmatrix.cpp.
|
inline |