MFEM  v3.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::Hybridization Class Reference

#include <hybridization.hpp>

Collaboration diagram for mfem::Hybridization:
[legend]

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 Finalize ()
 Finalize the construction of the hybridized matrix. More...
 
SparseMatrixGetMatrix ()
 Return the serial hybridized matrix. More...
 
HypreParMatrixGetParallelMatrix ()
 Return the parallel hybridized matrix. More...
 
void ReduceRHS (const Vector &b, Vector &b_r) const
 
void ComputeSolution (const Vector &b, const Vector &sol_r, Vector &sol) const
 

Protected Member Functions

void GetIBDofs (int el, Array< int > &i_dofs, Array< int > &b_dofs) const
 
void MultAfInv (const Vector &b, const Vector &lambda, Vector &bf, int mode) const
 

Protected Attributes

FiniteElementSpacefes
 
FiniteElementSpacec_fes
 
BilinearFormIntegratorc_bfi
 
SparseMatrixCt
 
SparseMatrixH
 
Array< int > hat_offsets
 
Array< int > hat_dofs_marker
 
Array< int > Af_offsets
 
Array< int > Af_f_offsets
 
double * Af_data
 
int * Af_ipiv
 
HypreParMatrixpH
 

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:

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

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

Constructor & Destructor Documentation

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

Constructor.

Definition at line 30 of file hybridization.cpp.

mfem::Hybridization::~Hybridization ( )

Destructor.

Definition at line 40 of file hybridization.cpp.

Member Function Documentation

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

Assemble the element matrix A into the hybridized system matrix.

Definition at line 344 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 567 of file hybridization.cpp.

void mfem::Hybridization::Finalize ( )

Finalize the construction of the hybridized matrix.

Definition at line 434 of file hybridization.cpp.

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

Definition at line 325 of file hybridization.cpp.

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

Return the serial hybridized matrix.

Definition at line 115 of file hybridization.hpp.

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

Return the parallel hybridized matrix.

Definition at line 119 of file hybridization.hpp.

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

Prepare the Hybridization object for assembly.

Definition at line 52 of file hybridization.cpp.

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

Definition at line 456 of file hybridization.cpp.

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

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 538 of file hybridization.cpp.

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

Member Data Documentation

double* mfem::Hybridization::Af_data
protected

Definition at line 74 of file hybridization.hpp.

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

Definition at line 73 of file hybridization.hpp.

int* mfem::Hybridization::Af_ipiv
protected

Definition at line 75 of file hybridization.hpp.

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

Definition at line 73 of file hybridization.hpp.

BilinearFormIntegrator* mfem::Hybridization::c_bfi
protected

Definition at line 68 of file hybridization.hpp.

FiniteElementSpace * mfem::Hybridization::c_fes
protected

Definition at line 67 of file hybridization.hpp.

SparseMatrix* mfem::Hybridization::Ct
protected

Definition at line 70 of file hybridization.hpp.

FiniteElementSpace* mfem::Hybridization::fes
protected

Definition at line 67 of file hybridization.hpp.

SparseMatrix * mfem::Hybridization::H
protected

Definition at line 70 of file hybridization.hpp.

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

Definition at line 72 of file hybridization.hpp.

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

Definition at line 72 of file hybridization.hpp.

HypreParMatrix* mfem::Hybridization::pH
protected

Definition at line 78 of file hybridization.hpp.


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