jacobi_oo_STL
Public Member Functions | List of all members
CRS_Matrix Class Reference

#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
 

Detailed Description

Square matrix in CRS format (compressed row storage; also named CSR), see an introduction.

Definition at line 42 of file getmatrix.h.

Constructor & Destructor Documentation

◆ CRS_Matrix()

CRS_Matrix::CRS_Matrix ( Mesh const &  mesh)
explicit

Intializes the CRS matrix structure from the given discetization in mesh.

The sparse matrix pattern is generated but the values are 0.

Parameters
[in]meshgiven discretization
Warning
A reference to the discretization mesh is stored inside this class. Therefore, changing mesh outside requires also to call method Derive_Matrix_Pattern explicitely.
See also
Derive_Matrix_Pattern

Definition at line 88 of file getmatrix.cpp.

◆ ~CRS_Matrix()

CRS_Matrix::~CRS_Matrix ( )
inline

Destructor.

Definition at line 63 of file getmatrix.h.

Member Function Documentation

◆ AddElem_3()

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).

Parameters
[in]ialnode indices of the three element vertices
[in]skeelement stiffness matrix
[in]feelement load vector
[in,out]fdistributed local vector storing the right hand side
Warning
Algorithm assumes linear triangular elements (ndof_e==3).

Definition at line 325 of file getmatrix.cpp.

◆ ApplyDirichletBC()

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.

Parameters
[in]u(global) vector with Dirichlet data
[in,out]fload vector

Definition at line 202 of file getmatrix.cpp.

◆ CalculateLaplace()

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.

Parameters
[in,out]f(preallocated) rhs/load vector

Definition at line 169 of file getmatrix.cpp.

◆ 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]nnoderow number of external matrix
[in]idstart indices of matrix rows of external matrix
[in]ikcolumn indices of external matrix
[in]sknon-zero values of external matrix
Returns
true iff all data are identical.

Definition at line 233 of file getmatrix.cpp.

◆ Debug()

void CRS_Matrix::Debug ( ) const

Show the matrix entries.

Definition at line 148 of file getmatrix.cpp.

◆ Defect()

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.

Parameters
[in,out]wresulting vector (preallocated)
[in]fload vector
[in]uvector

Definition at line 286 of file getmatrix.cpp.

◆ Derive_Matrix_Pattern()

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.

◆ fetch()

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.

Parameters
[in]rowrow index
[in]colcolumn 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.

Definition at line 304 of file getmatrix.cpp.

◆ GetDiag()

void CRS_Matrix::GetDiag ( std::vector< double > &  d) const

Extracts the diagonal elemenst of the sparse matrix.

Parameters
[in,out]d(prellocated) vector of diagonal elements

Definition at line 220 of file getmatrix.cpp.

◆ Mult()

void CRS_Matrix::Mult ( std::vector< double > &  w,
std::vector< double > const &  u 
) const

Performs the matrix-vector product w := K*u.

Parameters
[in,out]wresulting vector (preallocated)
[in]uvector

Definition at line 269 of file getmatrix.cpp.

◆ Nrows()

int CRS_Matrix::Nrows ( ) const
inline

Number rows in matrix.

Returns
number of rows.

Definition at line 120 of file getmatrix.h.


The documentation for this class was generated from the following files: