MFEM  v4.2.0
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
Public Member Functions | Protected Attributes | List of all members
mfem::TimeDependentAdjointOperator Class Referenceabstract

#include <operator.hpp>

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

Public Member Functions

 TimeDependentAdjointOperator (int dim, int adjdim, double t=0., Type type=EXPLICIT)
 The TimedependentAdjointOperator extends the TimeDependentOperator class to use features in SUNDIALS CVODESSolver for computing quadratures and solving adjoint problems. More...
 
virtual ~TimeDependentAdjointOperator ()
 Destructor. More...
 
virtual void QuadratureIntegration (const Vector &y, Vector &qdot) const
 Provide the operator integration of a quadrature equation. More...
 
virtual void AdjointRateMult (const Vector &y, Vector &yB, Vector &yBdot) const =0
 Perform the action of the operator: yBdot = k = f(y,@2 yB, t), where. More...
 
virtual void QuadratureSensitivityMult (const Vector &y, const Vector &yB, Vector &qBdot) const
 Provides the sensitivity of the quadrature w.r.t to primal and adjoint solutions. More...
 
virtual int SUNImplicitSetupB (const double t, const Vector &x, const Vector &xB, const Vector &fxB, int jokB, int *jcurB, double gammaB)
 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)} \). More...
 
virtual int SUNImplicitSolveB (Vector &x, const Vector &b, double tol)
 Solve the ODE linear system \( A(x,xB,t) xB = b \) as setup by the method SUNImplicitSetup(). More...
 
int GetAdjointHeight ()
 Returns the size of the adjoint problem state space. More...
 
- Public Member Functions inherited from mfem::TimeDependentOperator
 TimeDependentOperator (int n=0, double t_=0.0, Type type_=EXPLICIT)
 Construct a "square" TimeDependentOperator y = f(x,t), where x and y have the same dimension n. More...
 
 TimeDependentOperator (int h, int w, double t_=0.0, Type type_=EXPLICIT)
 Construct a TimeDependentOperator y = f(x,t), where x and y have dimensions w and h, respectively. More...
 
virtual double GetTime () const
 Read the currently set time. More...
 
virtual void SetTime (const double _t)
 Set the current time. More...
 
bool isExplicit () const
 True if type is EXPLICIT. More...
 
bool isImplicit () const
 True if type is IMPLICIT or HOMOGENEOUS. More...
 
bool isHomogeneous () const
 True if type is HOMOGENEOUS. More...
 
EvalMode GetEvalMode () const
 Return the current evaluation mode. See SetEvalMode() for details. More...
 
virtual void SetEvalMode (const EvalMode new_eval_mode)
 Set the evaluation mode of the time-dependent operator. More...
 
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. More...
 
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. More...
 
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. More...
 
virtual void ImplicitSolve (const double dt, const Vector &x, Vector &k)
 Solve the equation: k = f(x + dt k, t), for the unknown k at the current time t. More...
 
virtual OperatorGetImplicitGradient (const Vector &x, const Vector &k, double shift) const
 Return an Operator representing (dF/dk shift + dF/dx) at the given x, k, and the currently set time. More...
 
virtual OperatorGetExplicitGradient (const Vector &x) const
 Return an Operator representing dG/dx at the given point x and the currently set time. More...
 
virtual int SUNImplicitSetup (const Vector &x, const Vector &fx, int jok, int *jcur, double 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)} \). More...
 
virtual int SUNImplicitSolve (const Vector &b, Vector &x, double tol)
 Solve the ODE linear system \( A x = b \) as setup by the method SUNImplicitSetup(). More...
 
virtual int SUNMassSetup ()
 Setup the mass matrix in the ODE system \( M y' = f(y,t) \) . More...
 
virtual int SUNMassSolve (const Vector &b, Vector &x, double tol)
 Solve the mass matrix linear system \( M x = b \) as setup by the method SUNMassSetup(). More...
 
virtual int SUNMassMult (const Vector &x, Vector &v)
 Compute the mass matrix-vector product \( v = M x \) . More...
 
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. 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 OperatorGetGradient (const Vector &x) const
 Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate an error. More...
 
virtual const OperatorGetProlongation () const
 Prolongation operator from linear algebra (linear system) vectors, to input vectors for the operator. NULL means identity. More...
 
virtual const OperatorGetRestriction () const
 Restriction operator from input vectors for the operator to linear algebra (linear system) vectors. NULL means identity. 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 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=0, int m=0) const
 Prints operator with input size n and output size m in Matlab format. More...
 
virtual ~Operator ()
 Virtual destructor. More...
 
Type GetType () const
 Return the type ID of the Operator class. More...
 

Protected Attributes

int adjoint_height
 
- Protected Attributes inherited from mfem::TimeDependentOperator
double t
 Current time. More...
 
Type type
 Describes the form of the TimeDependentOperator. More...
 
EvalMode eval_mode
 Current evaluation mode. 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::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
}
 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() 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, taking in input/output Prolongation matrices. More...
 

Detailed Description

TimeDependentAdjointOperator is a TimeDependentOperator with Adjoint rate equations to be used with CVODESSolver.

Definition at line 456 of file operator.hpp.

Constructor & Destructor Documentation

mfem::TimeDependentAdjointOperator::TimeDependentAdjointOperator ( int  dim,
int  adjdim,
double  t = 0.,
Type  type = EXPLICIT 
)
inline

The TimedependentAdjointOperator extends the TimeDependentOperator class to use features in SUNDIALS CVODESSolver for computing quadratures and solving adjoint problems.

To solve adjoint problems one needs to implement the AdjointRateMult method to tell CVODES what the adjoint rate equation is.

QuadratureIntegration (optional) can be used to compute values over the forward problem

QuadratureSensitivityMult (optional) can be used to find the sensitivity of the quadrature using the adjoint solution in part.

SUNImplicitSetupB (optional) can be used to setup custom solvers for the newton solve for the adjoint problem.

SUNImplicitSolveB (optional) actually uses the solvers from SUNImplicitSetupB to solve the adjoint problem.

See SUNDIALS user manuals for specifics.

Parameters
[in]dimDimension of the forward operator
[in]adjdimDimension of the adjoint operator. Typically it is the same size as dim. However, SUNDIALS allows users to specify the size if one wants to perform custom operations.
[in]tStarting time to set
[in]typeThe TimeDependentOperator type

Definition at line 489 of file operator.hpp.

virtual mfem::TimeDependentAdjointOperator::~TimeDependentAdjointOperator ( )
inlinevirtual

Destructor.

Definition at line 496 of file operator.hpp.

Member Function Documentation

virtual void mfem::TimeDependentAdjointOperator::AdjointRateMult ( const Vector y,
Vector yB,
Vector yBdot 
) const
pure virtual

Perform the action of the operator: yBdot = k = f(y,@2 yB, t), where.

Parameters
[in]yThe primal solution at time t
[in]yBThe adjoint solution at time t
[out]yBdotthe rate at time t
int mfem::TimeDependentAdjointOperator::GetAdjointHeight ( )
inline

Returns the size of the adjoint problem state space.

Definition at line 572 of file operator.hpp.

virtual void mfem::TimeDependentAdjointOperator::QuadratureIntegration ( const Vector y,
Vector qdot 
) const
inlinevirtual

Provide the operator integration of a quadrature equation.

Parameters
[in]yThe current value at time t
[out]qdotThe current quadrature rate value at t

Definition at line 504 of file operator.hpp.

virtual void mfem::TimeDependentAdjointOperator::QuadratureSensitivityMult ( const Vector y,
const Vector yB,
Vector qBdot 
) const
inlinevirtual

Provides the sensitivity of the quadrature w.r.t to primal and adjoint solutions.

Parameters
[in]ythe value of the primal solution at time t
[in]yBthe value of the adjoint solution at time t
[out]qBdotthe value of the sensitivity of the quadrature rate at time t

Definition at line 525 of file operator.hpp.

virtual int mfem::TimeDependentAdjointOperator::SUNImplicitSetupB ( const double  t,
const Vector x,
const Vector xB,
const Vector fxB,
int  jokB,
int *  jcurB,
double  gammaB 
)
inlinevirtual

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]tThe current time
[in]xThe state at which \(A(x,xB,t)\) should be evaluated.
[in]xBThe state at which \(A(x,xB,t)\) should be evaluated.
[in]fxBThe current value of the ODE rhs function, \(f(x,t)\).
[in]jokBFlag indicating if the Jacobian should be updated.
[out]jcurBFlag to signal if the Jacobian was updated.
[in]gammaBThe 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 544 of file operator.hpp.

virtual int mfem::TimeDependentAdjointOperator::SUNImplicitSolveB ( Vector x,
const Vector b,
double  tol 
)
inlinevirtual

Solve the ODE linear system \( A(x,xB,t) xB = 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 564 of file operator.hpp.

Member Data Documentation

int mfem::TimeDependentAdjointOperator::adjoint_height
protected

Definition at line 575 of file operator.hpp.


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