14 #ifdef MFEM_USE_SUNDIALS
21 #include <sundials/sundials_config.h>
23 #ifndef SUNDIALS_VERSION_MAJOR
25 #define MFEM_SUNDIALS_VERSION 20700
27 #define SUNFALSE FALSE
29 #define MFEM_SUNDIALS_VERSION \
30 ((SUNDIALS_VERSION_MAJOR)*10000 + \
31 (SUNDIALS_VERSION_MINOR)*100 + \
32 (SUNDIALS_VERSION_PATCH))
35 #include <nvector/nvector_serial.h>
37 #include <nvector/nvector_parallel.h>
40 #include <cvode/cvode_impl.h>
49 #include <arkode/arkode_impl.h>
50 #include <kinsol/kinsol_impl.h>
53 #if MFEM_SUNDIALS_VERSION < 30000
55 #include <cvode/cvode_spgmr.h>
56 #include <arkode/arkode_spgmr.h>
57 #include <kinsol/kinsol_spgmr.h>
60 #include <sunlinsol/sunlinsol_spgmr.h>
61 #include <cvode/cvode_spils.h>
62 #include <arkode/arkode_spils.h>
63 #include <kinsol/kinsol_spils.h>
71 double SundialsODELinearSolver::GetTimeStep(
void *sundials_mem)
73 return (type == CVODE) ?
74 ((CVodeMem)sundials_mem)->cv_gamma :
75 ((ARKodeMem)sundials_mem)->ark_gamma;
79 SundialsODELinearSolver::GetTimeDependentOperator(
void *sundials_mem)
81 void *user_data = (type == CVODE) ?
82 ((CVodeMem)sundials_mem)->cv_user_data :
83 ((ARKodeMem)sundials_mem)->ark_user_data;
92 static int cvLinSysInit(CVodeMem cv_mem)
94 return to_solver(cv_mem->cv_lmem)->
InitSystem(cv_mem);
97 static int cvLinSysSetup(CVodeMem cv_mem,
int convfail,
98 N_Vector ypred, N_Vector fpred, booleantype *jcurPtr,
99 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
101 Vector yp(ypred), fp(fpred), vt1(vtemp1), vt2(vtemp2), vt3(vtemp3);
102 return to_solver(cv_mem->cv_lmem)->
SetupSystem(cv_mem, convfail, yp, fp,
103 *jcurPtr, vt1, vt2, vt3);
106 static int cvLinSysSolve(CVodeMem cv_mem, N_Vector b, N_Vector weight,
107 N_Vector ycur, N_Vector fcur)
109 Vector bb(b), w(weight), yc(ycur), fc(fcur);
110 return to_solver(cv_mem->cv_lmem)->
SolveSystem(cv_mem, bb, w, yc, fc);
113 static int cvLinSysFree(CVodeMem cv_mem)
115 return to_solver(cv_mem->cv_lmem)->
FreeSystem(cv_mem);
118 static int arkLinSysInit(ARKodeMem ark_mem)
120 return to_solver(ark_mem->ark_lmem)->
InitSystem(ark_mem);
123 static int arkLinSysSetup(ARKodeMem ark_mem,
int convfail,
124 N_Vector ypred, N_Vector fpred, booleantype *jcurPtr,
125 N_Vector vtemp1, N_Vector vtemp2, N_Vector vtemp3)
127 Vector yp(ypred), fp(fpred), vt1(vtemp1), vt2(vtemp2), vt3(vtemp3);
128 return to_solver(ark_mem->ark_lmem)->
SetupSystem(ark_mem, convfail, yp, fp,
129 *jcurPtr, vt1, vt2, vt3);
132 #if MFEM_SUNDIALS_VERSION < 30000
133 static int arkLinSysSolve(ARKodeMem ark_mem, N_Vector b, N_Vector weight,
134 N_Vector ycur, N_Vector fcur)
136 Vector bb(b), w(weight), yc(ycur), fc(fcur);
137 return to_solver(ark_mem->ark_lmem)->
SolveSystem(ark_mem, bb, w, yc, fc);
140 static int arkLinSysSolve(ARKodeMem ark_mem, N_Vector b,
141 N_Vector ycur, N_Vector fcur)
143 Vector bb(b), w(ark_mem->ark_rwt), yc(ycur), fc(fcur);
144 return to_solver(ark_mem->ark_lmem)->
SolveSystem(ark_mem, bb, w, yc, fc);
148 static int arkLinSysFree(ARKodeMem ark_mem)
150 return to_solver(ark_mem->ark_lmem)->
FreeSystem(ark_mem);
153 const double SundialsSolver::default_rel_tol = 1e-4;
154 const double SundialsSolver::default_abs_tol = 1e-9;
157 int SundialsSolver::ODEMult(realtype t,
const N_Vector y,
158 N_Vector ydot,
void *td_oper)
166 f->
Mult(mfem_y, mfem_ydot);
170 static inline CVodeMem Mem(
const CVODESolver *
self)
172 return CVodeMem(self->SundialsMem());
175 CVODESolver::CVODESolver(
int lmm,
int iter)
178 y = N_VNewEmpty_Serial(0);
179 MFEM_ASSERT(y,
"error in N_VNewEmpty_Serial()");
182 sundials_mem = CVodeCreate(lmm, iter);
183 MFEM_ASSERT(sundials_mem,
"error in CVodeCreate()");
185 SetStepMode(CV_NORMAL);
188 SetSStolerances(default_rel_tol, default_abs_tol);
195 CVODESolver::CVODESolver(MPI_Comm comm,
int lmm,
int iter)
197 if (comm == MPI_COMM_NULL)
200 y = N_VNewEmpty_Serial(0);
201 MFEM_ASSERT(y,
"error in N_VNewEmpty_Serial()");
206 y = N_VNewEmpty_Parallel(comm, 0, 0);
207 MFEM_ASSERT(y,
"error in N_VNewEmpty_Parallel()");
211 sundials_mem = CVodeCreate(lmm, iter);
212 MFEM_ASSERT(sundials_mem,
"error in CVodeCreate()");
214 SetStepMode(CV_NORMAL);
217 SetSStolerances(default_rel_tol, default_abs_tol);
222 #endif // MFEM_USE_MPI
224 void CVODESolver::SetSStolerances(
double reltol,
double abstol)
226 CVodeMem mem = Mem(
this);
228 mem->cv_reltol = reltol;
229 mem->cv_Sabstol = abstol;
235 CVodeMem mem = Mem(
this);
236 MFEM_ASSERT(mem->cv_iter == CV_NEWTON,
237 "The function is applicable only to CV_NEWTON iteration type.");
239 if (mem->cv_lfree != NULL) { (mem->cv_lfree)(mem); }
243 mem->cv_linit = cvLinSysInit;
244 mem->cv_lsetup = cvLinSysSetup;
245 mem->cv_lsolve = cvLinSysSolve;
246 mem->cv_lfree = cvLinSysFree;
247 mem->cv_lmem = &ls_spec;
248 #if MFEM_SUNDIALS_VERSION < 30000
249 mem->cv_setupNonNull = TRUE;
251 ls_spec.
type = SundialsODELinearSolver::CVODE;
254 void CVODESolver::SetStepMode(
int itask)
256 Mem(
this)->cv_taskc = itask;
259 void CVODESolver::SetMaxOrder(
int max_order)
261 flag = CVodeSetMaxOrd(sundials_mem, max_order);
262 if (flag == CV_ILL_INPUT)
264 MFEM_WARNING(
"CVodeSetMaxOrd() did not change the maximum order!");
269 static inline void cvCopyInit(CVodeMem src, CVodeMem dest)
271 dest->cv_lmm = src->cv_lmm;
272 dest->cv_iter = src->cv_iter;
274 dest->cv_linit = src->cv_linit;
275 dest->cv_lsetup = src->cv_lsetup;
276 dest->cv_lsolve = src->cv_lsolve;
277 dest->cv_lfree = src->cv_lfree;
278 dest->cv_lmem = src->cv_lmem;
279 #if MFEM_SUNDIALS_VERSION < 30000
280 dest->cv_setupNonNull = src->cv_setupNonNull;
283 dest->cv_reltol = src->cv_reltol;
284 dest->cv_Sabstol = src->cv_Sabstol;
286 dest->cv_taskc = src->cv_taskc;
287 dest->cv_qmax = src->cv_qmax;
294 CVodeMem mem = Mem(
this);
297 if (mem->cv_MallocDone == SUNTRUE)
300 cvCopyInit(mem, &backup);
301 CVodeFree(&sundials_mem);
302 sundials_mem = CVodeCreate(backup.cv_lmm, backup.cv_iter);
303 MFEM_ASSERT(sundials_mem,
"error in CVodeCreate()");
304 cvCopyInit(&backup, mem);
310 int loc_size = f_.
Height();
313 NV_LENGTH_S(y) = loc_size;
314 NV_DATA_S(y) =
new double[loc_size]();
319 long local_size = loc_size, global_size;
320 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
322 NV_LOCLENGTH_P(y) = local_size;
323 NV_GLOBLENGTH_P(y) = global_size;
324 NV_DATA_P(y) =
new double[loc_size]();
329 cvCopyInit(mem, &backup);
330 flag = CVodeInit(mem, ODEMult, f_.
GetTime(), y);
331 MFEM_ASSERT(flag >= 0,
"CVodeInit() failed!");
332 cvCopyInit(&backup, mem);
337 delete [] NV_DATA_S(y);
343 delete [] NV_DATA_P(y);
349 flag = CVodeSetUserData(sundials_mem, f);
350 MFEM_ASSERT(flag >= 0,
"CVodeSetUserData() failed!");
352 flag = CVodeSStolerances(mem, mem->cv_reltol, mem->cv_Sabstol);
353 MFEM_ASSERT(flag >= 0,
"CVodeSStolerances() failed!");
356 void CVODESolver::Step(
Vector &x,
double &t,
double &dt)
358 CVodeMem mem = Mem(
this);
363 MFEM_VERIFY(NV_LENGTH_S(y) == x.
Size(),
"");
369 MFEM_VERIFY(NV_LOCLENGTH_P(y) == x.
Size(),
"");
373 if (mem->cv_nst == 0)
376 if (mem->cv_iter == CV_NEWTON && mem->cv_lsolve == NULL)
378 #if MFEM_SUNDIALS_VERSION < 30000
379 flag = CVSpgmr(sundials_mem, PREC_NONE, 0);
382 LS = SUNSPGMR(y, PREC_NONE, 0);
383 flag = CVSpilsSetLinearSolver(sundials_mem, LS);
388 N_VScale(ONE, y, mem->cv_zn[0]);
391 double tout = t + dt;
393 flag = CVode(sundials_mem, tout, y, &t, mem->cv_taskc);
394 MFEM_ASSERT(flag >= 0,
"CVode() failed!");
400 void CVODESolver::PrintInfo()
const
402 CVodeMem mem = Mem(
this);
406 "num steps: " << mem->cv_nst <<
", "
407 "num evals: " << mem->cv_nfe <<
", "
408 "num lin setups: " << mem->cv_nsetups <<
", "
409 "num nonlin sol iters: " << mem->cv_nni <<
"\n "
410 "last order: " << mem->cv_qu <<
", "
411 "next order: " << mem->cv_next_q <<
", "
412 "last dt: " << mem->cv_hu <<
", "
413 "next dt: " << mem->cv_next_h
417 CVODESolver::~CVODESolver()
420 CVodeFree(&sundials_mem);
425 return ARKodeMem(self->SundialsMem());
428 ARKODESolver::ARKODESolver(
Type type)
429 : use_implicit(type == IMPLICIT), irk_table(-1), erk_table(-1)
432 y = N_VNewEmpty_Serial(0);
433 MFEM_ASSERT(
y,
"error in N_VNewEmpty_Serial()");
449 : use_implicit(type == IMPLICIT), irk_table(-1), erk_table(-1)
451 if (comm == MPI_COMM_NULL)
454 y = N_VNewEmpty_Serial(0);
455 MFEM_ASSERT(
y,
"error in N_VNew_Serial()");
460 y = N_VNewEmpty_Parallel(comm, 0, 0);
461 MFEM_ASSERT(
y,
"error in N_VNewEmpty_Parallel()");
479 ARKodeMem mem = Mem(
this);
481 mem->ark_reltol = reltol;
482 mem->ark_Sabstol = abstol;
488 ARKodeMem mem = Mem(
this);
490 "The function is applicable only to implicit time integration.");
492 if (mem->ark_lfree != NULL) { mem->ark_lfree(mem); }
495 mem->ark_lsolve_type = 4;
498 mem->ark_linit = arkLinSysInit;
499 mem->ark_lsetup = arkLinSysSetup;
500 mem->ark_lsolve = arkLinSysSolve;
501 mem->ark_lfree = arkLinSysFree;
502 mem->ark_lmem = &ls_spec;
503 #if MFEM_SUNDIALS_VERSION < 30000
504 mem->ark_setupNonNull = TRUE;
511 Mem(
this)->ark_taskc = itask;
516 ARKodeMem mem = Mem(
this);
537 MFEM_ASSERT(
flag >= 0,
"ARKodeSetFixedStep() failed!");
541 static inline void arkCopyInit(ARKodeMem src, ARKodeMem dest)
543 dest->ark_lsolve_type = src->ark_lsolve_type;
544 dest->ark_linit = src->ark_linit;
545 dest->ark_lsetup = src->ark_lsetup;
546 dest->ark_lsolve = src->ark_lsolve;
547 dest->ark_lfree = src->ark_lfree;
548 dest->ark_lmem = src->ark_lmem;
549 #if MFEM_SUNDIALS_VERSION < 30000
550 dest->ark_setupNonNull = src->ark_setupNonNull;
553 dest->ark_reltol = src->ark_reltol;
554 dest->ark_Sabstol = src->ark_Sabstol;
556 dest->ark_taskc = src->ark_taskc;
557 dest->ark_q = src->ark_q;
558 dest->ark_fixedstep = src->ark_fixedstep;
559 dest->ark_hin = src->ark_hin;
564 ARKodeMem mem = Mem(
this);
568 if (mem->ark_MallocDone == SUNTRUE)
571 arkCopyInit(mem, &backup);
575 arkCopyInit(&backup, mem);
581 int loc_size = f_.
Height();
584 NV_LENGTH_S(
y) = loc_size;
585 NV_DATA_S(
y) =
new double[loc_size]();
590 long local_size = loc_size, global_size;
591 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
593 NV_LOCLENGTH_P(
y) = local_size;
594 NV_GLOBLENGTH_P(
y) = global_size;
595 NV_DATA_P(
y) =
new double[loc_size]();
600 arkCopyInit(mem, &backup);
606 MFEM_ASSERT(
flag >= 0,
"CVodeInit() failed!");
607 arkCopyInit(&backup, mem);
612 delete [] NV_DATA_S(
y);
618 delete [] NV_DATA_P(
y);
625 MFEM_ASSERT(
flag >= 0,
"ARKodeSetUserData() failed!");
627 flag = ARKodeSStolerances(mem, mem->ark_reltol, mem->ark_Sabstol);
628 MFEM_ASSERT(
flag >= 0,
"CVodeSStolerances() failed!");
631 MFEM_ASSERT(
flag >= 0,
"ARKodeSetOrder() failed!");
636 MFEM_ASSERT(
flag >= 0,
"ARKodeSetIRKTableNum() failed!");
641 MFEM_ASSERT(
flag >= 0,
"ARKodeSetERKTableNum() failed!");
647 ARKodeMem mem = Mem(
this);
652 MFEM_VERIFY(NV_LENGTH_S(
y) == x.
Size(),
"");
658 MFEM_VERIFY(NV_LOCLENGTH_P(
y) == x.
Size(),
"");
662 if (mem->ark_nst == 0)
665 if (mem->ark_implicit && mem->ark_linit == NULL)
667 #if MFEM_SUNDIALS_VERSION < 30000
671 LS = SUNSPGMR(
y, PREC_NONE, 0);
679 N_VScale(ONE,
y, mem->ark_ycur);
682 double tout = t + dt;
685 MFEM_ASSERT(
flag >= 0,
"ARKode() failed!");
693 ARKodeMem mem = Mem(
this);
697 "num steps: " << mem->ark_nst <<
", "
698 "num evals: " << mem->ark_nfe <<
", "
699 "num lin setups: " << mem->ark_nsetups <<
", "
700 "num nonlin sol iters: " << mem->ark_nni <<
"\n "
701 "method order: " << mem->ark_q <<
", "
702 "last dt: " << mem->ark_h <<
", "
703 "next dt: " << mem->ark_next_h
714 static inline KINMem Mem(
const KinSolver *
self)
716 return KINMem(self->SundialsMem());
732 booleantype *new_u,
void *user_data)
740 self->jacobian = &
self->oper->GetGradient(mfem_u);
743 self->jacobian->Mult(mfem_v, mfem_Jv);
750 const Vector u(kin_mem->kin_uu);
755 self->prec->SetOperator(*self->jacobian);
762 realtype *sJpnorm, realtype *sFdotJp)
771 if ( (kin_mem->kin_globalstrategy == KIN_LINESEARCH) ||
772 (kin_mem->kin_globalstrategy != KIN_FP &&
773 kin_mem->kin_etaflag == KIN_ETACHOICE1) )
776 self->jacobian->Mult(mx, mb);
778 *sJpnorm = N_VWL2Norm(b, kin_mem->kin_fscale);
779 N_VProd(b, kin_mem->kin_fscale, b);
780 N_VProd(b, kin_mem->kin_fscale, b);
781 *sFdotJp = N_VDotProd(kin_mem->kin_fval, b);
789 : use_oper_grad(oper_grad), jacobian(NULL)
792 y = N_VNewEmpty_Serial(0);
793 y_scale = N_VNewEmpty_Serial(0);
794 f_scale = N_VNewEmpty_Serial(0);
795 MFEM_ASSERT(
y &&
y_scale &&
f_scale,
"Error in N_VNewEmpty_Serial().");
800 Mem(
this)->kin_globalstrategy = strategy;
802 abs_tol = Mem(
this)->kin_fnormtol;
811 : use_oper_grad(oper_grad), jacobian(NULL)
813 if (comm == MPI_COMM_NULL)
816 y = N_VNewEmpty_Serial(0);
817 y_scale = N_VNewEmpty_Serial(0);
818 f_scale = N_VNewEmpty_Serial(0);
819 MFEM_ASSERT(
y &&
y_scale &&
f_scale,
"Error in N_VNewEmpty_Serial().");
824 y = N_VNewEmpty_Parallel(comm, 0, 0);
825 y_scale = N_VNewEmpty_Parallel(comm, 0, 0);
826 f_scale = N_VNewEmpty_Parallel(comm, 0, 0);
827 MFEM_ASSERT(
y &&
y_scale &&
f_scale,
"Error in N_VNewEmpty_Parallel().");
833 Mem(
this)->kin_globalstrategy = strategy;
835 abs_tol = Mem(
this)->kin_fnormtol;
844 static inline void kinCopyInit(KINMem src, KINMem dest)
846 dest->kin_linit = src->kin_linit;
847 dest->kin_lsetup = src->kin_lsetup;
848 dest->kin_lsolve = src->kin_lsolve;
849 dest->kin_lfree = src->kin_lfree;
850 dest->kin_lmem = src->kin_lmem;
851 #if MFEM_SUNDIALS_VERSION < 30000
852 dest->kin_setupNonNull = src->kin_setupNonNull;
854 dest->kin_msbset = src->kin_msbset;
856 dest->kin_globalstrategy = src->kin_globalstrategy;
857 dest->kin_printfl = src->kin_printfl;
858 dest->kin_mxiter = src->kin_mxiter;
859 dest->kin_scsteptol = src->kin_scsteptol;
860 dest->kin_fnormtol = src->kin_fnormtol;
865 KINMem mem = Mem(
this);
869 if (mem->kin_MallocDone == SUNTRUE)
872 kinCopyInit(mem, &backup);
876 kinCopyInit(&backup, mem);
886 NV_DATA_S(
y) =
new double[
height]();
895 long local_size =
height, global_size;
896 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
898 NV_LOCLENGTH_P(
y) = local_size;
899 NV_GLOBLENGTH_P(
y) = global_size;
900 NV_DATA_P(
y) =
new double[local_size]();
901 NV_LOCLENGTH_P(
y_scale) = local_size;
902 NV_GLOBLENGTH_P(
y_scale) = global_size;
904 NV_LOCLENGTH_P(
f_scale) = local_size;
905 NV_GLOBLENGTH_P(
f_scale) = global_size;
910 kinCopyInit(mem, &backup);
915 N_VConst(ZERO, mem->kin_pp);
916 MFEM_ASSERT(
flag >= 0,
"KINInit() failed!");
917 kinCopyInit(&backup, mem);
922 delete [] NV_DATA_S(
y);
928 delete [] NV_DATA_P(
y);
935 MFEM_ASSERT(
flag >= 0,
"KINSetUserData() failed!");
940 #if MFEM_SUNDIALS_VERSION < 30000
942 MFEM_ASSERT(
flag >= 0,
"KINSpgmr() failed!");
944 SUNLinearSolver LS = NULL;
945 LS = SUNSPGMR(
y, PREC_NONE, 0);
947 MFEM_ASSERT(
flag >= 0,
"KINSpilsSetLinearSolver() failed!");
953 MFEM_ASSERT(
flag >= 0,
"KINSpilsSetJacTimesVecFn() failed!");
962 KINMem mem = Mem(
this);
964 mem->kin_linit = NULL;
967 mem->kin_lfree = NULL;
968 mem->kin_lmem =
this;
969 #if MFEM_SUNDIALS_VERSION < 30000
970 mem->kin_setupNonNull = TRUE;
977 Mem(
this)->kin_scsteptol = sstol;
982 Mem(
this)->kin_msbset = max_calls;
1003 MPI_Allreduce(&lnorm, &norm, 1, MPI_DOUBLE, MPI_MAX, NV_COMM_P(
y));
1014 Mem(
this)->kin_fnormtol =
abs_tol;
1019 Mem(
this)->kin_fnormtol =
rel_tol;
1024 Mem(
this)->kin_fnormtol =
abs_tol;
1034 KINMem mem = Mem(
this);
1037 MFEM_ASSERT(
flag >= 0,
"KINSetPrintLevel() failed!");
1040 MFEM_ASSERT(
flag >= 0,
"KINSetNumMaxIters() failed!");
1043 MFEM_ASSERT(
flag >= 0,
"KINSetScaledStepTol() failed!");
1046 MFEM_ASSERT(
flag >= 0,
"KINSetFuncNormTol() failed!");
1051 MFEM_VERIFY(NV_LENGTH_S(
y) == x.
Size(),
"");
1059 MFEM_VERIFY(NV_LOCLENGTH_P(
y) == x.
Size(),
"");
void SetFixedStep(double dt)
Use a fixed time step size, instead of performing any form of temporal adaptivity.
virtual double GetTime() const
Read the currently set time.
virtual int FreeSystem(void *sundials_mem)=0
virtual int SolveSystem(void *sundials_mem, Vector &b, const Vector &w, const Vector &y_cur, const Vector &f_cur)=0
static const double default_abs_tol
virtual void Init(TimeDependentOperator &f_)
Set the ODE right-hand-side operator.
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 Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
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 void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
virtual void SetTime(const double _t)
Set the current time.
int Size() const
Returns the size of the vector.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void * sundials_mem
Pointer to the SUNDIALS mem object.
double Normlinf() const
Returns the l_infinity norm of the vector.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
double * GetData() const
Return a pointer to the beginning of the Vector data.
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.
Wrapper for SUNDIALS' KINSOL library – Nonlinear solvers.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
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)
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.
enum mfem::SundialsODELinearSolver::@1 type
Is CVODE or ARKODE using this object?
void SetIRKTableNum(int table_num)
Choose a specific Butcher table for implicit RK method.
const Operator * jacobian
static int LinSysSolve(KINMemRec *kin_mem, N_Vector x, N_Vector b, realtype *sJpnorm, realtype *sFdotJp)
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...
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...
int height
Dimension of the output / number of rows in the matrix.
virtual ~ARKODESolver()
Destroy the associated ARKODE memory.
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
int flag
Flag returned by the last call to SUNDIALS.
void SetScaledStepTol(double sstol)
Set KINSOL's scaled step tolerance.
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
N_Vector y
Auxiliary N_Vector.
KinSolver(int strategy, bool oper_grad=true)
Construct a serial KinSolver, a wrapper for SUNDIALS' KINSOL solver.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
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.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.