MFEM v4.7.0
Finite element discretization library
|
Base abstract class for second order time dependent operators. More...
#include <operator.hpp>
Public Member Functions | |
SecondOrderTimeDependentOperator (int n=0, real_t t_=0.0, Type type_=EXPLICIT) | |
Construct a "square" SecondOrderTimeDependentOperator y = f(x,dxdt,t), where x, dxdt and y have the same dimension n. | |
SecondOrderTimeDependentOperator (int h, int w, real_t t_=0.0, Type type_=EXPLICIT) | |
Construct a SecondOrderTimeDependentOperator y = f(x,dxdt,t), where x, dxdt and y have the same dimension n. | |
virtual void | Mult (const Vector &x, const Vector &dxdt, Vector &y) const |
Perform the action of the operator: y = k = f(x,@ dxdt, t), where k solves the algebraic equation F(x,@ dxdt, k, t) = G(x,@ dxdt, t) and t is the current time. | |
virtual void | ImplicitSolve (const real_t fac0, const real_t fac1, const Vector &x, const Vector &dxdt, Vector &k) |
Solve the equation: k = f(x + fac0 k, dxdt + fac1 k, t), for the unknown k at the current time t. | |
virtual | ~SecondOrderTimeDependentOperator () |
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. | |
Public Member Functions inherited from mfem::TimeDependentOperator | |
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 Operator & | GetImplicitGradient (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 Operator & | GetExplicitGradient (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 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. | |
Additional Inherited Members | |
Public Types inherited from mfem::TimeDependentOperator | |
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... | |
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() | |
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". | |
Protected Attributes inherited from mfem::TimeDependentOperator | |
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. | |
Base abstract class for second order time dependent operators.
Operator of the form: (x,dxdt,t) -> f(x,dxdt,t), where k = f(x,dxdt,t) generally solves the algebraic equation F(x,dxdt,k,t) = G(x,dxdt,t). The functions F and G represent the_implicit_ and explicit parts of the operator, respectively. For explicit operators, F(x,dxdt,k,t) = k, so f(x,dxdt,t) = G(x,dxdt,t).
Definition at line 634 of file operator.hpp.
|
inlineexplicit |
Construct a "square" SecondOrderTimeDependentOperator y = f(x,dxdt,t), where x, dxdt and y have the same dimension n.
Definition at line 639 of file operator.hpp.
|
inline |
Construct a SecondOrderTimeDependentOperator y = f(x,dxdt,t), where x, dxdt and y have the same dimension n.
Definition at line 645 of file operator.hpp.
|
inlinevirtual |
Definition at line 677 of file operator.hpp.
|
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 from mfem::TimeDependentOperator.
Definition at line 421 of file operator.cpp.
|
virtual |
Solve the equation: k = f(x + fac0 k, dxdt + fac1 k, t), for the unknown k at the current time t.
For general F and G, the equation for k becomes: F(x + fac0 k, dxdt + fac1 k, t) = G(x + fac0 k, dxdt + fac1 k, t).
The input vectors x and dxdt corresponds to time index (or cycle) n, while the currently set time, t, and the result vector k correspond to time index n+1.
This method allows for the abstract implementation of some time integration methods.
If not re-implemented, this method simply generates an error.
Definition at line 359 of file operator.cpp.
|
virtual |
Perform the action of the operator: y = k = f(x,@ dxdt, t), where k solves the algebraic equation F(x,@ dxdt, k, t) = G(x,@ dxdt, t) and t is the current time.
Definition at line 352 of file operator.cpp.
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.
Reimplemented from mfem::TimeDependentOperator.
Definition at line 403 of file operator.cpp.