![]() |
MFEM v4.8.0
Finite element discretization library
|
Class representing the parallel weak formulation. (Convenient for DPG Equations) More...
#include <pcomplexweakform.hpp>
Public Member Functions | |
| ParComplexDPGWeakForm () | |
| ParComplexDPGWeakForm (Array< ParFiniteElementSpace * > &trial_pfes_, Array< FiniteElementCollection * > &fecol_) | |
| Creates bilinear form associated with FE spaces trial_pfes_. | |
| void | SetParSpaces (Array< ParFiniteElementSpace * > &trial_pfes_, Array< FiniteElementCollection * > &fecol_) |
| void | Assemble (int skip_zeros=1) |
| Assembles the form i.e. sums over all domain integrators. | |
| void | ParallelAssemble (BlockMatrix *mat_r, BlockMatrix *mat_i) |
| Returns the matrix assembled on the true dofs, i.e. P^t A P. | |
| void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, OperatorHandle &A, Vector &X, Vector &B, int copy_interior=0) |
| void | FormSystemMatrix (const Array< int > &ess_tdof_list, OperatorHandle &A) |
| virtual void | RecoverFEMSolution (const Vector &X, Vector &x) |
| virtual void | Update () |
| virtual | ~ParComplexDPGWeakForm () |
| Destroys bilinear form. | |
Public Member Functions inherited from mfem::ComplexDPGWeakForm | |
| 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. | |
| template<typename OpType > | |
| void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, OpType &A, Vector &X, Vector &B, int copy_interior=0) |
| 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) |
| void | SetDiagonalPolicy (Operator::DiagonalPolicy policy) |
| Sets diagonal policy used upon construction of the linear system. | |
| void | StoreMatrices (bool store_matrices_=true) |
| void | EnableStaticCondensation () |
| Vector & | ComputeResidual (const Vector &x) |
| virtual | ~ComplexDPGWeakForm () |
| Destroys bilinear form. | |
Protected Member Functions | |
| void | FillEssTdofLists (const Array< int > &ess_tdof_list) |
| void | BuildProlongation () |
Protected Member Functions inherited from mfem::ComplexDPGWeakForm | |
| void | Init () |
| void | ReleaseInitMemory () |
| void | AllocMat () |
| void | ConformingAssemble () |
| void | ComputeOffsets () |
Protected Attributes | |
| Array< ParFiniteElementSpace * > | trial_pfes |
| Trial FE spaces. | |
| Array< Array< int > * > | ess_tdofs |
| ess_tdof list for each space | |
| BlockOperator * | P = nullptr |
| BlockMatrix * | R = nullptr |
| ComplexOperator * | p_mat = nullptr |
| BlockOperator * | p_mat_r = nullptr |
| BlockOperator * | p_mat_i = nullptr |
| BlockOperator * | p_mat_e_r = nullptr |
| BlockOperator * | p_mat_e_i = nullptr |
Protected Attributes inherited from mfem::ComplexDPGWeakForm | |
| 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 parallel weak formulation. (Convenient for DPG Equations)
Definition at line 25 of file pcomplexweakform.hpp.
|
inline |
Definition at line 54 of file pcomplexweakform.hpp.
|
inline |
Creates bilinear form associated with FE spaces trial_pfes_.
Definition at line 57 of file pcomplexweakform.hpp.
|
virtual |
Destroys bilinear form.
Definition at line 265 of file pcomplexweakform.cpp.
| void mfem::ParComplexDPGWeakForm::Assemble | ( | int | skip_zeros = 1 | ) |
Assembles the form i.e. sums over all domain integrators.
Definition at line 34 of file pcomplexweakform.cpp.
|
protectedvirtual |
Reimplemented from mfem::ComplexDPGWeakForm.
Definition at line 105 of file pcomplexweakform.cpp.
|
protected |
split ess_tdof_list given in global tdof (for all spaces) to individual lists for each space
Definition at line 19 of file pcomplexweakform.cpp.
|
virtual |
Reimplemented from mfem::ComplexDPGWeakForm.
Definition at line 121 of file pcomplexweakform.cpp.
|
virtual |
Reimplemented from mfem::ComplexDPGWeakForm.
Definition at line 183 of file pcomplexweakform.cpp.
| void mfem::ParComplexDPGWeakForm::ParallelAssemble | ( | BlockMatrix * | mat_r, |
| BlockMatrix * | mat_i ) |
Returns the matrix assembled on the true dofs, i.e. P^t A P.
The returned matrix has to be deleted by the caller.
Definition at line 39 of file pcomplexweakform.cpp.
Call this method after solving a linear system constructed using the FormLinearSystem method to recover the solution as a ParGridFunction-size vector in x. Use the same arguments as in the FormLinearSystem call.
Reimplemented from mfem::ComplexDPGWeakForm.
Definition at line 218 of file pcomplexweakform.cpp.
|
inline |
Definition at line 64 of file pcomplexweakform.hpp.
|
virtual |
Reimplemented from mfem::ComplexDPGWeakForm.
Definition at line 241 of file pcomplexweakform.cpp.
ess_tdof list for each space
Definition at line 33 of file pcomplexweakform.hpp.
|
protected |
Definition at line 40 of file pcomplexweakform.hpp.
|
protected |
Definition at line 44 of file pcomplexweakform.hpp.
|
protected |
Definition at line 48 of file pcomplexweakform.hpp.
|
protected |
Definition at line 47 of file pcomplexweakform.hpp.
|
protected |
Definition at line 46 of file pcomplexweakform.hpp.
|
protected |
Definition at line 45 of file pcomplexweakform.hpp.
|
protected |
Definition at line 42 of file pcomplexweakform.hpp.
|
protected |
Trial FE spaces.
Definition at line 30 of file pcomplexweakform.hpp.