MFEM v4.7.0
Finite element discretization library
Loading...
Searching...
No Matches
mfem::ParDiscreteLinearOperator Class Reference

#include <pbilinearform.hpp>

Inheritance diagram for mfem::ParDiscreteLinearOperator:
[legend]
Collaboration diagram for mfem::ParDiscreteLinearOperator:
[legend]

Public Member Functions

 ParDiscreteLinearOperator (ParFiniteElementSpace *dfes, ParFiniteElementSpace *rfes)
 Construct a ParDiscreteLinearOperator on the given FiniteElementSpaces dfes (domain FE space) and rfes (range FE space).
 
HypreParMatrixParallelAssemble () const
 Returns the matrix "assembled" on the true dofs.
 
void ParallelAssemble (OperatorHandle &A)
 Returns the matrix assembled on the true dofs, i.e. A = R_test A_local P_trial, in the format (type id) specified by A.
 
void GetParBlocks (Array2D< HypreParMatrix * > &blocks) const
 
virtual void FormRectangularSystemMatrix (OperatorHandle &A)
 Return in A a parallel (on truedofs) version of this operator.
 
virtual ~ParDiscreteLinearOperator ()
 
virtual void FormRectangularSystemMatrix (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OperatorHandle &A)
 Return in A that is column-constrained.
 
template<typename OpType >
void FormRectangularSystemMatrix (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OpType &A)
 Form the column-constrained linear system matrix A.
 
- Public Member Functions inherited from mfem::DiscreteLinearOperator
 DiscreteLinearOperator (FiniteElementSpace *domain_fes, FiniteElementSpace *range_fes)
 Construct a DiscreteLinearOperator on the given FiniteElementSpaces domain_fes and range_fes.
 
void AddDomainInterpolator (DiscreteInterpolator *di)
 Adds a domain interpolator. Assumes ownership of di.
 
void AddDomainInterpolator (DiscreteInterpolator *di, Array< int > &elem_marker)
 
void AddTraceFaceInterpolator (DiscreteInterpolator *di)
 Adds a trace face interpolator. Assumes ownership of di.
 
Array< BilinearFormIntegrator * > * GetDI ()
 Access all interpolators added with AddDomainInterpolator().
 
Array< Array< int > * > * GetDI_Marker ()
 
void SetAssemblyLevel (AssemblyLevel assembly_level)
 Set the desired assembly level. The default is AssemblyLevel::FULL.
 
virtual void Assemble (int skip_zeros=1)
 Construct the internal matrix representation of the discrete linear operator.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Get the output finite element space restriction matrix in transposed form.
 
- Public Member Functions inherited from mfem::MixedBilinearForm
 MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes)
 Construct a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces.
 
 MixedBilinearForm (FiniteElementSpace *tr_fes, FiniteElementSpace *te_fes, MixedBilinearForm *mbf)
 Create a MixedBilinearForm on the given trial, tr_fes, and test, te_fes, FiniteElementSpaces, using the same integrators as the MixedBilinearForm mbf.
 
virtual real_tElem (int i, int j)
 Returns a reference to: \( M_{ij} \).
 
virtual const real_tElem (int i, int j) const
 Returns a reference to: \( M_{ij} \).
 
virtual void Mult (const Vector &x, Vector &y) const
 Matrix multiplication: \( y = M x \).
 
virtual void AddMult (const Vector &x, Vector &y, const real_t a=1.0) const
 Add the matrix vector multiple to a vector: \( y += a M x \).
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Matrix transpose vector multiplication: \( y = M^T x \).
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const
 Add the matrix transpose vector multiplication: \( y += a M^T x \).
 
virtual MatrixInverseInverse () const
 Returns a pointer to (approximation) of the matrix inverse: \( M^{-1} \) (currently unimplemented and returns NULL)
 
virtual void Finalize (int skip_zeros=1)
 Finalizes the matrix initialization if the AssemblyLevel is AssemblyLevel::LEGACY.
 
void GetBlocks (Array2D< SparseMatrix * > &blocks) const
 Extract the associated matrix as SparseMatrix blocks. The number of block rows and columns is given by the vector dimensions (vdim) of the test and trial spaces, respectively.
 
const SparseMatrixSpMat () const
 Returns a const reference to the sparse matrix: \( M \).
 
SparseMatrixSpMat ()
 Returns a reference to the sparse matrix: \( M \).
 
SparseMatrixLoseMat ()
 Nullifies the internal matrix \( M \) and returns a pointer to it. Used for transferring ownership.
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi)
 Adds a domain integrator. Assumes ownership of bfi.
 
void AddDomainIntegrator (BilinearFormIntegrator *bfi, Array< int > &elem_marker)
 Adds a domain integrator. Assumes ownership of bfi.
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi)
 Adds a boundary integrator. Assumes ownership of bfi.
 
void AddBoundaryIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds a boundary integrator. Assumes ownership of bfi.
 
void AddTraceFaceIntegrator (BilinearFormIntegrator *bfi)
 Add a trace face integrator. Assumes ownership of bfi.
 
void AddBdrTraceFaceIntegrator (BilinearFormIntegrator *bfi)
 Adds a boundary trace face integrator. Assumes ownership of bfi.
 
void AddBdrTraceFaceIntegrator (BilinearFormIntegrator *bfi, Array< int > &bdr_marker)
 Adds a boundary trace face integrator. Assumes ownership of bfi.
 
Array< BilinearFormIntegrator * > * GetDBFI ()
 Access all integrators added with AddDomainIntegrator().
 
Array< Array< int > * > * GetDBFI_Marker ()
 Access all domain markers added with AddDomainIntegrator(). If no marker was specified when the integrator was added, the corresponding pointer (to Array<int>) will be NULL.
 
Array< BilinearFormIntegrator * > * GetBBFI ()
 Access all integrators added with AddBoundaryIntegrator().
 
Array< Array< int > * > * GetBBFI_Marker ()
 Access all boundary markers added with AddBoundaryIntegrator().
 
Array< BilinearFormIntegrator * > * GetTFBFI ()
 Access all integrators added with AddTraceFaceIntegrator().
 
Array< BilinearFormIntegrator * > * GetBTFBFI ()
 Access all integrators added with AddBdrTraceFaceIntegrator().
 
Array< Array< int > * > * GetBTFBFI_Marker ()
 Access all boundary markers added with AddBdrTraceFaceIntegrator()
 
void operator= (const real_t a)
 Sets all sparse values of \( M \) to a.
 
void SetAssemblyLevel (AssemblyLevel assembly_level)
 Set the desired assembly level. The default is AssemblyLevel::LEGACY.
 
void Assemble (int skip_zeros=1)
 
void AssembleDiagonal_ADAt (const Vector &D, Vector &diag) const
 Assemble the diagonal of ADA^T into diag, where A is this mixed bilinear form and D is a diagonal.
 
virtual const OperatorGetProlongation () const
 Get the input finite element space prolongation matrix.
 
virtual const OperatorGetRestriction () const
 Get the input finite element space restriction matrix.
 
virtual const OperatorGetOutputProlongation () const
 Get the test finite element space prolongation matrix.
 
virtual const OperatorGetOutputRestriction () const
 Get the test finite element space restriction matrix.
 
void ConformingAssemble ()
 For partially conforming trial and/or test FE spaces, complete the assembly process by performing \( P2^t A P1 \) where \( A \) is the internal sparse matrix; \( P1 \) and \( P2 \) are the conforming prolongation matrices of the trial and test FE spaces, respectively. After this call the MixedBilinearForm becomes an operator on the conforming FE spaces.
 
void ComputeElementMatrix (int i, DenseMatrix &elmat)
 Compute the element matrix of the given element.
 
void ComputeBdrElementMatrix (int i, DenseMatrix &elmat)
 Compute the boundary element matrix of the given boundary element.
 
void AssembleElementMatrix (int i, const DenseMatrix &elmat, int skip_zeros=1)
 Assemble the given element matrix.
 
void AssembleElementMatrix (int i, const DenseMatrix &elmat, Array< int > &trial_vdofs, Array< int > &test_vdofs, int skip_zeros=1)
 Assemble the given element matrix.
 
void AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, int skip_zeros=1)
 Assemble the given boundary element matrix.
 
void AssembleBdrElementMatrix (int i, const DenseMatrix &elmat, Array< int > &trial_vdofs, Array< int > &test_vdofs, int skip_zeros=1)
 Assemble the given boundary element matrix.
 
void EliminateTrialDofs (const Array< int > &bdr_attr_is_ess, const Vector &sol, Vector &rhs)
 Eliminate essential boundary DOFs from the columns of the system.
 
void EliminateEssentialBCFromTrialDofs (const Array< int > &marked_vdofs, const Vector &sol, Vector &rhs)
 Eliminate the list of DOFs from the columns of the system.
 
virtual void EliminateTestDofs (const Array< int > &bdr_attr_is_ess)
 Eliminate essential boundary DOFs from the rows of the system.
 
template<typename OpType >
void FormRectangularSystemMatrix (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, OpType &A)
 Form the column-constrained linear system matrix A.
 
virtual void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OperatorHandle &A, Vector &X, Vector &B)
 Form the linear system A X = B, corresponding to this mixed bilinear form and the linear form b(.).
 
template<typename OpType >
void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, OpType &A, Vector &X, Vector &B)
 Form the linear system A X = B, corresponding to this bilinear form and the linear form b(.).
 
void Update ()
 Must be called after making changes to trial_fes or test_fes.
 
FiniteElementSpaceTrialFESpace ()
 Return the trial FE space associated with the BilinearForm.
 
const FiniteElementSpaceTrialFESpace () const
 Read-only access to the associated trial FiniteElementSpace.
 
FiniteElementSpaceTestFESpace ()
 Return the test FE space associated with the BilinearForm.
 
const FiniteElementSpaceTestFESpace () const
 Read-only access to the associated test FiniteElementSpace.
 
virtual ~MixedBilinearForm ()
 Deletes internal matrices, bilinear integrators, and the BilinearFormExtension.
 
- Public Member Functions inherited from mfem::Matrix
 Matrix (int s)
 Creates a square matrix of size s.
 
 Matrix (int h, int w)
 Creates a matrix of the given height and width.
 
bool IsSquare () const
 Returns whether the matrix is a square matrix.
 
virtual void Print (std::ostream &out=mfem::out, int width_=4) const
 Prints matrix to stream out.
 
virtual ~Matrix ()
 Destroys matrix.
 
- Public Member Functions inherited from mfem::Operator
void InitTVectors (const Operator *Po, const Operator *Ri, const Operator *Pi, Vector &x, Vector &b, Vector &X, Vector &B) const
 Initializes memory for true vectors of linear system.
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0).
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size).
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows().
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height().
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols().
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width().
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator.
 
virtual void ArrayMult (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Operator application on a matrix: Y=A(X).
 
virtual void ArrayMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X).
 
virtual void ArrayAddMult (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X).
 
virtual void ArrayAddMultTranspose (const Array< const Vector * > &X, Array< Vector * > &Y, const real_t a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X).
 
virtual OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error.
 
virtual void AssembleDiagonal (Vector &diag) const
 Computes the diagonal entries into diag. Typically, this operation only makes sense for linear Operators. In some cases, only an approximation of the diagonal is computed.
 
void FormLinearSystem (const Array< int > &ess_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B, int copy_interior=0)
 Form a constrained linear system using a matrix-free approach.
 
void FormRectangularLinearSystem (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Vector &x, Vector &b, Operator *&A, Vector &X, Vector &B)
 Form a column-constrained linear system using a matrix-free approach.
 
virtual void RecoverFEMSolution (const Vector &X, const Vector &b, Vector &x)
 Reconstruct a solution vector x (e.g. a GridFunction) from the solution X of a constrained linear system obtained from Operator::FormLinearSystem() or Operator::FormRectangularLinearSystem().
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator.
 
void FormRectangularSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator (including constraints).
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator.
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format.
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format.
 
virtual ~Operator ()
 Virtual destructor.
 
Type GetType () const
 Return the type ID of the Operator class.
 

Protected Attributes

ParFiniteElementSpacedomain_fes
 Points to the same object as trial_fes.
 
ParFiniteElementSpacerange_fes
 Points to the same object as test_fes.
 
- Protected Attributes inherited from mfem::MixedBilinearForm
SparseMatrixmat
 Owned.
 
SparseMatrixmat_e
 Owned.
 
FiniteElementSpacetrial_fes
 Not owned.
 
FiniteElementSpacetest_fes
 Not owned.
 
AssemblyLevel assembly
 The form assembly level (full, partial, etc.)
 
MixedBilinearFormExtensionext
 
int extern_bfs
 Indicates the BilinearFormIntegrators stored in MixedBilinearForm::domain_integs, MixedBilinearForm::boundary_integs, MixedBilinearForm::trace_face_integs and MixedBilinearForm::boundary_trace_face_integs are owned by another MixedBilinearForm.
 
Array< BilinearFormIntegrator * > domain_integs
 Domain integrators.
 
Array< Array< int > * > domain_integs_marker
 Entries are not owned.
 
Array< BilinearFormIntegrator * > boundary_integs
 Boundary integrators.
 
Array< Array< int > * > boundary_integs_marker
 Entries are not owned.
 
Array< BilinearFormIntegrator * > trace_face_integs
 Trace face (skeleton) integrators.
 
Array< BilinearFormIntegrator * > boundary_trace_face_integs
 Boundary trace face (skeleton) integrators.
 
Array< Array< int > * > boundary_trace_face_integs_marker
 Entries are not owned.
 
DenseMatrix elemmat
 
Array< int > trial_vdofs
 
Array< int > test_vdofs
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix.
 
int width
 Dimension of the input / number of columns in the matrix.
 

Additional Inherited Members

- Public Types inherited from mfem::Operator
enum  DiagonalPolicy { DIAG_ZERO , DIAG_ONE , DIAG_KEEP }
 Defines operator diagonal policy upon elimination of rows and/or columns. More...
 
enum  Type {
  ANY_TYPE , MFEM_SPARSEMAT , Hypre_ParCSR , PETSC_MATAIJ ,
  PETSC_MATIS , PETSC_MATSHELL , PETSC_MATNEST , PETSC_MATHYPRE ,
  PETSC_MATGENERIC , Complex_Operator , MFEM_ComplexSparseMat , Complex_Hypre_ParCSR ,
  Complex_DenseMat , MFEM_Block_Matrix , MFEM_Block_Operator
}
 Enumeration defining IDs for some classes derived from Operator. More...
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator()
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator()
 
OperatorSetupRAP (const Operator *Pi, const Operator *Po)
 Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt".
 

Detailed Description

The parallel matrix representation a linear operator between parallel finite element spaces

Definition at line 338 of file pbilinearform.hpp.

Constructor & Destructor Documentation

◆ ParDiscreteLinearOperator()

mfem::ParDiscreteLinearOperator::ParDiscreteLinearOperator ( ParFiniteElementSpace * dfes,
ParFiniteElementSpace * rfes )
inline

Construct a ParDiscreteLinearOperator on the given FiniteElementSpaces dfes (domain FE space) and rfes (range FE space).

The pointers dfes and rfes are not owned by the newly constructed object.

Definition at line 359 of file pbilinearform.hpp.

◆ ~ParDiscreteLinearOperator()

virtual mfem::ParDiscreteLinearOperator::~ParDiscreteLinearOperator ( )
inlinevirtual

Definition at line 380 of file pbilinearform.hpp.

Member Function Documentation

◆ FormRectangularSystemMatrix() [1/3]

void mfem::MixedBilinearForm::FormRectangularSystemMatrix ( const Array< int > & trial_tdof_list,
const Array< int > & test_tdof_list,
OperatorHandle & A )
virtual

Return in A that is column-constrained.

This returns the same operator as FormRectangularLinearSystem(), but does without the transformations of the right-hand side.

Reimplemented from mfem::MixedBilinearForm.

Definition at line 1015 of file bilinearform.cpp.

◆ FormRectangularSystemMatrix() [2/3]

template<typename OpType >
void mfem::MixedBilinearForm::FormRectangularSystemMatrix ( const Array< int > & trial_tdof_list,
const Array< int > & test_tdof_list,
OpType & A )
inline

Form the column-constrained linear system matrix A.

Version of the method FormRectangularSystemMatrix() where the system matrix is returned in the variable A, of type OpType, holding a reference* to the system matrix (created with the method OpType::MakeRef()). The reference will be invalidated when SetOperatorType(), Update(), or the destructor is called.

Definition at line 1027 of file bilinearform.hpp.

◆ FormRectangularSystemMatrix() [3/3]

void mfem::ParDiscreteLinearOperator::FormRectangularSystemMatrix ( OperatorHandle & A)
virtual

Return in A a parallel (on truedofs) version of this operator.

Definition at line 741 of file pbilinearform.cpp.

◆ GetParBlocks()

void mfem::ParDiscreteLinearOperator::GetParBlocks ( Array2D< HypreParMatrix * > & blocks) const

Extract the parallel blocks corresponding to the vector dimensions of the domain and range parallel finite element spaces

Definition at line 753 of file pbilinearform.cpp.

◆ ParallelAssemble() [1/2]

HypreParMatrix * mfem::ParDiscreteLinearOperator::ParallelAssemble ( ) const

Returns the matrix "assembled" on the true dofs.

Definition at line 702 of file pbilinearform.cpp.

◆ ParallelAssemble() [2/2]

void mfem::ParDiscreteLinearOperator::ParallelAssemble ( OperatorHandle & A)

Returns the matrix assembled on the true dofs, i.e. A = R_test A_local P_trial, in the format (type id) specified by A.

Definition at line 713 of file pbilinearform.cpp.

Member Data Documentation

◆ domain_fes

ParFiniteElementSpace* mfem::ParDiscreteLinearOperator::domain_fes
protected

Points to the same object as trial_fes.

Definition at line 342 of file pbilinearform.hpp.

◆ range_fes

ParFiniteElementSpace* mfem::ParDiscreteLinearOperator::range_fes
protected

Points to the same object as test_fes.

Definition at line 344 of file pbilinearform.hpp.


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