ug4
lib_algebra

Algebra Library. More...

Modules

 CPU Algebra
 
 CRS Algebra
 
 lib_algebra Common
 Utilities for libAlgebra.
 
 Parallel Algebra
 
 Small Algebra
 

Classes

class  ug::AlgebraType
 class describing the type of an algebra More...
 
class  ug::DefaultAlgebra
 Singleton, providing the current default algebra. More...
 
struct  ug::DenseMatrixFromSparseMatrix< TSparseMatrix >
 
class  ug::MapSparseMatrix< TValueType >
 sparse matrix for big, variable sparse matrices. More...
 
class  ug::MapVector< TValueType >
 
struct  ug::matrix_algebra_type_traits< MapSparseMatrix< T > >
 
struct  ug::sortStruct< T >
 

Macros

#define FORCE_CREATION   volatile size_t ___never_happens___ = 0; if(___never_happens___)
 
#define UNUSED_VARIABLE(var)   ((void) var);
 prevent unused variable-warnings More...
 

Typedefs

typedef TValueType ug::MapVector< TValueType >::value_type
 
typedef MapVector< TValueType > ug::MapVector< TValueType >::vector_type
 

Functions

template<typename ABC_type , typename A_type , typename B_type , typename C_type >
void ug::AddMultiplyOf (ABC_type &M, const A_type &A, const B_type &B, const C_type &C, double epsilonTruncation=0.0)
 Calculates M += A*B*C. More...
 
template<typename TSparseMatrix , typename vector_t >
bool ug::Axpy_transposedCommonSparseMatrix (const TSparseMatrix &A, vector_t &dest, const number &alpha1, const vector_t &v1, const number &beta1, const vector_t &w1)
 
template<typename TSparseMatrix , typename vector_t >
bool ug::AxpyCommonSparseMatrix (const TSparseMatrix &A, vector_t &dest, const number &alpha1, const vector_t &v1, const number &beta1, const vector_t &w1)
 
template<typename TSparseMatrix >
bool ug::CheckDiagonalInvertible (const ParallelMatrix< TSparseMatrix > &m)
 
template<typename TSparseMatrix >
bool ug::CheckDiagonalInvertible (const TSparseMatrix &A)
 
template<typename TSparseMatrix >
bool ug::CheckRowIterators (const TSparseMatrix &A)
 
template<typename TVector >
bool ug::CheckVectorInvertible (const ParallelVector< TVector > &v)
 
template<typename TVector >
bool ug::CheckVectorInvertible (const TVector &v)
 
template<typename AB_type , typename A_type , typename B_type >
void ug::CreateAsMultiplyOf (AB_type &M, const A_type &A, const B_type &B)
 Calculates M = A*B. More...
 
template<typename ABC_type , typename A_type , typename B_type , typename C_type >
void ug::CreateAsMultiplyOf (ABC_type &M, const A_type &A, const B_type &B, const C_type &C, double epsilonTruncation=0.0)
 Calculates M = A*B*C. More...
 
template<typename TSparseMatrix , class TIStream >
void ug::DeserializeMatrix (TIStream &buf, TSparseMatrix &A)
 
template<typename Matrix_type , typename Vector_type >
void ug::diag_step (const Matrix_type &A, Vector_type &c, const Vector_type &d, number damp)
 Performs a jacobi-step. More...
 
template<typename TDenseType , typename TSparseMatrix >
size_t ug::GetDenseDoubleFromSparse (TDenseType &A, const TSparseMatrix &S)
 
template<typename TSparseMatrix >
DenseMatrixFromSparseMatrix< TSparseMatrix >::type & ug::GetDenseFromSparse (typename DenseMatrixFromSparseMatrix< TSparseMatrix >::type &A, const TSparseMatrix &S)
 
template<typename TDoubleType , typename TSparseMatrix >
void ug::GetDoubleFromSparseBlock (TDoubleType &A, const TSparseMatrix &S)
 
template<typename TSparseMatrix >
size_t ug::GetDoubleSize (const TSparseMatrix &S)
 
template<typename TDoubleSparse , typename TSparseMatrix >
size_t ug::GetDoubleSparseFromBlockSparse (TDoubleSparse &A, const TSparseMatrix &S)
 
template<typename TSparseMatrix >
size_t ug::GetMaxConnections (const TSparseMatrix &A)
 returns max number of non-zero connections in rows More...
 
template<typename TSparseMatrix >
void ug::GetNeighborhood (const TSparseMatrix &A, size_t node, size_t depth, std::vector< size_t > &indices)
 
template<typename TSparseMatrix >
void ug::GetNeighborhood (const TSparseMatrix &A, size_t node, size_t depth, std::vector< size_t > &indices, std::vector< bool > &bVisited, bool bResetVisitedFlags=true)
 
template<typename TSparseMatrix >
void ug::GetNeighborhood_worker (const TSparseMatrix &A, size_t node, size_t depth, std::vector< size_t > &indices, std::vector< bool > &bVisited)
 
template<typename TSparseMatrix >
void ug::GetNeighborhoodHierachy (const TSparseMatrix &A, size_t node, size_t depth, std::vector< std::vector< size_t > > &indices)
 
template<typename TSparseMatrix >
void ug::GetNeighborhoodHierachy (const TSparseMatrix &A, size_t node, size_t depth, std::vector< std::vector< size_t > > &indices, std::vector< bool > &bVisited, bool bResetVisitedFlags=true)
 
template<typename TSparseMatrix >
void ug::GetNeighborhoodHierachy_worker (const TSparseMatrix &A, size_t node, size_t depth, size_t maxdepth, std::vector< std::vector< size_t > > &indices, std::vector< bool > &bVisited)
 
template<typename TSparseMatrix >
size_t ug::GetNNZs (const TSparseMatrix &A)
 returns the number of non-zeroes (!= number of connections) More...
 
template<typename Matrix_type , typename Vector_type >
void ug::gs_step_LL (const Matrix_type &A, Vector_type &c, const Vector_type &d, const number relaxFactor)
 Gauss-Seidel-Iterations. More...
 
template<typename Matrix_type , typename Vector_type >
bool ug::gs_step_LL (const Matrix_type &A, Vector_type &x, const Vector_type &b)
 Performs a forward gauss-seidel-step, that is, solve on the lower left of A. More...
 
template<typename Matrix_type , typename Vector_type >
void ug::gs_step_UR (const Matrix_type &A, Vector_type &c, const Vector_type &d, const number relaxFactor)
 Performs a backward gauss-seidel-step, that is, solve on the upper right of A. Using gs_step_UR within a preconditioner-scheme leads to the fact that we get the correction by successive inserting the already computed values of c in c = N * d, with c being the correction and d being the defect. N denotes the matrix of the second normal-form of a linear iteration scheme. More...
 
template<typename Matrix_type , typename Vector_type >
bool ug::gs_step_UR (const Matrix_type &A, Vector_type &x, const Vector_type &b)
 Performs a backward gauss-seidel-step, that is, solve on the upper right of A. More...
 
template<typename TSparseMatrix >
bool ug::IsCloseToBoundary (const TSparseMatrix &A, size_t node, size_t distance)
 gets the neighborhood of a node in the connectivity graph of a SparseMatrix. More...
 
template<typename TSparseMatrix >
bool ug::IsDirichletRow (const TSparseMatrix &A, size_t i, size_t alpha)
 Evaluates 'true', iff corresponding row is Dirichlet. More...
 
 ug::MapVector< TValueType >::MapVector ()
 constructor More...
 
 ug::MapVector< TValueType >::MapVector (size_t _length)
 constructor with length More...
 
template<typename TSparseMatrix >
void ug::MarkNeighbors (const TSparseMatrix &A, size_t node, size_t depth, std::vector< bool > &bVisited)
 
template<typename matrix_type >
void ug::MatAdd (matrix_type &M, number alpha1, const matrix_type &A, number alpha2, const matrix_type &B)
 Calculates M = A + B. More...
 
template<typename matrix_type >
void ug::MatAddNonDirichlet (matrix_type &M, number alpha1, const matrix_type &A, number alpha2, const matrix_type &B)
 Calculates M = A + B. More...
 
template<typename vector_t , typename matrix_t >
void ug::MatMultTransposedAdd (vector_t &dest, const number &alpha1, const vector_t &v1, const number &beta1, const MapSparseMatrix< matrix_t > &A1, const vector_t &w1)
 calculates dest = alpha1*v1 + beta1 * A1^T *w1; More...
 
std::ostream & ug::operator<< (std::ostream &out, const AlgebraType &v)
 writes the Identifier to the output stream More...
 
template<typename TSparseMatrix >
void ug::ScaleRow (TSparseMatrix &A, size_t i, number fac)
 
template<typename TSparseMatrix >
void ug::ScaleSparseMatrixCommon (TSparseMatrix &A, double d)
 
template<typename TSparseMatrix , class TOStream >
void ug::SerializeMatrix (TOStream &buf, const TSparseMatrix &A)
 
template<typename TSparseMatrix >
void ug::SetCol (TSparseMatrix &A, size_t i, size_t alpha, number val=0.0)
 
template<typename TSparseMatrix >
void ug::SetDirichletRow (TSparseMatrix &A, const std::vector< size_t > vIndex)
 
template<typename TSparseMatrix >
void ug::SetDirichletRow (TSparseMatrix &A, size_t i)
 
template<typename TSparseMatrix >
void ug::SetDirichletRow (TSparseMatrix &A, size_t i, size_t alpha)
 
template<typename TSparseMatrix >
void ug::SetRow (TSparseMatrix &A, size_t i, number val=0.0)
 
template<typename TSparseMatrix >
void ug::SetRow (TSparseMatrix &A, size_t i, size_t alpha, number val=0.0)
 
template<typename Matrix_type , typename Vector_type >
void ug::sgs_step (const Matrix_type &A, Vector_type &c, const Vector_type &d, const number relaxFactor)
 Performs a symmetric gauss-seidel step. Using sgs_step within a preconditioner-scheme leads to the fact that we get the correction by successive inserting the already computed values of c in c = N * d, with c being the correction and d being the defect. N denotes the matrix of the second normal-form of a linear iteration scheme. More...
 
template<typename Matrix_type , typename Vector_type >
bool ug::sgs_step (const Matrix_type &A, Vector_type &x, const Vector_type &b)
 Performs a symmetric gauss-seidel step. More...
 

Detailed Description

Algebra Library.

Macro Definition Documentation

◆ FORCE_CREATION

#define FORCE_CREATION   volatile size_t ___never_happens___ = 0; if(___never_happens___)

use this to force the creation of prsize_t routines or similar for use in gdb.

◆ UNUSED_VARIABLE

#define UNUSED_VARIABLE (   var)    ((void) var);

prevent unused variable-warnings

Typedef Documentation

◆ value_type

template<typename TValueType >
typedef TValueType ug::MapVector< TValueType >::value_type

◆ vector_type

template<typename TValueType >
typedef MapVector<TValueType> ug::MapVector< TValueType >::vector_type

Function Documentation

◆ AddMultiplyOf()

template<typename ABC_type , typename A_type , typename B_type , typename C_type >
void ug::AddMultiplyOf ( ABC_type M,
const A_type A,
const B_type B,
const C_type &  C,
double  epsilonTruncation = 0.0 
)

Calculates M += A*B*C.

Parameters
M(in/out) Matrix M, M += A*B*C$
A(in) Matrix A
B(in) Matrix B
C(in) Matrix C

Complete formula for calculating M=A*B*C:

\[ M_{ij} += \sum_{kl} A_{ik} * B_{kl} * C_{lj} \]

References ug::AddMult(), ug::AssignMult(), ug::BlockNorm(), PROFILE_FUNC_GROUP, UG_ASSERT, and UG_THROW.

Referenced by ug::AssembledMultiGridCycle< TDomain, TAlgebra >::init_rap_operator().

◆ Axpy_transposedCommonSparseMatrix()

template<typename TSparseMatrix , typename vector_t >
bool ug::Axpy_transposedCommonSparseMatrix ( const TSparseMatrix &  A,
vector_t &  dest,
const number alpha1,
const vector_t &  v1,
const number beta1,
const vector_t &  w1 
)

◆ AxpyCommonSparseMatrix()

template<typename TSparseMatrix , typename vector_t >
bool ug::AxpyCommonSparseMatrix ( const TSparseMatrix &  A,
vector_t &  dest,
const number alpha1,
const vector_t &  v1,
const number beta1,
const vector_t &  w1 
)

◆ CheckDiagonalInvertible() [1/2]

template<typename TSparseMatrix >
bool ug::CheckDiagonalInvertible ( const ParallelMatrix< TSparseMatrix > &  m)

◆ CheckDiagonalInvertible() [2/2]

template<typename TSparseMatrix >
bool ug::CheckDiagonalInvertible ( const TSparseMatrix &  A)

◆ CheckRowIterators()

template<typename TSparseMatrix >
bool ug::CheckRowIterators ( const TSparseMatrix &  A)

◆ CheckVectorInvertible() [1/2]

template<typename TVector >
bool ug::CheckVectorInvertible ( const ParallelVector< TVector > &  v)

◆ CheckVectorInvertible() [2/2]

template<typename TVector >
bool ug::CheckVectorInvertible ( const TVector &  v)

◆ CreateAsMultiplyOf() [1/2]

template<typename AB_type , typename A_type , typename B_type >
void ug::CreateAsMultiplyOf ( AB_type &  M,
const A_type A,
const B_type B 
)

Calculates M = A*B.

Parameters
M(out) Matrix M, M = A*B$
A(in) Matrix A
B(in) Matrix B \( M_{ij} = \sum_k A_{ik} * B_{kj} \) For implementation details, see also CreateAsMultiplyOf(M, A, B, C).

References ug::AddMult(), PROFILE_FUNC_GROUP, and UG_ASSERT.

◆ CreateAsMultiplyOf() [2/2]

template<typename ABC_type , typename A_type , typename B_type , typename C_type >
void ug::CreateAsMultiplyOf ( ABC_type M,
const A_type A,
const B_type B,
const C_type &  C,
double  epsilonTruncation = 0.0 
)

Calculates M = A*B*C.

Parameters
M(out) Matrix M, M = A*B*C$
A(in) Matrix A
B(in) Matrix B
C(in) Matrix C

Complete formula for calculating M=A*B*C:

\[ M_{ij} = \sum_{kl} A_{ik} * B_{kl} * C_{lj} \]

Calculation is done on row-basis without a temporary BC or AB matrix. This has shown to be much faster than implementations with temporary matrices due to cache effects. We also added an improved way of storing the results of the calculation: when we go through connections of B and C and want to add the connection (i, j) to M, we need to know if this connection already exists. For this we have an array posInConnections, needs n=A.num_rows() memory. posInConnections[i]: index in the connections for current row (if not in row: -1) tried this also with std::map, but took 1511.53 ms instead of 393.972 ms searching in the connections array is also slower

References ug::AddMult(), ug::AssignMult(), ug::BlockNorm(), PROFILE_FUNC_GROUP, and UG_ASSERT.

◆ DeserializeMatrix()

template<typename TSparseMatrix , class TIStream >
void ug::DeserializeMatrix ( TIStream &  buf,
TSparseMatrix &  A 
)

◆ diag_step()

template<typename Matrix_type , typename Vector_type >
void ug::diag_step ( const Matrix_type &  A,
Vector_type &  c,
const Vector_type &  d,
number  damp 
)

Performs a jacobi-step.

Parameters
AMatrix \(A = D - L - R\)
cwill be \(c = N * d = D^{-1} d \)
dthe vector d.
dampthe damping factor

References ug::InverseMatMult(), and UG_ASSERT.

◆ GetDenseDoubleFromSparse()

template<typename TDenseType , typename TSparseMatrix >
size_t ug::GetDenseDoubleFromSparse ( TDenseType &  A,
const TSparseMatrix &  S 
)

◆ GetDenseFromSparse()

template<typename TSparseMatrix >
DenseMatrixFromSparseMatrix<TSparseMatrix>::type& ug::GetDenseFromSparse ( typename DenseMatrixFromSparseMatrix< TSparseMatrix >::type &  A,
const TSparseMatrix &  S 
)

◆ GetDoubleFromSparseBlock()

template<typename TDoubleType , typename TSparseMatrix >
void ug::GetDoubleFromSparseBlock ( TDoubleType &  A,
const TSparseMatrix &  S 
)

◆ GetDoubleSize()

template<typename TSparseMatrix >
size_t ug::GetDoubleSize ( const TSparseMatrix &  S)

◆ GetDoubleSparseFromBlockSparse()

template<typename TDoubleSparse , typename TSparseMatrix >
size_t ug::GetDoubleSparseFromBlockSparse ( TDoubleSparse &  A,
const TSparseMatrix &  S 
)

◆ GetMaxConnections()

template<typename TSparseMatrix >
size_t ug::GetMaxConnections ( const TSparseMatrix &  A)

returns max number of non-zero connections in rows

◆ GetNeighborhood() [1/2]

template<typename TSparseMatrix >
void ug::GetNeighborhood ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
std::vector< size_t > &  indices 
)

◆ GetNeighborhood() [2/2]

template<typename TSparseMatrix >
void ug::GetNeighborhood ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
std::vector< size_t > &  indices,
std::vector< bool > &  bVisited,
bool  bResetVisitedFlags = true 
)

◆ GetNeighborhood_worker()

template<typename TSparseMatrix >
void ug::GetNeighborhood_worker ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
std::vector< size_t > &  indices,
std::vector< bool > &  bVisited 
)

Referenced by ug::GetNeighborhood().

◆ GetNeighborhoodHierachy() [1/2]

template<typename TSparseMatrix >
void ug::GetNeighborhoodHierachy ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
std::vector< std::vector< size_t > > &  indices 
)

◆ GetNeighborhoodHierachy() [2/2]

template<typename TSparseMatrix >
void ug::GetNeighborhoodHierachy ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
std::vector< std::vector< size_t > > &  indices,
std::vector< bool > &  bVisited,
bool  bResetVisitedFlags = true 
)

◆ GetNeighborhoodHierachy_worker()

template<typename TSparseMatrix >
void ug::GetNeighborhoodHierachy_worker ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
size_t  maxdepth,
std::vector< std::vector< size_t > > &  indices,
std::vector< bool > &  bVisited 
)

◆ GetNNZs()

template<typename TSparseMatrix >
size_t ug::GetNNZs ( const TSparseMatrix &  A)

returns the number of non-zeroes (!= number of connections)

◆ gs_step_LL() [1/2]

template<typename Matrix_type , typename Vector_type >
void ug::gs_step_LL ( const Matrix_type &  A,
Vector_type &  c,
const Vector_type &  d,
const number  relaxFactor 
)

Gauss-Seidel-Iterations.

Here, iteration schemes of gauss-seidel-type are described for solving the linear equation

\f$ A * x = b.      A \in R^{nxn}, x \in R^n, b \in R^n \f$.

Most of the common linear iteration-methods base on the decomposition of A into its diagonal (D) and strict-upper(-U) and strict-lower part (-L),

\( A = D - L - U \).

Among others, W. Hackbusch ('Iterative Loesung grosser Gleichungssysteme'), distinguishes three different forms for describing a linear iteration scheme. The general 'first normal-form' of a linear iteration scheme takes the form

\( x^{m+1} = M * x^m + N * b \),

with some Matrices \( M \) and \( N \in R^{nxn} \). m denotes the iteration index. The general 'second normal-form' of a linear iteration scheme takes the form

\( x^{m+1} = x^m - N * (A * x^m - b) \).

Those linear iteration schemes, which can be represented by the second normal-form are the linear, consistent iteration schemes.

Introducing the correction \( c{m+1} := x^{m+1} - x^m \) and the defect \( d^m := b - A * x^m \) the second normal-form can be rewritten as

\( c = N * d \).

Below, methods for the (forward) Gauss-Seidel step, the backward Gauss-Seidel step and the symmetric Gauss-Seidel step are implemented ('gs_step_LL', 'gs_step_UR' resp. 'sgs_step'). The matrices of the second normal-form are

\( N = (D-L)^{-1}\) for the (forward) Gauss-Seidel step, \( N = (D-U)^{-1}\) for the backward Gauss-Seidel step and \( N = (D-u)^{-1}D(D-U)^{-1} \) for the symmetric Gauss-Seidel step.

References:

  • W. Hackbusch. Iterative Loesung grosser Gleichungssysteme
See also
gs_step_LL, gs_step_UR, sgs_step

Performs a forward gauss-seidel-step, that is, solve on the lower left of A. Using gs_step_LL within a preconditioner-scheme leads to the fact that we get the correction by successive inserting the already computed values of c in c = N * d, with c being the correction and d being the defect. N denotes the matrix of the second normal-form of a linear iteration scheme.

Parameters
AMatrix \(A = D - L - U\)
cVector. \( c = N * d = (D-L)^{-1} * d \)
dVector d.
See also
gs_step_UR, sgs_step

References ug::InverseMatMult(), ug::MatMultAdd(), and s.

Referenced by ug::sgs_step(), and ug::GaussSeidel< TAlgebra >::step().

◆ gs_step_LL() [2/2]

template<typename Matrix_type , typename Vector_type >
bool ug::gs_step_LL ( const Matrix_type &  A,
Vector_type &  x,
const Vector_type &  b 
)

Performs a forward gauss-seidel-step, that is, solve on the lower left of A.

Parameters
AMatrix \(A = D - L - U\)
xwill be \(x = (D-L)^{-1}b \)
bthe vector b.
See also
gs_step_UR, sgs_step

References ug::InverseMatMult(), ug::MatMultAdd(), and s.

◆ gs_step_UR() [1/2]

template<typename Matrix_type , typename Vector_type >
void ug::gs_step_UR ( const Matrix_type &  A,
Vector_type &  c,
const Vector_type &  d,
const number  relaxFactor 
)

Performs a backward gauss-seidel-step, that is, solve on the upper right of A. Using gs_step_UR within a preconditioner-scheme leads to the fact that we get the correction by successive inserting the already computed values of c in c = N * d, with c being the correction and d being the defect. N denotes the matrix of the second normal-form of a linear iteration scheme.

Parameters
AMatrix \(A = D - L - U\)
cwill be \(c = N * d = (D-U)^{-1} * d \)
dthe vector d.
See also
gs_step_LL, sgs_step

References ug::InverseMatMult(), ug::MatMultAdd(), and s.

Referenced by ug::sgs_step(), and ug::BackwardGaussSeidel< TAlgebra >::step().

◆ gs_step_UR() [2/2]

template<typename Matrix_type , typename Vector_type >
bool ug::gs_step_UR ( const Matrix_type &  A,
Vector_type &  x,
const Vector_type &  b 
)

Performs a backward gauss-seidel-step, that is, solve on the upper right of A.

Parameters
AMatrix \(A = D - L - U\)
xwill be \(x = (D-U)^{-1}b \)
bthe vector b.
See also
gs_step_LL, sgs_step

References ug::InverseMatMult(), ug::MatMultAdd(), and s.

◆ IsCloseToBoundary()

template<typename TSparseMatrix >
bool ug::IsCloseToBoundary ( const TSparseMatrix &  A,
size_t  node,
size_t  distance 
)

gets the neighborhood of a node in the connectivity graph of a SparseMatrix.

Parameters
A(in) Matrix A
node(in) the node where to start
depth(in) the depth of neighborhood. 0 = empty.
indices(out) the indices of the neighbors
posInConnectionsarray to speed up computation. Has to be posInConnections[i] = 0 for all i=0..A.num_rows(). Can be NULL.
Note
the node itself is only included if there is a connection from node to node.

determines if a node is close to a unconnected node in the connectivity graph of a SparseMatrix.

Parameters
A(in) Matrix A
node(in) the node where to start
distance(in) up to which distance "close" is.
Returns
if there is a distance long path in graph(A) to an unconnected node, true. otherwise false.

◆ IsDirichletRow()

template<typename TSparseMatrix >
bool ug::IsDirichletRow ( const TSparseMatrix &  A,
size_t  i,
size_t  alpha 
)

Evaluates 'true', iff corresponding row is Dirichlet.

References alpha, ug::BlockRef(), and ug::GetCols().

Referenced by ug::MatAddNonDirichlet().

◆ MapVector() [1/2]

template<typename TValueType >
ug::MapVector< TValueType >::MapVector ( )
inline

constructor

◆ MapVector() [2/2]

template<typename TValueType >
ug::MapVector< TValueType >::MapVector ( size_t  _length)
inline

constructor with length

◆ MarkNeighbors()

template<typename TSparseMatrix >
void ug::MarkNeighbors ( const TSparseMatrix &  A,
size_t  node,
size_t  depth,
std::vector< bool > &  bVisited 
)

◆ MatAdd()

template<typename matrix_type >
void ug::MatAdd ( matrix_type M,
number  alpha1,
const matrix_type A,
number  alpha2,
const matrix_type B 
)

Calculates M = A + B.

Parameters
M(out) Matrix M, M = A + B
A(in) Matrix A
B(in) Matrix B note: A and/or B may be equal to M.

References PROFILE_FUNC_GROUP, and UG_ASSERT.

◆ MatAddNonDirichlet()

template<typename matrix_type >
void ug::MatAddNonDirichlet ( matrix_type M,
number  alpha1,
const matrix_type A,
number  alpha2,
const matrix_type B 
)

Calculates M = A + B.

Parameters
M(out) Matrix M, M = A + B
A(in) Matrix A
B(in) Matrix B note: A and/or B may be equal to M.

References alpha, ug::BlockRef(), ug::GetRows(), ug::IsDirichletRow(), PROFILE_FUNC_GROUP, and UG_ASSERT.

◆ MatMultTransposedAdd()

template<typename vector_t , typename matrix_t >
void ug::MatMultTransposedAdd ( vector_t &  dest,
const number alpha1,
const vector_t &  v1,
const number beta1,
const MapSparseMatrix< matrix_t > &  A1,
const vector_t &  w1 
)
inline

calculates dest = alpha1*v1 + beta1 * A1^T *w1;

References ug::MapSparseMatrix< TValueType >::axpy_transposed().

◆ operator<<()

std::ostream & ug::operator<< ( std::ostream &  out,
const AlgebraType v 
)
inline

◆ ScaleRow()

template<typename TSparseMatrix >
void ug::ScaleRow ( TSparseMatrix &  A,
size_t  i,
number  fac 
)

scales (block-)row i.

Parameters
A(in) Matrix A
i(in) row to scales
fac(in) Scaling factor

◆ ScaleSparseMatrixCommon()

template<typename TSparseMatrix >
void ug::ScaleSparseMatrixCommon ( TSparseMatrix &  A,
double  d 
)

◆ SerializeMatrix()

template<typename TSparseMatrix , class TOStream >
void ug::SerializeMatrix ( TOStream &  buf,
const TSparseMatrix &  A 
)

◆ SetCol()

template<typename TSparseMatrix >
void ug::SetCol ( TSparseMatrix &  A,
size_t  i,
size_t  alpha,
number  val = 0.0 
)

set value for col for entry (i,alpha).

Parameters
A(in) Matrix A
i(in) col to set
alphathe alpha index
valthe value to be set

References alpha, ug::BlockRef(), and ug::GetRows().

Referenced by ug::IObstacleConstraint< TDomain, TAlgebra >::adjust_restriction().

◆ SetDirichletRow() [1/3]

template<typename TSparseMatrix >
void ug::SetDirichletRow ( TSparseMatrix &  A,
const std::vector< size_t >  vIndex 
)

set Dirichlet row for block i.

Parameters
[in/out]Matrix A
[in]vIndexvector of row indices to set dirichlet, that is A(i,i) = 1.0, A(i,k) = 0.0 for all k != i.

References UG_ASSERT.

◆ SetDirichletRow() [2/3]

template<typename TSparseMatrix >
void ug::SetDirichletRow ( TSparseMatrix &  A,
size_t  i 
)

set Dirichlet row for block i.

Parameters
A(in) Matrix A
i(in) row to set dirichlet, that is A(i,i) = 1.0, A(i,k) = 0.0 for all k != i.

◆ SetDirichletRow() [3/3]

◆ SetRow() [1/2]

template<typename TSparseMatrix >
void ug::SetRow ( TSparseMatrix &  A,
size_t  i,
number  val = 0.0 
)

set value for (block-)row i.

Parameters
A(in) Matrix A
i(in) row to scales
val(in) value to be set

◆ SetRow() [2/2]

◆ sgs_step() [1/2]

template<typename Matrix_type , typename Vector_type >
void ug::sgs_step ( const Matrix_type &  A,
Vector_type &  c,
const Vector_type &  d,
const number  relaxFactor 
)

Performs a symmetric gauss-seidel step. Using sgs_step within a preconditioner-scheme leads to the fact that we get the correction by successive inserting the already computed values of c in c = N * d, with c being the correction and d being the defect. N denotes the matrix of the second normal-form of a linear iteration scheme.

Parameters
AMatrix \(A = D - L - R\)
cwill be \(c = N * d = (D-U)^{-1} D (D-L)^{-1} d \)
dthe vector d.
See also
gs_step_LL, gs_step_LL

References ug::gs_step_LL(), ug::gs_step_UR(), ug::MatMult(), and s.

Referenced by ug::SymmetricGaussSeidel< TAlgebra >::step().

◆ sgs_step() [2/2]

template<typename Matrix_type , typename Vector_type >
bool ug::sgs_step ( const Matrix_type &  A,
Vector_type &  x,
const Vector_type &  b 
)

Performs a symmetric gauss-seidel step.

Parameters
AMatrix \(A = D - L - R\)
xwill be \(x = (D-U)^{-1} D (D-L)^{-1} b \)
bthe vector b.
See also
gs_step_LL, gs_step_LL

References ug::gs_step_LL(), ug::gs_step_UR(), and ug::MatMult().