MFEM
v4.6.0
Finite element discretization library
|
#include <hybridization.hpp>
Public Member Functions | |
Hybridization (FiniteElementSpace *fespace, FiniteElementSpace *c_fespace) | |
Constructor. More... | |
~Hybridization () | |
Destructor. More... | |
void | SetConstraintIntegrator (BilinearFormIntegrator *c_integ) |
void | Init (const Array< int > &ess_tdof_list) |
Prepare the Hybridization object for assembly. More... | |
void | AssembleMatrix (int el, const DenseMatrix &A) |
Assemble the element matrix A into the hybridized system matrix. More... | |
void | AssembleBdrMatrix (int bdr_el, const DenseMatrix &A) |
Assemble the boundary element matrix A into the hybridized system matrix. More... | |
void | Finalize () |
Finalize the construction of the hybridized matrix. More... | |
SparseMatrix & | GetMatrix () |
Return the serial hybridized matrix. More... | |
HypreParMatrix & | GetParallelMatrix () |
Return the parallel hybridized matrix. More... | |
void | GetParallelMatrix (OperatorHandle &H_h) const |
Return the parallel hybridized matrix in the format specified by SetOperatorType(). More... | |
void | SetOperatorType (Operator::Type tid) |
Set the operator type id for the parallel hybridized matrix/operator. More... | |
void | ReduceRHS (const Vector &b, Vector &b_r) const |
void | ComputeSolution (const Vector &b, const Vector &sol_r, Vector &sol) const |
void | Reset () |
Destroy the current hybridization matrix while preserving the computed constraint matrix and the set of essential true dofs. After Reset(), a new hybridized matrix can be assembled via AssembleMatrix() and Finalize(). The Mesh and FiniteElementSpace objects are assumed to be un-modified. If that is not the case, a new Hybridization object must be created. More... | |
Protected Member Functions | |
void | ConstructC () |
void | GetIBDofs (int el, Array< int > &i_dofs, Array< int > &b_dofs) const |
void | GetBDofs (int el, int &num_idofs, Array< int > &b_dofs) const |
void | ComputeH () |
void | MultAfInv (const Vector &b, const Vector &lambda, Vector &bf, int mode) const |
Protected Attributes | |
FiniteElementSpace * | fes |
FiniteElementSpace * | c_fes |
BilinearFormIntegrator * | c_bfi |
SparseMatrix * | Ct |
SparseMatrix * | H |
Array< int > | hat_offsets |
Array< int > | hat_dofs_marker |
Array< int > | Af_offsets |
Array< int > | Af_f_offsets |
double * | Af_data |
int * | Af_ipiv |
HypreParMatrix * | pC |
HypreParMatrix * | P_pc |
OperatorHandle | pH |
Auxiliary class Hybridization, used to implement BilinearForm hybridization.
Hybridization can be viewed as a technique for solving linear systems obtained through finite element assembly. The assembled matrix \( A \) can be written as:
\[ A = P^T \hat{A} P, \]
where \( P \) is the matrix mapping the conforming finite element space to the purely local finite element space without any inter-element constraints imposed, and \( \hat{A} \) is the block-diagonal matrix of all element matrices.
We assume that:
Under these conditions, the linear system \( A x = b \) can be solved using the following procedure:
\[ (C \hat{A}^{-1} C^T) \lambda = C \hat{A}^{-1} R^T b \]
Hybridization is advantageous when the matrix \( H = (C \hat{A}^{-1} C^T) \) of the hybridized system is either smaller than the original system, or is simpler to invert with a known method.
In some cases, e.g. high-order elements, the matrix \( C \) can be written as
\[ C = \begin{pmatrix} 0 & C_b \end{pmatrix}, \]
and then the hybridized matrix \( H \) can be assembled using the identity
\[ H = C_b S_b^{-1} C_b^T, \]
where \( S_b \) is the Schur complement of \( \hat{A} \) with respect to the same decomposition as the columns of \( C \):
\[ S_b = \hat{A}_b - \hat{A}_{bf} \hat{A}_{f}^{-1} \hat{A}_{fb}. \]
Hybridization can also be viewed as a discretization method for imposing (weak) continuity constraints between neighboring elements.
Definition at line 60 of file hybridization.hpp.
mfem::Hybridization::Hybridization | ( | FiniteElementSpace * | fespace, |
FiniteElementSpace * | c_fespace | ||
) |
Constructor.
Definition at line 30 of file hybridization.cpp.
mfem::Hybridization::~Hybridization | ( | ) |
Destructor.
Definition at line 41 of file hybridization.cpp.
void mfem::Hybridization::AssembleBdrMatrix | ( | int | bdr_el, |
const DenseMatrix & | A | ||
) |
Assemble the boundary element matrix A into the hybridized system matrix.
Definition at line 482 of file hybridization.cpp.
void mfem::Hybridization::AssembleMatrix | ( | int | el, |
const DenseMatrix & | A | ||
) |
Assemble the element matrix A into the hybridized system matrix.
Definition at line 442 of file hybridization.cpp.
|
protected |
Definition at line 570 of file hybridization.cpp.
void mfem::Hybridization::ComputeSolution | ( | const Vector & | b, |
const Vector & | sol_r, | ||
Vector & | sol | ||
) | const |
Reconstruct the solution of the original system, sol, from solution of the hybridized system, sol_r, and the original r.h.s. vector, b. It is assumed that the vector sol has the right essential b.c.
Definition at line 841 of file hybridization.cpp.
|
protected |
Definition at line 54 of file hybridization.cpp.
void mfem::Hybridization::Finalize | ( | ) |
Finalize the construction of the hybridized matrix.
Definition at line 717 of file hybridization.cpp.
|
protected |
Definition at line 426 of file hybridization.cpp.
|
protected |
Definition at line 406 of file hybridization.cpp.
|
inline |
Return the serial hybridized matrix.
Definition at line 121 of file hybridization.hpp.
|
inline |
Return the parallel hybridized matrix.
Definition at line 125 of file hybridization.hpp.
|
inline |
Return the parallel hybridized matrix in the format specified by SetOperatorType().
Definition at line 129 of file hybridization.hpp.
void mfem::Hybridization::Init | ( | const Array< int > & | ess_tdof_list | ) |
Prepare the Hybridization object for assembly.
Definition at line 215 of file hybridization.cpp.
|
protected |
Definition at line 726 of file hybridization.cpp.
Perform the reduction of the given r.h.s. vector, b, to a r.h.s vector, b_r, for the hybridized system.
Definition at line 807 of file hybridization.cpp.
void mfem::Hybridization::Reset | ( | ) |
Destroy the current hybridization matrix while preserving the computed constraint matrix and the set of essential true dofs. After Reset(), a new hybridized matrix can be assembled via AssembleMatrix() and Finalize(). The Mesh and FiniteElementSpace objects are assumed to be un-modified. If that is not the case, a new Hybridization object must be created.
Definition at line 880 of file hybridization.cpp.
|
inline |
Set the integrator that will be used to construct the constraint matrix C. The Hybridization object assumes ownership of the integrator, i.e. it will delete the integrator when destroyed.
Definition at line 105 of file hybridization.hpp.
|
inline |
Set the operator type id for the parallel hybridized matrix/operator.
Definition at line 132 of file hybridization.hpp.
|
protected |
Definition at line 70 of file hybridization.hpp.
|
protected |
Definition at line 69 of file hybridization.hpp.
|
protected |
Definition at line 71 of file hybridization.hpp.
|
protected |
Definition at line 69 of file hybridization.hpp.
|
protected |
Definition at line 64 of file hybridization.hpp.
|
protected |
Definition at line 63 of file hybridization.hpp.
|
protected |
Definition at line 66 of file hybridization.hpp.
|
protected |
Definition at line 63 of file hybridization.hpp.
|
protected |
Definition at line 66 of file hybridization.hpp.
|
protected |
Definition at line 68 of file hybridization.hpp.
|
protected |
Definition at line 68 of file hybridization.hpp.
|
protected |
Definition at line 74 of file hybridization.hpp.
|
protected |
Definition at line 74 of file hybridization.hpp.
|
protected |
Definition at line 75 of file hybridization.hpp.