MFEM  v4.4.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | List of all members
mfem::StaticCondensation Class Reference

#include <staticcond.hpp>

Public Member Functions

 StaticCondensation (FiniteElementSpace *fespace)
 Construct a StaticCondensation object. More...
 
 ~StaticCondensation ()
 Destroy a StaticCondensation object. More...
 
int GetNPrDofs () const
 Return the number of vector private dofs. More...
 
int GetNExDofs () const
 Return the number of vector exposed/reduced dofs. More...
 
bool ReducesTrueVSize () const
 
void Init (bool symmetric, bool block_diagonal)
 
FiniteElementSpaceGetTraceFESpace ()
 Return a pointer to the reduced/trace FE space. More...
 
ParFiniteElementSpaceGetParTraceFESpace ()
 Return a pointer to the parallel reduced/trace FE space. More...
 
void AssembleMatrix (int el, const DenseMatrix &elmat)
 
void AssembleBdrMatrix (int el, const DenseMatrix &elmat)
 
void Finalize ()
 Finalize the construction of the Schur complement matrix. More...
 
void SetEssentialTrueDofs (const Array< int > &ess_tdof_list)
 Determine and save internally essential reduced true dofs. More...
 
void EliminateReducedTrueDofs (const Array< int > &ess_rtdof_list, Matrix::DiagonalPolicy dpolicy)
 Eliminate the given reduced true dofs from the Schur complement matrix S. More...
 
void EliminateReducedTrueDofs (Matrix::DiagonalPolicy dpolicy)
 Eliminate the internal reduced true dofs (set using SetEssentialTrueDofs()) from the Schur complement matrix S. More...
 
bool HasEliminatedBC () const
 Return true if essential boundary conditions have been eliminated from the Schur complement matrix. More...
 
SparseMatrixGetMatrix ()
 Return the serial Schur complement matrix. More...
 
SparseMatrixGetMatrixElim ()
 Return the eliminated part of the serial Schur complement matrix. More...
 
HypreParMatrixGetParallelMatrix ()
 Return the parallel Schur complement matrix. More...
 
HypreParMatrixGetParallelMatrixElim ()
 Return the eliminated part of the parallel Schur complement matrix. More...
 
void GetParallelMatrix (OperatorHandle &S_h) const
 Return the parallel Schur complement matrix in the format specified by SetOperatorType(). More...
 
void GetParallelMatrixElim (OperatorHandle &S_e_h) const
 Return the eliminated part of the parallel Schur complement matrix in the format specified by SetOperatorType(). More...
 
void SetOperatorType (Operator::Type tid)
 Set the operator type id for the parallel reduced matrix/operator. More...
 
void ReduceRHS (const Vector &b, Vector &sc_b) const
 
void ReduceSolution (const Vector &sol, Vector &sc_sol) const
 
void ReduceSystem (Vector &x, Vector &b, Vector &X, Vector &B, int copy_interior=0) const
 Set the reduced solution X and r.h.s B vectors from the full linear system solution x and r.h.s. b vectors. More...
 
void ConvertMarkerToReducedTrueDofs (const Array< int > &ess_tdof_marker, Array< int > &ess_rtdof_marker) const
 
void ConvertListToReducedTrueDofs (const Array< int > &ess_tdof_list_, Array< int > &ess_rtdof_list_) const
 
void ComputeSolution (const Vector &b, const Vector &sc_sol, Vector &sol) const
 

Detailed Description

Auxiliary class StaticCondensation, used to implement static condensation in class BilinearForm.

Static condensation is a technique for solving linear systems by eliminating groups/blocks of unknowns and reducing the original system to the remaining interfacial unknowns. The assumption is that unknowns in one group are connected (in the graph of the matrix) only to unknowns in the same group or to interfacial unknowns but not to other groups.

For finite element systems, the groups correspond to degrees of freedom (DOFs) associated with the interior of the elements. The rest of the DOFs (associated with the element boundaries) are interfacial.

In block form the matrix of the system can be written as

\[ A = \begin{pmatrix} A_{11} & A_{12} \\ A_{21} & A_{22} \end{pmatrix} \begin{array}{l} \text{-- groups: element interior/private DOFs} \\ \text{-- interface: element boundary/exposed DOFs} \end{array} \]

where the block \( A_1 \) is itself block diagonal with small local blocks and it is, therefore, easily invertible.

Starting with the block system

\[ \begin{pmatrix} A_{11} & A_{12} \\ A_{21} & A_{22} \end{pmatrix} \begin{pmatrix} X_1 \\ X_2 \end{pmatrix} = \begin{pmatrix} B_1 \\ B_2 \end{pmatrix} \]

the reduced, statically condensed system is given by

\[ S_{22} X_2 = B_2 - A_{21} A_{11}^{-1} B_1 \]

where the Schur complement matrix \( S_{22} \) is given by

\[ S_{22} = A_{22} - A_{21} A_{11}^{-1} A_{12}. \]

After solving the Schur complement system, the \( X_1 \) part of the solution can be recovered using the formula

\[ X_1 = A_{11}^{-1} ( B_1 - A_{12} X_2 ). \]

Definition at line 65 of file staticcond.hpp.

Constructor & Destructor Documentation

mfem::StaticCondensation::StaticCondensation ( FiniteElementSpace fespace)

Construct a StaticCondensation object.

Definition at line 17 of file staticcond.cpp.

mfem::StaticCondensation::~StaticCondensation ( )

Destroy a StaticCondensation object.

Definition at line 103 of file staticcond.cpp.

Member Function Documentation

void mfem::StaticCondensation::AssembleBdrMatrix ( int  el,
const DenseMatrix elmat 
)

Assemble the contribution to the Schur complement from the given boundary element matrix 'elmat'.

Definition at line 225 of file staticcond.cpp.

void mfem::StaticCondensation::AssembleMatrix ( int  el,
const DenseMatrix elmat 
)

Assemble the contribution to the Schur complement from the given element matrix 'elmat'; save the other blocks internally: A_pp_inv, A_pe, and A_ep.

Definition at line 187 of file staticcond.cpp.

void mfem::StaticCondensation::ComputeSolution ( const Vector b,
const Vector sc_sol,
Vector sol 
) const

Given a solution of the reduced system 'sc_sol' and the RHS 'b' for the full linear system, compute the solution of the full system 'sol'.

Definition at line 476 of file staticcond.cpp.

void mfem::StaticCondensation::ConvertListToReducedTrueDofs ( const Array< int > &  ess_tdof_list_,
Array< int > &  ess_rtdof_list_ 
) const
inline

Restrict a list of true FE space dofs to a list of reduced/trace true FE space dofs.

Definition at line 202 of file staticcond.hpp.

void mfem::StaticCondensation::ConvertMarkerToReducedTrueDofs ( const Array< int > &  ess_tdof_marker,
Array< int > &  ess_rtdof_marker 
) const

Restrict a marker Array on the true FE space dofs to a marker Array on the reduced/trace true FE space dofs.

Definition at line 439 of file staticcond.cpp.

void mfem::StaticCondensation::EliminateReducedTrueDofs ( const Array< int > &  ess_rtdof_list,
Matrix::DiagonalPolicy  dpolicy 
)

Eliminate the given reduced true dofs from the Schur complement matrix S.

Definition at line 286 of file staticcond.cpp.

void mfem::StaticCondensation::EliminateReducedTrueDofs ( Matrix::DiagonalPolicy  dpolicy)
inline

Eliminate the internal reduced true dofs (set using SetEssentialTrueDofs()) from the Schur complement matrix S.

Definition at line 137 of file staticcond.hpp.

void mfem::StaticCondensation::Finalize ( )

Finalize the construction of the Schur complement matrix.

Definition at line 233 of file staticcond.cpp.

SparseMatrix& mfem::StaticCondensation::GetMatrix ( )
inline

Return the serial Schur complement matrix.

Definition at line 152 of file staticcond.hpp.

SparseMatrix& mfem::StaticCondensation::GetMatrixElim ( )
inline

Return the eliminated part of the serial Schur complement matrix.

Definition at line 155 of file staticcond.hpp.

int mfem::StaticCondensation::GetNExDofs ( ) const
inline

Return the number of vector exposed/reduced dofs.

Definition at line 99 of file staticcond.hpp.

int mfem::StaticCondensation::GetNPrDofs ( ) const
inline

Return the number of vector private dofs.

Definition at line 97 of file staticcond.hpp.

HypreParMatrix& mfem::StaticCondensation::GetParallelMatrix ( )
inline

Return the parallel Schur complement matrix.

Definition at line 159 of file staticcond.hpp.

void mfem::StaticCondensation::GetParallelMatrix ( OperatorHandle S_h) const
inline

Return the parallel Schur complement matrix in the format specified by SetOperatorType().

Definition at line 167 of file staticcond.hpp.

HypreParMatrix& mfem::StaticCondensation::GetParallelMatrixElim ( )
inline

Return the eliminated part of the parallel Schur complement matrix.

Definition at line 162 of file staticcond.hpp.

void mfem::StaticCondensation::GetParallelMatrixElim ( OperatorHandle S_e_h) const
inline

Return the eliminated part of the parallel Schur complement matrix in the format specified by SetOperatorType().

Definition at line 171 of file staticcond.hpp.

ParFiniteElementSpace* mfem::StaticCondensation::GetParTraceFESpace ( )
inline

Return a pointer to the parallel reduced/trace FE space.

Definition at line 113 of file staticcond.hpp.

FiniteElementSpace* mfem::StaticCondensation::GetTraceFESpace ( )
inline

Return a pointer to the reduced/trace FE space.

Definition at line 109 of file staticcond.hpp.

bool mfem::StaticCondensation::HasEliminatedBC ( ) const
inline

Return true if essential boundary conditions have been eliminated from the Schur complement matrix.

Definition at line 142 of file staticcond.hpp.

void mfem::StaticCondensation::Init ( bool  symmetric,
bool  block_diagonal 
)

Prepare the StaticCondensation object to assembly: allocate the Schur complement matrix and the other element-wise blocks.

Definition at line 132 of file staticcond.cpp.

void mfem::StaticCondensation::ReduceRHS ( const Vector b,
Vector sc_b 
) const

Given a RHS vector for the full linear system, compute the RHS for the reduced linear system: sc_b = b_e - A_ep A_pp_inv b_p.

Definition at line 309 of file staticcond.cpp.

void mfem::StaticCondensation::ReduceSolution ( const Vector sol,
Vector sc_sol 
) const

Restrict a solution vector on the full FE space dofs to a vector on the reduced/trace true FE space dofs.

Definition at line 389 of file staticcond.cpp.

bool mfem::StaticCondensation::ReducesTrueVSize ( ) const

Return true if applying the static condensation actually reduces the (global) number of true vector dofs.

Definition at line 116 of file staticcond.cpp.

void mfem::StaticCondensation::ReduceSystem ( Vector x,
Vector b,
Vector X,
Vector B,
int  copy_interior = 0 
) const

Set the reduced solution X and r.h.s B vectors from the full linear system solution x and r.h.s. b vectors.

This method should be called after the internal reduced essential dofs have been set using SetEssentialTrueDofs() and both the Schur complement and its eliminated part have been finalized.

Definition at line 416 of file staticcond.cpp.

void mfem::StaticCondensation::SetEssentialTrueDofs ( const Array< int > &  ess_tdof_list)
inline

Determine and save internally essential reduced true dofs.

Definition at line 128 of file staticcond.hpp.

void mfem::StaticCondensation::SetOperatorType ( Operator::Type  tid)
inline

Set the operator type id for the parallel reduced matrix/operator.

Definition at line 174 of file staticcond.hpp.


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