MFEM v4.7.0
Finite element discretization library
|
Class representing the DPG weak formulation for complex valued systems (see the class DPGWeakForm). More...
#include <complexweakform.hpp>
Public Member Functions | |
ComplexDPGWeakForm () | |
ComplexDPGWeakForm (Array< FiniteElementSpace * > &fes_, Array< FiniteElementCollection * > &fecol_) | |
Creates bilinear form associated with FE spaces fes_. | |
void | SetTestFECollVdim (int test_fec, int vdim) |
void | SetSpaces (Array< FiniteElementSpace * > &fes_, Array< FiniteElementCollection * > &fecol_) |
int | Size () const |
void | AllocateMatrix () |
void | Finalize (int skip_zeros=1) |
Finalizes the matrix initialization. | |
BlockMatrix & | BlockMat_r () |
Returns a reference to the BlockMatrix: \( M_r \). | |
BlockMatrix & | BlockMat_i () |
Returns a reference to the BlockMatrix: \( M_i \). | |
BlockMatrix & | BlockMatElim_r () |
Returns a reference to the BlockMatrix of eliminated b.c.: \( M_e_r \). | |
BlockMatrix & | BlockMatElim_i () |
Returns a reference to the BlockMatrix of eliminated b.c.: \( M_e_i \). | |
void | AddTrialIntegrator (BilinearFormIntegrator *bfi_r, BilinearFormIntegrator *bfi_i, int n, int m) |
Adds new Domain BF Integrator. Assumes ownership of bfi. | |
void | AddTestIntegrator (BilinearFormIntegrator *bfi_r, BilinearFormIntegrator *bfi_i, int n, int m) |
Adds new Test Integrator. Assumes ownership of bfi_r and bfi_i. | |
void | AddDomainLFIntegrator (LinearFormIntegrator *lfi_r, LinearFormIntegrator *lfi_i, int n) |
Adds new Domain LF Integrator. Assumes ownership of lfi_r and lfi_i. | |
void | Assemble (int skip_zeros=1) |
Assembles the form i.e. sums over all integrators. | |
virtual void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0) |
template<typename OpType > | |
void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, OpType &A, Vector &X, Vector &B, int copy_interior=0) |
virtual void | FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A) |
template<typename OpType > | |
void | FormSystemMatrix (const Array< int > &ess_tdof_list, OpType &A) |
void | EliminateVDofs (const Array< int > &vdofs, Operator::DiagonalPolicy dpolicy=Operator::DIAG_ONE) |
void | EliminateVDofsInRHS (const Array< int > &vdofs, const Vector &x_r, const Vector &x_i, Vector &b_r, Vector &b_i) |
virtual void | RecoverFEMSolution (const Vector &X, Vector &x) |
void | SetDiagonalPolicy (Operator::DiagonalPolicy policy) |
Sets diagonal policy used upon construction of the linear system. | |
virtual void | Update () |
void | StoreMatrices (bool store_matrices_=true) |
void | EnableStaticCondensation () |
Vector & | ComputeResidual (const Vector &x) |
virtual | ~ComplexDPGWeakForm () |
Destroys bilinear form. | |
Protected Member Functions | |
void | Init () |
void | ReleaseInitMemory () |
void | AllocMat () |
void | ConformingAssemble () |
void | ComputeOffsets () |
virtual void | BuildProlongation () |
Protected Attributes | |
ComplexBlockStaticCondensation * | static_cond |
Owned. | |
bool | initialized = false |
Mesh * | mesh = nullptr |
int | height |
int | width |
int | nblocks |
Array< int > | dof_offsets |
Array< int > | tdof_offsets |
BlockMatrix * | mat_r = nullptr |
Block matrix \( M \) to be associated with the real/imag Block bilinear form. Owned. | |
BlockMatrix * | mat_i = nullptr |
ComplexOperator * | mat = nullptr |
BlockVector * | y_r = nullptr |
BlockVectors to be associated with the real/imag Block linear form. | |
BlockVector * | y_i = nullptr |
Vector * | y = nullptr |
BlockMatrix * | mat_e_r = nullptr |
Block Matrix \( M_e \) used to store the eliminations from the b.c. Owned. \( M + M_e = M_{original} \). | |
BlockMatrix * | mat_e_i = nullptr |
Array< FiniteElementSpace * > | trial_fes |
Trial FE spaces. | |
Array< int > | IsTraceFes |
Flags to determine if a FiniteElementSpace is Trace. | |
Array< FiniteElementCollection * > | test_fecols |
Test FE Collections (Broken) | |
Array< int > | test_fecols_vdims |
Array2D< Array< BilinearFormIntegrator * > * > | trial_integs_r |
Set of Trial Integrators to be applied for matrix B. | |
Array2D< Array< BilinearFormIntegrator * > * > | trial_integs_i |
Array2D< Array< BilinearFormIntegrator * > * > | test_integs_r |
Set of Test Space (broken) Integrators to be applied for matrix G. | |
Array2D< Array< BilinearFormIntegrator * > * > | test_integs_i |
Array< Array< LinearFormIntegrator * > * > | lfis_r |
Set of LinearForm Integrators to be applied. | |
Array< Array< LinearFormIntegrator * > * > | lfis_i |
BlockMatrix * | P = nullptr |
Block Prolongation. | |
BlockMatrix * | R = nullptr |
Block Restriction. | |
mfem::Operator::DiagonalPolicy | diag_policy |
bool | store_matrices = false |
Array< ComplexDenseMatrix * > | Bmat |
Array< Vector * > | fvec |
Vector | residuals |
Class representing the DPG weak formulation for complex valued systems (see the class DPGWeakForm).
Definition at line 23 of file complexweakform.hpp.
|
inline |
Definition at line 107 of file complexweakform.hpp.
|
inline |
Creates bilinear form associated with FE spaces fes_.
Definition at line 114 of file complexweakform.hpp.
|
virtual |
Destroys bilinear form.
Definition at line 1024 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::AddDomainLFIntegrator | ( | LinearFormIntegrator * | lfi_r, |
LinearFormIntegrator * | lfi_i, | ||
int | n ) |
Adds new Domain LF Integrator. Assumes ownership of lfi_r and lfi_i.
Adds new Domain LF Integrator. Assumes ownership of bfi.
Definition at line 162 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::AddTestIntegrator | ( | BilinearFormIntegrator * | bfi_r, |
BilinearFormIntegrator * | bfi_i, | ||
int | n, | ||
int | m ) |
Adds new Test Integrator. Assumes ownership of bfi_r and bfi_i.
Adds new Domain BF Integrator. Assumes ownership of bfi.
Definition at line 144 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::AddTrialIntegrator | ( | BilinearFormIntegrator * | bfi_r, |
BilinearFormIntegrator * | bfi_i, | ||
int | n, | ||
int | m ) |
Adds new Domain BF Integrator. Assumes ownership of bfi.
Adds new Trial Integrator. Assumes ownership of bfi_r and bfi_i. n and m correspond to the trial FESpace and test FEColl respectively
Definition at line 124 of file complexweakform.cpp.
|
inline |
Definition at line 152 of file complexweakform.hpp.
|
protected |
Definition at line 83 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::Assemble | ( | int | skip_zeros = 1 | ) |
Assembles the form i.e. sums over all integrators.
Assembles the form i.e. sums over all domain integrators.
Definition at line 318 of file complexweakform.cpp.
|
inline |
Returns a reference to the BlockMatrix: \( M_i \).
Definition at line 164 of file complexweakform.hpp.
|
inline |
Returns a reference to the BlockMatrix: \( M_r \).
Definition at line 158 of file complexweakform.hpp.
|
inline |
Returns a reference to the BlockMatrix of eliminated b.c.: \( M_e_i \).
Definition at line 178 of file complexweakform.hpp.
|
inline |
Returns a reference to the BlockMatrix of eliminated b.c.: \( M_e_r \).
Definition at line 171 of file complexweakform.hpp.
|
protectedvirtual |
Reimplemented in mfem::ParComplexDPGWeakForm.
Definition at line 178 of file complexweakform.cpp.
|
protected |
Definition at line 67 of file complexweakform.cpp.
Definition at line 924 of file complexweakform.cpp.
|
protected |
Definition at line 196 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::EliminateVDofs | ( | const Array< int > & | vdofs, |
Operator::DiagonalPolicy | dpolicy = Operator::DIAG_ONE ) |
Definition at line 740 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::EliminateVDofsInRHS | ( | const Array< int > & | vdofs, |
const Vector & | x_r, | ||
const Vector & | x_i, | ||
Vector & | b_r, | ||
Vector & | b_i ) |
Definition at line 727 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::EnableStaticCondensation | ( | ) |
Definition at line 918 of file complexweakform.cpp.
void mfem::ComplexDPGWeakForm::Finalize | ( | int | skip_zeros = 1 | ) |
Finalizes the matrix initialization.
Definition at line 108 of file complexweakform.cpp.
|
virtual |
Reimplemented in mfem::ParComplexDPGWeakForm.
Definition at line 613 of file complexweakform.cpp.
|
inline |
Definition at line 210 of file complexweakform.hpp.
|
virtual |
Reimplemented in mfem::ParComplexDPGWeakForm.
Definition at line 690 of file complexweakform.cpp.
|
inline |
Definition at line 226 of file complexweakform.hpp.
|
protected |
Definition at line 17 of file complexweakform.cpp.
Reimplemented in mfem::ParComplexDPGWeakForm.
Definition at line 768 of file complexweakform.cpp.
|
protected |
Definition at line 807 of file complexweakform.cpp.
|
inline |
Sets diagonal policy used upon construction of the linear system.
Policies include:
Definition at line 250 of file complexweakform.hpp.
|
inline |
Definition at line 125 of file complexweakform.hpp.
|
inline |
Definition at line 120 of file complexweakform.hpp.
|
inline |
Definition at line 149 of file complexweakform.hpp.
|
inline |
Definition at line 257 of file complexweakform.hpp.
|
virtual |
Reimplemented in mfem::ParComplexDPGWeakForm.
Definition at line 867 of file complexweakform.cpp.
|
protected |
Store the matrix L^-1 B and Vector L^-1 l where G = L L^t
Definition at line 99 of file complexweakform.hpp.
|
protected |
Definition at line 81 of file complexweakform.hpp.
|
protected |
Definition at line 35 of file complexweakform.hpp.
Definition at line 100 of file complexweakform.hpp.
|
protected |
Definition at line 33 of file complexweakform.hpp.
|
protected |
Definition at line 30 of file complexweakform.hpp.
|
protected |
Flags to determine if a FiniteElementSpace is Trace.
Definition at line 58 of file complexweakform.hpp.
|
protected |
Definition at line 74 of file complexweakform.hpp.
|
protected |
Set of LinearForm Integrators to be applied.
Definition at line 73 of file complexweakform.hpp.
|
protected |
Definition at line 41 of file complexweakform.hpp.
|
protected |
Definition at line 52 of file complexweakform.hpp.
|
protected |
Block Matrix \( M_e \) used to store the eliminations from the b.c. Owned. \( M + M_e = M_{original} \).
Definition at line 51 of file complexweakform.hpp.
|
protected |
Definition at line 40 of file complexweakform.hpp.
|
protected |
Block matrix \( M \) to be associated with the real/imag Block bilinear form. Owned.
Definition at line 39 of file complexweakform.hpp.
|
protected |
Definition at line 32 of file complexweakform.hpp.
|
protected |
Definition at line 34 of file complexweakform.hpp.
|
protected |
Block Prolongation.
Definition at line 77 of file complexweakform.hpp.
|
protected |
Block Restriction.
Definition at line 79 of file complexweakform.hpp.
|
protected |
Definition at line 101 of file complexweakform.hpp.
|
protected |
Owned.
Definition at line 28 of file complexweakform.hpp.
|
protected |
Definition at line 95 of file complexweakform.hpp.
|
protected |
Definition at line 36 of file complexweakform.hpp.
|
protected |
Test FE Collections (Broken)
Definition at line 61 of file complexweakform.hpp.
|
protected |
Definition at line 62 of file complexweakform.hpp.
|
protected |
Definition at line 70 of file complexweakform.hpp.
|
protected |
Set of Test Space (broken) Integrators to be applied for matrix G.
Definition at line 69 of file complexweakform.hpp.
|
protected |
Trial FE spaces.
Definition at line 55 of file complexweakform.hpp.
|
protected |
Definition at line 66 of file complexweakform.hpp.
|
protected |
Set of Trial Integrators to be applied for matrix B.
Definition at line 65 of file complexweakform.hpp.
|
protected |
Definition at line 33 of file complexweakform.hpp.
|
protected |
Definition at line 46 of file complexweakform.hpp.
|
protected |
Definition at line 45 of file complexweakform.hpp.
|
protected |
BlockVectors to be associated with the real/imag Block linear form.
Definition at line 44 of file complexweakform.hpp.