15 #include "../config/config.hpp"
17 #ifdef MFEM_USE_SUNDIALS
26 #include <cvode/cvode.h>
27 #include <arkode/arkode.h>
28 #include <kinsol/kinsol.h>
72 virtual int InitSystem(
void *sundials_mem) = 0;
73 virtual int SetupSystem(
void *sundials_mem,
int conv_fail,
75 int &jac_cur,
Vector &v_temp1,
79 virtual int FreeSystem(
void *sundials_mem) = 0;
96 {
return (
y->ops->nvgetvectorid != N_VGetVectorID_Serial); }
106 static int ODEMult(realtype t,
const N_Vector
y,
107 N_Vector ydot,
void *td_oper);
201 virtual void Step(
Vector &x,
double &t,
double &dt);
294 virtual void Step(
Vector &x,
double &t,
double &dt);
323 static int Mult(
const N_Vector u, N_Vector fu,
void *user_data);
326 static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u,
327 booleantype *new_u,
void *user_data);
331 static int LinSysSolve(KINMemRec *kin_mem, N_Vector x, N_Vector b,
332 realtype *sJpnorm, realtype *sFdotJp);
342 KinSolver(
int strategy,
bool oper_grad =
true);
352 KinSolver(MPI_Comm comm,
int strategy,
bool oper_grad =
true);
409 #endif // MFEM_USE_SUNDIALS
411 #endif // MFEM_SUNDIALS
void SetFixedStep(double dt)
Use a fixed time step size, instead of performing any form of temporal adaptivity.
virtual void Init(TimeDependentOperator &f_)
Set the ODE right-hand-side operator.
virtual int FreeSystem(void *sundials_mem)=0
virtual ~SundialsODELinearSolver()
static const double default_abs_tol
virtual void Init(TimeDependentOperator &f_)
Set the ODE right-hand-side operator.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Base abstract class for time dependent operators.
virtual void Step(Vector &x, double &t, double &dt)
Use CVODE to integrate over [t, t + dt], with the specified step mode.
SundialsODELinearSolver()
static int ODEMult(realtype t, const N_Vector y, N_Vector ydot, void *td_oper)
Callback function used in CVODESolver and ARKODESolver.
virtual ~KinSolver()
Destroy the associated KINSOL memory.
void SetStepMode(int itask)
ARKode supports two modes, specified by itask: ARK_NORMAL (default) and ARK_ONE_STEP.
virtual ~CVODESolver()
Destroy the associated CVODE memory.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
int GetFlag() const
Return the flag returned by the last call to a SUNDIALS function.
void SetLinearSolver(SundialsODELinearSolver &ls_spec)
void * sundials_mem
Pointer to the SUNDIALS mem object.
TimeDependentOperator * GetTimeDependentOperator(void *sundials_mem)
Get the TimeDependentOperator associated with sundials_mem.
Type
Types of ARKODE solvers.
void SetSStolerances(double reltol, double abstol)
Specify the scalar relative and scalar absolute tolerances.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for explicit RK method.
virtual int InitSystem(void *sundials_mem)=0
virtual int SetupSystem(void *sundials_mem, int conv_fail, const Vector &y_pred, const Vector &f_pred, int &jac_cur, Vector &v_temp1, Vector &v_temp2, Vector &v_temp3)=0
ARKODESolver(Type type=EXPLICIT)
Construct a serial ARKODESolver, a wrapper for SUNDIALS' ARKODE solver.
enum mfem::SundialsODELinearSolver::@0 type
Is CVODE or ARKODE using this object?
void * SundialsMem() const
Access the underlying SUNDIALS object.
Wrapper for SUNDIALS' KINSOL library – Nonlinear solvers.
Wrapper for SUNDIALS' CVODE library – Multi-step time integration.
static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u, booleantype *new_u, void *user_data)
void SetMaxStep(double dt_max)
Set the maximum time step of the linear multistep method.
Newton's method for solving F(x)=b for a given operator F.
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system. This method calls KINInit().
Wrapper for SUNDIALS' ARKODE library – Runge-Kutta time integration.
void SetLinearSolver(SundialsODELinearSolver &ls_spec)
Set a custom Jacobian system solver for implicit methods.
void SetIRKTableNum(int table_num)
Choose a specific Butcher table for implicit RK method.
const Operator * jacobian
virtual int SolveSystem(void *sundials_mem, Vector &b, const Vector &weight, const Vector &y_cur, const Vector &f_cur)=0
A base class for the MFEM classes wrapping SUNDIALS' solvers.
static int LinSysSolve(KINMemRec *kin_mem, N_Vector x, N_Vector b, realtype *sJpnorm, realtype *sFdotJp)
SundialsSolver(void *mem)
void SetFuncNormTol(double ftol)
Set KINSOL's functional norm tolerance.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
Abstract base class, wrapping the custom linear solvers interface in SUNDIALS' CVODE and ARKODE solve...
virtual ~ARKODESolver()
Destroy the associated ARKODE memory.
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
virtual void SetPreconditioner(Solver &solver)
Equivalent to SetSolver(Solver).
double GetTimeStep(void *sundials_mem)
Get the current scaled time step, gamma, from sundials_mem.
int flag
Flag returned by the last call to SUNDIALS.
CVODESolver(int lmm, int iter)
Construct a serial CVODESolver, a wrapper for SUNDIALS' CVODE solver.
void SetScaledStepTol(double sstol)
Set KINSOL's scaled step tolerance.
void PrintInfo() const
Print CVODE statistics.
void SetMaxOrder(int max_order)
Set the maximum order of the linear multistep method.
N_Vector y
Auxiliary N_Vector.
KinSolver(int strategy, bool oper_grad=true)
Construct a serial KinSolver, a wrapper for SUNDIALS' KINSOL solver.
void SetMaxStep(double dt_max)
Set the maximum time step of the Runge-Kutta method.
static int LinSysSetup(KINMemRec *kin_mem)
void PrintInfo() const
Print ARKODE statistics.
static const double default_rel_tol
virtual void Step(Vector &x, double &t, double &dt)
Use ARKODE to integrate over [t, t + dt], with the specified step mode.
void SetStepMode(int itask)
CVode supports two modes, specified by itask: CV_NORMAL (default) and CV_ONE_STEP.