15 #include "../config/config.hpp"
17 #ifdef MFEM_USE_SUNDIALS
27 #include <sundials/sundials_config.h>
29 #if !defined(SUNDIALS_VERSION_MAJOR) || (SUNDIALS_VERSION_MAJOR < 5)
30 #error MFEM requires SUNDIALS version 5.0.0 or newer!
32 #if defined(MFEM_USE_CUDA) && ((SUNDIALS_VERSION_MAJOR == 5) && (SUNDIALS_VERSION_MINOR < 4))
33 #error MFEM requires SUNDIALS version 5.4.0 or newer when MFEM_USE_CUDA=TRUE!
35 #include <sundials/sundials_matrix.h>
36 #include <sundials/sundials_linearsolver.h>
37 #include <arkode/arkode_arkstep.h>
38 #include <cvodes/cvodes.h>
39 #include <kinsol/kinsol.h>
90 SundialsNVector(MPI_Comm comm,
double *_data,
int loc_size,
long glob_size);
104 inline N_Vector_ID
GetNVectorID(N_Vector _x)
const {
return N_VGetVectorID(_x); }
108 inline MPI_Comm
GetComm()
const {
return *
static_cast<MPI_Comm*
>(N_VGetCommunicator(
x)); }
115 void SetSize(
int s,
long glob_size = 0);
147 operator N_Vector()
const {
return x; }
161 using Vector::operator=;
178 static N_Vector
MakeNVector(MPI_Comm comm,
bool use_device);
260 static int RHS(realtype t,
const N_Vector y, N_Vector ydot,
void *user_data);
263 static int LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix
A,
264 booleantype jok, booleantype *jcur,
265 realtype gamma,
void *user_data, N_Vector tmp1,
266 N_Vector tmp2, N_Vector tmp3);
269 static int LinSysSolve(SUNLinearSolver LS, SUNMatrix
A, N_Vector x,
270 N_Vector
b, realtype tol);
273 static int root(realtype t, N_Vector y, realtype *gout,
void *user_data);
276 typedef std::function<int(realtype t, Vector y, Vector gout, CVODESolver *)>
283 typedef std::function<int(Vector y, Vector w, CVODESolver*)>
EWTFunction;
330 virtual void Step(
Vector &x,
double &t,
double &dt);
395 static int RHSQ(realtype t,
const N_Vector y, N_Vector qdot,
void *user_data);
398 static int RHSB(realtype t, N_Vector y,
399 N_Vector
yB, N_Vector yBdot,
void *user_dataB);
402 static int RHSQB(realtype t, N_Vector y, N_Vector
yB,
403 N_Vector qBdot,
void *user_dataB);
406 static int ewt(N_Vector y, N_Vector w,
void *user_data);
456 virtual void Step(
Vector &x,
double &t,
double &dt);
459 virtual void StepB(
Vector &w,
double &t,
double &dt);
466 double reltolQ = 1e-3,
467 double abstolQ = 1e-8);
471 double abstolQB = 1e-8);
518 static int LinSysSetupB(realtype t, N_Vector y, N_Vector
yB, N_Vector fyB,
520 booleantype jok, booleantype *jcur,
521 realtype gamma,
void *user_data, N_Vector tmp1,
522 N_Vector tmp2, N_Vector tmp3);
525 static int LinSysSolveB(SUNLinearSolver LS, SUNMatrix
A, N_Vector x,
526 N_Vector
b, realtype tol);
559 static int RHS1(realtype t,
const N_Vector y, N_Vector ydot,
void *user_data);
560 static int RHS2(realtype t,
const N_Vector y, N_Vector ydot,
void *user_data);
564 static int LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix
A,
565 SUNMatrix
M, booleantype jok, booleantype *jcur,
566 realtype gamma,
void *user_data, N_Vector tmp1,
567 N_Vector tmp2, N_Vector tmp3);
570 static int LinSysSolve(SUNLinearSolver LS, SUNMatrix
A, N_Vector x,
571 N_Vector
b, realtype tol);
574 static int MassSysSetup(realtype t, SUNMatrix
M,
void *user_data,
575 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
578 static int MassSysSolve(SUNLinearSolver LS, SUNMatrix
M, N_Vector x,
579 N_Vector
b, realtype tol);
582 static int MassMult1(SUNMatrix
M, N_Vector x, N_Vector v);
585 static int MassMult2(N_Vector x, N_Vector v, realtype t,
633 virtual void Step(
Vector &x,
double &t,
double &dt);
727 static int Mult(
const N_Vector u, N_Vector fu,
void *user_data);
730 static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u,
731 booleantype *new_u,
void *user_data);
734 static int LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J,
735 void *user_data, N_Vector tmp1, N_Vector tmp2);
738 static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u,
739 N_Vector
b, realtype tol);
766 KINSolver(
int strategy,
bool oper_grad =
true);
776 KINSolver(MPI_Comm comm,
int strategy,
bool oper_grad =
true);
865 #endif // MFEM_USE_SUNDIALS
867 #endif // MFEM_SUNDIALS
static int RHSQ(realtype t, const N_Vector y, N_Vector qdot, void *user_data)
Wrapper to compute the ODE RHS Quadrature function.
void SetSize(int s, long glob_size=0)
Resize the vector to size s.
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
void SetJFNKSolver(Solver &solver)
int GetOwnership() const
Gets ownership of the internal N_Vector.
SundialsNVector * q
Quadrature vector.
void GetForwardSolution(double tB, mfem::Vector &yy)
Get Interpolated Forward solution y at backward integration time tB.
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
std::function< int(Vector y, Vector w, CVODESolver *)> EWTFunction
Typedef declaration for error weight functions.
static int LinSysSolveB(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system A x = b.
~SundialsNVector()
Calls SUNDIALS N_VDestroy function if the N_Vector is owned by 'this'.
void SetScaledStepTol(double sstol)
Set KINSOL's scaled step tolerance.
RootFunction root_func
A class member to facilitate pointing to a user-specified root function.
SundialsNVector * Y
State vector.
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Base abstract class for first order time dependent operators.
virtual void Step(Vector &x, double &t, double &dt)
Integrate the ODE with CVODE using the specified step mode.
static constexpr double default_abs_tolB
Default scalar backward absolute tolerance.
SundialsNVector * f_scale
scaling vectors
void SetSVtolerances(double reltol, Vector abstol)
Set the scalar relative and vector of absolute tolerances.
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.
void MakeRef(Vector &base, int offset, int s)
Reset the Vector to be a reference to a sub-vector of base.
void * GetMem() const
Access the SUNDIALS memory structure.
int indexB
backward problem index
virtual ~CVODESolver()
Destroy the associated CVODE memory and SUNDIALS objects.
bool use_oper_grad
use the Jv prod function
void SetIRKTableNum(int table_num)
Choose a specific Butcher table for a diagonally implicit RK method.
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
void SetIMEXTableNum(int etable_num, int itable_num)
Choose a specific Butcher table for an IMEX RK method.
static int LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix A, SUNMatrix M, booleantype jok, booleantype *jcur, realtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
void SetMaxStep(double dt_max)
Set the maximum time step.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
int GetFlag() const
Returns the last flag returned by a call to a SUNDIALS function.
void * sundials_mem
SUNDIALS mem structure.
static int PrecSolve(N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, N_Vector vv, void *user_data)
Solve the preconditioner equation .
long saved_global_size
Global vector length on last initialization.
virtual void Step(Vector &x, double &t, double &dt)
SundialsNVector()
Creates an empty SundialsNVector.
void UseSundialsLinearSolverB()
Use built in SUNDIALS Newton solver.
static int MassSysSetup(realtype t, SUNMatrix M, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
void SetSStolerancesB(double reltol, double abstol)
Tolerance specification functions for the adjoint problem.
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
void _SetDataAndSize_()
Set data and length from the internal N_Vector x.
ARKStepSolver(Type type=EXPLICIT)
Construct a serial wrapper to SUNDIALS' ARKode integrator.
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
int global_strategy
KINSOL solution strategy.
SUNLinearSolver LSA
Linear solver for A.
static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
Compute the matrix-vector product .
SundialsSolver()
Protected constructor: objects of this type should be constructed only as part of a derived class...
SUNMatrix M
Mass matrix M.
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
void SetSVtolerancesB(double reltol, Vector abstol)
Tolerance specification functions for the adjoint problem.
int lmm_type
Linear multistep method type.
Interface to the CVODE library – linear multi-step methods.
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)
Interface to the KINSOL library – nonlinear solver methods.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
static constexpr double default_rel_tolB
Default scalar backward relative tolerance.
SundialsNVector * yy
State vector.
Vector interface for SUNDIALS N_Vectors.
static constexpr double default_abs_tolQB
Default scalar backward absolute quadrature tolerance.
void SetMaxStep(double dt_max)
Set the maximum time step.
Type
Types of ARKODE solvers.
Type rk_type
Runge-Kutta type.
virtual ~ARKStepSolver()
Destroy the associated ARKode memory and SUNDIALS objects.
Newton's method for solving F(x)=b for a given operator F.
SundialsNVector * qB
State vector.
bool use_implicit
True for implicit or imex integration.
void SetJFNK(bool use_jfnk)
Set the Jacobian Free Newton Krylov flag. The default is false.
SUNNonlinearSolver NLS
Nonlinear solver.
CVODESolver(int lmm)
Construct a serial wrapper to SUNDIALS' CVODE integrator.
N_Vector StealNVector()
Changes the ownership of the the vector.
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Wrapper to compute the nonlinear residual .
int step_mode
ARKStep step mode (ARK_NORMAL or ARK_ONE_STEP).
SundialsNVector * y_scale
Wrapper for hypre's parallel vector class.
static int LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2)
Setup the linear system .
SUNLinearSolver LSM
Linear solver for M.
int ncheck
number of checkpoints used so far
static int RHS1(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
N_Vector x
The actual SUNDIALS object.
Base class for interfacing with SUNDIALS packages.
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to CVODE.
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
virtual ~KINSolver()
Destroy the associated KINSOL memory.
void MakeRef(Vector &base, int offset)
Reset the Vector to be a reference to a sub-vector of base without changing its current size...
void SetRootFinder(int components, RootFunction func)
Initialize Root Finder.
Vector wrk
Work vector needed for the JFNK PC.
EWTFunction ewt_func
A class member to facilitate pointing to a user-specified error weight function.
void SetFixedStep(double dt)
Use a fixed time step size (disable temporal adaptivity).
static int root(realtype t, N_Vector y, realtype *gout, void *user_data)
Prototype to define root finding for CVODE.
void SetLSMaxIter(int m)
Set the maximum number of linear solver iterations.
bool reinit
Flag to signal memory reinitialization is need.
int maa
number of acceleration vectors
static N_Vector MakeNVector(bool use_device)
Create a N_Vector.
void UseMFEMLinearSolverB()
Set Linear Solver for the backward problem.
std::function< int(realtype t, Vector y, Vector gout, CVODESolver *)> RootFunction
Typedef for root finding functions.
void _SetNvecDataAndSize_(long glob_size=0)
Set data and length of internal N_Vector x from 'this'.
void SetMAA(int maa)
Set the number of acceleration vectors to use with KIN_FP or KIN_PICARD.
int maxli
Maximum linear iterations.
void SetOwnership(int own)
Sets ownership of the internal N_Vector.
long GlobalSize() const
Returns the MPI global length for the internal N_Vector x.
long GetNumSteps()
Get the number of internal steps taken so far.
MPI_Comm GetComm() const
Returns the MPI communicator for the internal N_Vector x.
const double * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
void SetDataAndSize(double *d, int s, long glob_size=0)
Set the vector data and size.
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 .
void PrintInfo() const
Print various ARKStep statistics.
N_Vector_ID GetNVectorID() const
Returns the N_Vector_ID for the internal N_Vector.
static int RHS2(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u, booleantype *new_u, void *user_data)
Wrapper to compute the Jacobian-vector product .
void AllocateEmptyNVector(N_Vector &y)
const Operator * jacobian
stores oper->GetGradient()
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
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.
SUNLinearSolver LSB
Linear solver for A.
SundialsNVector * yB
State vector.
static constexpr double default_rel_tol
Default scalar relative tolerance.
void SetWFTolerances(EWTFunction func)
Set multiplicative error weights.
SUNMatrix AB
Linear system A = I - gamma J, M - gamma J, or J.
static bool UseManagedMemory()
int step_mode
CVODE step mode (CV_NORMAL or CV_ONE_STEP).
static int MassMult2(N_Vector x, N_Vector v, realtype t, void *mtimes_data)
Compute the matrix-vector product at time t.
KINSolver(int strategy, bool oper_grad=true)
Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.
Implicit-explicit ARK method.
N_Vector_ID GetNVectorID(N_Vector _x) const
Returns the N_Vector_ID for the N_Vector _x.
static constexpr double default_abs_tol
Default scalar absolute tolerance.
int flag
Last flag returned from a call to SUNDIALS.
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to ARKode.
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
void PrintInfo() const
Print various CVODE statistics.
void Init(TimeDependentAdjointOperator &f_)
static int PrecSetup(N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, void *user_data)
Setup the preconditioner.
void SetMaxOrder(int max_order)
Set the maximum method order.
virtual ~CVODESSolver()
Destroy the associated CVODES memory and SUNDIALS objects.
static int ewt(N_Vector y, N_Vector w, void *user_data)
Error control function.
void SetLSMaxRestarts(int m)
Set the maximum number of linear solver restarts.
void UseMFEMMassLinearSolver(int tdep)
Attach mass matrix linear system setup, solve, and matrix-vector product methods from the TimeDepende...
void SetMaxNStepsB(int mxstepsB)
Set the maximum number of backward steps.
void SetERKTableNum(int table_num)
Choose a specific Butcher table for an explicit RK method.
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
virtual void SetPreconditioner(Solver &solver)
Equivalent to SetSolver(solver).
static int RHS(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Number of components in gout.
int maxlrs
Maximum linear solver restarts.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
long GetNumSteps()
Get Number of Steps for ForwardSolve.
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, realtype tol)
Solve the linear system .
static int MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
virtual void Step(Vector &x, double &t, double &dt)
Integrate the ODE with ARKode using the specified step mode.
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.
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
void SetMaxNSteps(int steps)
Set the maximum number of time steps.
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
void UseSundialsMassLinearSolver(int tdep)
Attach the SUNDIALS GMRES linear solver and the mass matrix matrix-vector product method from the Tim...