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

Base abstract class for first order time dependent operators. More...

#include <operator.hpp>

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

Public Types

enum  Type { EXPLICIT , IMPLICIT , HOMOGENEOUS }
 
enum  EvalMode { NORMAL , ADDITIVE_TERM_1 , ADDITIVE_TERM_2 }
 Evaluation mode. See SetEvalMode() for details. More...
 
- 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...
 

Public Member Functions

 TimeDependentOperator (int n=0, real_t t_=0.0, Type type_=EXPLICIT)
 Construct a "square" TimeDependentOperator y = f(x,t), where x and y have the same dimension n.
 
 TimeDependentOperator (int h, int w, real_t t_=0.0, Type type_=EXPLICIT)
 Construct a TimeDependentOperator y = f(x,t), where x and y have dimensions w and h, respectively.
 
virtual real_t GetTime () const
 Read the currently set time.
 
virtual void SetTime (const real_t t_)
 Set the current time.
 
bool isExplicit () const
 True if type is EXPLICIT.
 
bool isImplicit () const
 True if type is IMPLICIT or HOMOGENEOUS.
 
bool isHomogeneous () const
 True if type is HOMOGENEOUS.
 
EvalMode GetEvalMode () const
 Return the current evaluation mode. See SetEvalMode() for details.
 
virtual void SetEvalMode (const EvalMode new_eval_mode)
 Set the evaluation mode of the time-dependent operator.
 
virtual void ExplicitMult (const Vector &x, Vector &y) const
 Perform the action of the explicit part of the operator, G: y = G(x, t) where t is the current time.
 
virtual void ImplicitMult (const Vector &x, const Vector &k, Vector &y) const
 Perform the action of the implicit part of the operator, F: y = F(x, k, t) where t is the current time.
 
virtual void Mult (const Vector &x, Vector &y) const
 Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x, k, t) = G(x, t) and t is the current time.
 
virtual void ImplicitSolve (const real_t dt, const Vector &x, Vector &k)
 Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.
 
virtual OperatorGetImplicitGradient (const Vector &x, const Vector &k, real_t shift) const
 Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time.
 
virtual OperatorGetExplicitGradient (const Vector &x) const
 Return an Operator representing dG/dx at the given point x and the currently set time.
 
virtual int SUNImplicitSetup (const Vector &x, const Vector &fx, int jok, int *jcur, real_t gamma)
 Setup the ODE linear system \( A(x,t) = (I - gamma J) \) or \( A = (M - gamma J) \), where \( J(x,t) = \frac{df}{dt(x,t)} \).
 
virtual int SUNImplicitSolve (const Vector &b, Vector &x, real_t tol)
 Solve the ODE linear system \( A x = b \) as setup by the method SUNImplicitSetup().
 
virtual int SUNMassSetup ()
 Setup the mass matrix in the ODE system \( M y' = f(y,t) \) .
 
virtual int SUNMassSolve (const Vector &b, Vector &x, real_t tol)
 Solve the mass matrix linear system \( M x = b \) as setup by the method SUNMassSetup().
 
virtual int SUNMassMult (const Vector &x, Vector &v)
 Compute the mass matrix-vector product \( v = M x \) .
 
virtual ~TimeDependentOperator ()
 
- 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 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.
 
virtual void AddMult (const Vector &x, Vector &y, const real_t a=1.0) const
 Operator application: y+=A(x) (default) or y+=a*A(x).
 
virtual void AddMultTranspose (const Vector &x, Vector &y, const real_t a=1.0) const
 Operator transpose application: y+=A^t(x) (default) or y+=a*A^t(x).
 
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.
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity.
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
 
virtual const OperatorGetOutputProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity.
 
virtual const OperatorGetOutputRestrictionTranspose () const
 Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators.
 
virtual const OperatorGetOutputRestriction () const
 Restriction operator from output vectors for the operator to linear algebra (linear system) vectors. NULL means identity.
 
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

real_t t
 Current time.
 
Type type
 Describes the form of the TimeDependentOperator.
 
EvalMode eval_mode
 Current evaluation mode.
 
- 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

- 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

Base abstract class for first order time dependent operators.

Operator of the form: (x,t) -> f(x,t), where k = f(x,t) generally solves the algebraic equation F(x,k,t) = G(x,t). The functions F and G represent the implicit and explicit parts of the operator, respectively. For explicit operators, F(x,k,t) = k, so f(x,t) = G(x,t).

Definition at line 316 of file operator.hpp.

Member Enumeration Documentation

◆ EvalMode

Evaluation mode. See SetEvalMode() for details.

Enumerator
NORMAL 

Normal evaluation.

ADDITIVE_TERM_1 

Assuming additive split, f(x,t) = f1(x,t) + f2(x,t), evaluate the first term, f1.

ADDITIVE_TERM_2 

Assuming additive split, f(x,t) = f1(x,t) + f2(x,t), evaluate the second term, f2.

Definition at line 327 of file operator.hpp.

◆ Type

Enumerator
EXPLICIT 

This type assumes F(x,k,t) = k, i.e. k = f(x,t) = G(x,t).

IMPLICIT 

This is the most general type, no assumptions on F and G.

HOMOGENEOUS 

This type assumes that G(x,t) = 0.

Definition at line 319 of file operator.hpp.

Constructor & Destructor Documentation

◆ TimeDependentOperator() [1/2]

mfem::TimeDependentOperator::TimeDependentOperator ( int n = 0,
real_t t_ = 0.0,
Type type_ = EXPLICIT )
inlineexplicit

Construct a "square" TimeDependentOperator y = f(x,t), where x and y have the same dimension n.

Definition at line 347 of file operator.hpp.

◆ TimeDependentOperator() [2/2]

mfem::TimeDependentOperator::TimeDependentOperator ( int h,
int w,
real_t t_ = 0.0,
Type type_ = EXPLICIT )
inline

Construct a TimeDependentOperator y = f(x,t), where x and y have dimensions w and h, respectively.

Definition at line 353 of file operator.hpp.

◆ ~TimeDependentOperator()

virtual mfem::TimeDependentOperator::~TimeDependentOperator ( )
inlinevirtual

Definition at line 499 of file operator.hpp.

Member Function Documentation

◆ ExplicitMult()

void mfem::TimeDependentOperator::ExplicitMult ( const Vector & x,
Vector & y ) const
virtual

Perform the action of the explicit part of the operator, G: y = G(x, t) where t is the current time.

Presently, this method is used by some PETSc ODE solvers, for more details, see the PETSc Manual.

Definition at line 282 of file operator.cpp.

◆ GetEvalMode()

EvalMode mfem::TimeDependentOperator::GetEvalMode ( ) const
inline

Return the current evaluation mode. See SetEvalMode() for details.

Definition at line 370 of file operator.hpp.

◆ GetExplicitGradient()

Operator & mfem::TimeDependentOperator::GetExplicitGradient ( const Vector & x) const
virtual

Return an Operator representing dG/dx at the given point x and the currently set time.

Presently, this method is used by some PETSc ODE solvers, for more details, see the PETSc Manual.

Definition at line 312 of file operator.cpp.

◆ GetImplicitGradient()

Operator & mfem::TimeDependentOperator::GetImplicitGradient ( const Vector & x,
const Vector & k,
real_t shift ) const
virtual

Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time.

Presently, this method is used by some PETSc ODE solvers, for more details, see the PETSc Manual.

Definition at line 304 of file operator.cpp.

◆ GetTime()

virtual real_t mfem::TimeDependentOperator::GetTime ( ) const
inlinevirtual

Read the currently set time.

Definition at line 357 of file operator.hpp.

◆ ImplicitMult()

void mfem::TimeDependentOperator::ImplicitMult ( const Vector & x,
const Vector & k,
Vector & y ) const
virtual

Perform the action of the implicit part of the operator, F: y = F(x, k, t) where t is the current time.

Presently, this method is used by some PETSc ODE solvers, for more details, see the PETSc Manual.

Definition at line 287 of file operator.cpp.

◆ ImplicitSolve()

void mfem::TimeDependentOperator::ImplicitSolve ( const real_t dt,
const Vector & x,
Vector & k )
virtual

Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t.

For general F and G, the equation for k becomes: F(x + dt k, k, t) = G(x + dt k, t).

The input vector x corresponds to time index (or cycle) n, while the currently set time, t, and the result vector k correspond to time index n+1. The time step dt corresponds to the time interval between cycles n and n+1.

This method allows for the abstract implementation of some time integration methods, including diagonal implicit Runge-Kutta (DIRK) methods and the backward Euler method in particular.

If not re-implemented, this method simply generates an error.

Reimplemented in mfem::electromagnetics::MagneticDiffusionEOperator, mfem::electromagnetics::MaxwellSolver, and mfem::SecondOrderTimeDependentOperator.

Definition at line 298 of file operator.cpp.

◆ isExplicit()

bool mfem::TimeDependentOperator::isExplicit ( ) const
inline

True if type is EXPLICIT.

Definition at line 363 of file operator.hpp.

◆ isHomogeneous()

bool mfem::TimeDependentOperator::isHomogeneous ( ) const
inline

True if type is HOMOGENEOUS.

Definition at line 367 of file operator.hpp.

◆ isImplicit()

bool mfem::TimeDependentOperator::isImplicit ( ) const
inline

True if type is IMPLICIT or HOMOGENEOUS.

Definition at line 365 of file operator.hpp.

◆ Mult()

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

Perform the action of the operator: y = k = f(x, t), where k solves the algebraic equation F(x, k, t) = G(x, t) and t is the current time.

Implements mfem::Operator.

Reimplemented in mfem::AdvectionOper, mfem::DGHyperbolicConservationLaws, mfem::electromagnetics::MagneticDiffusionEOperator, mfem::electromagnetics::MaxwellSolver, mfem::ParAdvectorCGOper, mfem::SecondOrderTimeDependentOperator, and mfem::SerialAdvectorCGOper.

Definition at line 293 of file operator.cpp.

◆ SetEvalMode()

virtual void mfem::TimeDependentOperator::SetEvalMode ( const EvalMode new_eval_mode)
inlinevirtual

Set the evaluation mode of the time-dependent operator.

The evaluation mode is a switch that allows time-stepping methods to request evaluation of separate components/terms of the time-dependent operator. For example, IMEX methods typically assume additive split of the operator: f(x,t) = f1(x,t) + f2(x,t) and they rely on the ability to evaluate the two terms separately.

Generally, setting the evaluation mode should affect the behavior of all evaluation-related methods in the class, such as Mult(), ImplicitSolve(), etc. However, the exact list of methods that need to support a specific mode will depend on the used time-stepping method.

Definition at line 383 of file operator.hpp.

◆ SetTime()

virtual void mfem::TimeDependentOperator::SetTime ( const real_t t_)
inlinevirtual

Set the current time.

Reimplemented in mfem::electromagnetics::MagneticDiffusionEOperator.

Definition at line 360 of file operator.hpp.

◆ SUNImplicitSetup()

int mfem::TimeDependentOperator::SUNImplicitSetup ( const Vector & x,
const Vector & fx,
int jok,
int * jcur,
real_t gamma )
virtual

Setup the ODE linear system \( A(x,t) = (I - gamma J) \) or \( A = (M - gamma J) \), where \( J(x,t) = \frac{df}{dt(x,t)} \).

Parameters
[in]xThe state at which \(A(x,t)\) should be evaluated.
[in]fxThe current value of the ODE rhs function, \(f(x,t)\).
[in]jokFlag indicating if the Jacobian should be updated.
[out]jcurFlag to signal if the Jacobian was updated.
[in]gammaThe scaled time step value.

If not re-implemented, this method simply generates an error.

Presently, this method is used by SUNDIALS ODE solvers, for more details, see the SUNDIALS User Guides.

Definition at line 319 of file operator.cpp.

◆ SUNImplicitSolve()

int mfem::TimeDependentOperator::SUNImplicitSolve ( const Vector & b,
Vector & x,
real_t tol )
virtual

Solve the ODE linear system \( A x = b \) as setup by the method SUNImplicitSetup().

Parameters
[in]bThe linear system right-hand side.
[in,out]xOn input, the initial guess. On output, the solution.
[in]tolLinear solve tolerance.

If not re-implemented, this method simply generates an error.

Presently, this method is used by SUNDIALS ODE solvers, for more details, see the SUNDIALS User Guides.

Definition at line 327 of file operator.cpp.

◆ SUNMassMult()

int mfem::TimeDependentOperator::SUNMassMult ( const Vector & x,
Vector & v )
virtual

Compute the mass matrix-vector product \( v = M x \) .

Parameters
[in]xThe vector to multiply.
[out]vThe result of the matrix-vector product.

If not re-implemented, this method simply generates an error.

Presently, this method is used by SUNDIALS ARKStep integrator, for more details, see the ARKode User Guide.

Definition at line 345 of file operator.cpp.

◆ SUNMassSetup()

int mfem::TimeDependentOperator::SUNMassSetup ( )
virtual

Setup the mass matrix in the ODE system \( M y' = f(y,t) \) .

If not re-implemented, this method simply generates an error.

Presently, this method is used by SUNDIALS ARKStep integrator, for more details, see the ARKode User Guide.

Definition at line 333 of file operator.cpp.

◆ SUNMassSolve()

int mfem::TimeDependentOperator::SUNMassSolve ( const Vector & b,
Vector & x,
real_t tol )
virtual

Solve the mass matrix linear system \( M x = b \) as setup by the method SUNMassSetup().

Parameters
[in]bThe linear system right-hand side.
[in,out]xOn input, the initial guess. On output, the solution.
[in]tolLinear solve tolerance.

If not re-implemented, this method simply generates an error.

Presently, this method is used by SUNDIALS ARKStep integrator, for more details, see the ARKode User Guide.

Definition at line 339 of file operator.cpp.

Member Data Documentation

◆ eval_mode

EvalMode mfem::TimeDependentOperator::eval_mode
protected

Current evaluation mode.

Definition at line 342 of file operator.hpp.

◆ t

real_t mfem::TimeDependentOperator::t
protected

Current time.

Definition at line 340 of file operator.hpp.

◆ type

Type mfem::TimeDependentOperator::type
protected

Describes the form of the TimeDependentOperator.

Definition at line 341 of file operator.hpp.


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