14#ifdef MFEM_USE_SUNDIALS
22#include <nvector/nvector_serial.h>
23#if defined(MFEM_USE_CUDA)
24#include <nvector/nvector_cuda.h>
25#elif defined(MFEM_USE_HIP)
26#include <nvector/nvector_hip.h>
29#include <nvector/nvector_mpiplusx.h>
30#include <nvector/nvector_parallel.h>
34#include <sunlinsol/sunlinsol_spgmr.h>
35#include <sunlinsol/sunlinsol_spfgmr.h>
38#define GET_CONTENT(X) ( X->content )
40#if defined(MFEM_USE_CUDA)
41#define SUN_Hip_OR_Cuda(X) X##_Cuda
42#define SUN_HIP_OR_CUDA(X) X##_CUDA
43#elif defined(MFEM_USE_HIP)
44#define SUN_Hip_OR_Cuda(X) X##_Hip
45#define SUN_HIP_OR_CUDA(X) X##_HIP
50#if (SUNDIALS_VERSION_MAJOR < 6)
116 sunindextype local_length,
117 sunindextype global_length,
125#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
131 SUNMemoryHelper helper,
146#if defined(MFEM_USE_MPI) && (defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP))
160#if MFEM_SUNDIALS_VERSION < 70100
161#define MFEM_ARKode(FUNC) ARKStep##FUNC
163#define MFEM_ARKode(FUNC) ARKode##FUNC
168#define STR(s) STR1(s)
176 Sundials::Instance();
187 return Sundials::Instance().context;
192 return Sundials::Instance().memHelper;
195#if (SUNDIALS_VERSION_MAJOR >= 6)
200 int mpi_initialized = 0;
201 MPI_Initialized(&mpi_initialized);
202 MPI_Comm communicator = mpi_initialized ? MPI_COMM_WORLD : MPI_COMM_NULL;
203#if SUNDIALS_VERSION_MAJOR < 7
204 int return_val = SUNContext_Create((
void*) &communicator, &context);
206 int return_val = SUNContext_Create(communicator, &context);
209#if SUNDIALS_VERSION_MAJOR < 7
210 int return_val = SUNContext_Create(
nullptr, &context);
212 int return_val = SUNContext_Create((SUNComm)(0), &context);
215 MFEM_VERIFY(return_val == 0,
"Call to SUNContext_Create failed");
217 memHelper = std::move(actual_helper);
218 isInitialized =
true;
229 if (sundials.isInitialized)
232 SUNContext_Free(&context);
233 sundials.isInitialized =
false;
256#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
265 h->ops->copy = SUN_Hip_OR_Cuda(SUNMemoryHelper_Copy);
266 h->ops->copyasync = SUN_Hip_OR_Cuda(SUNMemoryHelper_CopyAsync);
271 this->h = that_helper.h;
272 that_helper.h =
nullptr;
283 SUNMemory* memptr,
size_t memsize,
284 SUNMemoryType mem_type
285#
if (SUNDIALS_VERSION_MAJOR >= 6)
290#if (SUNDIALS_VERSION_MAJOR < 7)
291 SUNMemory sunmem = SUNMemoryNewEmpty();
293 SUNMemory sunmem = SUNMemoryNewEmpty(helper->sunctx);
297 sunmem->own = SUNTRUE;
300 if (mem_type == SUNMEMTYPE_HOST)
305 sunmem->type = SUNMEMTYPE_HOST;
308 else if (mem_type == SUNMEMTYPE_DEVICE || mem_type == SUNMEMTYPE_UVM)
313 sunmem->type = mem_type;
328#
if (SUNDIALS_VERSION_MAJOR >= 6)
333 if (sunmem->ptr && sunmem->own && !
mm.
IsKnown(sunmem->ptr))
335 if (sunmem->type == SUNMEMTYPE_HOST)
341 else if (sunmem->type == SUNMEMTYPE_DEVICE || sunmem->type == SUNMEMTYPE_UVM)
349 MFEM_ABORT(
"Invalid SUNMEMTYPE");
367 N_Vector local_x =
MPIPlusX() ? N_VGetLocalVector_MPIPlusX(
x) :
x;
369 N_Vector local_x =
x;
371 N_Vector_ID
id = N_VGetVectorID(local_x);
376 case SUNDIALS_NVEC_SERIAL:
378 MFEM_ASSERT(NV_OWN_DATA_S(local_x) == SUNFALSE,
"invalid serial N_Vector");
380 NV_LENGTH_S(local_x) =
size;
383#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
384 case SUN_HIP_OR_CUDA(SUNDIALS_NVEC):
386 SUN_Hip_OR_Cuda(N_VSetHostArrayPointer)(
HostReadWrite(), local_x);
387 SUN_Hip_OR_Cuda(N_VSetDeviceArrayPointer)(
ReadWrite(), local_x);
388 static_cast<SUN_Hip_OR_Cuda(N_VectorContent)
>(GET_CONTENT(
389 local_x))->length =
size;
394 case SUNDIALS_NVEC_PARALLEL:
396 MFEM_ASSERT(NV_OWN_DATA_P(
x) == SUNFALSE,
"invalid parallel N_Vector");
398 NV_LOCLENGTH_P(
x) =
size;
403 if (glob_size == 0 && glob_size !=
size)
405 long local_size =
size;
406 MPI_Allreduce(&local_size, &glob_size, 1, MPI_LONG,
410 NV_GLOBLENGTH_P(
x) = glob_size;
415 MFEM_ABORT(
"N_Vector type " <<
id <<
" is not supported");
425 if (glob_size == 0 && glob_size !=
size)
427 long local_size =
size;
428 MPI_Allreduce(&local_size, &glob_size, 1, MPI_LONG,
432 static_cast<N_VectorContent_MPIManyVector
>(GET_CONTENT(
x))->global_length =
441 N_Vector local_x =
MPIPlusX() ? N_VGetLocalVector_MPIPlusX(
x) :
x;
443 N_Vector local_x =
x;
445 N_Vector_ID
id = N_VGetVectorID(local_x);
450 case SUNDIALS_NVEC_SERIAL:
452 const bool known =
mm.
IsKnown(NV_DATA_S(local_x));
453 size = NV_LENGTH_S(local_x);
458#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
459 case SUN_HIP_OR_CUDA(SUNDIALS_NVEC):
461 double *h_ptr = SUN_Hip_OR_Cuda(N_VGetHostArrayPointer)(local_x);
462 double *d_ptr = SUN_Hip_OR_Cuda(N_VGetDeviceArrayPointer)(local_x);
464 size = SUN_Hip_OR_Cuda(N_VGetLength)(local_x);
472 case SUNDIALS_NVEC_PARALLEL:
475 size = NV_LENGTH_S(
x);
476 data.
Wrap(NV_DATA_P(
x), NV_LOCLENGTH_P(
x),
false);
482 MFEM_ABORT(
"N_Vector type " <<
id <<
" is not supported");
541 :
SundialsNVector(vec.GetComm(), vec.GetData(), vec.Size(), vec.GlobalSize())
552 N_VDestroy(N_VGetLocalVector_MPIPlusX(
x));
580#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
595 MFEM_VERIFY(
x,
"Error in SundialsNVector::MakeNVector.");
605 if (comm == MPI_COMM_NULL)
611#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
629 MFEM_VERIFY(
x,
"Error in SundialsNVector::MakeNVector.");
641static SUNMatrix_ID MatGetID(SUNMatrix)
643 return (SUNMATRIX_CUSTOM);
646static void MatDestroy(SUNMatrix A)
648 if (A->content) { A->content = NULL; }
649 if (A->ops) { free(A->ops); A->ops = NULL; }
659static SUNLinearSolver_Type LSGetType(SUNLinearSolver)
661 return (SUNLINEARSOLVER_MATRIX_ITERATIVE);
664static int LSFree(SUNLinearSolver LS)
666 if (LS->content) { LS->content = NULL; }
667 if (LS->ops) { free(LS->ops); LS->ops = NULL; }
686 self->
f->
Mult(mfem_y, mfem_ydot);
697 if (!self->
root_func) {
return CV_RTFUNC_FAIL; }
702 return self->
root_func(t, mfem_y, mfem_gout, self);
710 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in SetRootFinder()");
716 void*, N_Vector, N_Vector, N_Vector)
739 : lmm_type(lmm), step_mode(CV_NORMAL)
746 : lmm_type(lmm), step_mode(CV_NORMAL)
758 long local_size = f_.
Height();
761 long global_size = 0;
764 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
778 resize = (
Y->
Size() != local_size);
783 int l_resize = (
Y->
Size() != local_size) ||
785 MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
811 Y->
SetSize(local_size, global_size);
822 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeInit()");
826 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetUserData()");
830 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetSStolerances()");
843 MFEM_VERIFY(
Y->
Size() == x.
Size(),
"size mismatch");
849 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeReInit()");
855 double tout = t + dt;
857 MFEM_VERIFY(
flag >= 0,
"error in CVode()");
864 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetLastStep()");
870 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
871 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
875 MFEM_VERIFY(
LSA,
"error in SUNLinSolNewEmpty()");
878 LSA->ops->gettype = LSGetType;
880 LSA->ops->free = LSFree;
883 MFEM_VERIFY(
A,
"error in SUNMatNewEmpty()");
886 A->ops->getid = MatGetID;
887 A->ops->destroy = MatDestroy;
891 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolver()");
895 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinSysFn()");
901 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
902 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
906 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPGMR()");
910 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolver()");
921 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSStolerances()");
927 "abs tolerance is not the same size.");
933 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSVtolerances()");
939 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetMaxStep()");
945 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetMaxNumSteps()");
952 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetNumSteps()");
959 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetMaxOrd()");
964 long int nsteps, nfevals, nlinsetups, netfails;
966 double hinused, hlast, hcur, tcur;
967 long int nniters, nncfails;
981 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetIntegratorStats()");
987 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetNonlinSolvStats()");
991 "num steps: " << nsteps <<
"\n"
992 "num rhs evals: " << nfevals <<
"\n"
993 "num lin setups: " << nlinsetups <<
"\n"
994 "num nonlin sol iters: " << nniters <<
"\n"
995 "num nonlin conv fail: " << nncfails <<
"\n"
996 "num error test fails: " << netfails <<
"\n"
997 "last order: " << qlast <<
"\n"
998 "current order: " << qcur <<
"\n"
999 "initial dt: " << hinused <<
"\n"
1000 "last dt: " << hlast <<
"\n"
1001 "current dt: " << hcur <<
"\n"
1002 "current t: " << tcur <<
"\n" << endl;
1012 SUNNonlinSolFree(
NLS);
1050 MFEM_VERIFY(t <= f->GetTime(),
"t > current forward solver time");
1053 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetQuad()");
1060 MFEM_VERIFY(t <= f->GetTime(),
"t > current forward solver time");
1063 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetQuadB()");
1065 dG_dp.
Set(-1., *
qB);
1073 MFEM_VERIFY(
flag >= 0,
"error in CVodeGetAdjY()");
1093 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeCreateB()");
1097 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeInit()");
1101 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetUserDataB()");
1106 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetSStolerancesB()");
1118 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeAdjInit");
1124 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeSetMaxNumStepsB()");
1133 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadInit()");
1136 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeSetQuadErrCon");
1139 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadSStolerances");
1148 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadInitB()");
1151 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeSetQuadErrConB");
1154 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadSStolerancesB");
1160 if (
AB != NULL) { SUNMatDestroy(
AB);
AB = NULL; }
1161 if (
LSB != NULL) { SUNLinSolFree(
LSB);
LSB = NULL; }
1165 MFEM_VERIFY(
LSB,
"error in SUNLinSolNewEmpty()");
1167 LSB->content =
this;
1168 LSB->ops->gettype = LSGetType;
1170 LSB->ops->free = LSFree;
1173 MFEM_VERIFY(
AB,
"error in SUNMatNewEmpty()");
1176 AB->ops->getid = MatGetID;
1177 AB->ops->destroy = MatDestroy;
1181 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolverB()");
1186 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinSysFn()");
1192 if (
AB != NULL) { SUNMatDestroy(
AB);
AB = NULL; }
1193 if (
LSB != NULL) { SUNLinSolFree(
LSB);
LSB = NULL; }
1197 MFEM_VERIFY(
LSB,
"error in SUNLinSol_SPGMR()");
1201 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolverB()");
1205 N_Vector fyB, SUNMatrix AB,
1208 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
1219 return (
f->SUNImplicitSetupB(t, mfem_y, mfem_yB, mfem_fyB, jokB, jcurB,
1232 int ret =
f->SUNImplicitSolveB(mfem_yB, mfem_Rb, tol);
1239 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSStolerancesB()");
1245 "abs tolerance is not the same size.");
1251 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSVtolerancesB()");
1271 f->QuadratureIntegration(mfem_y, mfem_qdot);
1285 f->QuadratureSensitivityMult(mfem_y, mfem_yB, mfem_qBdot);
1301 f->AdjointRateMult(mfem_y, mfem_yB, mfem_yBdot);
1312 return self->
ewt_func(mfem_y, mfem_w, self);
1319 MFEM_VERIFY(
Y->
Size() == x.
Size(),
"size mismatch");
1325 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeReInit()");
1332 double tout = t + dt;
1334 MFEM_VERIFY(
flag >= 0,
"error in CVodeF()");
1341 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetLastStep()");
1353 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeReInit()");
1360 double tout = tB - dtB;
1362 MFEM_VERIFY(
flag >= 0,
"error in CVodeB()");
1366 MFEM_VERIFY(
flag >= 0,
"error in CVodeGetB()");
1408 self->
f->
Mult(mfem_y, mfem_result);
1434 self->
f->
Mult(mfem_y, mfem_result);
1448 void*, N_Vector, N_Vector, N_Vector)
1480 void*, N_Vector, N_Vector, N_Vector)
1523 : rk_type(type), step_mode(ARK_NORMAL),
1524 use_implicit(type == IMPLICIT || type == IMEX)
1531 : rk_type(type), step_mode(ARK_NORMAL),
1532 use_implicit(type == IMPLICIT || type == IMEX)
1544 long local_size = f_.
Height();
1552 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
1566 resize = (
Y->
Size() != local_size);
1571 int l_resize = (
Y->
Size() != local_size) ||
1573 MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
1595 Y->
SetSize(local_size, global_size);
1620 MFEM_VERIFY(
flag == ARK_SUCCESS,
1621 "error in " STR(MFEM_ARKode(SetUserData))
"()");
1626 MFEM_VERIFY(
flag == ARK_SUCCESS,
1627 "error in " STR(MFEM_ARKode(SStolerances))
"()");
1640 MFEM_VERIFY(
Y->
Size() == x.
Size(),
"size mismatch");
1658 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepReInit()");
1665 double tout = t + dt;
1667 MFEM_VERIFY(
flag >= 0,
"error in " STR(MFEM_ARKode(Evolve))
"()");
1674 MFEM_VERIFY(
flag == ARK_SUCCESS,
1675 "error in " STR(MFEM_ARKode(GetLastStep))
"()");
1681 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
1682 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1686 MFEM_VERIFY(
LSA,
"error in SUNLinSolNewEmpty()");
1688 LSA->content =
this;
1689 LSA->ops->gettype = LSGetType;
1691 LSA->ops->free = LSFree;
1694 MFEM_VERIFY(
A,
"error in SUNMatNewEmpty()");
1697 A->ops->getid = MatGetID;
1698 A->ops->destroy = MatDestroy;
1702 MFEM_VERIFY(
flag == ARK_SUCCESS,
1703 "error in " STR(MFEM_ARKode(SetLinearSolver))
"()");
1707 MFEM_VERIFY(
flag == ARK_SUCCESS,
1708 "error in " STR(MFEM_ARKode(SetLinSysFn))
"()");
1714 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
1715 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1719 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPGMR()");
1723 MFEM_VERIFY(
flag == ARK_SUCCESS,
1724 "error in " STR(MFEM_ARKode(SetLinearSolver))
"()");
1730 if (
M != NULL) { SUNMatDestroy(
M);
M = NULL; }
1731 if (
LSM != NULL) { SUNLinSolFree(
LSM);
LSM = NULL; }
1735 MFEM_VERIFY(
LSM,
"error in SUNLinSolNewEmpty()");
1737 LSM->content =
this;
1738 LSM->ops->gettype = LSGetType;
1740 LSM->ops->free = LSFree;
1743 MFEM_VERIFY(
M,
"error in SUNMatNewEmpty()");
1746 M->ops->getid = MatGetID;
1748 M->ops->destroy = MatDestroy;
1752 MFEM_VERIFY(
flag == ARK_SUCCESS,
1753 "error in " STR(MFEM_ARKode(SetMassLinearSolver))
"()");
1757 MFEM_VERIFY(
flag == ARK_SUCCESS,
1758 "error in " STR(MFEM_ARKode(SetMassFn))
"()");
1761 MFEM_VERIFY(!
f->
isExplicit(),
"ODE operator is expressed in EXPLICIT form")
1767 if (
M != NULL) { SUNMatDestroy(
A);
M = NULL; }
1768 if (
LSM != NULL) { SUNLinSolFree(
LSM);
LSM = NULL; }
1772 MFEM_VERIFY(
LSM,
"error in SUNLinSol_SPGMR()");
1776 MFEM_VERIFY(
flag == ARK_SUCCESS,
1777 "error in " STR(MFEM_ARKode(SetMassLinearSolver))
"()");
1782 MFEM_VERIFY(
flag == ARK_SUCCESS,
1783 "error in " STR(MFEM_ARKode(SetMassTimes))
"()");
1786 MFEM_VERIFY(!
f->
isExplicit(),
"ODE operator is expressed in EXPLICIT form")
1797 MFEM_VERIFY(
flag == ARK_SUCCESS,
1798 "error in " STR(MFEM_ARKode(SStolerances))
"()");
1804 MFEM_VERIFY(
flag == ARK_SUCCESS,
1805 "error in " STR(MFEM_ARKode(
SetMaxStep))
"()");
1811 MFEM_VERIFY(
flag == ARK_SUCCESS,
1812 "error in " STR(MFEM_ARKode(
SetOrder))
"()");
1818 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetTableNum()");
1824 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetTableNum()");
1831 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetTableNum()");
1837 MFEM_VERIFY(
flag == ARK_SUCCESS,
1843 long int nsteps, expsteps, accsteps, step_attempts;
1844 long int nfe_evals, nfi_evals;
1845 long int nlinsetups, netfails;
1846 double hinused, hlast, hcur, tcur;
1847 long int nniters, nncfails;
1858 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepGetTimestepperStats()");
1871 MFEM_VERIFY(
flag == ARK_SUCCESS,
1872 "error in " STR(MFEM_ARKode(GetNonlinSolvStats))
"()");
1876 "num steps: " << nsteps <<
"\n"
1877 "num exp rhs evals: " << nfe_evals <<
"\n"
1878 "num imp rhs evals: " << nfi_evals <<
"\n"
1879 "num lin setups: " << nlinsetups <<
"\n"
1880 "num nonlin sol iters: " << nniters <<
"\n"
1881 "num nonlin conv fail: " << nncfails <<
"\n"
1882 "num steps attempted: " << step_attempts <<
"\n"
1883 "num acc limited steps: " << accsteps <<
"\n"
1884 "num exp limited stepfails: " << expsteps <<
"\n"
1885 "num error test fails: " << netfails <<
"\n"
1886 "initial dt: " << hinused <<
"\n"
1887 "last dt: " << hlast <<
"\n"
1888 "current dt: " << hcur <<
"\n"
1889 "current t: " << tcur <<
"\n" << endl;
1899 SUNNonlinSolFree(
NLS);
1946 void *, N_Vector, N_Vector )
2014 : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
2015 f_scale(NULL), jacobian(NULL)
2022#if MFEM_SUNDIALS_VERSION < 70000
2023 abs_tol = pow(UNIT_ROUNDOFF, 1.0/3.0);
2025 abs_tol = pow(SUN_UNIT_ROUNDOFF, 1.0/3.0);
2032 : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
2033 f_scale(NULL), jacobian(NULL)
2040#if MFEM_SUNDIALS_VERSION < 70000
2041 abs_tol = pow(UNIT_ROUNDOFF, 1.0/3.0);
2043 abs_tol = pow(SUN_UNIT_ROUNDOFF, 1.0/3.0);
2057 long local_size =
height;
2065 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
2076 resize = (
Y->
Size() != local_size);
2081 int l_resize = (
Y->
Size() != local_size) ||
2083 MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
2105 Y->
SetSize(local_size, global_size);
2120 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetMAA()");
2123 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetDelayAA()");
2126 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetDampingAA()");
2128#if SUNDIALS_VERSION_MAJOR >= 6
2130 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetOrthAA()");
2136 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINInit()");
2140 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetUserData()");
2143 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetDamping()");
2153 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
2154 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
2157 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPGMR()");
2160 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetLinearSolver()");
2166 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetJacTimesVecFn()");
2184 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
2185 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
2189 MFEM_VERIFY(
LSA,
"error in SUNLinSolNewEmpty()");
2191 LSA->content =
this;
2192 LSA->ops->gettype = LSGetType;
2194 LSA->ops->free = LSFree;
2197 MFEM_VERIFY(
A,
"error in SUNMatNewEmpty()");
2200 A->ops->getid = MatGetID;
2201 A->ops->destroy = MatDestroy;
2205 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetLinearSolver()");
2209 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetJacFn()");
2221 if (
A != NULL) { SUNMatDestroy(
A);
A = NULL; }
2222 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
2227 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPFGMR()");
2233 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetLinearSolver()");
2240 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetPreconditioner()");
2247 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetScaledStepTol()");
2253 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetMaxSetupCalls()");
2262 MFEM_ABORT(
"Subsequent calls to EnableAndersonAcc() must set"
2263 " the subspace size to less or equal to the initially requested size."
2264 " If SetOperator() has already been called, the subspace size can't be"
2269 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetMAA()");
2272 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetDelayAA()");
2275 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetDampingAA()");
2277#if SUNDIALS_VERSION_MAJOR >= 6
2279 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetOrthAA()");
2281 if (orth != KIN_ORTH_MGS)
2283 MFEM_WARNING(
"SUNDIALS < v6 does not support setting the Anderson"
2284 " acceleration orthogonalization routine!");
2301 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetDamping()");
2307 MFEM_ABORT(
"this method is not supported! Use SetPrintLevel(int) instead.");
2332 double lnorm =
norm;
2356 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetFuncNormTol()");
2367 MFEM_ASSERT(
flag == KIN_SUCCESS,
"KINSetNumMaxIters() failed!");
2381 MPI_Comm_rank(
Y->
GetComm(), &rank);
2387#if MFEM_SUNDIALS_VERSION < 70000
2389 MFEM_VERIFY(
flag == KIN_SUCCESS,
"KINSetPrintLevel() failed!");
2393#ifdef SUNDIALS_BUILD_WITH_MONITORING
2396 flag = SUNLinSolSetInfoFile_SPFGMR(
LSA, stdout);
2398 "error in SUNLinSolSetInfoFile_SPFGMR()");
2400 flag = SUNLinSolSetPrintLevel_SPFGMR(
LSA, 1);
2402 "error in SUNLinSolSetPrintLevel_SPFGMR()");
2419 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINGetNumNonlinSolvIters()");
2424 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINGetFuncNorm()");
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
static int RHS2(sunrealtype t, const N_Vector y, N_Vector ydot, void *user_data)
void SetMaxStep(double dt_max)
Set the maximum time step.
void PrintInfo() const
Print various ARKStep statistics.
Type rk_type
Runge-Kutta type.
ARKStepSolver(Type type=EXPLICIT)
Construct a serial wrapper to SUNDIALS' ARKode integrator.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Type
Types of ARKODE solvers.
@ IMPLICIT
Implicit RK method.
@ IMEX
Implicit-explicit ARK method.
@ EXPLICIT
Explicit RK method.
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system .
static int LinSysSetup(sunrealtype t, N_Vector y, N_Vector fy, SUNMatrix A, SUNMatrix M, sunbooleantype jok, sunbooleantype *jcur, sunrealtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
void Init(TimeDependentOperator &f_) override
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
static int MassMult2(N_Vector x, N_Vector v, sunrealtype t, void *mtimes_data)
Compute the matrix-vector product at time t.
void SetIRKTableNum(ARKODE_DIRKTableID table_id)
Choose a specific Butcher table for a diagonally implicit RK method.
void SetFixedStep(double dt)
Use a fixed time step size (disable temporal adaptivity).
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
int step_mode
ARKStep step mode (ARK_NORMAL or ARK_ONE_STEP).
static int RHS1(sunrealtype t, const N_Vector y, N_Vector ydot, void *user_data)
void UseMFEMMassLinearSolver(int tdep)
Attach mass matrix linear system setup, solve, and matrix-vector product methods from the TimeDepende...
bool use_implicit
True for implicit or imex integration.
void UseSundialsMassLinearSolver(int tdep)
Attach the SUNDIALS GMRES linear solver and the mass matrix matrix-vector product method from the Tim...
virtual ~ARKStepSolver()
Destroy the associated ARKode memory and SUNDIALS objects.
void SetIMEXTableNum(ARKODE_ERKTableID etable_id, ARKODE_DIRKTableID itable_id)
Choose a specific Butcher table for an IMEX RK method.
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i....
static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
Compute the matrix-vector product .
static int MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system .
static int MassSysSetup(sunrealtype t, SUNMatrix M, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
void Step(Vector &x, real_t &t, real_t &dt) override
Integrate the ODE with ARKode using the specified step mode.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
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.
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
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 constexpr double default_abs_tolB
Default scalar backward absolute tolerance.
void SetMaxNStepsB(int mxstepsB)
Set the maximum number of backward steps.
static int RHSQ(sunrealtype t, const N_Vector y, N_Vector qdot, void *user_data)
Wrapper to compute the ODE RHS Quadrature function.
void Step(Vector &x, double &t, double &dt) override
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
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.
SundialsNVector * q
Quadrature vector.
int indexB
backward problem index
void GetForwardSolution(double tB, mfem::Vector &yy)
Get Interpolated Forward solution y at backward integration time tB.
SUNLinearSolver LSB
Linear solver for A.
void SetSVtolerancesB(double reltol, Vector abstol)
Tolerance specification functions for the adjoint problem.
void UseSundialsLinearSolverB()
Use built in SUNDIALS Newton solver.
void SetWFTolerances(EWTFunction func)
Set multiplicative error weights.
SundialsNVector * yy
State vector.
void Init(TimeDependentAdjointOperator &f_)
SUNMatrix AB
Linear system A = I - gamma J, M - gamma J, or J.
int ncheck
number of checkpoints used so far
SundialsNVector * qB
State vector.
static int ewt(N_Vector y, N_Vector w, void *user_data)
Error control function.
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
Initialize Quadrature Integration (Adjoint)
static constexpr double default_rel_tolB
Default scalar backward relative tolerance.
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
SundialsNVector * yB
State vector.
void InitQuadIntegration(mfem::Vector &q0, double reltolQ=1e-3, double abstolQ=1e-8)
virtual ~CVODESSolver()
Destroy the associated CVODES memory and SUNDIALS objects.
void SetSStolerancesB(double reltol, double abstol)
Tolerance specification functions for the adjoint problem.
void UseMFEMLinearSolverB()
Set Linear Solver for the backward problem.
Interface to the CVODE library – linear multi-step methods.
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
void SetRootFinder(int components, RootFunction func)
Initialize Root Finder.
std::function< int(sunrealtype t, Vector y, Vector gout, CVODESolver *)> RootFunction
Typedef for root finding functions.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
void Init(TimeDependentOperator &f_) override
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system .
virtual ~CVODESolver()
Destroy the associated CVODE memory and SUNDIALS objects.
static int root(sunrealtype t, N_Vector y, sunrealtype *gout, void *user_data)
Prototype to define root finding for CVODE.
void SetMaxNSteps(int steps)
Set the maximum number of time steps.
CVODESolver(int lmm)
Construct a serial wrapper to SUNDIALS' CVODE integrator.
long GetNumSteps()
Get the number of internal steps taken so far.
void Step(Vector &x, double &t, double &dt) override
Integrate the ODE with CVODE using the specified step mode.
EWTFunction ewt_func
A class member to facilitate pointing to a user-specified error weight function.
static int RHS(sunrealtype t, const N_Vector y, N_Vector ydot, void *user_data)
Number of components in gout.
void SetMaxStep(double dt_max)
Set the maximum time step.
int lmm_type
Linear multistep method type.
void PrintInfo() const
Print various CVODE statistics.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i....
RootFunction root_func
A class member to facilitate pointing to a user-specified root function.
std::function< int(Vector y, Vector w, CVODESolver *)> EWTFunction
Typedef declaration for error weight functions.
int step_mode
CVODE step mode (CV_NORMAL or CV_ONE_STEP).
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 .
void SetSVtolerances(double reltol, Vector abstol)
Set the scalar relative and vector of absolute tolerances.
void SetMaxOrder(int max_order)
Set the maximum method order.
static MemoryType GetHostMemoryType()
Get the current Host MemoryType. This is the MemoryType used by most MFEM classes when allocating mem...
static bool IsAvailable()
Return true if an actual device (e.g. GPU) has been configured.
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Wrapper for hypre's parallel vector class.
real_t abs_tol
Absolute tolerance.
real_t rel_tol
Relative tolerance.
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
int max_iter
Limit for the number of iterations the solver is allowed to do.
Interface to the KINSOL library – nonlinear solver methods.
SundialsNVector * f_scale
scaling vectors
KINSolver(int strategy, bool oper_grad=true)
Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.
double aa_damping
Anderson Acceleration damping.
void SetJFNKSolver(Solver &solver)
virtual ~KINSolver()
Destroy the associated KINSOL memory.
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Wrapper to compute the nonlinear residual .
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, sunrealtype tol)
Solve the linear system .
int aa_delay
Anderson Acceleration delay.
int global_strategy
KINSOL solution strategy.
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 .
int maxlrs
Maximum linear solver restarts.
const Operator * jacobian
stores oper->GetGradient()
int maxli
Maximum linear iterations.
Vector wrk
Work vector needed for the JFNK PC.
void SetScaledStepTol(double sstol)
Set KINSOL's scaled step tolerance.
void SetSolver(Solver &solver) override
Set the linear solver for inverting the Jacobian.
SundialsNVector * y_scale
void SetOperator(const Operator &op) override
Set the nonlinear Operator of the system and initialize KINSOL.
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
static int LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2)
Setup the linear system .
int aa_orth
Anderson Acceleration orthogonalization routine.
static int PrecSetup(N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, void *user_data)
Setup the preconditioner.
void SetDamping(double damping)
double fp_damping
Fixed Point or Picard damping parameter.
void EnableAndersonAcc(int n, int orth=KIN_ORTH_MGS, int delay=0, double damping=1.0)
Enable Anderson Acceleration for KIN_FP or KIN_PICARD.
void SetPrintLevel(int print_lvl) override
Set the print level for the KINSetPrintLevel function.
bool use_oper_grad
use the Jv prod function
int aa_n
number of acceleration vectors
static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u, sunbooleantype *new_u, void *user_data)
Wrapper to compute the Jacobian-vector product .
bool IsKnown(const void *h_ptr)
Return true if the pointer is known by the memory manager.
Class used by MFEM to store pointers to host and/or device memory.
void SetHostPtrOwner(bool own) const
Set/clear the ownership flag for the host pointer. Ownership indicates whether the pointer will be de...
void SetDevicePtrOwner(bool own) const
Set/clear the ownership flag for the device pointer. Ownership indicates whether the pointer will be ...
void Wrap(T *ptr, int size, bool own)
Wrap an externally allocated host pointer, ptr with the current host memory type returned by MemoryMa...
void Delete()
Delete the owned pointers and reset the Memory object.
void ClearOwnerFlags() const
Clear the ownership flags for the host and device pointers, as well as any internal data allocated by...
void SetOperator(const Operator &op) override
Also calls SetOperator for the preconditioner if there is one.
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
int height
Dimension of the output / number of rows in the matrix.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
virtual void SetOperator(const Operator &op)=0
Set/update the solver for the given operator.
SundialsMemHelper & operator=(const SundialsMemHelper &)=delete
Disable copy assignment.
SundialsMemHelper()=default
Default constructor – object must be moved to.
static int SundialsMemHelper_Alloc(SUNMemoryHelper helper, SUNMemory *memptr, size_t memsize, SUNMemoryType mem_type #if(SUNDIALS_VERSION_MAJOR >=6), void *queue #endif)
static int SundialsMemHelper_Dealloc(SUNMemoryHelper helper, SUNMemory sunmem #if(SUNDIALS_VERSION_MAJOR >=6), void *queue #endif)
Vector interface for SUNDIALS N_Vectors.
long GlobalSize() const
Returns the MPI global length for the internal N_Vector x.
MPI_Comm GetComm() const
Returns the MPI communicator for the internal N_Vector x.
void MakeRef(Vector &base, int offset, int s)
Reset the Vector to be a reference to a sub-vector of base.
void SetSize(int s, long glob_size=0)
Resize the vector to size s.
static bool UseManagedMemory()
static N_Vector MakeNVector(bool use_device)
Create a N_Vector.
~SundialsNVector()
Calls SUNDIALS N_VDestroy function if the N_Vector is owned by 'this'.
void SetDataAndSize(double *d, int s, long glob_size=0)
Set the vector data and size.
void _SetNvecDataAndSize_(long glob_size=0)
Set data and length of internal N_Vector x from 'this'.
void _SetDataAndSize_()
Set data and length from the internal N_Vector x.
N_Vector x
The actual SUNDIALS object.
SundialsNVector()
Creates an empty SundialsNVector.
static constexpr double default_abs_tol
Default scalar absolute tolerance.
SUNMatrix M
Mass matrix M.
int flag
Last flag returned from a call to SUNDIALS.
long saved_global_size
Global vector length on last initialization.
static constexpr double default_rel_tol
Default scalar relative tolerance.
bool reinit
Flag to signal memory reinitialization is need.
SundialsNVector * Y
State vector.
void * sundials_mem
SUNDIALS mem structure.
SUNLinearSolver LSA
Linear solver for A.
SUNLinearSolver LSM
Linear solver for M.
SUNNonlinearSolver NLS
Nonlinear solver.
Singleton class for SUNContext and SundialsMemHelper objects.
static SUNContext & GetContext()
Provides access to the SUNContext object.
static SundialsMemHelper & GetMemHelper()
Provides access to the SundialsMemHelper object.
static void Finalize()
Finalize sundials (called automatically at program exit if Sundials::Init() has been called).
static void Init()
Initializes SUNContext and SundialsMemHelper objects. Should be called at the beginning of the callin...
int GetAdjointHeight()
Returns the size of the adjoint problem state space.
Base abstract class for first order time dependent operators.
bool isExplicit() const
True if type is EXPLICIT.
virtual int SUNImplicitSolve(const Vector &r, Vector &dk, real_t tol)
Solve the ODE linear system A dk = r , where A and r are defined by the method SUNImplicitSetup().
virtual int SUNMassMult(const Vector &x, Vector &v)
Compute the mass matrix-vector product v = M x, where M is defined by the method SUNMassSetup().
virtual int SUNMassSetup()
Setup the mass matrix in the ODE system .
void Mult(const Vector &u, Vector &k) const override
Perform the action of the operator (u,t) -> k(u,t) where t is the current time set by SetTime() and k...
virtual int SUNMassSolve(const Vector &b, Vector &x, real_t tol)
Solve the mass matrix linear system M x = b, where M is defined by the method SUNMassSetup().
virtual void ExplicitMult(const Vector &u, Vector &v) const
Perform the action of the explicit part of the operator, G: v = G(u, t) where t is the current time.
virtual void SetEvalMode(const EvalMode new_eval_mode)
Set the evaluation mode of the time-dependent operator.
virtual void SetTime(const real_t t_)
Set the current time.
virtual int SUNImplicitSetup(const Vector &y, const Vector &v, int jok, int *jcur, real_t gamma)
Setup a linear system as needed by some SUNDIALS ODE solvers to perform a similar action to ImplicitS...
virtual real_t GetTime() const
Read the currently set time.
virtual const real_t * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
virtual real_t * ReadWrite(bool on_dev=true)
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), on_dev).
void SetDataAndSize(real_t *d, int s)
Set the Vector data and size.
real_t Normlinf() const
Returns the l_infinity norm of the vector.
Vector & Set(const real_t a, const Vector &x)
(*this) = a * x
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
int Size() const
Returns the size of the vector.
void SetSize(int s)
Resize the vector to size s.
virtual real_t * HostReadWrite()
Shortcut for mfem::ReadWrite(vec.GetMemory(), vec.Size(), false).
T * HostReadWrite(Memory< T > &mem, int size)
Shortcut to ReadWrite(Memory<T> &mem, int size, false)
real_t u(const Vector &xvec)
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
T * ReadWrite(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read+write access to mem with the mfem::Device's DeviceMemoryClass,...
MemoryManager mm
The (single) global memory manager object.
MFEM_HOST_DEVICE real_t norm(const Complex &z)
Settings for the output behavior of the IterativeSolver.
Helper struct to convert a C++ type to an MPI type.
MFEM_DEPRECATED SUNLinearSolver SUNLinSolNewEmpty(SUNContext)
MFEM_DEPRECATED N_Vector N_VMake_MPIPlusX(MPI_Comm comm, N_Vector local_vector, SUNContext)
MFEM_DEPRECATED void * KINCreate(SUNContext)
MFEM_DEPRECATED SUNLinearSolver SUNLinSol_SPFGMR(N_Vector y, int pretype, int maxl, SUNContext)
MFEM_DEPRECATED void * CVodeCreate(int lmm, SUNContext)
MFEM_DEPRECATED N_Vector N_VNewEmpty_Parallel(MPI_Comm comm, sunindextype local_length, sunindextype global_length, SUNContext)
MFEM_DEPRECATED SUNMatrix SUNMatNewEmpty(SUNContext)
MFEM_DEPRECATED SUNMemoryHelper SUNMemoryHelper_NewEmpty(SUNContext)
MFEM_DEPRECATED N_Vector SUN_Hip_OR_Cuda N_VNewWithMemHelp(sunindextype length, sunbooleantype use_managed_mem, SUNMemoryHelper helper, SUNContext)
MFEM_DEPRECATED void * ARKStepCreate(ARKRhsFn fe, ARKRhsFn fi, sunrealtype t0, N_Vector y0, SUNContext)
MFEM_DEPRECATED N_Vector N_VNewEmpty_Serial(sunindextype vec_length, SUNContext)
MFEM_DEPRECATED SUNLinearSolver SUNLinSol_SPGMR(N_Vector y, int pretype, int maxl, SUNContext)
realtype sunrealtype
'sunrealtype' was first introduced in v6.0.0
booleantype sunbooleantype
'sunbooleantype' was first introduced in v6.0.0
constexpr ARKODE_ERKTableID ARKODE_ERK_NONE
constexpr ARKODE_DIRKTableID ARKODE_DIRK_NONE