MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
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 fes, using the same integrators as the BilinearForms bfr and bfi .
 
ComplexOperator::Convention GetConvention () const
 
void SetConvention (const ComplexOperator::Convention &convention)
 
void SetAssemblyLevel (AssemblyLevel assembly_level)
 Set the desired assembly level.
 
BilinearFormreal ()
 
BilinearFormimag ()
 
const BilinearFormreal () const
 
const BilinearFormimag () 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 the given 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.
 
ComplexSparseMatrixAssembleComplexSparseMatrix ()
 Returns the matrix assembled on the true dofs, i.e. P^t A P.
 
FiniteElementSpaceFESpace () 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)
 
void SetDiagonalPolicy (mfem::Matrix::DiagonalPolicy dpolicy)
 Sets diagonal policy used upon construction of the linear system.
 
Matrix::DiagonalPolicy GetDiagonalPolicy () const
 Returns the diagonal policy of the sesquilinear form.
 
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 289 of file complex_fem.hpp.

Constructor & Destructor Documentation

◆ SesquilinearForm() [1/2]

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

Definition at line 413 of file complex_fem.cpp.

◆ SesquilinearForm() [2/2]

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

Create a SesquilinearForm on the FiniteElementSpace fes, 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 420 of file complex_fem.cpp.

◆ ~SesquilinearForm()

mfem::SesquilinearForm::~SesquilinearForm ( )
virtual

Definition at line 433 of file complex_fem.cpp.

Member Function Documentation

◆ AddBdrFaceIntegrator() [1/2]

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

Adds new boundary Face Integrator. Assumes ownership of bfi.

Definition at line 479 of file complex_fem.cpp.

◆ AddBdrFaceIntegrator() [2/2]

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

◆ AddBoundaryIntegrator() [1/2]

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

Adds new Boundary Integrator.

Definition at line 455 of file complex_fem.cpp.

◆ AddBoundaryIntegrator() [2/2]

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

◆ AddDomainIntegrator() [1/2]

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

Adds new Domain Integrator.

Definition at line 439 of file complex_fem.cpp.

◆ AddDomainIntegrator() [2/2]

void mfem::SesquilinearForm::AddDomainIntegrator ( BilinearFormIntegrator * bfi_real,
BilinearFormIntegrator * bfi_imag,
Array< int > & elem_marker )

Adds new Domain Integrator, restricted to the given attributes.

Definition at line 446 of file complex_fem.cpp.

◆ AddInteriorFaceIntegrator()

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

Adds new interior Face Integrator. Assumes ownership of bfi.

Definition at line 472 of file complex_fem.cpp.

◆ Assemble()

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

Assemble the local matrix.

Definition at line 495 of file complex_fem.cpp.

◆ AssembleComplexSparseMatrix()

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

◆ FESpace()

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

Return the parallel FE space associated with the ParBilinearForm.

Definition at line 395 of file complex_fem.hpp.

◆ Finalize()

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

Finalizes the matrix initialization.

Definition at line 502 of file complex_fem.cpp.

◆ FormLinearSystem()

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

◆ FormSystemMatrix()

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

Definition at line 672 of file complex_fem.cpp.

◆ GetConvention()

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

Definition at line 322 of file complex_fem.hpp.

◆ GetDiagonalPolicy()

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

Returns the diagonal policy of the sesquilinear form.

Definition at line 415 of file complex_fem.hpp.

◆ imag() [1/2]

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

Definition at line 343 of file complex_fem.hpp.

◆ imag() [2/2]

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

Definition at line 345 of file complex_fem.hpp.

◆ real() [1/2]

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

Definition at line 342 of file complex_fem.hpp.

◆ real() [2/2]

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

Definition at line 344 of file complex_fem.hpp.

◆ RecoverFEMSolution()

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

◆ SetAssemblyLevel()

void mfem::SesquilinearForm::SetAssemblyLevel ( AssemblyLevel assembly_level)
inline

Set the desired assembly level.

Valid choices are:

This method must be called before assembly.

Definition at line 336 of file complex_fem.hpp.

◆ SetConvention()

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

Definition at line 323 of file complex_fem.hpp.

◆ SetDiagonalPolicy()

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

Sets diagonal policy used upon construction of the linear system.

Definition at line 428 of file complex_fem.cpp.

◆ Update()

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

Definition at line 765 of file complex_fem.cpp.


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