![]() |
MFEM v4.9.0
Finite element discretization library
|
#include <complex_fem.hpp>
Public Member Functions | |
| ParSesquilinearForm (ParFiniteElementSpace *pf, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN) | |
| ParSesquilinearForm (ParFiniteElementSpace *pf, ParBilinearForm *pbfr, ParBilinearForm *pbfi, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN) | |
| Create a ParSesquilinearForm on the ParFiniteElementSpace pf, using the same integrators as the ParBilinearForms pbfr and pbfi . | |
| ComplexOperator::Convention | GetConvention () const |
| void | SetConvention (const ComplexOperator::Convention &convention) |
| void | SetAssemblyLevel (AssemblyLevel assembly_level) |
| Set the desired assembly level. | |
| ParBilinearForm & | real () |
| ParBilinearForm & | imag () |
| const ParBilinearForm & | real () const |
| const ParBilinearForm & | imag () const |
| void | AddDomainIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag) |
| Adds new Domain Integrator. | |
| void | AddDomainIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag, Array< int > &elem_marker) |
| Adds new Domain Integrator, restricted to specific attributes. | |
| void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag) |
| Adds new Boundary Integrator. | |
| void | AddBoundaryIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag, Array< int > &bdr_marker) |
| Adds new boundary Integrator, restricted to specific boundary attributes. | |
| void | AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag) |
| Adds new interior Face Integrator. Assumes ownership of bfi. | |
| void | AddBdrFaceIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag) |
| Adds new boundary Face Integrator. Assumes ownership of bfi. | |
| void | AddBdrFaceIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag, Array< int > &bdr_marker) |
| Adds new boundary Face Integrator, restricted to specific boundary attributes. | |
| void | Assemble (int skip_zeros=1) |
| Assemble the local matrix. | |
| void | Finalize (int skip_zeros=1) |
| Finalizes the matrix initialization. | |
| ComplexHypreParMatrix * | ParallelAssemble () |
| Returns the matrix assembled on the true dofs, i.e. P^t A P. | |
| ParFiniteElementSpace * | ParFESpace () const |
| Return the parallel FE space associated with the ParBilinearForm. | |
| void | FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, 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, const Vector &b, Vector &x) |
| virtual void | Update (FiniteElementSpace *nfes=NULL) |
| virtual | ~ParSesquilinearForm () |
Class for a parallel sesquilinear form
A sesquilinear form is a generalization of a bilinear form to complex-valued fields. Sesquilinear forms are linear in the second argument but the first argument involves a complex conjugate in the sense that:
a(alpha u, beta v) = conj(alpha) beta a(u, v)
The convention argument in the class's constructor is documented in the mfem::ComplexOperator class found in linalg/complex_operator.hpp.
When supplying integrators to the ParSesquilinearForm either the real or imaginary integrator can be NULL. This indicates that the corresponding portion of the complex-valued material coefficient is equal to zero.
Definition at line 708 of file complex_fem.hpp.
| mfem::ParSesquilinearForm::ParSesquilinearForm | ( | ParFiniteElementSpace * | pf, |
| ComplexOperator::Convention | convention = ComplexOperator::HERMITIAN ) |
Definition at line 1312 of file complex_fem.cpp.
| mfem::ParSesquilinearForm::ParSesquilinearForm | ( | ParFiniteElementSpace * | pf, |
| ParBilinearForm * | pbfr, | ||
| ParBilinearForm * | pbfi, | ||
| ComplexOperator::Convention | convention = ComplexOperator::HERMITIAN ) |
Create a ParSesquilinearForm on the ParFiniteElementSpace pf, using the same integrators as the ParBilinearForms pbfr and pbfi .
The pointer pf is not owned by the newly constructed object.
The integrators are copied as pointers and they are not owned by the newly constructed ParSesquilinearForm.
Definition at line 1320 of file complex_fem.cpp.
|
virtual |
Definition at line 1329 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddBdrFaceIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag ) |
Adds new boundary Face Integrator. Assumes ownership of bfi.
Definition at line 1376 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddBdrFaceIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag, | ||
| Array< int > & | bdr_marker ) |
Adds new boundary Face Integrator, restricted to specific boundary attributes.
Assumes ownership of bfi.
The array bdr_marker is stored internally as a pointer to the given Array<int> object.
Definition at line 1384 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddBoundaryIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag ) |
Adds new Boundary Integrator.
Definition at line 1351 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddBoundaryIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag, | ||
| Array< int > & | bdr_marker ) |
Adds new boundary Integrator, restricted to specific boundary attributes.
Assumes ownership of bfi.
The array bdr_marker is stored internally as a pointer to the given Array<int> object.
Definition at line 1359 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddDomainIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag ) |
Adds new Domain Integrator.
Definition at line 1335 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddDomainIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag, | ||
| Array< int > & | elem_marker ) |
Adds new Domain Integrator, restricted to specific attributes.
Definition at line 1342 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::AddInteriorFaceIntegrator | ( | BilinearFormIntegrator * | bfi_real, |
| BilinearFormIntegrator * | bfi_imag ) |
Adds new interior Face Integrator. Assumes ownership of bfi.
Definition at line 1368 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::Assemble | ( | int | skip_zeros = 1 | ) |
Assemble the local matrix.
Definition at line 1393 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::Finalize | ( | int | skip_zeros = 1 | ) |
Finalizes the matrix initialization.
Definition at line 1400 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::FormLinearSystem | ( | const Array< int > & | ess_tdof_list, |
| Vector & | x, | ||
| Vector & | b, | ||
| OperatorHandle & | A, | ||
| Vector & | X, | ||
| Vector & | B, | ||
| int | copy_interior = 0 ) |
Definition at line 1415 of file complex_fem.cpp.
| void mfem::ParSesquilinearForm::FormSystemMatrix | ( | const Array< int > & | ess_tdof_list, |
| OperatorHandle & | A ) |
Definition at line 1581 of file complex_fem.cpp.
|
inline |
Definition at line 738 of file complex_fem.hpp.
|
inline |
Definition at line 759 of file complex_fem.hpp.
|
inline |
Definition at line 761 of file complex_fem.hpp.
| ComplexHypreParMatrix * mfem::ParSesquilinearForm::ParallelAssemble | ( | ) |
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 1407 of file complex_fem.cpp.
|
inline |
Return the parallel FE space associated with the ParBilinearForm.
Definition at line 817 of file complex_fem.hpp.
|
inline |
Definition at line 758 of file complex_fem.hpp.
|
inline |
Definition at line 760 of file complex_fem.hpp.
|
virtual |
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.
Definition at line 1649 of file complex_fem.cpp.
|
inline |
Set the desired assembly level.
Valid choices are:
This method must be called before assembly.
Definition at line 752 of file complex_fem.hpp.
|
inline |
Definition at line 739 of file complex_fem.hpp.
|
virtual |
Definition at line 1676 of file complex_fem.cpp.