MFEM v4.8.0
Finite element discretization library
Loading...
Searching...
No Matches
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.
 
void Step (Vector &x, double &t, double &dt) override
 
virtual void StepB (Vector &w, double &t, double &dt)
 Solve one adjoint time step.
 
void SetWFTolerances (EWTFunction func)
 Set multiplicative error weights.
 
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)
 
void InitAdjointSolve (int steps, int interpolation)
 Initialize Adjoint.
 
void SetMaxNStepsB (int mxstepsB)
 Set the maximum number of backward steps.
 
long GetNumSteps ()
 Get Number of Steps for ForwardSolve.
 
void EvalQuadIntegration (double t, Vector &q)
 Evaluate Quadrature.
 
void EvalQuadIntegrationB (double t, Vector &dG_dp)
 Evaluate Quadrature solution.
 
void GetForwardSolution (double tB, mfem::Vector &yy)
 Get Interpolated Forward solution y at backward integration time tB.
 
void UseMFEMLinearSolverB ()
 Set Linear Solver for the backward problem.
 
void UseSundialsLinearSolverB ()
 Use built in SUNDIALS Newton solver.
 
void SetSStolerancesB (double reltol, double abstol)
 Tolerance specification functions for the adjoint problem.
 
void SetSVtolerancesB (double reltol, Vector abstol)
 Tolerance specification functions for the adjoint problem.
 
virtual ~CVODESSolver ()
 Destroy the associated CVODES memory and SUNDIALS objects.
 
- Public Member Functions inherited from mfem::CVODESolver
 CVODESolver (int lmm)
 Construct a serial wrapper to SUNDIALS' CVODE integrator.
 
 CVODESolver (MPI_Comm comm, int lmm)
 Construct a parallel wrapper to SUNDIALS' CVODE integrator.
 
void Step (Vector &x, double &t, double &dt) override
 Integrate the ODE with CVODE using the specified step mode.
 
void UseMFEMLinearSolver ()
 Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to CVODE.
 
void UseSundialsLinearSolver ()
 Attach SUNDIALS GMRES linear solver to CVODE.
 
void SetStepMode (int itask)
 Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
 
void SetSStolerances (double reltol, double abstol)
 Set the scalar relative and scalar absolute tolerances.
 
void SetSVtolerances (double reltol, Vector abstol)
 Set the scalar relative and vector of absolute tolerances.
 
void SetRootFinder (int components, RootFunction func)
 Initialize Root Finder.
 
void SetMaxStep (double dt_max)
 Set the maximum time step.
 
void SetMaxNSteps (int steps)
 Set the maximum number of time steps.
 
long GetNumSteps ()
 Get the number of internal steps taken so far.
 
void SetMaxOrder (int max_order)
 Set the maximum method order.
 
void PrintInfo () const
 Print various CVODE statistics.
 
virtual ~CVODESolver ()
 Destroy the associated CVODE memory and SUNDIALS objects.
 
- Public Member Functions inherited from mfem::ODESolver
 ODESolver ()
 
virtual void Step (Vector &x, real_t &t, real_t &dt)=0
 Perform a time step from time t [in] to time t [out] based on the requested step size dt [in].
 
virtual void Run (Vector &x, real_t &t, real_t &dt, real_t tf)
 Perform time integration from time t [in] to time tf [in].
 
virtual int GetStateSize ()
 Returns how many State vectors the ODE requires.
 
virtual ~ODESolver ()
 
- Public Member Functions inherited from mfem::SundialsSolver
void * GetMem () const
 Access the SUNDIALS memory structure.
 
int GetFlag () const
 Returns the last flag returned by a call to a SUNDIALS function.
 

Static Public Member Functions

static int LinSysSetupB (sunrealtype t, N_Vector y, N_Vector yB, N_Vector fyB, SUNMatrix A, sunbooleantype jok, sunbooleantype *jcur, sunrealtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
 Setup the linear system A x = b.
 
static int LinSysSolveB (SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, sunrealtype tol)
 Solve the linear system A x = b.
 
- Static Public Member Functions inherited from mfem::ODESolver
static MFEM_EXPORT std::unique_ptr< ODESolverSelect (const int ode_solver_type)
 
static MFEM_EXPORT std::unique_ptr< ODESolverSelectExplicit (const int ode_solver_type)
 
static MFEM_EXPORT std::unique_ptr< ODESolverSelectImplicit (const int ode_solver_type)
 

Static Protected Member Functions

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

Protected Attributes

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

Static Protected Attributes

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

Additional Inherited Members

- Static Public Attributes inherited from mfem::ODESolver
static MFEM_EXPORT std::string ExplicitTypes
 
static MFEM_EXPORT std::string ImplicitTypes
 
static MFEM_EXPORT std::string Types
 
- Protected Types inherited from mfem::CVODESolver
typedef std::function< int(sunrealtype t, Vector y, Vector gout, CVODESolver *)> RootFunction
 Typedef for root finding functions.
 
typedef std::function< int(Vector y, Vector w, CVODESolver *)> EWTFunction
 Typedef declaration for error weight functions.
 
- 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.
 
void AllocateEmptyNVector (N_Vector &y)
 
void AllocateEmptyNVector (N_Vector &y, MPI_Comm comm)
 

Detailed Description

Definition at line 555 of file sundials.hpp.

Constructor & Destructor Documentation

◆ CVODESSolver() [1/2]

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 1003 of file sundials.cpp.

◆ CVODESSolver() [2/2]

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 1017 of file sundials.cpp.

◆ ~CVODESSolver()

mfem::CVODESSolver::~CVODESSolver ( )
virtual

Destroy the associated CVODES memory and SUNDIALS objects.

Definition at line 1355 of file sundials.cpp.

Member Function Documentation

◆ EvalQuadIntegration()

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

Evaluate Quadrature.

Definition at line 1031 of file sundials.cpp.

◆ EvalQuadIntegrationB()

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

Evaluate Quadrature solution.

Definition at line 1041 of file sundials.cpp.

◆ ewt()

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

Error control function.

Definition at line 1288 of file sundials.cpp.

◆ GetForwardSolution()

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

Get Interpolated Forward solution y at backward integration time tB.

Definition at line 1051 of file sundials.cpp.

◆ GetNumSteps()

long mfem::CVODESSolver::GetNumSteps ( )

Get Number of Steps for ForwardSolve.

◆ Init()

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 1060 of file sundials.cpp.

◆ InitAdjointSolve()

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

Initialize Adjoint.

Definition at line 1098 of file sundials.cpp.

◆ InitB()

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

Initialize the adjoint problem.

Definition at line 1065 of file sundials.cpp.

◆ InitQuadIntegration()

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

Definition at line 1110 of file sundials.cpp.

◆ InitQuadIntegrationB()

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

Initialize Quadrature Integration (Adjoint)

Definition at line 1125 of file sundials.cpp.

◆ LinSysSetupB()

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

Setup the linear system A x = b.

Definition at line 1187 of file sundials.cpp.

◆ LinSysSolveB()

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

Solve the linear system A x = b.

Definition at line 1206 of file sundials.cpp.

◆ RHSB()

int mfem::CVODESSolver::RHSB ( sunrealtype 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 1272 of file sundials.cpp.

◆ RHSQ()

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

Wrapper to compute the ODE RHS Quadrature function.

Definition at line 1245 of file sundials.cpp.

◆ RHSQB()

int mfem::CVODESSolver::RHSQB ( sunrealtype 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 1258 of file sundials.cpp.

◆ SetMaxNStepsB()

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

Set the maximum number of backward steps.

Definition at line 1104 of file sundials.cpp.

◆ SetSStolerancesB()

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 1219 of file sundials.cpp.

◆ SetSVtolerancesB()

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 1225 of file sundials.cpp.

◆ SetWFTolerances()

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

Set multiplicative error weights.

Definition at line 1237 of file sundials.cpp.

◆ Step()

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

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.

Definition at line 1299 of file sundials.cpp.

◆ StepB()

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

Solve one adjoint time step.

Definition at line 1327 of file sundials.cpp.

◆ UseMFEMLinearSolverB()

void mfem::CVODESSolver::UseMFEMLinearSolverB ( )

Set Linear Solver for the backward problem.

Definition at line 1140 of file sundials.cpp.

◆ UseSundialsLinearSolverB()

void mfem::CVODESSolver::UseSundialsLinearSolverB ( )

Use built in SUNDIALS Newton solver.

Definition at line 1172 of file sundials.cpp.

Member Data Documentation

◆ AB

SUNMatrix mfem::CVODESSolver::AB
protected

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

Definition at line 579 of file sundials.hpp.

◆ default_abs_tolB

double mfem::CVODESSolver::default_abs_tolB = 1e-9
staticconstexprprotected

Default scalar backward absolute tolerance.

Definition at line 589 of file sundials.hpp.

◆ default_abs_tolQB

double mfem::CVODESSolver::default_abs_tolQB = 1e-9
staticconstexprprotected

Default scalar backward absolute quadrature tolerance.

Definition at line 591 of file sundials.hpp.

◆ default_rel_tolB

double mfem::CVODESSolver::default_rel_tolB = 1e-4
staticconstexprprotected

Default scalar backward relative tolerance.

Definition at line 587 of file sundials.hpp.

◆ indexB

int mfem::CVODESSolver::indexB
protected

backward problem index

Definition at line 562 of file sundials.hpp.

◆ LSB

SUNLinearSolver mfem::CVODESSolver::LSB
protected

Linear solver for A.

Definition at line 580 of file sundials.hpp.

◆ ncheck

int mfem::CVODESSolver::ncheck
protected

number of checkpoints used so far

Definition at line 561 of file sundials.hpp.

◆ q

SundialsNVector* mfem::CVODESSolver::q
protected

Quadrature vector.

Definition at line 581 of file sundials.hpp.

◆ qB

SundialsNVector* mfem::CVODESSolver::qB
protected

State vector.

Definition at line 584 of file sundials.hpp.

◆ yB

SundialsNVector* mfem::CVODESSolver::yB
protected

State vector.

Definition at line 582 of file sundials.hpp.

◆ yy

SundialsNVector* mfem::CVODESSolver::yy
protected

State vector.

Definition at line 583 of file sundials.hpp.


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