MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::Hybridization Class Reference

Auxiliary class Hybridization, used to implement BilinearForm hybridization. More...

#include <hybridization.hpp>

Collaboration diagram for mfem::Hybridization:
[legend]

Public Member Functions

 Hybridization (FiniteElementSpace *fespace, FiniteElementSpace *c_fespace)
 Constructor.
 
 ~Hybridization ()
 Destructor.
 
void EnableDeviceExecution ()
 Turns on device execution.
 
void SetConstraintIntegrator (BilinearFormIntegrator *c_integ)
 Set the integrator that will be used to construct the constraint matrix C.
 
void AddBdrConstraintIntegrator (BilinearFormIntegrator *c_integ)
 
void AddBdrConstraintIntegrator (BilinearFormIntegrator *c_integ, Array< int > &bdr_marker)
 
BilinearFormIntegratorGetBdrConstraintIntegrator (int i)
 Access all integrators added with AddBdrConstraintIntegrator().
 
Array< int > * GetBdrConstraintIntegratorMarker (int i)
 Access all boundary markers added with AddBdrConstraintIntegrator().
 
void UseExternalBdrConstraintIntegrators ()
 Indicate that boundary constraint integrators are not owned.
 
void Init (const Array< int > &ess_tdof_list)
 Prepare the Hybridization object for assembly.
 
void AssembleMatrix (int el, const DenseMatrix &A)
 Assemble the element matrix A into the hybridized system matrix.
 
void AssembleElementMatrices (const class DenseTensor &el_mats)
 Assemble all of the element matrices given in the form of a DenseTensor.
 
void AssembleBdrMatrix (int bdr_el, const DenseMatrix &A)
 Assemble the boundary element matrix A into the hybridized system matrix.
 
void Finalize ()
 Finalize the construction of the hybridized matrix.
 
SparseMatrixGetMatrix ()
 Return the serial hybridized matrix.
 
HypreParMatrixGetParallelMatrix ()
 Return the parallel hybridized matrix.
 
void GetParallelMatrix (OperatorHandle &H_h) const
 Return the parallel hybridized matrix in the format specified by SetOperatorType().
 
void SetOperatorType (Operator::Type tid)
 Set the operator type id for the parallel hybridized matrix/operator.
 
void ReduceRHS (const Vector &b, Vector &b_r) const
 Perform the reduction of the given right-hand side b to a right-hand side vector b_r for the hybridized system.
 
void 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 right-hand side b.
 
void Reset ()
 Destroy the current hybridization matrix while preserving the computed constraint matrix and the set of essential true dofs.
 

Protected Member Functions

void ConstructC ()
 Construct the constraint matrix.
 
void GetIBDofs (int el, Array< int > &i_dofs, Array< int > &b_dofs) const
 Returns the local indices of the i-dofs and b-dofs of element el.
 
void GetBDofs (int el, int &num_idofs, Array< int > &b_dofs) const
 Returns global indices of the b-dofs of element el.
 
void ComputeH ()
 Construct the Schur complement system.
 
void MultAfInv (const Vector &b, const Vector &lambda, Vector &bf, int mode) const
 

Protected Attributes

FiniteElementSpacefes
 The finite element space.
 
FiniteElementSpacec_fes
 
std::unique_ptr< class HybridizationExtensionext
 Extension for device execution.
 
std::unique_ptr< BilinearFormIntegratorc_bfi
 The constraint integrator.
 
std::vector< BilinearFormIntegrator * > boundary_constraint_integs
 The constraint boundary face integrators.
 
std::vector< Array< int > * > boundary_constraint_integs_marker
 Boundary markers for constraint face integrators.
 
bool extern_bdr_constr_integs {false}
 Indicates if the boundary_constraint_integs integrators are owned externally.
 
std::unique_ptr< SparseMatrixCt
 The constraint matrix.
 
std::unique_ptr< SparseMatrixH
 The Schur complement system for the Lagrange multiplier.
 
Array< int > hat_offsets
 
Array< int > hat_dofs_marker
 
Array< int > Af_offsets
 
Array< int > Af_f_offsets
 
Array< real_tAf_data
 
Array< int > Af_ipiv
 
std::unique_ptr< HypreParMatrixpC
 
std::unique_ptr< HypreParMatrixP_pc
 
OperatorHandle pH
 

Friends

class HybridizationExtension
 

Detailed Description

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:

  • \( \hat{A} \) is invertible,
  • \( P \) has a left inverse \( R \), such that \( R P = I \),
  • a constraint matrix \( C \) can be constructed, such that \( \operatorname{Ker}(C) = \operatorname{Im}(P) \).

Under these conditions, the linear system \( A x = b \) can be solved using the following procedure:

  • solve for \( \lambda \) in the linear system:

    \[ (C \hat{A}^{-1} C^T) \lambda = C \hat{A}^{-1} R^T b \]

  • compute \( x = R \hat{A}^{-1} (R^T b - C^T \lambda) \)

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 62 of file hybridization.hpp.

Constructor & Destructor Documentation

◆ Hybridization()

mfem::Hybridization::Hybridization ( FiniteElementSpace * fespace,
FiniteElementSpace * c_fespace )

Constructor.

Definition at line 31 of file hybridization.cpp.

◆ ~Hybridization()

mfem::Hybridization::~Hybridization ( )

Destructor.

Definition at line 40 of file hybridization.cpp.

Member Function Documentation

◆ AddBdrConstraintIntegrator() [1/2]

void mfem::Hybridization::AddBdrConstraintIntegrator ( BilinearFormIntegrator * c_integ)
inline

Add the boundary face 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 137 of file hybridization.hpp.

◆ AddBdrConstraintIntegrator() [2/2]

void mfem::Hybridization::AddBdrConstraintIntegrator ( BilinearFormIntegrator * c_integ,
Array< int > & bdr_marker )
inline

Definition at line 142 of file hybridization.hpp.

◆ AssembleBdrMatrix()

void mfem::Hybridization::AssembleBdrMatrix ( int bdr_el,
const DenseMatrix & A )

Assemble the boundary element matrix A into the hybridized system matrix.

Definition at line 581 of file hybridization.cpp.

◆ AssembleElementMatrices()

void mfem::Hybridization::AssembleElementMatrices ( const class DenseTensor & el_mats)

Assemble all of the element matrices given in the form of a DenseTensor.

Definition at line 567 of file hybridization.cpp.

◆ AssembleMatrix()

void mfem::Hybridization::AssembleMatrix ( int el,
const DenseMatrix & A )

Assemble the element matrix A into the hybridized system matrix.

Definition at line 521 of file hybridization.cpp.

◆ ComputeH()

void mfem::Hybridization::ComputeH ( )
protected

Construct the Schur complement system.

Definition at line 675 of file hybridization.cpp.

◆ ComputeSolution()

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 right-hand side b.

It is assumed that the vector sol has the correct essential boundary conditions.

Definition at line 962 of file hybridization.cpp.

◆ ConstructC()

void mfem::Hybridization::ConstructC ( )
protected

Construct the constraint matrix.

Definition at line 54 of file hybridization.cpp.

◆ EnableDeviceExecution()

void mfem::Hybridization::EnableDeviceExecution ( )

Turns on device execution.

Definition at line 49 of file hybridization.cpp.

◆ Finalize()

void mfem::Hybridization::Finalize ( )

Finalize the construction of the hybridized matrix.

Definition at line 830 of file hybridization.cpp.

◆ GetBDofs()

void mfem::Hybridization::GetBDofs ( int el,
int & num_idofs,
Array< int > & b_dofs ) const
protected

Returns global indices of the b-dofs of element el.

Definition at line 505 of file hybridization.cpp.

◆ GetBdrConstraintIntegrator()

BilinearFormIntegrator & mfem::Hybridization::GetBdrConstraintIntegrator ( int i)
inline

Access all integrators added with AddBdrConstraintIntegrator().

Definition at line 150 of file hybridization.hpp.

◆ GetBdrConstraintIntegratorMarker()

Array< int > * mfem::Hybridization::GetBdrConstraintIntegratorMarker ( int i)
inline

Access all boundary markers added with AddBdrConstraintIntegrator().

If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.

Definition at line 156 of file hybridization.hpp.

◆ GetIBDofs()

void mfem::Hybridization::GetIBDofs ( int el,
Array< int > & i_dofs,
Array< int > & b_dofs ) const
protected

Returns the local indices of the i-dofs and b-dofs of element el.

Definition at line 485 of file hybridization.cpp.

◆ GetMatrix()

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

Return the serial hybridized matrix.

Definition at line 178 of file hybridization.hpp.

◆ GetParallelMatrix() [1/2]

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

Return the parallel hybridized matrix.

Definition at line 182 of file hybridization.hpp.

◆ GetParallelMatrix() [2/2]

void mfem::Hybridization::GetParallelMatrix ( OperatorHandle & H_h) const
inline

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

Definition at line 186 of file hybridization.hpp.

◆ Init()

void mfem::Hybridization::Init ( const Array< int > & ess_tdof_list)

Prepare the Hybridization object for assembly.

Definition at line 288 of file hybridization.cpp.

◆ MultAfInv()

void mfem::Hybridization::MultAfInv ( const Vector & b,
const Vector & lambda,
Vector & bf,
int mode ) const
protected

Definition at line 839 of file hybridization.cpp.

◆ ReduceRHS()

void mfem::Hybridization::ReduceRHS ( const Vector & b,
Vector & b_r ) const

Perform the reduction of the given right-hand side b to a right-hand side vector b_r for the hybridized system.

Definition at line 922 of file hybridization.cpp.

◆ Reset()

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 unmodified. If that is not the case, a new Hybridization object must be created.

Definition at line 1007 of file hybridization.cpp.

◆ SetConstraintIntegrator()

void mfem::Hybridization::SetConstraintIntegrator ( BilinearFormIntegrator * c_integ)
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 131 of file hybridization.hpp.

◆ SetOperatorType()

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

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

Definition at line 189 of file hybridization.hpp.

◆ UseExternalBdrConstraintIntegrators()

void mfem::Hybridization::UseExternalBdrConstraintIntegrators ( )
inline

Indicate that boundary constraint integrators are not owned.

Definition at line 160 of file hybridization.hpp.

Friends And Related Symbol Documentation

◆ HybridizationExtension

friend class HybridizationExtension
friend

Definition at line 64 of file hybridization.hpp.

Member Data Documentation

◆ Af_data

Array<real_t> mfem::Hybridization::Af_data
protected

Definition at line 86 of file hybridization.hpp.

◆ Af_f_offsets

Array<int> mfem::Hybridization::Af_f_offsets
protected

Definition at line 85 of file hybridization.hpp.

◆ Af_ipiv

Array<int> mfem::Hybridization::Af_ipiv
protected

Definition at line 87 of file hybridization.hpp.

◆ Af_offsets

Array<int> mfem::Hybridization::Af_offsets
protected

Definition at line 85 of file hybridization.hpp.

◆ boundary_constraint_integs

std::vector<BilinearFormIntegrator*> mfem::Hybridization::boundary_constraint_integs
protected

The constraint boundary face integrators.

Definition at line 73 of file hybridization.hpp.

◆ boundary_constraint_integs_marker

std::vector<Array<int>*> mfem::Hybridization::boundary_constraint_integs_marker
protected

Boundary markers for constraint face integrators.

Definition at line 75 of file hybridization.hpp.

◆ c_bfi

std::unique_ptr<BilinearFormIntegrator> mfem::Hybridization::c_bfi
protected

The constraint integrator.

Definition at line 71 of file hybridization.hpp.

◆ c_fes

FiniteElementSpace& mfem::Hybridization::c_fes
protected

The constraint finite element space.

Definition at line 67 of file hybridization.hpp.

◆ Ct

std::unique_ptr<SparseMatrix> mfem::Hybridization::Ct
protected

The constraint matrix.

Definition at line 80 of file hybridization.hpp.

◆ ext

std::unique_ptr<class HybridizationExtension> mfem::Hybridization::ext
protected

Extension for device execution.

Definition at line 69 of file hybridization.hpp.

◆ extern_bdr_constr_integs

bool mfem::Hybridization::extern_bdr_constr_integs {false}
protected

Indicates if the boundary_constraint_integs integrators are owned externally.

Definition at line 77 of file hybridization.hpp.

◆ fes

FiniteElementSpace& mfem::Hybridization::fes
protected

The finite element space.

Definition at line 66 of file hybridization.hpp.

◆ H

std::unique_ptr<SparseMatrix> mfem::Hybridization::H
protected

The Schur complement system for the Lagrange multiplier.

Definition at line 82 of file hybridization.hpp.

◆ hat_dofs_marker

Array<int> mfem::Hybridization::hat_dofs_marker
protected

Definition at line 84 of file hybridization.hpp.

◆ hat_offsets

Array<int> mfem::Hybridization::hat_offsets
protected

Definition at line 84 of file hybridization.hpp.

◆ P_pc

std::unique_ptr<HypreParMatrix> mfem::Hybridization::P_pc
protected

Definition at line 90 of file hybridization.hpp.

◆ pC

std::unique_ptr<HypreParMatrix> mfem::Hybridization::pC
protected

Definition at line 90 of file hybridization.hpp.

◆ pH

OperatorHandle mfem::Hybridization::pH
protected

Definition at line 91 of file hybridization.hpp.


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