14 #ifdef MFEM_USE_SUNDIALS
22 #include <nvector/nvector_serial.h>
24 #include <nvector/nvector_cuda.h>
25 #include <sunmemory/sunmemory_cuda.h>
28 #include <nvector/nvector_mpiplusx.h>
29 #include <nvector/nvector_parallel.h>
33 #include <sunlinsol/sunlinsol_spgmr.h>
34 #include <sunlinsol/sunlinsol_spfgmr.h>
37 #define GET_CONTENT(X) ( X->content )
49 class SundialsMemHelper
55 friend class SundialsNVector;
61 h = SUNMemoryHelper_NewEmpty();
64 h->ops->alloc = SundialsMemHelper_Alloc;
65 h->ops->dealloc = SundialsMemHelper_Dealloc;
67 h->ops->copy = SUNMemoryHelper_Copy_Cuda;
68 h->ops->copyasync = SUNMemoryHelper_CopyAsync_Cuda;
74 SUNMemoryHelper_Destroy(h);
78 operator SUNMemoryHelper()
const {
return h; }
80 static int SundialsMemHelper_Alloc(SUNMemoryHelper helper,
83 SUNMemoryType mem_type)
85 int length = memsize/
sizeof(double);
86 SUNMemory sunmem = SUNMemoryNewEmpty();
89 sunmem->own = SUNTRUE;
91 if (mem_type == SUNMEMTYPE_HOST)
93 Memory<double> mem(length, Device::GetHostMemoryType());
94 mem.SetHostPtrOwner(
false);
96 sunmem->type = SUNMEMTYPE_HOST;
99 else if (mem_type == SUNMEMTYPE_DEVICE || mem_type == SUNMEMTYPE_UVM)
101 Memory<double> mem(length, Device::GetDeviceMemoryType());
102 mem.SetDevicePtrOwner(
false);
104 sunmem->type = mem_type;
117 static int SundialsMemHelper_Dealloc(SUNMemoryHelper helper, SUNMemory sunmem)
119 if (sunmem->ptr && sunmem->own && !
mm.
IsKnown(sunmem->ptr))
121 if (sunmem->type == SUNMEMTYPE_HOST)
123 Memory<double> mem(static_cast<double*>(sunmem->ptr), 1,
124 Device::GetHostMemoryType(),
true);
127 else if (sunmem->type == SUNMEMTYPE_DEVICE || sunmem->type == SUNMEMTYPE_UVM)
129 Memory<double> mem(static_cast<double*>(sunmem->ptr), 1,
130 Device::GetDeviceMemoryType(),
true);
135 MFEM_ABORT(
"Invalid SUNMEMTYPE");
152 void SundialsNVector::_SetNvecDataAndSize_(
long glob_size)
155 N_Vector local_x = MPIPlusX() ? N_VGetLocalVector_MPIPlusX(x) : x;
157 N_Vector local_x = x;
159 N_Vector_ID
id = N_VGetVectorID(local_x);
164 case SUNDIALS_NVEC_SERIAL:
166 MFEM_ASSERT(NV_OWN_DATA_S(local_x) == SUNFALSE,
"invalid serial N_Vector");
168 NV_LENGTH_S(local_x) = size;
172 case SUNDIALS_NVEC_CUDA:
175 N_VSetDeviceArrayPointer_Cuda(
ReadWrite(), local_x);
176 static_cast<N_VectorContent_Cuda
>(GET_CONTENT(local_x))->length = size;
181 case SUNDIALS_NVEC_PARALLEL:
183 MFEM_ASSERT(NV_OWN_DATA_P(x) == SUNFALSE,
"invalid parallel N_Vector");
185 NV_LOCLENGTH_P(x) = size;
188 glob_size = GlobalSize();
190 if (glob_size == 0 && glob_size != size)
192 long local_size = size;
193 MPI_Allreduce(&local_size, &glob_size, 1, MPI_LONG,
197 NV_GLOBLENGTH_P(x) = glob_size;
202 MFEM_ABORT(
"N_Vector type " <<
id <<
" is not supported");
210 glob_size = GlobalSize();
212 if (glob_size == 0 && glob_size != size)
214 long local_size = size;
215 MPI_Allreduce(&local_size, &glob_size, 1, MPI_LONG,
219 static_cast<N_VectorContent_MPIManyVector
>(GET_CONTENT(x))->global_length =
225 void SundialsNVector::_SetDataAndSize_()
228 N_Vector local_x = MPIPlusX() ? N_VGetLocalVector_MPIPlusX(x) : x;
230 N_Vector local_x = x;
232 N_Vector_ID
id = N_VGetVectorID(local_x);
237 case SUNDIALS_NVEC_SERIAL:
239 const bool known =
mm.
IsKnown(NV_DATA_S(local_x));
240 size = NV_LENGTH_S(local_x);
241 data.Wrap(NV_DATA_S(local_x), size,
false);
242 if (known) { data.ClearOwnerFlags(); }
246 case SUNDIALS_NVEC_CUDA:
248 double *h_ptr = N_VGetHostArrayPointer_Cuda(local_x);
249 double *d_ptr = N_VGetDeviceArrayPointer_Cuda(local_x);
251 size = N_VGetLength_Cuda(local_x);
252 data.Wrap(h_ptr, d_ptr, size, Device::GetHostMemoryType(),
false);
253 if (known) { data.ClearOwnerFlags(); }
259 case SUNDIALS_NVEC_PARALLEL:
261 const bool known =
mm.
IsKnown(NV_DATA_P(x));
262 size = NV_LENGTH_S(x);
263 data.Wrap(NV_DATA_P(x), NV_LOCLENGTH_P(x),
false);
264 if (known) { data.ClearOwnerFlags(); }
269 MFEM_ABORT(
"N_Vector type " <<
id <<
" is not supported");
273 SundialsNVector::SundialsNVector()
328 :
SundialsNVector(vec.GetComm(), vec.GetData(), vec.Size(), vec.GlobalSize())
339 N_VDestroy(N_VGetLocalVector_MPIPlusX(
x));
374 x = N_VNewEmpty_Serial(0);
377 x = N_VNewEmpty_Serial(0);
380 MFEM_VERIFY(x,
"Error in SundialsNVector::MakeNVector.");
390 if (comm == MPI_COMM_NULL)
404 x = N_VNewEmpty_Parallel(comm, 0, 0);
407 x = N_VNewEmpty_Parallel(comm, 0, 0);
408 #endif // MFEM_USE_CUDA
411 MFEM_VERIFY(x,
"Error in SundialsNVector::MakeNVector.");
415 #endif // MFEM_USE_MPI
423 static SUNMatrix_ID MatGetID(SUNMatrix)
425 return (SUNMATRIX_CUSTOM);
428 static void MatDestroy(SUNMatrix A)
430 if (A->content) { A->content = NULL; }
431 if (A->ops) { free(A->ops); A->ops = NULL; }
441 static SUNLinearSolver_Type LSGetType(SUNLinearSolver)
443 return (SUNLINEARSOLVER_MATRIX_ITERATIVE);
446 static int LSFree(SUNLinearSolver LS)
448 if (LS->content) { LS->content = NULL; }
449 if (LS->ops) { free(LS->ops); LS->ops = NULL; }
468 self->f->Mult(mfem_y, mfem_ydot);
478 if (!self->root_func) {
return CV_RTFUNC_FAIL; }
483 return self->root_func(t, mfem_y, mfem_gout,
self);
491 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in SetRootFinder()");
495 booleantype jok, booleantype *jcur, realtype gamma,
496 void*, N_Vector, N_Vector, N_Vector)
505 return (self->f->SUNImplicitSetup(mfem_y, mfem_fy, jok, jcur, gamma));
509 N_Vector
b, realtype tol)
515 return (self->f->SUNImplicitSolve(mfem_b, mfem_x, tol));
519 : lmm_type(lmm), step_mode(CV_NORMAL)
526 : lmm_type(lmm), step_mode(CV_NORMAL)
538 long local_size = f_.
Height();
541 long global_size = 0;
544 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
558 resize = (
Y->
Size() != local_size);
563 int l_resize = (
Y->
Size() != local_size) ||
565 MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
591 Y->
SetSize(local_size, global_size);
602 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeInit()");
606 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetUserData()");
610 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetSStolerances()");
623 MFEM_VERIFY(
Y->
Size() == x.
Size(),
"size mismatch");
629 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeReInit()");
635 double tout = t + dt;
637 MFEM_VERIFY(
flag >= 0,
"error in CVode()");
644 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetLastStep()");
650 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
651 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
654 LSA = SUNLinSolNewEmpty();
655 MFEM_VERIFY(
LSA,
"error in SUNLinSolNewEmpty()");
658 LSA->ops->gettype = LSGetType;
660 LSA->ops->free = LSFree;
662 A = SUNMatNewEmpty();
663 MFEM_VERIFY(A,
"error in SUNMatNewEmpty()");
666 A->ops->getid = MatGetID;
667 A->ops->destroy = MatDestroy;
671 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolver()");
675 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinSysFn()");
681 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
682 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
685 LSA = SUNLinSol_SPGMR(*
Y, PREC_NONE, 0);
686 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPGMR()");
690 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolver()");
701 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSStolerances()");
707 "abs tolerance is not the same size.");
713 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSVtolerances()");
719 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetMaxStep()");
725 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetMaxNumSteps()");
732 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetNumSteps()");
739 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetMaxOrd()");
744 long int nsteps, nfevals, nlinsetups, netfails;
746 double hinused, hlast, hcur, tcur;
747 long int nniters, nncfails;
761 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetIntegratorStats()");
767 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetNonlinSolvStats()");
771 "num steps: " << nsteps <<
"\n"
772 "num rhs evals: " << nfevals <<
"\n"
773 "num lin setups: " << nlinsetups <<
"\n"
774 "num nonlin sol iters: " << nniters <<
"\n"
775 "num nonlin conv fail: " << nncfails <<
"\n"
776 "num error test fails: " << netfails <<
"\n"
777 "last order: " << qlast <<
"\n"
778 "current order: " << qcur <<
"\n"
779 "initial dt: " << hinused <<
"\n"
780 "last dt: " << hlast <<
"\n"
781 "current dt: " << hcur <<
"\n"
782 "current t: " << tcur <<
"\n" << endl;
792 SUNNonlinSolFree(
NLS);
830 MFEM_VERIFY(t <= f->GetTime(),
"t > current forward solver time");
833 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetQuad()");
840 MFEM_VERIFY(t <= f->GetTime(),
"t > current forward solver time");
843 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetQuadB()");
853 MFEM_VERIFY(
flag >= 0,
"error in CVodeGetAdjY()");
873 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeCreateB()");
877 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeInit()");
881 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetUserDataB()");
886 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetSStolerancesB()");
898 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeAdjInit");
904 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeSetMaxNumStepsB()");
913 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadInit()");
916 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeSetQuadErrCon");
919 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadSStolerances");
928 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadInitB()");
931 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeSetQuadErrConB");
934 MFEM_VERIFY(
flag == CV_SUCCESS,
"Error in CVodeQuadSStolerancesB");
940 if (
AB != NULL) { SUNMatDestroy(
AB);
AB = NULL; }
941 if (
LSB != NULL) { SUNLinSolFree(
LSB);
LSB = NULL; }
944 LSB = SUNLinSolNewEmpty();
945 MFEM_VERIFY(
LSB,
"error in SUNLinSolNewEmpty()");
948 LSB->ops->gettype = LSGetType;
950 LSB->ops->free = LSFree;
952 AB = SUNMatNewEmpty();
953 MFEM_VERIFY(
AB,
"error in SUNMatNewEmpty()");
956 AB->ops->getid = MatGetID;
957 AB->ops->destroy = MatDestroy;
961 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolverB()");
966 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinSysFn()");
972 if (
AB != NULL) { SUNMatDestroy(
AB);
AB = NULL; }
973 if (
LSB != NULL) { SUNLinSolFree(
LSB);
LSB = NULL; }
976 LSB = SUNLinSol_SPGMR(*
yB, PREC_NONE, 0);
977 MFEM_VERIFY(
LSB,
"error in SUNLinSol_SPGMR()");
981 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSetLinearSolverB()");
985 N_Vector fyB, SUNMatrix AB,
986 booleantype jokB, booleantype *jcurB,
987 realtype gammaB,
void *user_data, N_Vector tmp1,
988 N_Vector tmp2, N_Vector tmp3)
1004 N_Vector Rb, realtype tol)
1019 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSStolerancesB()");
1025 "abs tolerance is not the same size.");
1031 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeSVtolerancesB()");
1092 return self->ewt_func(mfem_y, mfem_w,
self);
1099 MFEM_VERIFY(
Y->
Size() == x.
Size(),
"size mismatch");
1105 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeReInit()");
1112 double tout = t + dt;
1114 MFEM_VERIFY(
flag >= 0,
"error in CVodeF()");
1121 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeGetLastStep()");
1133 MFEM_VERIFY(
flag == CV_SUCCESS,
"error in CVodeReInit()");
1140 double tout = tB - dtB;
1142 MFEM_VERIFY(
flag >= 0,
"error in CVodeB()");
1146 MFEM_VERIFY(
flag >= 0,
"error in CVodeGetB()");
1177 if (self->rk_type ==
IMEX)
1181 self->f->Mult(mfem_y, mfem_ydot);
1198 self->f->Mult(mfem_y, mfem_ydot);
1205 SUNMatrix, booleantype jok, booleantype *jcur,
1207 void*, N_Vector, N_Vector, N_Vector)
1216 if (self->rk_type ==
IMEX)
1220 return (self->f->SUNImplicitSetup(mfem_y, mfem_fy, jok, jcur, gamma));
1224 N_Vector
b, realtype tol)
1231 if (self->rk_type ==
IMEX)
1235 return (self->f->SUNImplicitSolve(mfem_b, mfem_x, tol));
1239 void*, N_Vector, N_Vector, N_Vector)
1245 return (self->f->SUNMassSetup());
1249 N_Vector
b, realtype tol)
1256 return (self->f->SUNMassSolve(mfem_b, mfem_x, tol));
1266 return (self->f->SUNMassMult(mfem_x, mfem_v));
1278 return (self->f->SUNMassMult(mfem_x, mfem_v));
1282 : rk_type(type), step_mode(ARK_NORMAL),
1283 use_implicit(type == IMPLICIT || type == IMEX)
1290 : rk_type(type), step_mode(ARK_NORMAL),
1291 use_implicit(type == IMPLICIT || type == IMEX)
1303 long local_size = f_.
Height();
1311 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
1325 resize = (
Y->
Size() != local_size);
1330 int l_resize = (
Y->
Size() != local_size) ||
1332 MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
1354 Y->
SetSize(local_size, global_size);
1377 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetUserData()");
1381 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetSStolerances()");
1394 MFEM_VERIFY(
Y->
Size() == x.
Size(),
"size mismatch");
1412 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepReInit()");
1419 double tout = t + dt;
1421 MFEM_VERIFY(
flag >= 0,
"error in ARKStepEvolve()");
1428 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepGetLastStep()");
1434 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1435 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1438 LSA = SUNLinSolNewEmpty();
1439 MFEM_VERIFY(
LSA,
"error in SUNLinSolNewEmpty()");
1441 LSA->content =
this;
1442 LSA->ops->gettype = LSGetType;
1444 LSA->ops->free = LSFree;
1446 A = SUNMatNewEmpty();
1447 MFEM_VERIFY(A,
"error in SUNMatNewEmpty()");
1450 A->ops->getid = MatGetID;
1451 A->ops->destroy = MatDestroy;
1455 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetLinearSolver()");
1459 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetLinSysFn()");
1465 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1466 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1469 LSA = SUNLinSol_SPGMR(*
Y, PREC_NONE, 0);
1470 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPGMR()");
1474 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetLinearSolver()");
1480 if (
M != NULL) { SUNMatDestroy(
M);
M = NULL; }
1481 if (
LSM != NULL) { SUNLinSolFree(
LSM);
LSM = NULL; }
1484 LSM = SUNLinSolNewEmpty();
1485 MFEM_VERIFY(
LSM,
"error in SUNLinSolNewEmpty()");
1487 LSM->content =
this;
1488 LSM->ops->gettype = LSGetType;
1490 LSA->ops->free = LSFree;
1492 M = SUNMatNewEmpty();
1493 MFEM_VERIFY(
M,
"error in SUNMatNewEmpty()");
1496 M->ops->getid = SUNMatGetID;
1498 M->ops->destroy = MatDestroy;
1502 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetLinearSolver()");
1506 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetMassFn()");
1512 if (
M != NULL) { SUNMatDestroy(A);
M = NULL; }
1513 if (
LSM != NULL) { SUNLinSolFree(
LSM);
LSM = NULL; }
1516 LSM = SUNLinSol_SPGMR(*
Y, PREC_NONE, 0);
1517 MFEM_VERIFY(
LSM,
"error in SUNLinSol_SPGMR()");
1521 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetMassLinearSolver()");
1526 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetMassTimes()");
1537 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSStolerances()");
1543 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetMaxStep()");
1549 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetOrder()");
1555 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetTableNum()");
1561 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetTableNum()");
1567 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetTableNum()");
1573 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepSetFixedStep()");
1578 long int nsteps, expsteps, accsteps, step_attempts;
1579 long int nfe_evals, nfi_evals;
1580 long int nlinsetups, netfails;
1581 double hinused, hlast, hcur, tcur;
1582 long int nniters, nncfails;
1593 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepGetTimestepperStats()");
1606 MFEM_VERIFY(
flag == ARK_SUCCESS,
"error in ARKStepGetNonlinSolvStats()");
1610 "num steps: " << nsteps <<
"\n"
1611 "num exp rhs evals: " << nfe_evals <<
"\n"
1612 "num imp rhs evals: " << nfi_evals <<
"\n"
1613 "num lin setups: " << nlinsetups <<
"\n"
1614 "num nonlin sol iters: " << nniters <<
"\n"
1615 "num nonlin conv fail: " << nncfails <<
"\n"
1616 "num steps attempted: " << step_attempts <<
"\n"
1617 "num acc limited steps: " << accsteps <<
"\n"
1618 "num exp limited stepfails: " << expsteps <<
"\n"
1619 "num error test fails: " << netfails <<
"\n"
1620 "initial dt: " << hinused <<
"\n"
1621 "last dt: " << hlast <<
"\n"
1622 "current dt: " << hcur <<
"\n"
1623 "current t: " << tcur <<
"\n" << endl;
1633 SUNNonlinSolFree(
NLS);
1657 booleantype *new_u,
void *user_data)
1667 self->jacobian = &
self->oper->GetGradient(mfem_u);
1672 self->jacobian->Mult(mfem_v, mfem_Jv);
1680 void *, N_Vector, N_Vector )
1689 self->prec->SetOperator(*self->jacobian);
1697 N_Vector
b, realtype)
1722 self->prec->SetOperator(*self->jacobian);
1740 self->prec->Mult(mfem_v, self->wrk);
1748 : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
1749 f_scale(NULL), jacobian(NULL), maa(0)
1762 : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
1763 f_scale(NULL), jacobian(NULL), maa(0)
1783 long local_size =
height;
1791 MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
1802 resize = (
Y->
Size() != local_size);
1807 int l_resize = (
Y->
Size() != local_size) ||
1809 MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
1831 Y->
SetSize(local_size, global_size);
1846 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetMAA()");
1851 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINInit()");
1855 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetUserData()");
1865 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1866 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1868 LSA = SUNLinSol_SPGMR(*
Y, PREC_NONE, 0);
1869 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPGMR()");
1872 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetLinearSolver()");
1878 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetJacTimesVecFn()");
1896 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1897 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1900 LSA = SUNLinSolNewEmpty();
1901 MFEM_VERIFY(
LSA,
"error in SUNLinSolNewEmpty()");
1903 LSA->content =
this;
1904 LSA->ops->gettype = LSGetType;
1906 LSA->ops->free = LSFree;
1908 A = SUNMatNewEmpty();
1909 MFEM_VERIFY(A,
"error in SUNMatNewEmpty()");
1912 A->ops->getid = MatGetID;
1913 A->ops->destroy = MatDestroy;
1917 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetLinearSolver()");
1921 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetJacFn()");
1933 if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1934 if (
LSA != NULL) { SUNLinSolFree(
LSA);
LSA = NULL; }
1937 LSA = SUNLinSol_SPFGMR(*
Y,
prec ? PREC_RIGHT : PREC_NONE,
maxli);
1938 MFEM_VERIFY(
LSA,
"error in SUNLinSol_SPFGMR()");
1941 MFEM_VERIFY(
flag == SUNLS_SUCCESS,
"error in SUNLinSol_SPFGMR()");
1944 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetLinearSolver()");
1951 MFEM_VERIFY(
flag == KIN_SUCCESS,
"error in KINSetPreconditioner()");
1958 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetScaledStepTol()");
1964 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetMaxSetupCalls()");
1975 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetMAA()");
1981 MFEM_ABORT(
"this method is not supported! Use SetPrintLevel(int) instead.");
2006 double lnorm = norm;
2007 MPI_Allreduce(&lnorm, &norm, 1, MPI_DOUBLE, MPI_MAX,
2030 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINSetFuncNormTol()");
2041 MFEM_ASSERT(
flag == KIN_SUCCESS,
"KINSetNumMaxIters() failed!");
2055 MPI_Comm_rank(
Y->
GetComm(), &rank);
2062 MFEM_VERIFY(
flag == KIN_SUCCESS,
"KINSetPrintLevel() failed!");
2064 #ifdef SUNDIALS_BUILD_WITH_MONITORING
2067 flag = SUNLinSolSetInfoFile_SPFGMR(
LSA, stdout);
2068 MFEM_VERIFY(
flag == SUNLS_SUCCESS,
2069 "error in SUNLinSolSetInfoFile_SPFGMR()");
2071 flag = SUNLinSolSetPrintLevel_SPFGMR(
LSA, 1);
2072 MFEM_VERIFY(
flag == SUNLS_SUCCESS,
2073 "error in SUNLinSolSetPrintLevel_SPFGMR()");
2090 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINGetNumNonlinSolvIters()");
2095 MFEM_ASSERT(
flag == KIN_SUCCESS,
"error in KINGetFuncNorm()");
2110 #endif // MFEM_USE_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.
virtual double GetTime() const
Read the currently set time.
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
static bool IsAvailable()
Return true if an actual device (e.g. GPU) has been configured.
void SetJFNKSolver(Solver &solver)
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.
virtual int SUNImplicitSolveB(Vector &x, const Vector &b, double tol)
Solve the ODE linear system as setup by the method SUNImplicitSetup().
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.
virtual void QuadratureSensitivityMult(const Vector &y, const Vector &yB, Vector &qBdot) const
Provides the sensitivity of the quadrature w.r.t to primal and adjoint solutions. ...
~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.
void SetSize(int s)
Resize the vector to size s.
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 Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
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.
virtual void AdjointRateMult(const Vector &y, Vector &yB, Vector &yBdot) const =0
Perform the action of the operator: yBdot = k = f(y,@2 yB, t), where.
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.
int indexB
backward problem index
virtual ~CVODESolver()
Destroy the associated CVODE memory and SUNDIALS objects.
bool use_oper_grad
use the Jv prod function
int Size() const
Returns the size of the vector.
void SetIRKTableNum(int table_num)
Choose a specific Butcher table for a diagonally implicit RK method.
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
void SetIMEXTableNum(int etable_num, int itable_num)
Choose a specific Butcher table for an IMEX RK method.
virtual void SetPrintLevel(int print_lvl)
Set the print level for the KINSetPrintLevel function.
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.
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
virtual void SetTime(const double t_)
Set the current time.
void * sundials_mem
SUNDIALS mem structure.
double Normlinf() const
Returns the l_infinity norm of the vector.
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.
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
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.
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
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.
SundialsMemHelper sunmemHelper
SUNLinearSolver LSA
Linear solver for A.
static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
Compute the matrix-vector product .
virtual void SetEvalMode(const EvalMode new_eval_mode)
Set the evaluation mode of the time-dependent operator.
SUNMatrix M
Mass matrix M.
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
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 constexpr double default_rel_tolB
Default scalar backward relative tolerance.
SundialsNVector * yy
State vector.
Vector interface for SUNDIALS N_Vectors.
int GetAdjointHeight()
Returns the size of the adjoint problem state space.
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.
SundialsNVector * qB
State vector.
bool use_implicit
True for implicit or imex integration.
SUNNonlinearSolver NLS
Nonlinear solver.
CVODESolver(int lmm)
Construct a serial wrapper to SUNDIALS' CVODE integrator.
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
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).
FDualNumber< tbase > pow(const FDualNumber< tbase > &a, const FDualNumber< tbase > &b)
pow([dual number],[dual number])
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.
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.
int max_iter
Limit for the number of iterations the solver is allowed to do.
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.
bool IsKnown(const void *h_ptr)
Return true if the pointer is known by the memory manager.
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.
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
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'.
Vector & Set(const double a, const Vector &x)
(*this) = a * x
int height
Dimension of the output / number of rows in the matrix.
void SetMAA(int maa)
Set the number of acceleration vectors to use with KIN_FP or KIN_PICARD.
virtual void QuadratureIntegration(const Vector &y, Vector &qdot) const
Provide the operator integration of a quadrature equation.
int maxli
Maximum linear iterations.
MemoryManager mm
The (single) global memory manager object.
long GetNumSteps()
Get the number of internal steps taken so far.
double rel_tol
Relative tolerance.
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, if on_dev = true, or the mfem::Device's HostMemoryClass, otherwise.
MPI_Comm GetComm() const
Returns the MPI communicator for the internal N_Vector x.
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.
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 .
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.
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.
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
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_)
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
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.
double u(const Vector &xvec)
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.
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
T * HostReadWrite(Memory< T > &mem, int size)
Shortcut to ReadWrite(Memory<T> &mem, int size, false)
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.
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
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 int SUNImplicitSetupB(const double t, const Vector &x, const Vector &xB, const Vector &fxB, int jokB, int *jcurB, double gammaB)
Setup the ODE linear system or , where .
Settings for the output behavior of the IterativeSolver.
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.
double abs_tol
Absolute tolerance.
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
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...