![]() |
MFEM v4.8.0
Finite element discretization library
|
Base abstract class for first order time dependent operators. More...
#include <operator.hpp>
Public Types | |
enum | Type { EXPLICIT , IMPLICIT , HOMOGENEOUS } |
Enum used to describe the form of the time-dependent operator. More... | |
enum | EvalMode { NORMAL , ADDITIVE_TERM_1 , ADDITIVE_TERM_2 } |
Evaluation mode. See SetEvalMode() for details. More... | |
![]() | |
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 (u,t) -> k(u,t), where u and k have the same dimension n. | |
TimeDependentOperator (int h, int w, double t_=0.0, Type type_=EXPLICIT) | |
Construct a TimeDependentOperator (u,t) -> k(u,t), where u and k 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 &u, Vector &v) const |
Perform the action of the explicit part of the operator, G: v = G(u, t) where t is the current time. | |
virtual void | ImplicitMult (const Vector &u, const Vector &k, Vector &v) const |
Perform the action of the implicit part of the operator, F: v = F(u, k, t) where t is the current time. | |
void | Mult (const Vector &u, Vector &k) const override |
Perform the action of the operator (u,t) -> k(u,t) where t is the current time set by SetTime() and k satisfies F(u, k, t) = G(u, t). | |
virtual void | ImplicitSolve (const real_t gamma, const Vector &u, Vector &k) |
Solve for the unknown k, at the current time t, the following equation: F(u + gamma k, k, t) = G(u + gamma k, t). | |
virtual Operator & | GetImplicitGradient (const Vector &u, const Vector &k, real_t shift) const |
Return an Operator representing (dF/dk shift + dF/du) at the given u, k, and the currently set time. | |
virtual Operator & | GetExplicitGradient (const Vector &u) const |
Return an Operator representing dG/du at the given point u and the currently set time. | |
virtual int | SUNImplicitSetup (const Vector &y, const Vector &v, int jok, int *jcur, real_t gamma) |
Setup a linear system as needed by some SUNDIALS ODE solvers to perform a similar action to ImplicitSolve, i.e., solve for k, at the current time t, in F(u + gamma k, k, t) = G(u + gamma k, t). | |
virtual int | SUNImplicitSolve (const Vector &r, Vector &dk, real_t tol) |
Solve the ODE linear system A dk = r , where A and r are defined by the method SUNImplicitSetup(). | |
virtual int | SUNMassSetup () |
Setup the mass matrix in the ODE system \( M \frac{dy}{dt} = g(y,t) \) . | |
virtual int | SUNMassSolve (const Vector &b, Vector &x, real_t tol) |
Solve the mass matrix linear system M x = b, where M is defined by the method SUNMassSetup(). | |
virtual int | SUNMassMult (const Vector &x, Vector &v) |
Compute the mass matrix-vector product v = M x, where M is defined by the method SUNMassSetup(). | |
virtual | ~TimeDependentOperator () |
![]() | |
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 Operator & | GetGradient (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 Operator * | GetProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. | |
virtual const Operator * | GetRestriction () const |
Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. | |
virtual const Operator * | GetOutputProlongation () const |
Prolongation operator from linear algebra (linear system) vectors, to output vectors for the operator. NULL means identity. | |
virtual const Operator * | GetOutputRestrictionTranspose () const |
Transpose of GetOutputRestriction, directly available in this form to facilitate matrix-free RAP-type operators. | |
virtual const Operator * | GetOutputRestriction () 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, see the documentation of Type. | |
EvalMode | eval_mode |
Current evaluation mode. | |
![]() | |
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 | |
![]() | |
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() | |
Operator * | SetupRAP (const Operator *Pi, const Operator *Po) |
Returns RAP Operator of this, using input/output Prolongation matrices Pi corresponds to "P", Po corresponds to "Rt". | |
Base abstract class for first order time dependent operators.
Operator of the form: (u,t) -> k(u,t), where k generally solves the algebraic equation F(u,k,t) = G(u,t). The functions F and G represent the implicit and explicit parts of the operator, respectively.
A common use for this class is representing a differential algebraic equation of the form \( F(y,\frac{dy}{dt},t) = G(y,t) \).
For example, consider an ordinary differential equation of the form \( M \frac{dy}{dt} = g(y,t) \). There are various ways of expressing this ODE as a TimeDependentOperator depending on the choices for F and G. Here are some common choices:
Note that depending on the ODE solver, some of the above choices may be preferable to the others.
Definition at line 331 of file operator.hpp.
Evaluation mode. See SetEvalMode() for details.
Definition at line 360 of file operator.hpp.
Enum used to describe the form of the time-dependent operator.
The type should be set by classes derived from TimeDependentOperator to describe the form, in terms of the functions F and G, used by the specific derived class. This information can be queried by classes or functions (like time stepping algorithms) to make choices about the algorithm to use, or to ensure that the TimeDependentOperator uses the form expected by the class/function.
For example, assume that a derived class is implementing the ODE \(M \frac{dy}{dt} = g(y,t)\) and chooses to define \(F(u,k,t) = M k\) and \(G(u,t) = g(u,t)\). Then it cannot use type EXPLICIT, unless \(M = I\), or type HOMOGENEOUS, unless \(g(u,t) = 0\). If, on the other hand, the derived class chooses to define \(F(u,k,t) = k\) and \(G(u,t) = M^{-1} g(y,t)\), then the natural choice is to set the type to EXPLICIT, even though setting it to IMPLICIT is also not wrong – doing so will simply fail to inform methods that query this information that it uses a more specific implementation, EXPLICIT, that may allow the use of algorithms that support only the EXPLICIT type.
Enumerator | |
---|---|
EXPLICIT | This type assumes F(u,k,t) = k. |
IMPLICIT | This is the most general type, no assumptions on F and G. |
HOMOGENEOUS | This type assumes that G(u,t) = 0. |
Definition at line 352 of file operator.hpp.
|
inlineexplicit |
Construct a "square" TimeDependentOperator (u,t) -> k(u,t), where u and k have the same dimension n.
Definition at line 381 of file operator.hpp.
|
inline |
Construct a TimeDependentOperator (u,t) -> k(u,t), where u and k have dimensions w and h, respectively.
Definition at line 387 of file operator.hpp.
|
inlinevirtual |
Definition at line 596 of file operator.hpp.
Perform the action of the explicit part of the operator, G: v = G(u, t) where t is the current time.
Presently, this method is used by some PETSc ODE solvers and the SUNDIALS ARKStep integrator, for more details, see either the PETSc Manual or the ARKode User Guide, respectively.
Definition at line 282 of file operator.cpp.
|
inline |
Return the current evaluation mode. See SetEvalMode() for details.
Definition at line 404 of file operator.hpp.
Return an Operator representing dG/du at the given point u 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.
|
virtual |
Return an Operator representing (dF/dk shift + dF/du) at the given u, 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.
|
inlinevirtual |
Read the currently set time.
Definition at line 391 of file operator.hpp.
|
virtual |
Perform the action of the implicit part of the operator, F: v = F(u, 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.
|
virtual |
Solve for the unknown k, at the current time t, the following equation: F(u + gamma k, k, t) = G(u + gamma k, t).
For solving an ordinary differential equation of the form \( M \frac{dy}{dt} = g(y,t) \), recall that F and G can be defined in various ways, e.g.:
Regardless of the choice of F and G, this function should solve for k in M k = g(u + gamma k, t).
To see how k can be useful, consider the backward Euler method defined by \( y(t + \Delta t) = y(t) + \Delta t k_0 \) where \( M k_0 = g \big( y(t) + \Delta t k_0, t + \Delta t \big) \). A backward Euler integrator can use k from this function for \(k_0\), with the call using u set to \( y(t) \), gamma set to \( \Delta t\), and time set to \(t + \Delta t\). See class BackwardEulerSolver.
Generalizing further, consider a diagonally implicit Runge-Kutta (DIRK) method defined by \( y(t + \Delta t) = y(t) + \Delta t \sum_{i=1}^s b_i k_i \) where \( M k_i = g \big( y(t) + \Delta t \sum_{j=1}^i a_{ij} k_j, t + c_i \Delta t \big) \). A DIRK integrator can use k from this function, with u set to \( y(t) + \Delta t \sum_{j=1}^{i-1} a_{ij} k_j \) and gamma set to \( a_{ii} \Delta t \), for \( k_i \). For example, see class SDIRK33Solver.
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.
|
inline |
Definition at line 397 of file operator.hpp.
|
inline |
True if type is HOMOGENEOUS.
Definition at line 401 of file operator.hpp.
|
inline |
True if type is IMPLICIT or HOMOGENEOUS.
Definition at line 399 of file operator.hpp.
Perform the action of the operator (u,t) -> k(u,t) where t is the current time set by SetTime() and k satisfies F(u, k, t) = G(u, t).
For solving an ordinary differential equation of the form \( M \frac{dy}{dt} = g(y,t) \), recall that F and G can be defined in various ways, e.g.:
Regardless of the choice of F and G, this function should always compute k = inv(M) g(u, t).
Implements mfem::Operator.
Definition at line 293 of file operator.cpp.
|
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: k(u,t) = k1(u,t) + k2(u,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 417 of file operator.hpp.
|
inlinevirtual |
Set the current time.
Reimplemented in mfem::electromagnetics::MagneticDiffusionEOperator.
Definition at line 394 of file operator.hpp.
|
virtual |
Setup a linear system as needed by some SUNDIALS ODE solvers to perform a similar action to ImplicitSolve, i.e., solve for k, at the current time t, in F(u + gamma k, k, t) = G(u + gamma k, t).
The SUNDIALS ODE solvers iteratively solve for k, as knew = kold + dk. The linear system here is for dk, obtained by linearizing the nonlinear system F(u + gamma knew, knew, t) = G(u + gamma knew, t) about dk = 0: F(u + gamma (kold + dk), kold + dk, t) = G(u + gamma (kold + dk), t) => [dF/dk + gamma (dF/du - dG/du)] dk = G - F + O(dk^2) In other words, the linear system to be setup here is A dk = r, where A = [dF/dk + gamma (dF/du - dG/du)] and r = G - F.
For solving an ordinary differential equation of the form \( M \frac{dy}{dt} = g(y,t) \), recall that F and G can be defined as one of the following:
This function performs setup to solve \( A dk = r \) where A is either
with J = dg/dy (or a reasonable approximation thereof).
[in] | y | The state at which A(y,t) should be evaluated. |
[in] | v | The value of inv(M) g(y,t) for 1 or g(y,t) for 2 & 3. |
[in] | jok | Flag indicating if the Jacobian should be updated. |
[out] | jcur | Flag to signal if the Jacobian was updated. |
[in] | gamma | The 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.
|
virtual |
Solve the ODE linear system A dk = r , where A and r are defined by the method SUNImplicitSetup().
For solving an ordinary differential equation of the form \( M \frac{dy}{dt} = g(y,t) \), recall that F and G can be defined as one of the following:
[in] | r | inv(M) g(y,t) - k for 1 or g(y,t) - M k for 2 & 3. |
[in,out] | dk | On input, the initial guess. On output, the solution. |
[in] | tol | Linear 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.
Compute the mass matrix-vector product v = M x, where M is defined by the method SUNMassSetup().
[in] | x | The vector to multiply. |
[out] | v | The 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.
|
virtual |
Setup the mass matrix in the ODE system \( M \frac{dy}{dt} = g(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.
Solve the mass matrix linear system M x = b, where M is defined by the method SUNMassSetup().
[in] | b | The linear system right-hand side. |
[in,out] | x | On input, the initial guess. On output, the solution. |
[in] | tol | Linear 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.
|
protected |
Current evaluation mode.
Definition at line 376 of file operator.hpp.
|
protected |
Current time.
Definition at line 373 of file operator.hpp.
|
protected |
Describes the form of the TimeDependentOperator, see the documentation of Type.
Definition at line 374 of file operator.hpp.