MFEM
v3.1
Finite element discretization library
|
#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) |
FiniteElementSpace * | GetTraceFESpace () |
Return a pointer to the reduced/trace FE space. More... | |
ParFiniteElementSpace * | GetParTraceFESpace () |
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 | EliminateReducedTrueDofs (const Array< int > &ess_rtdof_list, int keep_diagonal) |
Eliminate the given reduced true dofs from the Schur complement matrix S. More... | |
bool | HasEliminatedBC () const |
SparseMatrix & | GetMatrix () |
Return the serial Schur complement matrix. More... | |
SparseMatrix & | GetMatrixElim () |
Return the eliminated part of the serial Schur complement matrix. More... | |
HypreParMatrix & | GetParallelMatrix () |
Return the parallel Schur complement matrix. More... | |
HypreParMatrix & | GetParallelMatrixElim () |
Return the eliminated part of the parallel Schur complement matrix. More... | |
void | ReduceRHS (const Vector &b, Vector &sc_b) const |
void | ReduceSolution (const Vector &sol, Vector &sc_sol) const |
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 |
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
where the block is itself block diagonal with small local blocks and it is, therefore, easily invertible.
Starting with the block system
the reduced, statically condensed system is given by
where the Schur complement matrix is given by
After solving the Schur complement system, the part of the solution can be recovered using the formula
Definition at line 65 of file staticcond.hpp.
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 104 of file staticcond.cpp.
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 227 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 189 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 452 of file staticcond.cpp.
|
inline |
Restrict a list of true FE space dofs to a list of reduced/trace true FE space dofs.
Definition at line 169 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 415 of file staticcond.cpp.
void mfem::StaticCondensation::EliminateReducedTrueDofs | ( | const Array< int > & | ess_rtdof_list, |
int | keep_diagonal | ||
) |
Eliminate the given reduced true dofs from the Schur complement matrix S.
Definition at line 286 of file staticcond.cpp.
void mfem::StaticCondensation::Finalize | ( | ) |
Finalize the construction of the Schur complement matrix.
Definition at line 235 of file staticcond.cpp.
|
inline |
Return the serial Schur complement matrix.
Definition at line 141 of file staticcond.hpp.
|
inline |
Return the eliminated part of the serial Schur complement matrix.
Definition at line 144 of file staticcond.hpp.
|
inline |
Return the number of vector exposed/reduced dofs.
Definition at line 97 of file staticcond.hpp.
|
inline |
Return the number of vector private dofs.
Definition at line 95 of file staticcond.hpp.
|
inline |
Return the parallel Schur complement matrix.
Definition at line 148 of file staticcond.hpp.
|
inline |
Return the eliminated part of the parallel Schur complement matrix.
Definition at line 151 of file staticcond.hpp.
|
inline |
Return a pointer to the parallel reduced/trace FE space.
Definition at line 111 of file staticcond.hpp.
|
inline |
Return a pointer to the reduced/trace FE space.
Definition at line 107 of file staticcond.hpp.
|
inline |
Return true if essential boundary conditions have been eliminated from the Schur complement matrix.
Definition at line 131 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 134 of file staticcond.cpp.
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.
Restrict a solution vector on the full FE space dofs to a vector on the reduced/trace true FE space dofs.
Definition at line 388 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 118 of file staticcond.cpp.