MFEM  v4.1.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | List of all members
mfem::ParSesquilinearForm Class Reference

#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 . More...
 
ComplexOperator::Convention GetConvention () const
 
void SetConvention (const ComplexOperator::Convention &convention)
 
ParBilinearFormreal ()
 
ParBilinearFormimag ()
 
const ParBilinearFormreal () const
 
const ParBilinearFormimag () const
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
 Adds new Domain Integrator. More...
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
 Adds new Boundary Integrator. More...
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag, Array< int > &bdr_marker)
 Adds new boundary Integrator, restricted to specific boundary attributes. More...
 
void AddInteriorFaceIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
 Adds new interior Face Integrator. Assumes ownership of bfi. More...
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag)
 Adds new boundary Face Integrator. Assumes ownership of bfi. More...
 
void AddBdrFaceIntegrator (BilinearFormIntegrator *bfi_real, BilinearFormIntegrator *bfi_imag, Array< int > &bdr_marker)
 Adds new boundary Face Integrator, restricted to specific boundary attributes. More...
 
void Assemble (int skip_zeros=1)
 Assemble the local matrix. More...
 
void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization. More...
 
ComplexHypreParMatrixParallelAssemble ()
 Returns the matrix assembled on the true dofs, i.e. P^t A P. More...
 
ParFiniteElementSpaceParFESpace () const
 Return the parallel FE space associated with the ParBilinearForm. More...
 
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 ()
 

Detailed Description

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 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 493 of file complex_fem.hpp.

Constructor & Destructor Documentation

mfem::ParSesquilinearForm::ParSesquilinearForm ( ParFiniteElementSpace pf,
ComplexOperator::Convention  convention = ComplexOperator::HERMITIAN 
)

Definition at line 811 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 819 of file complex_fem.cpp.

mfem::ParSesquilinearForm::~ParSesquilinearForm ( )
virtual

Definition at line 827 of file complex_fem.cpp.

Member Function Documentation

void mfem::ParSesquilinearForm::AddBdrFaceIntegrator ( BilinearFormIntegrator bfi_real,
BilinearFormIntegrator bfi_imag 
)

Adds new boundary Face Integrator. Assumes ownership of bfi.

Definition at line 866 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 874 of file complex_fem.cpp.

void mfem::ParSesquilinearForm::AddBoundaryIntegrator ( BilinearFormIntegrator bfi_real,
BilinearFormIntegrator bfi_imag 
)

Adds new Boundary Integrator.

Definition at line 841 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 849 of file complex_fem.cpp.

void mfem::ParSesquilinearForm::AddDomainIntegrator ( BilinearFormIntegrator bfi_real,
BilinearFormIntegrator bfi_imag 
)

Adds new Domain Integrator.

Definition at line 833 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 858 of file complex_fem.cpp.

void mfem::ParSesquilinearForm::Assemble ( int  skip_zeros = 1)

Assemble the local matrix.

Definition at line 883 of file complex_fem.cpp.

void mfem::ParSesquilinearForm::Finalize ( int  skip_zeros = 1)

Finalizes the matrix initialization.

Definition at line 890 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 906 of file complex_fem.cpp.

void mfem::ParSesquilinearForm::FormSystemMatrix ( const Array< int > &  ess_tdof_list,
OperatorHandle A 
)

Definition at line 1029 of file complex_fem.cpp.

ComplexOperator::Convention mfem::ParSesquilinearForm::GetConvention ( ) const
inline

Definition at line 523 of file complex_fem.hpp.

ParBilinearForm& mfem::ParSesquilinearForm::imag ( )
inline

Definition at line 528 of file complex_fem.hpp.

const ParBilinearForm& mfem::ParSesquilinearForm::imag ( ) const
inline

Definition at line 530 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 897 of file complex_fem.cpp.

ParFiniteElementSpace* mfem::ParSesquilinearForm::ParFESpace ( ) const
inline

Return the parallel FE space associated with the ParBilinearForm.

Definition at line 581 of file complex_fem.hpp.

ParBilinearForm& mfem::ParSesquilinearForm::real ( )
inline

Definition at line 527 of file complex_fem.hpp.

const ParBilinearForm& mfem::ParSesquilinearForm::real ( ) const
inline

Definition at line 529 of file complex_fem.hpp.

void mfem::ParSesquilinearForm::RecoverFEMSolution ( const Vector X,
const Vector b,
Vector x 
)
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 1093 of file complex_fem.cpp.

void mfem::ParSesquilinearForm::SetConvention ( const ComplexOperator::Convention convention)
inline

Definition at line 524 of file complex_fem.hpp.

void mfem::ParSesquilinearForm::Update ( FiniteElementSpace nfes = NULL)
virtual

Definition at line 1115 of file complex_fem.cpp.


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