MFEM  v4.5.2
Finite element discretization library
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
mfem::NonlinearForm Class Reference

#include <nonlinearform.hpp>

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

Public Member Functions

 NonlinearForm (FiniteElementSpace *f)
 Construct a NonlinearForm on the given FiniteElementSpace, f. More...
 
void SetAssemblyLevel (AssemblyLevel assembly_level)
 Set the desired assembly level. The default is AssemblyLevel::LEGACY. More...
 
FiniteElementSpaceFESpace ()
 
const FiniteElementSpaceFESpace () const
 
void AddDomainIntegrator (NonlinearFormIntegrator *nlfi)
 Adds new Domain Integrator. More...
 
Array< NonlinearFormIntegrator * > * GetDNFI ()
 Access all integrators added with AddDomainIntegrator(). More...
 
const Array< NonlinearFormIntegrator * > * GetDNFI () const
 
void AddInteriorFaceIntegrator (NonlinearFormIntegrator *nlfi)
 Adds new Interior Face Integrator. More...
 
const Array< NonlinearFormIntegrator * > & GetInteriorFaceIntegrators () const
 Access all interior face integrators added with AddInteriorFaceIntegrator(). More...
 
void AddBdrFaceIntegrator (NonlinearFormIntegrator *nlfi)
 Adds new Boundary Face Integrator. More...
 
void AddBdrFaceIntegrator (NonlinearFormIntegrator *nfi, Array< int > &bdr_marker)
 Adds new Boundary Face Integrator, restricted to specific boundary attributes. More...
 
const Array< NonlinearFormIntegrator * > & GetBdrFaceIntegrators () const
 Access all boundary face integrators added with AddBdrFaceIntegrator(). More...
 
void SetEssentialBC (const Array< int > &bdr_attr_is_ess, Vector *rhs=NULL)
 Specify essential boundary conditions. More...
 
void SetEssentialVDofs (const Array< int > &ess_vdofs_list)
 Specify essential boundary conditions. More...
 
void SetEssentialTrueDofs (const Array< int > &ess_tdof_list_)
 Specify essential boundary conditions. More...
 
const Array< int > & GetEssentialTrueDofs () const
 Return a (read-only) list of all essential true dofs. More...
 
double GetGridFunctionEnergy (const Vector &x) const
 Compute the energy corresponding to the state x. More...
 
virtual double GetEnergy (const Vector &x) const
 Compute the energy corresponding to the state x. More...
 
virtual void Mult (const Vector &x, Vector &y) const
 Evaluate the action of the NonlinearForm. More...
 
virtual OperatorGetGradient (const Vector &x) const
 Compute the gradient Operator of the NonlinearForm corresponding to the state x. More...
 
virtual void Update ()
 Update the NonlinearForm to propagate updates of the associated FE space. More...
 
virtual void Setup ()
 Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally, precompute and store data that will be reused in subsequent call to Mult(). More...
 
virtual const OperatorGetProlongation () const
 Get the finite element space prolongation matrix. More...
 
virtual const OperatorGetRestriction () const
 Get the finite element space restriction matrix. More...
 
virtual ~NonlinearForm ()
 Destroy the NonlinearForm including the owned NonlinearFormIntegrators and gradient Operator. More...
 
- 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. More...
 
 Operator (int s=0)
 Construct a square Operator with given size s (default 0). More...
 
 Operator (int h, int w)
 Construct an Operator with the given height (output size) and width (input size). More...
 
int Height () const
 Get the height (size of output) of the Operator. Synonym with NumRows(). More...
 
int NumRows () const
 Get the number of rows (size of output) of the Operator. Synonym with Height(). More...
 
int Width () const
 Get the width (size of input) of the Operator. Synonym with NumCols(). More...
 
int NumCols () const
 Get the number of columns (size of input) of the Operator. Synonym with Width(). More...
 
virtual MemoryClass GetMemoryClass () const
 Return the MemoryClass preferred by the Operator. More...
 
virtual void MultTranspose (const Vector &x, Vector &y) const
 Action of the transpose operator: y=A^t(x). The default behavior in class Operator is to generate an error. More...
 
virtual void AddMult (const Vector &x, Vector &y, const double a=1.0) const
 Operator application: y+=A(x) (default) or y+=a*A(x). More...
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const double a=1.0) const
 Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x). More...
 
virtual void ArrayMult (const Array< const Vector *> &X, Array< Vector *> &Y) const
 Operator application on a matrix: Y=A(X). More...
 
virtual void ArrayMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y) const
 Action of the transpose operator on a matrix: Y=A^t(X). More...
 
virtual void ArrayAddMult (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
 Operator application on a matrix: Y+=A(X) (default) or Y+=a*A(X). More...
 
virtual void ArrayAddMultTranspose (const Array< const Vector *> &X, Array< Vector *> &Y, const double a=1.0) const
 Operator transpose application on a matrix: Y+=A^t(X) (default) or Y+=a*A^t(X). More...
 
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. More...
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. More...
 
virtual const OperatorGetOutputRestriction () const
 Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity. More...
 
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. More...
 
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. More...
 
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(). More...
 
void FormSystemOperator (const Array< int > &ess_tdof_list, Operator *&A)
 Return in A a parallel (on truedofs) version of this square operator. More...
 
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). More...
 
void FormDiscreteOperator (Operator *&A)
 Return in A a parallel (on truedofs) version of this rectangular operator. More...
 
void PrintMatlab (std::ostream &out, int n, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual void PrintMatlab (std::ostream &out) const
 Prints operator in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 

Protected Member Functions

bool Serial () const
 
const VectorProlongate (const Vector &x) const
 
- Protected Member Functions inherited from mfem::Operator
void FormConstrainedSystemOperator (const Array< int > &ess_tdof_list, ConstrainedOperator *&Aout)
 see FormSystemOperator() More...
 
void FormRectangularConstrainedSystemOperator (const Array< int > &trial_tdof_list, const Array< int > &test_tdof_list, RectangularConstrainedOperator *&Aout)
 see FormRectangularSystemOperator() More...
 
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". More...
 

Protected Attributes

AssemblyLevel assembly
 The assembly level. More...
 
NonlinearFormExtensionext
 Extension for supporting different AssemblyLevels. More...
 
FiniteElementSpacefes
 FE space on which the form lives. More...
 
Array< NonlinearFormIntegrator * > dnfi
 Set of Domain Integrators to be assembled (added). More...
 
Array< NonlinearFormIntegrator * > fnfi
 Set of interior face Integrators to be assembled (added). More...
 
Array< NonlinearFormIntegrator * > bfnfi
 Set of boundary face Integrators to be assembled (added). More...
 
Array< Array< int > * > bfnfi_marker
 
SparseMatrixGrad
 
SparseMatrixcGrad
 
OperatorHandle hGrad
 Gradient Operator when not assembled as a matrix. More...
 
Array< int > ess_tdof_list
 A list of all essential true dofs. More...
 
long sequence
 Counter for updates propagated from the FiniteElementSpace. More...
 
Vector aux1
 Auxiliary Vectors. More...
 
Vector aux2
 
const OperatorP
 Pointer to the prolongation matrix of fes, may be NULL. More...
 
const SparseMatrixcP
 The result of dynamic-casting P to SparseMatrix pointer. More...
 
- Protected Attributes inherited from mfem::Operator
int height
 Dimension of the output / number of rows in the matrix. More...
 
int width
 Dimension of the input / number of columns in the matrix. More...
 

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...
 

Detailed Description

Definition at line 24 of file nonlinearform.hpp.

Constructor & Destructor Documentation

◆ NonlinearForm()

mfem::NonlinearForm::NonlinearForm ( FiniteElementSpace f)
inline

Construct a NonlinearForm on the given FiniteElementSpace, f.

As an Operator, the NonlinearForm has input and output size equal to the number of true degrees of freedom, i.e. f->GetTrueVSize().

Definition at line 73 of file nonlinearform.hpp.

◆ ~NonlinearForm()

mfem::NonlinearForm::~NonlinearForm ( )
virtual

Destroy the NonlinearForm including the owned NonlinearFormIntegrators and gradient Operator.

Definition at line 470 of file nonlinearform.cpp.

Member Function Documentation

◆ AddBdrFaceIntegrator() [1/2]

void mfem::NonlinearForm::AddBdrFaceIntegrator ( NonlinearFormIntegrator nlfi)
inline

Adds new Boundary Face Integrator.

Definition at line 127 of file nonlinearform.hpp.

◆ AddBdrFaceIntegrator() [2/2]

void mfem::NonlinearForm::AddBdrFaceIntegrator ( NonlinearFormIntegrator nfi,
Array< int > &  bdr_marker 
)
inline

Adds new Boundary Face Integrator, restricted to specific boundary attributes.

Definition at line 132 of file nonlinearform.hpp.

◆ AddDomainIntegrator()

void mfem::NonlinearForm::AddDomainIntegrator ( NonlinearFormIntegrator nlfi)
inline

Adds new Domain Integrator.

Definition at line 110 of file nonlinearform.hpp.

◆ AddInteriorFaceIntegrator()

void mfem::NonlinearForm::AddInteriorFaceIntegrator ( NonlinearFormIntegrator nlfi)
inline

Adds new Interior Face Integrator.

Definition at line 118 of file nonlinearform.hpp.

◆ FESpace() [1/2]

FiniteElementSpace* mfem::NonlinearForm::FESpace ( )
inline

Definition at line 106 of file nonlinearform.hpp.

◆ FESpace() [2/2]

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

Definition at line 107 of file nonlinearform.hpp.

◆ GetBdrFaceIntegrators()

const Array<NonlinearFormIntegrator*>& mfem::NonlinearForm::GetBdrFaceIntegrators ( ) const
inline

Access all boundary face integrators added with AddBdrFaceIntegrator().

Definition at line 138 of file nonlinearform.hpp.

◆ GetDNFI() [1/2]

Array<NonlinearFormIntegrator*>* mfem::NonlinearForm::GetDNFI ( )
inline

Access all integrators added with AddDomainIntegrator().

Definition at line 114 of file nonlinearform.hpp.

◆ GetDNFI() [2/2]

const Array<NonlinearFormIntegrator*>* mfem::NonlinearForm::GetDNFI ( ) const
inline

Definition at line 115 of file nonlinearform.hpp.

◆ GetEnergy()

virtual double mfem::NonlinearForm::GetEnergy ( const Vector x) const
inlinevirtual

Compute the energy corresponding to the state x.

In general, x may have non-homogeneous essential boundary values.

The state x must be a true-dof vector.

Reimplemented in mfem::ParNonlinearForm.

Definition at line 171 of file nonlinearform.hpp.

◆ GetEssentialTrueDofs()

const Array<int>& mfem::NonlinearForm::GetEssentialTrueDofs ( ) const
inline

Return a (read-only) list of all essential true dofs.

Definition at line 158 of file nonlinearform.hpp.

◆ GetGradient()

Operator & mfem::NonlinearForm::GetGradient ( const Vector x) const
virtual

Compute the gradient Operator of the NonlinearForm corresponding to the state x.

Any previously specified essential boundary conditions will be automatically imposed on the gradient operator.

The returned object is valid until the next call to this method or the destruction of this object.

In general, x may have non-homogeneous essential boundary values.

The state x must be a true-dof vector.

Reimplemented from mfem::Operator.

Reimplemented in mfem::ParNonlinearForm.

Definition at line 290 of file nonlinearform.cpp.

◆ GetGridFunctionEnergy()

double mfem::NonlinearForm::GetGridFunctionEnergy ( const Vector x) const

Compute the energy corresponding to the state x.

In general, x may have non-homogeneous essential boundary values.

The state x must be a "GridFunction size" vector, i.e. its size must be fes->GetVSize().

Definition at line 86 of file nonlinearform.cpp.

◆ GetInteriorFaceIntegrators()

const Array<NonlinearFormIntegrator*>& mfem::NonlinearForm::GetInteriorFaceIntegrators ( ) const
inline

Access all interior face integrators added with AddInteriorFaceIntegrator().

Definition at line 123 of file nonlinearform.hpp.

◆ GetProlongation()

virtual const Operator* mfem::NonlinearForm::GetProlongation ( ) const
inlinevirtual

Get the finite element space prolongation matrix.

Reimplemented from mfem::Operator.

Definition at line 209 of file nonlinearform.hpp.

◆ GetRestriction()

virtual const Operator* mfem::NonlinearForm::GetRestriction ( ) const
inlinevirtual

Get the finite element space restriction matrix.

Reimplemented from mfem::Operator.

Definition at line 211 of file nonlinearform.hpp.

◆ Mult()

void mfem::NonlinearForm::Mult ( const Vector x,
Vector y 
) const
virtual

Evaluate the action of the NonlinearForm.

The input essential dofs in x will, generally, be non-zero. However, the output essential dofs in y will always be set to zero.

Both the input and the output vectors, x and y, must be true-dof vectors, i.e. their size must be fes->GetTrueVSize().

Implements mfem::Operator.

Reimplemented in mfem::ParNonlinearForm.

Definition at line 143 of file nonlinearform.cpp.

◆ Prolongate()

const Vector & mfem::NonlinearForm::Prolongate ( const Vector x) const
protected

Definition at line 131 of file nonlinearform.cpp.

◆ Serial()

bool mfem::NonlinearForm::Serial ( ) const
inlineprotected

Definition at line 66 of file nonlinearform.hpp.

◆ SetAssemblyLevel()

void mfem::NonlinearForm::SetAssemblyLevel ( AssemblyLevel  assembly_level)

Set the desired assembly level. The default is AssemblyLevel::LEGACY.

For nonlinear operators, the "matrix" assembly levels usually do not make sense, so only LEGACY, NONE (matrix-free) and PARTIAL are supported.

Currently, AssemblyLevel::LEGACY uses the standard nonlinear action methods like AssembleElementVector of the NonlinearFormIntegrator class which work only on CPU and do not utilize features such as fast tensor-product basis evaluations. In this mode, the gradient operator is constructed as a SparseMatrix (or, in parallel, format such as HypreParMatrix).

When using AssemblyLevel::PARTIAL, the action is performed using methods like AddMultPA of the NonlinearFormIntegrator class which typically support both CPU and GPU backends and utilize features such as fast tensor-product basis evaluations. In this mode, the gradient operator also uses partial assembly with support for CPU and GPU backends.

When using AssemblyLevel::NONE, the action is performed using methods like AddMultMF of the NonlinearFormIntegrator class which typically support both CPU and GPU backends and utilize features such as fast tensor-product basis evaluations. In this mode, the gradient operator is currently not supported.

This method must be called before "assembly" with Setup().

Definition at line 18 of file nonlinearform.cpp.

◆ SetEssentialBC()

void mfem::NonlinearForm::SetEssentialBC ( const Array< int > &  bdr_attr_is_ess,
Vector rhs = NULL 
)

Specify essential boundary conditions.

This method calls FiniteElementSpace::GetEssentialTrueDofs() and stores the result internally for use by other methods. If the rhs pointer is not NULL, its essential true dofs will be set to zero. This makes it "compatible" with the output vectors from the Mult() method which also have zero entries at the essential true dofs.

Definition at line 41 of file nonlinearform.cpp.

◆ SetEssentialTrueDofs()

void mfem::NonlinearForm::SetEssentialTrueDofs ( const Array< int > &  ess_tdof_list_)
inline

Specify essential boundary conditions.

Definition at line 154 of file nonlinearform.hpp.

◆ SetEssentialVDofs()

void mfem::NonlinearForm::SetEssentialVDofs ( const Array< int > &  ess_vdofs_list)

Specify essential boundary conditions.

Use either SetEssentialBC() or SetEssentialTrueDofs() if possible.

Definition at line 56 of file nonlinearform.cpp.

◆ Setup()

void mfem::NonlinearForm::Setup ( )
virtual

Setup the NonlinearForm: based on the current AssemblyLevel and the current mesh, optionally, precompute and store data that will be reused in subsequent call to Mult().

Typically, this method has to be called before Mult() when using AssemblyLevel::PARTIAL, after calling Update(), or after modifying the mesh coordinates.

Definition at line 465 of file nonlinearform.cpp.

◆ Update()

void mfem::NonlinearForm::Update ( )
virtual

Update the NonlinearForm to propagate updates of the associated FE space.

After calling this method, the essential boundary conditions need to be set again.

Reimplemented in mfem::ParNonlinearForm.

Definition at line 448 of file nonlinearform.cpp.

Member Data Documentation

◆ assembly

AssemblyLevel mfem::NonlinearForm::assembly
protected

The assembly level.

Definition at line 28 of file nonlinearform.hpp.

◆ aux1

Vector mfem::NonlinearForm::aux1
mutableprotected

Auxiliary Vectors.

Definition at line 59 of file nonlinearform.hpp.

◆ aux2

Vector mfem::NonlinearForm::aux2
mutableprotected

Definition at line 59 of file nonlinearform.hpp.

◆ bfnfi

Array<NonlinearFormIntegrator*> mfem::NonlinearForm::bfnfi
protected

Set of boundary face Integrators to be assembled (added).

Definition at line 45 of file nonlinearform.hpp.

◆ bfnfi_marker

Array<Array<int>*> mfem::NonlinearForm::bfnfi_marker
protected

Definition at line 46 of file nonlinearform.hpp.

◆ cGrad

SparseMatrix * mfem::NonlinearForm::cGrad
mutableprotected

Definition at line 48 of file nonlinearform.hpp.

◆ cP

const SparseMatrix* mfem::NonlinearForm::cP
protected

The result of dynamic-casting P to SparseMatrix pointer.

Definition at line 64 of file nonlinearform.hpp.

◆ dnfi

Array<NonlinearFormIntegrator*> mfem::NonlinearForm::dnfi
protected

Set of Domain Integrators to be assembled (added).

Definition at line 39 of file nonlinearform.hpp.

◆ ess_tdof_list

Array<int> mfem::NonlinearForm::ess_tdof_list
protected

A list of all essential true dofs.

Definition at line 53 of file nonlinearform.hpp.

◆ ext

NonlinearFormExtension* mfem::NonlinearForm::ext
protected

Extension for supporting different AssemblyLevels.

For nonlinear operators, the "matrix" assembly levels usually do not make sense, so only PARTIAL and NONE (matrix-free) are supported.

Definition at line 33 of file nonlinearform.hpp.

◆ fes

FiniteElementSpace* mfem::NonlinearForm::fes
protected

FE space on which the form lives.

Definition at line 36 of file nonlinearform.hpp.

◆ fnfi

Array<NonlinearFormIntegrator*> mfem::NonlinearForm::fnfi
protected

Set of interior face Integrators to be assembled (added).

Definition at line 42 of file nonlinearform.hpp.

◆ Grad

SparseMatrix* mfem::NonlinearForm::Grad
mutableprotected

Definition at line 48 of file nonlinearform.hpp.

◆ hGrad

OperatorHandle mfem::NonlinearForm::hGrad
mutableprotected

Gradient Operator when not assembled as a matrix.

Definition at line 50 of file nonlinearform.hpp.

◆ P

const Operator* mfem::NonlinearForm::P
protected

Pointer to the prolongation matrix of fes, may be NULL.

Definition at line 62 of file nonlinearform.hpp.

◆ sequence

long mfem::NonlinearForm::sequence
protected

Counter for updates propagated from the FiniteElementSpace.

Definition at line 56 of file nonlinearform.hpp.


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