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::SesquilinearForm Class Reference

#include <complex_fem.hpp>

Public Member Functions

 SesquilinearForm (FiniteElementSpace *fes, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
 
 SesquilinearForm (FiniteElementSpace *fes, BilinearForm *bfr, BilinearForm *bfi, ComplexOperator::Convention convention=ComplexOperator::HERMITIAN)
 Create a SesquilinearForm on the FiniteElementSpace f, using the same integrators as the BilinearForms bfr and bfi . More...
 
ComplexOperator::Convention GetConvention () const
 
void SetConvention (const ComplexOperator::Convention &convention)
 
BilinearFormreal ()
 
BilinearFormimag ()
 
const BilinearFormreal () const
 
const BilinearFormimag () 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...
 
ComplexSparseMatrixAssembleComplexSparseMatrix ()
 Returns the matrix assembled on the true dofs, i.e. P^t A P. More...
 
FiniteElementSpaceFESpace () 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)
 
void SetDiagonalPolicy (mfem::Matrix::DiagonalPolicy dpolicy)
 Sets diagonal policy used upon construction of the linear system. More...
 
Matrix::DiagonalPolicy GetDiagonalPolicy () const
 Returns the diagonal policy of the sesquilinear form. More...
 
virtual ~SesquilinearForm ()
 

Detailed Description

Class for 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 SesquilinearForm 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 185 of file complex_fem.hpp.

Constructor & Destructor Documentation

mfem::SesquilinearForm::SesquilinearForm ( FiniteElementSpace fes,
ComplexOperator::Convention  convention = ComplexOperator::HERMITIAN 
)

Definition at line 242 of file complex_fem.cpp.

mfem::SesquilinearForm::SesquilinearForm ( FiniteElementSpace fes,
BilinearForm bfr,
BilinearForm bfi,
ComplexOperator::Convention  convention = ComplexOperator::HERMITIAN 
)

Create a SesquilinearForm on the FiniteElementSpace f, using the same integrators as the BilinearForms bfr and bfi .

The pointer fes is not owned by the newly constructed object.

The integrators are copied as pointers and they are not owned by the newly constructed SesquilinearForm.

Definition at line 249 of file complex_fem.cpp.

mfem::SesquilinearForm::~SesquilinearForm ( )
virtual

Definition at line 262 of file complex_fem.cpp.

Member Function Documentation

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

Adds new boundary Face Integrator. Assumes ownership of bfi.

Definition at line 300 of file complex_fem.cpp.

void mfem::SesquilinearForm::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 307 of file complex_fem.cpp.

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

Adds new Boundary Integrator.

Definition at line 276 of file complex_fem.cpp.

void mfem::SesquilinearForm::AddBoundaryIntegrator ( BilinearFormIntegrator bfi_real,
BilinearFormIntegrator bfi_imag,
Array< int > &  bdr_marker 
)

Adds new Boundary Integrator, restricted to specific boundary attributes.

Definition at line 284 of file complex_fem.cpp.

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

Adds new Domain Integrator.

Definition at line 268 of file complex_fem.cpp.

void mfem::SesquilinearForm::AddInteriorFaceIntegrator ( BilinearFormIntegrator bfi_real,
BilinearFormIntegrator bfi_imag 
)

Adds new interior Face Integrator. Assumes ownership of bfi.

Definition at line 293 of file complex_fem.cpp.

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

Assemble the local matrix.

Definition at line 316 of file complex_fem.cpp.

ComplexSparseMatrix * mfem::SesquilinearForm::AssembleComplexSparseMatrix ( )

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 330 of file complex_fem.cpp.

FiniteElementSpace* mfem::SesquilinearForm::FESpace ( ) const
inline

Return the parallel FE space associated with the ParBilinearForm.

Definition at line 270 of file complex_fem.hpp.

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

Finalizes the matrix initialization.

Definition at line 323 of file complex_fem.cpp.

void mfem::SesquilinearForm::FormLinearSystem ( const Array< int > &  ess_tdof_list,
Vector x,
Vector b,
OperatorHandle A,
Vector X,
Vector B,
int  copy_interior = 0 
)

Definition at line 338 of file complex_fem.cpp.

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

Definition at line 433 of file complex_fem.cpp.

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

Definition at line 218 of file complex_fem.hpp.

Matrix::DiagonalPolicy mfem::SesquilinearForm::GetDiagonalPolicy ( ) const
inline

Returns the diagonal policy of the sesquilinear form.

Definition at line 290 of file complex_fem.hpp.

BilinearForm& mfem::SesquilinearForm::imag ( )
inline

Definition at line 223 of file complex_fem.hpp.

const BilinearForm& mfem::SesquilinearForm::imag ( ) const
inline

Definition at line 225 of file complex_fem.hpp.

BilinearForm& mfem::SesquilinearForm::real ( )
inline

Definition at line 222 of file complex_fem.hpp.

const BilinearForm& mfem::SesquilinearForm::real ( ) const
inline

Definition at line 224 of file complex_fem.hpp.

void mfem::SesquilinearForm::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 465 of file complex_fem.cpp.

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

Definition at line 219 of file complex_fem.hpp.

void mfem::SesquilinearForm::SetDiagonalPolicy ( mfem::Matrix::DiagonalPolicy  dpolicy)

Sets diagonal policy used upon construction of the linear system.

Definition at line 257 of file complex_fem.cpp.

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

Definition at line 494 of file complex_fem.cpp.


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