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

#include <sundials.hpp>

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

Public Member Functions

 CVODESSolver (int lmm)
 
 CVODESSolver (MPI_Comm comm, int lmm)
 
void Init (TimeDependentAdjointOperator &f_)
 
void InitB (TimeDependentAdjointOperator &f_)
 Initialize the adjoint problem. More...
 
virtual void Step (Vector &x, double &t, double &dt)
 
virtual void StepB (Vector &w, double &t, double &dt)
 Solve one adjoint time step. More...
 
void SetWFTolerances (EWTFunction func)
 Set multiplicative error weights. More...
 
void InitQuadIntegration (mfem::Vector &q0, double reltolQ=1e-3, double abstolQ=1e-8)
 
void InitQuadIntegrationB (mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
 Initialize Quadrature Integration (Adjoint) More...
 
void InitAdjointSolve (int steps, int interpolation)
 Initialize Adjoint. More...
 
void SetMaxNStepsB (int mxstepsB)
 Set the maximum number of backward steps. More...
 
long GetNumSteps ()
 Get Number of Steps for ForwardSolve. More...
 
void EvalQuadIntegration (double t, Vector &q)
 Evaluate Quadrature. More...
 
void EvalQuadIntegrationB (double t, Vector &dG_dp)
 Evaluate Quadrature solution. More...
 
void GetForwardSolution (double tB, mfem::Vector &yy)
 Get Interpolated Forward solution y at backward integration time tB. More...
 
void UseMFEMLinearSolverB ()
 Set Linear Solver for the backward problem. More...
 
void UseSundialsLinearSolverB ()
 Use built in SUNDIALS Newton solver. More...
 
void SetSStolerancesB (double reltol, double abstol)
 Tolerance specification functions for the adjoint problem. More...
 
void SetSVtolerancesB (double reltol, Vector abstol)
 Tolerance specification functions for the adjoint problem. More...
 
virtual ~CVODESSolver ()
 Destroy the associated CVODES memory and SUNDIALS objects. More...
 
- Public Member Functions inherited from mfem::CVODESolver
 CVODESolver (int lmm)
 Construct a serial wrapper to SUNDIALS' CVODE integrator. More...
 
 CVODESolver (MPI_Comm comm, int lmm)
 Construct a parallel wrapper to SUNDIALS' CVODE integrator. More...
 
void Init (TimeDependentOperator &f_)
 Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults. More...
 
void UseMFEMLinearSolver ()
 Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to CVODE. More...
 
void UseSundialsLinearSolver ()
 Attach SUNDIALS GMRES linear solver to CVODE. More...
 
void SetStepMode (int itask)
 Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP. More...
 
void SetSStolerances (double reltol, double abstol)
 Set the scalar relative and scalar absolute tolerances. More...
 
void SetSVtolerances (double reltol, Vector abstol)
 Set the scalar relative and vector of absolute tolerances. More...
 
void SetRootFinder (int components, RootFunction func)
 Initialize Root Finder. More...
 
void SetMaxStep (double dt_max)
 Set the maximum time step. More...
 
void SetMaxNSteps (int steps)
 Set the maximum number of time steps. More...
 
long GetNumSteps ()
 Get the number of internal steps taken so far. More...
 
void SetMaxOrder (int max_order)
 Set the maximum method order. More...
 
void PrintInfo () const
 Print various CVODE statistics. More...
 
virtual ~CVODESolver ()
 Destroy the associated CVODE memory and SUNDIALS objects. More...
 
- Public Member Functions inherited from mfem::ODESolver
 ODESolver ()
 
virtual void Run (Vector &x, double &t, double &dt, double tf)
 Perform time integration from time t [in] to time tf [in]. More...
 
virtual int GetMaxStateSize ()
 Function for getting and setting the state vectors. More...
 
virtual int GetStateSize ()
 
virtual const VectorGetStateVector (int i)
 
virtual void GetStateVector (int i, Vector &state)
 
virtual void SetStateVector (int i, Vector &state)
 
virtual ~ODESolver ()
 
- Public Member Functions inherited from mfem::SundialsSolver
void * GetMem () const
 Access the SUNDIALS memory structure. More...
 
int GetFlag () const
 Returns the last flag returned by a call to a SUNDIALS function. More...
 

Static Public Member Functions

static int LinSysSetupB (realtype t, N_Vector y, N_Vector yB, N_Vector fyB, SUNMatrix A, booleantype jok, booleantype *jcur, realtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
 Setup the linear system A x = b. More...
 
static int LinSysSolveB (SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
 Solve the linear system A x = b. More...
 

Static Protected Member Functions

static int RHSQ (realtype t, const N_Vector y, N_Vector qdot, void *user_data)
 Wrapper to compute the ODE RHS Quadrature function. More...
 
static int RHSB (realtype t, N_Vector y, N_Vector yB, N_Vector yBdot, void *user_dataB)
 Wrapper to compute the ODE RHS backward function. More...
 
static int RHSQB (realtype t, N_Vector y, N_Vector yB, N_Vector qBdot, void *user_dataB)
 Wrapper to compute the ODE RHS Backwards Quadrature function. More...
 
static int ewt (N_Vector y, N_Vector w, void *user_data)
 Error control function. More...
 
- Static Protected Member Functions inherited from mfem::CVODESolver
static int RHS (realtype t, const N_Vector y, N_Vector ydot, void *user_data)
 Number of components in gout. More...
 
static int LinSysSetup (realtype t, N_Vector y, N_Vector fy, SUNMatrix A, booleantype jok, booleantype *jcur, realtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
 Setup the linear system \( A x = b \). More...
 
static int LinSysSolve (SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
 Solve the linear system \( A x = b \). More...
 
static int root (realtype t, N_Vector y, realtype *gout, void *user_data)
 Prototype to define root finding for CVODE. More...
 

Protected Attributes

int ncheck
 number of checkpoints used so far More...
 
int indexB
 backward problem index More...
 
SUNMatrix AB
 Linear system A = I - gamma J, M - gamma J, or J. More...
 
SUNLinearSolver LSB
 Linear solver for A. More...
 
SundialsNVectorq
 Quadrature vector. More...
 
SundialsNVectoryB
 State vector. More...
 
SundialsNVectoryy
 State vector. More...
 
SundialsNVectorqB
 State vector. More...
 
- Protected Attributes inherited from mfem::CVODESolver
int lmm_type
 Linear multistep method type. More...
 
int step_mode
 CVODE step mode (CV_NORMAL or CV_ONE_STEP). More...
 
int root_components
 
RootFunction root_func
 A class member to facilitate pointing to a user-specified root function. More...
 
EWTFunction ewt_func
 A class member to facilitate pointing to a user-specified error weight function. More...
 
- Protected Attributes inherited from mfem::ODESolver
TimeDependentOperatorf
 Pointer to the associated TimeDependentOperator. More...
 
MemoryType mem_type
 
- Protected Attributes inherited from mfem::SundialsSolver
void * sundials_mem
 SUNDIALS mem structure. More...
 
int flag
 Last flag returned from a call to SUNDIALS. More...
 
bool reinit
 Flag to signal memory reinitialization is need. More...
 
long saved_global_size
 Global vector length on last initialization. More...
 
SundialsNVectorY
 State vector. More...
 
SUNMatrix A
 
SUNMatrix M
 Mass matrix M. More...
 
SUNLinearSolver LSA
 Linear solver for A. More...
 
SUNLinearSolver LSM
 Linear solver for M. More...
 
SUNNonlinearSolver NLS
 Nonlinear solver. More...
 

Static Protected Attributes

static constexpr double default_rel_tolB = 1e-4
 Default scalar backward relative tolerance. More...
 
static constexpr double default_abs_tolB = 1e-9
 Default scalar backward absolute tolerance. More...
 
static constexpr double default_abs_tolQB = 1e-9
 Default scalar backward absolute quadrature tolerance. More...
 
- Static Protected Attributes inherited from mfem::SundialsSolver
static constexpr double default_rel_tol = 1e-4
 Default scalar relative tolerance. More...
 
static constexpr double default_abs_tol = 1e-9
 Default scalar absolute tolerance. More...
 

Additional Inherited Members

- Protected Types inherited from mfem::CVODESolver
typedef std::function< int(realtype
t, Vector y, Vector gout,
CVODESolver *)> 
RootFunction
 Typedef for root finding functions. More...
 
typedef std::function< int(Vector
y, Vector w, CVODESolver *)> 
EWTFunction
 Typedef declaration for error weight functions. More...
 
- Protected Member Functions inherited from mfem::SundialsSolver
bool Parallel () const
 
bool Parallel () const
 
 SundialsSolver ()
 Protected constructor: objects of this type should be constructed only as part of a derived class. More...
 
void AllocateEmptyNVector (N_Vector &y)
 
void AllocateEmptyNVector (N_Vector &y, MPI_Comm comm)
 

Detailed Description

Definition at line 385 of file sundials.hpp.

Constructor & Destructor Documentation

mfem::CVODESSolver::CVODESSolver ( int  lmm)

Construct a serial wrapper to SUNDIALS' CVODE integrator.

Parameters
[in]lmmSpecifies the linear multistep method, the options are: CV_ADAMS - implicit methods for non-stiff systems CV_BDF - implicit methods for stiff systems

Definition at line 800 of file sundials.cpp.

mfem::CVODESSolver::CVODESSolver ( MPI_Comm  comm,
int  lmm 
)

Construct a parallel wrapper to SUNDIALS' CVODE integrator.

Parameters
[in]commThe MPI communicator used to partition the ODE system
[in]lmmSpecifies the linear multistep method, the options are: CV_ADAMS - implicit methods for non-stiff systems CV_BDF - implicit methods for stiff systems

Definition at line 814 of file sundials.cpp.

mfem::CVODESSolver::~CVODESSolver ( )
virtual

Destroy the associated CVODES memory and SUNDIALS objects.

Definition at line 1152 of file sundials.cpp.

Member Function Documentation

void mfem::CVODESSolver::EvalQuadIntegration ( double  t,
Vector q 
)

Evaluate Quadrature.

Definition at line 828 of file sundials.cpp.

void mfem::CVODESSolver::EvalQuadIntegrationB ( double  t,
Vector dG_dp 
)

Evaluate Quadrature solution.

Definition at line 838 of file sundials.cpp.

int mfem::CVODESSolver::ewt ( N_Vector  y,
N_Vector  w,
void *  user_data 
)
staticprotected

Error control function.

Definition at line 1085 of file sundials.cpp.

void mfem::CVODESSolver::GetForwardSolution ( double  tB,
mfem::Vector yy 
)

Get Interpolated Forward solution y at backward integration time tB.

Definition at line 848 of file sundials.cpp.

long mfem::CVODESSolver::GetNumSteps ( )

Get Number of Steps for ForwardSolve.

void mfem::CVODESSolver::Init ( TimeDependentAdjointOperator f_)

Initialize CVODE: Calls CVodeInit() and sets some defaults. We define this to force the time dependent operator to be a TimeDependenAdjointOperator.

Parameters
[in]f_the TimeDependentAdjointOperator that defines the ODE system
Note
All other methods must be called after Init().

Definition at line 857 of file sundials.cpp.

void mfem::CVODESSolver::InitAdjointSolve ( int  steps,
int  interpolation 
)

Initialize Adjoint.

Definition at line 895 of file sundials.cpp.

void mfem::CVODESSolver::InitB ( TimeDependentAdjointOperator f_)

Initialize the adjoint problem.

Definition at line 862 of file sundials.cpp.

void mfem::CVODESSolver::InitQuadIntegration ( mfem::Vector q0,
double  reltolQ = 1e-3,
double  abstolQ = 1e-8 
)

Definition at line 907 of file sundials.cpp.

void mfem::CVODESSolver::InitQuadIntegrationB ( mfem::Vector qB0,
double  reltolQB = 1e-3,
double  abstolQB = 1e-8 
)

Initialize Quadrature Integration (Adjoint)

Definition at line 922 of file sundials.cpp.

int mfem::CVODESSolver::LinSysSetupB ( realtype  t,
N_Vector  y,
N_Vector  yB,
N_Vector  fyB,
SUNMatrix  A,
booleantype  jok,
booleantype *  jcur,
realtype  gamma,
void *  user_data,
N_Vector  tmp1,
N_Vector  tmp2,
N_Vector  tmp3 
)
static

Setup the linear system A x = b.

Definition at line 984 of file sundials.cpp.

int mfem::CVODESSolver::LinSysSolveB ( SUNLinearSolver  LS,
SUNMatrix  A,
N_Vector  x,
N_Vector  b,
realtype  tol 
)
static

Solve the linear system A x = b.

Definition at line 1003 of file sundials.cpp.

int mfem::CVODESSolver::RHSB ( realtype  t,
N_Vector  y,
N_Vector  yB,
N_Vector  yBdot,
void *  user_dataB 
)
staticprotected

Wrapper to compute the ODE RHS backward function.

Definition at line 1069 of file sundials.cpp.

int mfem::CVODESSolver::RHSQ ( realtype  t,
const N_Vector  y,
N_Vector  qdot,
void *  user_data 
)
staticprotected

Wrapper to compute the ODE RHS Quadrature function.

Definition at line 1042 of file sundials.cpp.

int mfem::CVODESSolver::RHSQB ( realtype  t,
N_Vector  y,
N_Vector  yB,
N_Vector  qBdot,
void *  user_dataB 
)
staticprotected

Wrapper to compute the ODE RHS Backwards Quadrature function.

Definition at line 1055 of file sundials.cpp.

void mfem::CVODESSolver::SetMaxNStepsB ( int  mxstepsB)

Set the maximum number of backward steps.

Definition at line 901 of file sundials.cpp.

void mfem::CVODESSolver::SetSStolerancesB ( double  reltol,
double  abstol 
)

Tolerance specification functions for the adjoint problem.

It should be called after InitB() is called.

Parameters
[in]reltolthe scalar relative error tolerance.
[in]abstolthe scalar absolute error tolerance.

Definition at line 1016 of file sundials.cpp.

void mfem::CVODESSolver::SetSVtolerancesB ( double  reltol,
Vector  abstol 
)

Tolerance specification functions for the adjoint problem.

It should be called after InitB() is called.

Parameters
[in]reltolthe scalar relative error tolerance
[in]abstolthe vector of absolute error tolerances

Definition at line 1022 of file sundials.cpp.

void mfem::CVODESSolver::SetWFTolerances ( EWTFunction  func)

Set multiplicative error weights.

Definition at line 1034 of file sundials.cpp.

void mfem::CVODESSolver::Step ( Vector x,
double &  t,
double &  dt 
)
virtual

Integrate the ODE with CVODE using the specified step mode.

Parameters
[out]xSolution vector at the requested output time x=x(t).
[in,out]tOn output, the output time reached.
[in,out]dtOn output, the last time step taken.
Note
On input, the values of t and dt are used to compute desired output time for the integration, tout = t + dt.

Reimplemented from mfem::CVODESolver.

Definition at line 1096 of file sundials.cpp.

void mfem::CVODESSolver::StepB ( Vector w,
double &  t,
double &  dt 
)
virtual

Solve one adjoint time step.

Definition at line 1124 of file sundials.cpp.

void mfem::CVODESSolver::UseMFEMLinearSolverB ( )

Set Linear Solver for the backward problem.

Definition at line 937 of file sundials.cpp.

void mfem::CVODESSolver::UseSundialsLinearSolverB ( )

Use built in SUNDIALS Newton solver.

Definition at line 969 of file sundials.cpp.

Member Data Documentation

SUNMatrix mfem::CVODESSolver::AB
protected

Linear system A = I - gamma J, M - gamma J, or J.

Definition at line 408 of file sundials.hpp.

constexpr double mfem::CVODESSolver::default_abs_tolB = 1e-9
staticprotected

Default scalar backward absolute tolerance.

Definition at line 418 of file sundials.hpp.

constexpr double mfem::CVODESSolver::default_abs_tolQB = 1e-9
staticprotected

Default scalar backward absolute quadrature tolerance.

Definition at line 420 of file sundials.hpp.

constexpr double mfem::CVODESSolver::default_rel_tolB = 1e-4
staticprotected

Default scalar backward relative tolerance.

Definition at line 416 of file sundials.hpp.

int mfem::CVODESSolver::indexB
protected

backward problem index

Definition at line 392 of file sundials.hpp.

SUNLinearSolver mfem::CVODESSolver::LSB
protected

Linear solver for A.

Definition at line 409 of file sundials.hpp.

int mfem::CVODESSolver::ncheck
protected

number of checkpoints used so far

Definition at line 391 of file sundials.hpp.

SundialsNVector* mfem::CVODESSolver::q
protected

Quadrature vector.

Definition at line 410 of file sundials.hpp.

SundialsNVector* mfem::CVODESSolver::qB
protected

State vector.

Definition at line 413 of file sundials.hpp.

SundialsNVector* mfem::CVODESSolver::yB
protected

State vector.

Definition at line 411 of file sundials.hpp.

SundialsNVector* mfem::CVODESSolver::yy
protected

State vector.

Definition at line 412 of file sundials.hpp.


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