MFEM  v4.5.1
Finite element discretization library
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Pages
sundials.cpp
Go to the documentation of this file.
1 // Copyright (c) 2010-2022, Lawrence Livermore National Security, LLC. Produced
2 // at the Lawrence Livermore National Laboratory. All Rights reserved. See files
3 // LICENSE and NOTICE for details. LLNL-CODE-806117.
4 //
5 // This file is part of the MFEM library. For more information and source code
6 // availability visit https://mfem.org.
7 //
8 // MFEM is free software; you can redistribute it and/or modify it under the
9 // terms of the BSD-3 license. We welcome feedback and contributions, see file
10 // CONTRIBUTING.md for details.
11 
12 #include "sundials.hpp"
13 
14 #ifdef MFEM_USE_SUNDIALS
15 
16 #include "solvers.hpp"
17 #ifdef MFEM_USE_MPI
18 #include "hypre.hpp"
19 #endif
20 
21 // SUNDIALS vectors
22 #include <nvector/nvector_serial.h>
23 #ifdef MFEM_USE_CUDA
24 #include <nvector/nvector_cuda.h>
25 #endif
26 #ifdef MFEM_USE_MPI
27 #include <nvector/nvector_mpiplusx.h>
28 #include <nvector/nvector_parallel.h>
29 #endif
30 
31 // SUNDIALS linear solvers
32 #include <sunlinsol/sunlinsol_spgmr.h>
33 #include <sunlinsol/sunlinsol_spfgmr.h>
34 
35 // Access SUNDIALS object's content pointer
36 #define GET_CONTENT(X) ( X->content )
37 
38 using namespace std;
39 
40 #if (SUNDIALS_VERSION_MAJOR < 6)
41 
42 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
43 /// version < 6
44 MFEM_DEPRECATED N_Vector N_VNewEmpty_Serial(sunindextype vec_length, SUNContext)
45 {
46  return N_VNewEmpty_Serial(vec_length);
47 }
48 
49 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
50 /// version < 6
51 MFEM_DEPRECATED SUNMatrix SUNMatNewEmpty(SUNContext)
52 {
53  return SUNMatNewEmpty();
54 }
55 
56 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
57 /// version < 6
58 MFEM_DEPRECATED SUNLinearSolver SUNLinSolNewEmpty(SUNContext)
59 {
60  return SUNLinSolNewEmpty();
61 }
62 
63 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
64 /// version < 6
65 MFEM_DEPRECATED SUNLinearSolver SUNLinSol_SPGMR(N_Vector y, int pretype,
66  int maxl, SUNContext)
67 {
68  return SUNLinSol_SPGMR(y, pretype, maxl);
69 }
70 
71 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
72 /// version < 6
73 MFEM_DEPRECATED SUNLinearSolver SUNLinSol_SPFGMR(N_Vector y, int pretype,
74  int maxl, SUNContext)
75 {
76  return SUNLinSol_SPFGMR(y, pretype, maxl);
77 }
78 
79 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
80 /// version < 6
81 MFEM_DEPRECATED void* CVodeCreate(int lmm, SUNContext)
82 {
83  return CVodeCreate(lmm);
84 }
85 
86 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
87 /// version < 6
88 MFEM_DEPRECATED void* ARKStepCreate(ARKRhsFn fe, ARKRhsFn fi, realtype t0,
89  N_Vector y0, SUNContext)
90 {
91  return ARKStepCreate(fe, fi, t0, y0);
92 }
93 
94 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
95 /// version < 6
96 MFEM_DEPRECATED void* KINCreate(SUNContext)
97 {
98  return KINCreate();
99 }
100 
101 #ifdef MFEM_USE_MPI
102 
103 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
104 /// version < 6
105 MFEM_DEPRECATED N_Vector N_VNewEmpty_Parallel(MPI_Comm comm,
106  sunindextype local_length,
107  sunindextype global_length,
108  SUNContext)
109 {
110  return N_VNewEmpty_Parallel(comm, local_length, global_length);
111 }
112 
113 #endif // MFEM_USE_MPI
114 
115 #ifdef MFEM_USE_CUDA
116 
117 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
118 /// version < 6
119 MFEM_DEPRECATED N_Vector N_VNewWithMemHelp_Cuda(sunindextype length,
120  booleantype use_managed_mem,
121  SUNMemoryHelper helper,
122  SUNContext)
123 {
124  return N_VNewWithMemHelp_Cuda(length, use_managed_mem, helper);
125 }
126 
127 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
128 /// version < 6
129 MFEM_DEPRECATED SUNMemoryHelper SUNMemoryHelper_NewEmpty(SUNContext)
130 {
131  return SUNMemoryHelper_NewEmpty();
132 }
133 
134 #endif // MFEM_USE_CUDA
135 
136 #if defined(MFEM_USE_MPI) && defined(MFEM_USE_CUDA)
137 
138 /// (DEPRECATED) Wrapper function for backwards compatibility with SUNDIALS
139 /// version < 6
140 MFEM_DEPRECATED N_Vector N_VMake_MPIPlusX(MPI_Comm comm, N_Vector local_vector,
141  SUNContext)
142 {
143  return N_VMake_MPIPlusX(comm, local_vector);
144 }
145 
146 #endif // MFEM_USE_MPI && MFEM_USE_CUDA
147 
148 #endif // SUNDIALS_VERSION_MAJOR < 6
149 
150 
151 namespace mfem
152 {
153 
154 void Sundials::Init()
155 {
156  Sundials::Instance();
157 }
158 
159 Sundials &Sundials::Instance()
160 {
161  static Sundials sundials;
162  return sundials;
163 }
164 
165 SUNContext &Sundials::GetContext()
166 {
167  return Sundials::Instance().context;
168 }
169 
170 SundialsMemHelper &Sundials::GetMemHelper()
171 {
172  return Sundials::Instance().memHelper;
173 }
174 
175 #if (SUNDIALS_VERSION_MAJOR >= 6)
176 
177 Sundials::Sundials()
178 {
179 #ifdef MFEM_USE_MPI
180  MPI_Comm communicator = MPI_COMM_WORLD;
181  int return_val = SUNContext_Create((void*) &communicator, &context);
182 #else
183  int return_val = SUNContext_Create(nullptr, &context);
184 #endif
185  MFEM_VERIFY(return_val == 0, "Call to SUNContext_Create failed");
186  SundialsMemHelper actual_helper(context);
187  memHelper = std::move(actual_helper);
188 }
189 
190 Sundials::~Sundials()
191 {
192  SUNContext_Free(&context);
193 }
194 
195 #else // SUNDIALS_VERSION_MAJOR >= 6
196 
197 Sundials::Sundials()
198 {
199  // Do nothing
200 }
201 
202 Sundials::~Sundials()
203 {
204  // Do nothing
205 }
206 
207 #endif // SUNDIALS_VERSION_MAJOR >= 6
208 
209 #ifdef MFEM_USE_CUDA
210 SundialsMemHelper::SundialsMemHelper(SUNContext context)
211 {
212  /* Allocate helper */
213  h = SUNMemoryHelper_NewEmpty(context);
214 
215  /* Set the ops */
216  h->ops->alloc = SundialsMemHelper_Alloc;
217  h->ops->dealloc = SundialsMemHelper_Dealloc;
218  h->ops->copy = SUNMemoryHelper_Copy_Cuda;
219  h->ops->copyasync = SUNMemoryHelper_CopyAsync_Cuda;
220 }
221 
222 SundialsMemHelper::SundialsMemHelper(SundialsMemHelper&& that_helper)
223 {
224  this->h = that_helper.h;
225  that_helper.h = nullptr;
226 }
227 
228 SundialsMemHelper& SundialsMemHelper::operator=(SundialsMemHelper&& rhs)
229 {
230  this->h = rhs.h;
231  rhs.h = nullptr;
232  return *this;
233 }
234 
235 int SundialsMemHelper::SundialsMemHelper_Alloc(SUNMemoryHelper helper,
236  SUNMemory* memptr, size_t memsize,
237  SUNMemoryType mem_type
238 #if (SUNDIALS_VERSION_MAJOR >= 6)
239  , void*
240 #endif
241  )
242 {
243  int length = memsize/sizeof(double);
244  SUNMemory sunmem = SUNMemoryNewEmpty();
245 
246  sunmem->ptr = NULL;
247  sunmem->own = SUNTRUE;
248 
249  if (mem_type == SUNMEMTYPE_HOST)
250  {
251  Memory<double> mem(length, Device::GetHostMemoryType());
252  mem.SetHostPtrOwner(false);
253  sunmem->ptr = mfem::HostReadWrite(mem, length);
254  sunmem->type = SUNMEMTYPE_HOST;
255  mem.Delete();
256  }
257  else if (mem_type == SUNMEMTYPE_DEVICE || mem_type == SUNMEMTYPE_UVM)
258  {
259  Memory<double> mem(length, Device::GetDeviceMemoryType());
260  mem.SetDevicePtrOwner(false);
261  sunmem->ptr = mfem::ReadWrite(mem, length);
262  sunmem->type = mem_type;
263  mem.Delete();
264  }
265  else
266  {
267  free(sunmem);
268  return -1;
269  }
270 
271  *memptr = sunmem;
272  return 0;
273 }
274 
275 int SundialsMemHelper::SundialsMemHelper_Dealloc(SUNMemoryHelper helper,
276  SUNMemory sunmem
277 #if (SUNDIALS_VERSION_MAJOR >= 6)
278  , void*
279 #endif
280  )
281 {
282  if (sunmem->ptr && sunmem->own && !mm.IsKnown(sunmem->ptr))
283  {
284  if (sunmem->type == SUNMEMTYPE_HOST)
285  {
286  Memory<double> mem(static_cast<double*>(sunmem->ptr), 1,
287  Device::GetHostMemoryType(), true);
288  mem.Delete();
289  }
290  else if (sunmem->type == SUNMEMTYPE_DEVICE || sunmem->type == SUNMEMTYPE_UVM)
291  {
292  Memory<double> mem(static_cast<double*>(sunmem->ptr), 1,
293  Device::GetDeviceMemoryType(), true);
294  mem.Delete();
295  }
296  else
297  {
298  MFEM_ABORT("Invalid SUNMEMTYPE");
299  return -1;
300  }
301  }
302  free(sunmem);
303  return 0;
304 }
305 
306 #endif // MFEM_USE_CUDA
307 
308 
309 // ---------------------------------------------------------------------------
310 // SUNDIALS N_Vector interface functions
311 // ---------------------------------------------------------------------------
312 
313 void SundialsNVector::_SetNvecDataAndSize_(long glob_size)
314 {
315 #ifdef MFEM_USE_MPI
316  N_Vector local_x = MPIPlusX() ? N_VGetLocalVector_MPIPlusX(x) : x;
317 #else
318  N_Vector local_x = x;
319 #endif
320  N_Vector_ID id = N_VGetVectorID(local_x);
321 
322  // Set the N_Vector data and length from the Vector data and size.
323  switch (id)
324  {
325  case SUNDIALS_NVEC_SERIAL:
326  {
327  MFEM_ASSERT(NV_OWN_DATA_S(local_x) == SUNFALSE, "invalid serial N_Vector");
328  NV_DATA_S(local_x) = HostReadWrite();
329  NV_LENGTH_S(local_x) = size;
330  break;
331  }
332 #ifdef MFEM_USE_CUDA
333  case SUNDIALS_NVEC_CUDA:
334  {
335  N_VSetHostArrayPointer_Cuda(HostReadWrite(), local_x);
336  N_VSetDeviceArrayPointer_Cuda(ReadWrite(), local_x);
337  static_cast<N_VectorContent_Cuda>(GET_CONTENT(local_x))->length = size;
338  break;
339  }
340 #endif
341 #ifdef MFEM_USE_MPI
342  case SUNDIALS_NVEC_PARALLEL:
343  {
344  MFEM_ASSERT(NV_OWN_DATA_P(x) == SUNFALSE, "invalid parallel N_Vector");
345  NV_DATA_P(x) = HostReadWrite();
346  NV_LOCLENGTH_P(x) = size;
347  if (glob_size == 0)
348  {
349  glob_size = GlobalSize();
350 
351  if (glob_size == 0 && glob_size != size)
352  {
353  long local_size = size;
354  MPI_Allreduce(&local_size, &glob_size, 1, MPI_LONG,
355  MPI_SUM, GetComm());
356  }
357  }
358  NV_GLOBLENGTH_P(x) = glob_size;
359  break;
360  }
361 #endif
362  default:
363  MFEM_ABORT("N_Vector type " << id << " is not supported");
364  }
365 
366 #ifdef MFEM_USE_MPI
367  if (MPIPlusX())
368  {
369  if (glob_size == 0)
370  {
371  glob_size = GlobalSize();
372 
373  if (glob_size == 0 && glob_size != size)
374  {
375  long local_size = size;
376  MPI_Allreduce(&local_size, &glob_size, 1, MPI_LONG,
377  MPI_SUM, GetComm());
378  }
379  }
380  static_cast<N_VectorContent_MPIManyVector>(GET_CONTENT(x))->global_length =
381  glob_size;
382  }
383 #endif
384 }
385 
386 void SundialsNVector::_SetDataAndSize_()
387 {
388 #ifdef MFEM_USE_MPI
389  N_Vector local_x = MPIPlusX() ? N_VGetLocalVector_MPIPlusX(x) : x;
390 #else
391  N_Vector local_x = x;
392 #endif
393  N_Vector_ID id = N_VGetVectorID(local_x);
394 
395  // The SUNDIALS NVector owns the data if it created it.
396  switch (id)
397  {
398  case SUNDIALS_NVEC_SERIAL:
399  {
400  const bool known = mm.IsKnown(NV_DATA_S(local_x));
401  size = NV_LENGTH_S(local_x);
402  data.Wrap(NV_DATA_S(local_x), size, false);
403  if (known) { data.ClearOwnerFlags(); }
404  break;
405  }
406 #ifdef MFEM_USE_CUDA
407  case SUNDIALS_NVEC_CUDA:
408  {
409  double *h_ptr = N_VGetHostArrayPointer_Cuda(local_x);
410  double *d_ptr = N_VGetDeviceArrayPointer_Cuda(local_x);
411  const bool known = mm.IsKnown(h_ptr);
412  size = N_VGetLength_Cuda(local_x);
413  data.Wrap(h_ptr, d_ptr, size, Device::GetHostMemoryType(), false);
414  if (known) { data.ClearOwnerFlags(); }
415  UseDevice(true);
416  break;
417  }
418 #endif
419 #ifdef MFEM_USE_MPI
420  case SUNDIALS_NVEC_PARALLEL:
421  {
422  const bool known = mm.IsKnown(NV_DATA_P(x));
423  size = NV_LENGTH_S(x);
424  data.Wrap(NV_DATA_P(x), NV_LOCLENGTH_P(x), false);
425  if (known) { data.ClearOwnerFlags(); }
426  break;
427  }
428 #endif
429  default:
430  MFEM_ABORT("N_Vector type " << id << " is not supported");
431  }
432 }
433 
434 SundialsNVector::SundialsNVector()
435  : Vector()
436 {
437  // MFEM creates and owns the data,
438  // and provides it to the SUNDIALS NVector.
440  x = MakeNVector(UseDevice());
441  own_NVector = 1;
442 }
443 
444 SundialsNVector::SundialsNVector(double *data_, int size_)
445  : Vector(data_, size_)
446 {
448  x = MakeNVector(UseDevice());
449  own_NVector = 1;
451 }
452 
454  : x(nv)
455 {
457  own_NVector = 0;
458 }
459 
460 #ifdef MFEM_USE_MPI
462  : Vector()
463 {
465  x = MakeNVector(comm, UseDevice());
466  own_NVector = 1;
467 }
468 
469 SundialsNVector::SundialsNVector(MPI_Comm comm, int loc_size, long glob_size)
470  : Vector(loc_size)
471 {
473  x = MakeNVector(comm, UseDevice());
474  own_NVector = 1;
475  _SetNvecDataAndSize_(glob_size);
476 }
477 
478 SundialsNVector::SundialsNVector(MPI_Comm comm, double *data_, int loc_size,
479  long glob_size)
480  : Vector(data_, loc_size)
481 {
483  x = MakeNVector(comm, UseDevice());
484  own_NVector = 1;
485  _SetNvecDataAndSize_(glob_size);
486 }
487 
489  : SundialsNVector(vec.GetComm(), vec.GetData(), vec.Size(), vec.GlobalSize())
490 {}
491 #endif
492 
494 {
495  if (own_NVector)
496  {
497 #ifdef MFEM_USE_MPI
498  if (MPIPlusX())
499  {
500  N_VDestroy(N_VGetLocalVector_MPIPlusX(x));
501  }
502 #endif
503  N_VDestroy(x);
504  }
505 }
506 
507 void SundialsNVector::SetSize(int s, long glob_size)
508 {
509  Vector::SetSize(s);
510  _SetNvecDataAndSize_(glob_size);
511 }
512 
514 {
515  Vector::SetData(d);
517 }
518 
519 void SundialsNVector::SetDataAndSize(double *d, int s, long glob_size)
520 {
522  _SetNvecDataAndSize_(glob_size);
523 }
524 
525 N_Vector SundialsNVector::MakeNVector(bool use_device)
526 {
527  N_Vector x;
528 #ifdef MFEM_USE_CUDA
529  if (use_device)
530  {
533  }
534  else
535  {
537  }
538 #else
540 #endif
541 
542  MFEM_VERIFY(x, "Error in SundialsNVector::MakeNVector.");
543 
544  return x;
545 }
546 
547 #ifdef MFEM_USE_MPI
548 N_Vector SundialsNVector::MakeNVector(MPI_Comm comm, bool use_device)
549 {
550  N_Vector x;
551 
552  if (comm == MPI_COMM_NULL)
553  {
554  x = MakeNVector(use_device);
555  }
556  else
557  {
558 #ifdef MFEM_USE_CUDA
559  if (use_device)
560  {
565  }
566  else
567  {
568  x = N_VNewEmpty_Parallel(comm, 0, 0, Sundials::GetContext());
569  }
570 #else
571  x = N_VNewEmpty_Parallel(comm, 0, 0, Sundials::GetContext());
572 #endif // MFEM_USE_CUDA
573  }
574 
575  MFEM_VERIFY(x, "Error in SundialsNVector::MakeNVector.");
576 
577  return x;
578 }
579 #endif // MFEM_USE_MPI
580 
581 
582 // ---------------------------------------------------------------------------
583 // SUNMatrix interface functions
584 // ---------------------------------------------------------------------------
585 
586 // Return the matrix ID
587 static SUNMatrix_ID MatGetID(SUNMatrix)
588 {
589  return (SUNMATRIX_CUSTOM);
590 }
591 
592 static void MatDestroy(SUNMatrix A)
593 {
594  if (A->content) { A->content = NULL; }
595  if (A->ops) { free(A->ops); A->ops = NULL; }
596  free(A); A = NULL;
597  return;
598 }
599 
600 // ---------------------------------------------------------------------------
601 // SUNLinearSolver interface functions
602 // ---------------------------------------------------------------------------
603 
604 // Return the linear solver type
605 static SUNLinearSolver_Type LSGetType(SUNLinearSolver)
606 {
607  return (SUNLINEARSOLVER_MATRIX_ITERATIVE);
608 }
609 
610 static int LSFree(SUNLinearSolver LS)
611 {
612  if (LS->content) { LS->content = NULL; }
613  if (LS->ops) { free(LS->ops); LS->ops = NULL; }
614  free(LS); LS = NULL;
615  return (0);
616 }
617 
618 // ---------------------------------------------------------------------------
619 // CVODE interface
620 // ---------------------------------------------------------------------------
621 int CVODESolver::RHS(realtype t, const N_Vector y, N_Vector ydot,
622  void *user_data)
623 {
624  // At this point the up-to-date data for N_Vector y and ydot is on the device.
625  const SundialsNVector mfem_y(y);
626  SundialsNVector mfem_ydot(ydot);
627 
628  CVODESolver *self = static_cast<CVODESolver*>(user_data);
629 
630  // Compute y' = f(t, y)
631  self->f->SetTime(t);
632  self->f->Mult(mfem_y, mfem_ydot);
633 
634  // Return success
635  return (0);
636 }
637 
638 int CVODESolver::root(realtype t, N_Vector y, realtype *gout, void *user_data)
639 {
640  CVODESolver *self = static_cast<CVODESolver*>(user_data);
641 
642  if (!self->root_func) { return CV_RTFUNC_FAIL; }
643 
644  SundialsNVector mfem_y(y);
645  SundialsNVector mfem_gout(gout, self->root_components);
646 
647  return self->root_func(t, mfem_y, mfem_gout, self);
648 }
649 
650 void CVODESolver::SetRootFinder(int components, RootFunction func)
651 {
652  root_func = func;
653 
654  flag = CVodeRootInit(sundials_mem, components, root);
655  MFEM_VERIFY(flag == CV_SUCCESS, "error in SetRootFinder()");
656 }
657 
658 int CVODESolver::LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix A,
659  booleantype jok, booleantype *jcur, realtype gamma,
660  void*, N_Vector, N_Vector, N_Vector)
661 {
662  // Get data from N_Vectors
663  const SundialsNVector mfem_y(y);
664  const SundialsNVector mfem_fy(fy);
665  CVODESolver *self = static_cast<CVODESolver*>(GET_CONTENT(A));
666 
667  // Compute the linear system
668  self->f->SetTime(t);
669  return (self->f->SUNImplicitSetup(mfem_y, mfem_fy, jok, jcur, gamma));
670 }
671 
672 int CVODESolver::LinSysSolve(SUNLinearSolver LS, SUNMatrix, N_Vector x,
673  N_Vector b, realtype tol)
674 {
675  SundialsNVector mfem_x(x);
676  const SundialsNVector mfem_b(b);
677  CVODESolver *self = static_cast<CVODESolver*>(GET_CONTENT(LS));
678  // Solve the linear system
679  return (self->f->SUNImplicitSolve(mfem_b, mfem_x, tol));
680 }
681 
683  : lmm_type(lmm), step_mode(CV_NORMAL)
684 {
685  Y = new SundialsNVector();
686 }
687 
688 #ifdef MFEM_USE_MPI
689 CVODESolver::CVODESolver(MPI_Comm comm, int lmm)
690  : lmm_type(lmm), step_mode(CV_NORMAL)
691 {
692  Y = new SundialsNVector(comm);
693 }
694 #endif
695 
697 {
698  // Initialize the base class
699  ODESolver::Init(f_);
700 
701  // Get the vector length
702  long local_size = f_.Height();
703 
704 #ifdef MFEM_USE_MPI
705  long global_size = 0;
706  if (Parallel())
707  {
708  MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
709  Y->GetComm());
710  }
711 #endif
712 
713  // Get current time
714  double t = f_.GetTime();
715 
716  if (sundials_mem)
717  {
718  // Check if the problem size has changed since the last Init() call
719  int resize = 0;
720  if (!Parallel())
721  {
722  resize = (Y->Size() != local_size);
723  }
724  else
725  {
726 #ifdef MFEM_USE_MPI
727  int l_resize = (Y->Size() != local_size) ||
728  (saved_global_size != global_size);
729  MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
730  Y->GetComm());
731 #endif
732  }
733 
734  // Free existing solver memory and re-create with new vector size
735  if (resize)
736  {
737  CVodeFree(&sundials_mem);
738  sundials_mem = NULL;
739  }
740  }
741 
742  if (!sundials_mem)
743  {
744  // Temporarily set N_Vector wrapper data to create CVODE. The correct
745  // initial condition will be set using CVodeReInit() when Step() is
746  // called.
747 
748  if (!Parallel())
749  {
750  Y->SetSize(local_size);
751  }
752 #ifdef MFEM_USE_MPI
753  else
754  {
755  Y->SetSize(local_size, global_size);
756  saved_global_size = global_size;
757  }
758 #endif
759 
760  // Create CVODE
762  MFEM_VERIFY(sundials_mem, "error in CVodeCreate()");
763 
764  // Initialize CVODE
765  flag = CVodeInit(sundials_mem, CVODESolver::RHS, t, *Y);
766  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeInit()");
767 
768  // Attach the CVODESolver as user-defined data
769  flag = CVodeSetUserData(sundials_mem, this);
770  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetUserData()");
771 
772  // Set default tolerances
773  flag = CVodeSStolerances(sundials_mem, default_rel_tol, default_abs_tol);
774  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetSStolerances()");
775 
776  // Attach MFEM linear solver by default
778  }
779 
780  // Set the reinit flag to call CVodeReInit() in the next Step() call.
781  reinit = true;
782 }
783 
784 void CVODESolver::Step(Vector &x, double &t, double &dt)
785 {
786  Y->MakeRef(x, 0, x.Size());
787  MFEM_VERIFY(Y->Size() == x.Size(), "size mismatch");
788 
789  // Reinitialize CVODE memory if needed
790  if (reinit)
791  {
792  flag = CVodeReInit(sundials_mem, t, *Y);
793  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeReInit()");
794  // reset flag
795  reinit = false;
796  }
797 
798  // Integrate the system
799  double tout = t + dt;
800  flag = CVode(sundials_mem, tout, *Y, &t, step_mode);
801  MFEM_VERIFY(flag >= 0, "error in CVode()");
802 
803  // Make sure host is up to date
804  Y->HostRead();
805 
806  // Return the last incremental step size
807  flag = CVodeGetLastStep(sundials_mem, &dt);
808  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetLastStep()");
809 }
810 
812 {
813  // Free any existing matrix and linear solver
814  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
815  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
816 
817  // Wrap linear solver as SUNLinearSolver and SUNMatrix
819  MFEM_VERIFY(LSA, "error in SUNLinSolNewEmpty()");
820 
821  LSA->content = this;
822  LSA->ops->gettype = LSGetType;
823  LSA->ops->solve = CVODESolver::LinSysSolve;
824  LSA->ops->free = LSFree;
825 
827  MFEM_VERIFY(A, "error in SUNMatNewEmpty()");
828 
829  A->content = this;
830  A->ops->getid = MatGetID;
831  A->ops->destroy = MatDestroy;
832 
833  // Attach the linear solver and matrix
834  flag = CVodeSetLinearSolver(sundials_mem, LSA, A);
835  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinearSolver()");
836 
837  // Set the linear system evaluation function
838  flag = CVodeSetLinSysFn(sundials_mem, CVODESolver::LinSysSetup);
839  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinSysFn()");
840 }
841 
843 {
844  // Free any existing matrix and linear solver
845  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
846  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
847 
848  // Create linear solver
849  LSA = SUNLinSol_SPGMR(*Y, PREC_NONE, 0, Sundials::GetContext());
850  MFEM_VERIFY(LSA, "error in SUNLinSol_SPGMR()");
851 
852  // Attach linear solver
853  flag = CVodeSetLinearSolver(sundials_mem, LSA, NULL);
854  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinearSolver()");
855 }
856 
858 {
859  step_mode = itask;
860 }
861 
862 void CVODESolver::SetSStolerances(double reltol, double abstol)
863 {
864  flag = CVodeSStolerances(sundials_mem, reltol, abstol);
865  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSStolerances()");
866 }
867 
868 void CVODESolver::SetSVtolerances(double reltol, Vector abstol)
869 {
870  MFEM_VERIFY(abstol.Size() == f->Height(),
871  "abs tolerance is not the same size.");
872 
873  SundialsNVector mfem_abstol;
874  mfem_abstol.MakeRef(abstol, 0, abstol.Size());
875 
876  flag = CVodeSVtolerances(sundials_mem, reltol, mfem_abstol);
877  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSVtolerances()");
878 }
879 
880 void CVODESolver::SetMaxStep(double dt_max)
881 {
882  flag = CVodeSetMaxStep(sundials_mem, dt_max);
883  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetMaxStep()");
884 }
885 
886 void CVODESolver::SetMaxNSteps(int mxsteps)
887 {
888  flag = CVodeSetMaxNumSteps(sundials_mem, mxsteps);
889  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetMaxNumSteps()");
890 }
891 
893 {
894  long nsteps;
895  flag = CVodeGetNumSteps(sundials_mem, &nsteps);
896  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetNumSteps()");
897  return nsteps;
898 }
899 
900 void CVODESolver::SetMaxOrder(int max_order)
901 {
902  flag = CVodeSetMaxOrd(sundials_mem, max_order);
903  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetMaxOrd()");
904 }
905 
907 {
908  long int nsteps, nfevals, nlinsetups, netfails;
909  int qlast, qcur;
910  double hinused, hlast, hcur, tcur;
911  long int nniters, nncfails;
912 
913  // Get integrator stats
914  flag = CVodeGetIntegratorStats(sundials_mem,
915  &nsteps,
916  &nfevals,
917  &nlinsetups,
918  &netfails,
919  &qlast,
920  &qcur,
921  &hinused,
922  &hlast,
923  &hcur,
924  &tcur);
925  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetIntegratorStats()");
926 
927  // Get nonlinear solver stats
928  flag = CVodeGetNonlinSolvStats(sundials_mem,
929  &nniters,
930  &nncfails);
931  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetNonlinSolvStats()");
932 
933  mfem::out <<
934  "CVODE:\n"
935  "num steps: " << nsteps << "\n"
936  "num rhs evals: " << nfevals << "\n"
937  "num lin setups: " << nlinsetups << "\n"
938  "num nonlin sol iters: " << nniters << "\n"
939  "num nonlin conv fail: " << nncfails << "\n"
940  "num error test fails: " << netfails << "\n"
941  "last order: " << qlast << "\n"
942  "current order: " << qcur << "\n"
943  "initial dt: " << hinused << "\n"
944  "last dt: " << hlast << "\n"
945  "current dt: " << hcur << "\n"
946  "current t: " << tcur << "\n" << endl;
947 
948  return;
949 }
950 
952 {
953  delete Y;
954  SUNMatDestroy(A);
955  SUNLinSolFree(LSA);
956  SUNNonlinSolFree(NLS);
957  CVodeFree(&sundials_mem);
958 }
959 
960 // ---------------------------------------------------------------------------
961 // CVODESSolver interface
962 // ---------------------------------------------------------------------------
963 
965  CVODESolver(lmm),
966  ncheck(0),
967  indexB(0),
968  AB(nullptr),
969  LSB(nullptr)
970 {
971  q = new SundialsNVector();
972  qB = new SundialsNVector();
973  yB = new SundialsNVector();
974  yy = new SundialsNVector();
975 }
976 
977 #ifdef MFEM_USE_MPI
978 CVODESSolver::CVODESSolver(MPI_Comm comm, int lmm) :
979  CVODESolver(comm, lmm),
980  ncheck(0),
981  indexB(0),
982  AB(nullptr),
983  LSB(nullptr)
984 {
985  q = new SundialsNVector(comm);
986  qB = new SundialsNVector(comm);
987  yB = new SundialsNVector(comm);
988  yy = new SundialsNVector(comm);
989 }
990 #endif
991 
993 {
994  MFEM_VERIFY(t <= f->GetTime(), "t > current forward solver time");
995 
996  flag = CVodeGetQuad(sundials_mem, &t, *q);
997  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetQuad()");
998 
999  Q.Set(1., *q);
1000 }
1001 
1003 {
1004  MFEM_VERIFY(t <= f->GetTime(), "t > current forward solver time");
1005 
1006  flag = CVodeGetQuadB(sundials_mem, indexB, &t, *qB);
1007  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetQuadB()");
1008 
1009  dG_dp.Set(-1., *qB);
1010 }
1011 
1013 {
1014  yy->MakeRef(yyy, 0, yyy.Size());
1015 
1016  flag = CVodeGetAdjY(sundials_mem, tB, *yy);
1017  MFEM_VERIFY(flag >= 0, "error in CVodeGetAdjY()");
1018 }
1019 
1020 // Implemented to enforce type checking for TimeDependentAdjointOperator
1022 {
1023  CVODESolver::Init(f_);
1024 }
1025 
1027 {
1028  long local_size = f_.GetAdjointHeight();
1029 
1030  // Get current time
1031  double tB = f_.GetTime();
1032 
1033  yB->SetSize(local_size);
1034 
1035  // Create the solver memory
1036  flag = CVodeCreateB(sundials_mem, CV_BDF, &indexB);
1037  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeCreateB()");
1038 
1039  // Initialize
1040  flag = CVodeInitB(sundials_mem, indexB, RHSB, tB, *yB);
1041  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeInit()");
1042 
1043  // Attach the CVODESSolver as user-defined data
1044  flag = CVodeSetUserDataB(sundials_mem, indexB, this);
1045  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetUserDataB()");
1046 
1047  // Set default tolerances
1048  flag = CVodeSStolerancesB(sundials_mem, indexB, default_rel_tolB,
1050  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetSStolerancesB()");
1051 
1052  // Attach MFEM linear solver by default
1054 
1055  // Set the reinit flag to call CVodeReInit() in the next Step() call.
1056  reinit = true;
1057 }
1058 
1059 void CVODESSolver::InitAdjointSolve(int steps, int interpolation)
1060 {
1061  flag = CVodeAdjInit(sundials_mem, steps, interpolation);
1062  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeAdjInit");
1063 }
1064 
1066 {
1067  flag = CVodeSetMaxNumStepsB(sundials_mem, indexB, mxstepsB);
1068  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeSetMaxNumStepsB()");
1069 }
1070 
1072  double abstolQ)
1073 {
1074  q->MakeRef(q0, 0, q0.Size());
1075 
1076  flag = CVodeQuadInit(sundials_mem, RHSQ, *q);
1077  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeQuadInit()");
1078 
1079  flag = CVodeSetQuadErrCon(sundials_mem, SUNTRUE);
1080  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeSetQuadErrCon");
1081 
1082  flag = CVodeQuadSStolerances(sundials_mem, reltolQ, abstolQ);
1083  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeQuadSStolerances");
1084 }
1085 
1087  double abstolQB)
1088 {
1089  qB->MakeRef(qB0, 0, qB0.Size());
1090 
1091  flag = CVodeQuadInitB(sundials_mem, indexB, RHSQB, *qB);
1092  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeQuadInitB()");
1093 
1094  flag = CVodeSetQuadErrConB(sundials_mem, indexB, SUNTRUE);
1095  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeSetQuadErrConB");
1096 
1097  flag = CVodeQuadSStolerancesB(sundials_mem, indexB, reltolQB, abstolQB);
1098  MFEM_VERIFY(flag == CV_SUCCESS, "Error in CVodeQuadSStolerancesB");
1099 }
1100 
1102 {
1103  // Free any existing linear solver
1104  if (AB != NULL) { SUNMatDestroy(AB); AB = NULL; }
1105  if (LSB != NULL) { SUNLinSolFree(LSB); LSB = NULL; }
1106 
1107  // Wrap linear solver as SUNLinearSolver and SUNMatrix
1109  MFEM_VERIFY(LSB, "error in SUNLinSolNewEmpty()");
1110 
1111  LSB->content = this;
1112  LSB->ops->gettype = LSGetType;
1113  LSB->ops->solve = CVODESSolver::LinSysSolveB; // JW change
1114  LSB->ops->free = LSFree;
1115 
1117  MFEM_VERIFY(AB, "error in SUNMatNewEmpty()");
1118 
1119  AB->content = this;
1120  AB->ops->getid = MatGetID;
1121  AB->ops->destroy = MatDestroy;
1122 
1123  // Attach the linear solver and matrix
1124  flag = CVodeSetLinearSolverB(sundials_mem, indexB, LSB, AB);
1125  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinearSolverB()");
1126 
1127  // Set the linear system evaluation function
1128  flag = CVodeSetLinSysFnB(sundials_mem, indexB,
1129  CVODESSolver::LinSysSetupB); // JW change
1130  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinSysFn()");
1131 }
1132 
1134 {
1135  // Free any existing matrix and linear solver
1136  if (AB != NULL) { SUNMatDestroy(AB); AB = NULL; }
1137  if (LSB != NULL) { SUNLinSolFree(LSB); LSB = NULL; }
1138 
1139  // Set default linear solver (Newton is the default Nonlinear Solver)
1140  LSB = SUNLinSol_SPGMR(*yB, PREC_NONE, 0, Sundials::GetContext());
1141  MFEM_VERIFY(LSB, "error in SUNLinSol_SPGMR()");
1142 
1143  /* Attach the matrix and linear solver */
1144  flag = CVodeSetLinearSolverB(sundials_mem, indexB, LSB, NULL);
1145  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinearSolverB()");
1146 }
1147 
1148 int CVODESSolver::LinSysSetupB(realtype t, N_Vector y, N_Vector yB,
1149  N_Vector fyB, SUNMatrix AB,
1150  booleantype jokB, booleantype *jcurB,
1151  realtype gammaB, void *user_data, N_Vector tmp1,
1152  N_Vector tmp2, N_Vector tmp3)
1153 {
1154  // Get data from N_Vectors
1155  const SundialsNVector mfem_y(y);
1156  const SundialsNVector mfem_yB(yB);
1157  SundialsNVector mfem_fyB(fyB);
1158  CVODESSolver *self = static_cast<CVODESSolver*>(GET_CONTENT(AB));
1160  (self->f);
1161  f->SetTime(t);
1162  // Compute the linear system
1163  return (f->SUNImplicitSetupB(t, mfem_y, mfem_yB, mfem_fyB, jokB, jcurB,
1164  gammaB));
1165 }
1166 
1167 int CVODESSolver::LinSysSolveB(SUNLinearSolver LS, SUNMatrix AB, N_Vector yB,
1168  N_Vector Rb, realtype tol)
1169 {
1170  SundialsNVector mfem_yB(yB);
1171  const SundialsNVector mfem_Rb(Rb);
1172  CVODESSolver *self = static_cast<CVODESSolver*>(GET_CONTENT(LS));
1174  (self->f);
1175  // Solve the linear system
1176  int ret = f->SUNImplicitSolveB(mfem_yB, mfem_Rb, tol);
1177  return (ret);
1178 }
1179 
1180 void CVODESSolver::SetSStolerancesB(double reltol, double abstol)
1181 {
1182  flag = CVodeSStolerancesB(sundials_mem, indexB, reltol, abstol);
1183  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSStolerancesB()");
1184 }
1185 
1186 void CVODESSolver::SetSVtolerancesB(double reltol, Vector abstol)
1187 {
1188  MFEM_VERIFY(abstol.Size() == f->Height(),
1189  "abs tolerance is not the same size.");
1190 
1191  SundialsNVector mfem_abstol;
1192  mfem_abstol.MakeRef(abstol, 0, abstol.Size());
1193 
1194  flag = CVodeSVtolerancesB(sundials_mem, indexB, reltol, mfem_abstol);
1195  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSVtolerancesB()");
1196 }
1197 
1199 {
1200  ewt_func = func;
1201  CVodeWFtolerances(sundials_mem, ewt);
1202 }
1203 
1204 // CVODESSolver static functions
1205 
1206 int CVODESSolver::RHSQ(realtype t, const N_Vector y, N_Vector qdot,
1207  void *user_data)
1208 {
1209  CVODESSolver *self = static_cast<CVODESSolver*>(user_data);
1210  const SundialsNVector mfem_y(y);
1211  SundialsNVector mfem_qdot(qdot);
1213  (self->f);
1214  f->SetTime(t);
1215  f->QuadratureIntegration(mfem_y, mfem_qdot);
1216  return 0;
1217 }
1218 
1219 int CVODESSolver::RHSQB(realtype t, N_Vector y, N_Vector yB, N_Vector qBdot,
1220  void *user_dataB)
1221 {
1222  CVODESSolver *self = static_cast<CVODESSolver*>(user_dataB);
1223  SundialsNVector mfem_y(y);
1224  SundialsNVector mfem_yB(yB);
1225  SundialsNVector mfem_qBdot(qBdot);
1227  (self->f);
1228  f->SetTime(t);
1229  f->QuadratureSensitivityMult(mfem_y, mfem_yB, mfem_qBdot);
1230  return 0;
1231 }
1232 
1233 int CVODESSolver::RHSB(realtype t, N_Vector y, N_Vector yB, N_Vector yBdot,
1234  void *user_dataB)
1235 {
1236  CVODESSolver *self = static_cast<CVODESSolver*>(user_dataB);
1237  SundialsNVector mfem_y(y);
1238  SundialsNVector mfem_yB(yB);
1239  SundialsNVector mfem_yBdot(yBdot);
1240 
1241  mfem_yBdot = 0.;
1243  (self->f);
1244  f->SetTime(t);
1245  f->AdjointRateMult(mfem_y, mfem_yB, mfem_yBdot);
1246  return 0;
1247 }
1248 
1249 int CVODESSolver::ewt(N_Vector y, N_Vector w, void *user_data)
1250 {
1251  CVODESSolver *self = static_cast<CVODESSolver*>(user_data);
1252 
1253  SundialsNVector mfem_y(y);
1254  SundialsNVector mfem_w(w);
1255 
1256  return self->ewt_func(mfem_y, mfem_w, self);
1257 }
1258 
1259 // Pretty much a copy of CVODESolver::Step except we use CVodeF instead of CVode
1260 void CVODESSolver::Step(Vector &x, double &t, double &dt)
1261 {
1262  Y->MakeRef(x, 0, x.Size());
1263  MFEM_VERIFY(Y->Size() == x.Size(), "size mismatch");
1264 
1265  // Reinitialize CVODE memory if needed, initializes the N_Vector y with x
1266  if (reinit)
1267  {
1268  flag = CVodeReInit(sundials_mem, t, *Y);
1269  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeReInit()");
1270 
1271  // reset flag
1272  reinit = false;
1273  }
1274 
1275  // Integrate the system
1276  double tout = t + dt;
1277  flag = CVodeF(sundials_mem, tout, *Y, &t, step_mode, &ncheck);
1278  MFEM_VERIFY(flag >= 0, "error in CVodeF()");
1279 
1280  // Make sure host is up to date
1281  Y->HostRead();
1282 
1283  // Return the last incremental step size
1284  flag = CVodeGetLastStep(sundials_mem, &dt);
1285  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetLastStep()");
1286 }
1287 
1288 void CVODESSolver::StepB(Vector &xB, double &tB, double &dtB)
1289 {
1290  yB->MakeRef(xB, 0, xB.Size());
1291  MFEM_VERIFY(yB->Size() == xB.Size(), "");
1292 
1293  // Reinitialize CVODE memory if needed
1294  if (reinit)
1295  {
1296  flag = CVodeReInitB(sundials_mem, indexB, tB, *yB);
1297  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeReInit()");
1298 
1299  // reset flag
1300  reinit = false;
1301  }
1302 
1303  // Integrate the system
1304  double tout = tB - dtB;
1305  flag = CVodeB(sundials_mem, tout, step_mode);
1306  MFEM_VERIFY(flag >= 0, "error in CVodeB()");
1307 
1308  // Call CVodeGetB to get yB of the backward ODE problem.
1309  flag = CVodeGetB(sundials_mem, indexB, &tB, *yB);
1310  MFEM_VERIFY(flag >= 0, "error in CVodeGetB()");
1311 
1312  // Make sure host is up to date
1313  yB->HostRead();
1314 }
1315 
1317 {
1318  delete yB;
1319  delete yy;
1320  delete qB;
1321  delete q;
1322  SUNMatDestroy(AB);
1323  SUNLinSolFree(LSB);
1324 }
1325 
1326 
1327 // ---------------------------------------------------------------------------
1328 // ARKStep interface
1329 // ---------------------------------------------------------------------------
1330 
1331 int ARKStepSolver::RHS1(realtype t, const N_Vector y, N_Vector ydot,
1332  void *user_data)
1333 {
1334  // Get data from N_Vectors
1335  const SundialsNVector mfem_y(y);
1336  SundialsNVector mfem_ydot(ydot);
1337  ARKStepSolver *self = static_cast<ARKStepSolver*>(user_data);
1338 
1339  // Compute f(t, y) in y' = f(t, y) or fe(t, y) in y' = fe(t, y) + fi(t, y)
1340  self->f->SetTime(t);
1341  if (self->rk_type == IMEX)
1342  {
1343  self->f->SetEvalMode(TimeDependentOperator::ADDITIVE_TERM_1);
1344  }
1345  self->f->Mult(mfem_y, mfem_ydot);
1346 
1347  // Return success
1348  return (0);
1349 }
1350 
1351 int ARKStepSolver::RHS2(realtype t, const N_Vector y, N_Vector ydot,
1352  void *user_data)
1353 {
1354  // Get data from N_Vectors
1355  const SundialsNVector mfem_y(y);
1356  SundialsNVector mfem_ydot(ydot);
1357  ARKStepSolver *self = static_cast<ARKStepSolver*>(user_data);
1358 
1359  // Compute fi(t, y) in y' = fe(t, y) + fi(t, y)
1360  self->f->SetTime(t);
1361  self->f->SetEvalMode(TimeDependentOperator::ADDITIVE_TERM_2);
1362  self->f->Mult(mfem_y, mfem_ydot);
1363 
1364  // Return success
1365  return (0);
1366 }
1367 
1368 int ARKStepSolver::LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix A,
1369  SUNMatrix, booleantype jok, booleantype *jcur,
1370  realtype gamma,
1371  void*, N_Vector, N_Vector, N_Vector)
1372 {
1373  // Get data from N_Vectors
1374  const SundialsNVector mfem_y(y);
1375  const SundialsNVector mfem_fy(fy);
1376  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(A));
1377 
1378  // Compute the linear system
1379  self->f->SetTime(t);
1380  if (self->rk_type == IMEX)
1381  {
1382  self->f->SetEvalMode(TimeDependentOperator::ADDITIVE_TERM_2);
1383  }
1384  return (self->f->SUNImplicitSetup(mfem_y, mfem_fy, jok, jcur, gamma));
1385 }
1386 
1387 int ARKStepSolver::LinSysSolve(SUNLinearSolver LS, SUNMatrix, N_Vector x,
1388  N_Vector b, realtype tol)
1389 {
1390  SundialsNVector mfem_x(x);
1391  const SundialsNVector mfem_b(b);
1392  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(LS));
1393 
1394  // Solve the linear system
1395  if (self->rk_type == IMEX)
1396  {
1398  }
1399  return (self->f->SUNImplicitSolve(mfem_b, mfem_x, tol));
1400 }
1401 
1402 int ARKStepSolver::MassSysSetup(realtype t, SUNMatrix M,
1403  void*, N_Vector, N_Vector, N_Vector)
1404 {
1405  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(M));
1406 
1407  // Compute the mass matrix system
1408  self->f->SetTime(t);
1409  return (self->f->SUNMassSetup());
1410 }
1411 
1412 int ARKStepSolver::MassSysSolve(SUNLinearSolver LS, SUNMatrix, N_Vector x,
1413  N_Vector b, realtype tol)
1414 {
1415  SundialsNVector mfem_x(x);
1416  const SundialsNVector mfem_b(b);
1417  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(LS));
1418 
1419  // Solve the mass matrix system
1420  return (self->f->SUNMassSolve(mfem_b, mfem_x, tol));
1421 }
1422 
1423 int ARKStepSolver::MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
1424 {
1425  const SundialsNVector mfem_x(x);
1426  SundialsNVector mfem_v(v);
1427  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(M));
1428 
1429  // Compute the mass matrix-vector product
1430  return (self->f->SUNMassMult(mfem_x, mfem_v));
1431 }
1432 
1433 int ARKStepSolver::MassMult2(N_Vector x, N_Vector v, realtype t,
1434  void* mtimes_data)
1435 {
1436  const SundialsNVector mfem_x(x);
1437  SundialsNVector mfem_v(v);
1438  ARKStepSolver *self = static_cast<ARKStepSolver*>(mtimes_data);
1439 
1440  // Compute the mass matrix-vector product
1441  self->f->SetTime(t);
1442  return (self->f->SUNMassMult(mfem_x, mfem_v));
1443 }
1444 
1446  : rk_type(type), step_mode(ARK_NORMAL),
1447  use_implicit(type == IMPLICIT || type == IMEX)
1448 {
1449  Y = new SundialsNVector();
1450 }
1451 
1452 #ifdef MFEM_USE_MPI
1454  : rk_type(type), step_mode(ARK_NORMAL),
1455  use_implicit(type == IMPLICIT || type == IMEX)
1456 {
1457  Y = new SundialsNVector(comm);
1458 }
1459 #endif
1460 
1462 {
1463  // Initialize the base class
1464  ODESolver::Init(f_);
1465 
1466  // Get the vector length
1467  long local_size = f_.Height();
1468 #ifdef MFEM_USE_MPI
1469  long global_size;
1470 #endif
1471 
1472  if (Parallel())
1473  {
1474 #ifdef MFEM_USE_MPI
1475  MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
1476  Y->GetComm());
1477 #endif
1478  }
1479 
1480  // Get current time
1481  double t = f_.GetTime();
1482 
1483  if (sundials_mem)
1484  {
1485  // Check if the problem size has changed since the last Init() call
1486  int resize = 0;
1487  if (!Parallel())
1488  {
1489  resize = (Y->Size() != local_size);
1490  }
1491  else
1492  {
1493 #ifdef MFEM_USE_MPI
1494  int l_resize = (Y->Size() != local_size) ||
1495  (saved_global_size != global_size);
1496  MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
1497  Y->GetComm());
1498 #endif
1499  }
1500 
1501  // Free existing solver memory and re-create with new vector size
1502  if (resize)
1503  {
1504  ARKStepFree(&sundials_mem);
1505  sundials_mem = NULL;
1506  }
1507  }
1508 
1509  if (!sundials_mem)
1510  {
1511  if (!Parallel())
1512  {
1513  Y->SetSize(local_size);
1514  }
1515 #ifdef MFEM_USE_MPI
1516  else
1517  {
1518  Y->SetSize(local_size, global_size);
1519  saved_global_size = global_size;
1520  }
1521 #endif
1522 
1523  // Create ARKStep memory
1524  if (rk_type == IMPLICIT)
1525  {
1528  }
1529  else if (rk_type == EXPLICIT)
1530  {
1533  }
1534  else
1535  {
1537  t, *Y, Sundials::GetContext());
1538  }
1539  MFEM_VERIFY(sundials_mem, "error in ARKStepCreate()");
1540 
1541  // Attach the ARKStepSolver as user-defined data
1542  flag = ARKStepSetUserData(sundials_mem, this);
1543  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetUserData()");
1544 
1545  // Set default tolerances
1546  flag = ARKStepSStolerances(sundials_mem, default_rel_tol, default_abs_tol);
1547  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetSStolerances()");
1548 
1549  // If implicit, attach MFEM linear solver by default
1550  if (use_implicit) { UseMFEMLinearSolver(); }
1551  }
1552 
1553  // Set the reinit flag to call ARKStepReInit() in the next Step() call.
1554  reinit = true;
1555 }
1556 
1557 void ARKStepSolver::Step(Vector &x, double &t, double &dt)
1558 {
1559  Y->MakeRef(x, 0, x.Size());
1560  MFEM_VERIFY(Y->Size() == x.Size(), "size mismatch");
1561 
1562  // Reinitialize ARKStep memory if needed
1563  if (reinit)
1564  {
1565  if (rk_type == IMPLICIT)
1566  {
1567  flag = ARKStepReInit(sundials_mem, NULL, ARKStepSolver::RHS1, t, *Y);
1568  }
1569  else if (rk_type == EXPLICIT)
1570  {
1571  flag = ARKStepReInit(sundials_mem, ARKStepSolver::RHS1, NULL, t, *Y);
1572  }
1573  else
1574  {
1575  flag = ARKStepReInit(sundials_mem,
1577  }
1578  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepReInit()");
1579 
1580  // reset flag
1581  reinit = false;
1582  }
1583 
1584  // Integrate the system
1585  double tout = t + dt;
1586  flag = ARKStepEvolve(sundials_mem, tout, *Y, &t, step_mode);
1587  MFEM_VERIFY(flag >= 0, "error in ARKStepEvolve()");
1588 
1589  // Make sure host is up to date
1590  Y->HostRead();
1591 
1592  // Return the last incremental step size
1593  flag = ARKStepGetLastStep(sundials_mem, &dt);
1594  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepGetLastStep()");
1595 }
1596 
1598 {
1599  // Free any existing matrix and linear solver
1600  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1601  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
1602 
1603  // Wrap linear solver as SUNLinearSolver and SUNMatrix
1605  MFEM_VERIFY(LSA, "error in SUNLinSolNewEmpty()");
1606 
1607  LSA->content = this;
1608  LSA->ops->gettype = LSGetType;
1609  LSA->ops->solve = ARKStepSolver::LinSysSolve;
1610  LSA->ops->free = LSFree;
1611 
1613  MFEM_VERIFY(A, "error in SUNMatNewEmpty()");
1614 
1615  A->content = this;
1616  A->ops->getid = MatGetID;
1617  A->ops->destroy = MatDestroy;
1618 
1619  // Attach the linear solver and matrix
1620  flag = ARKStepSetLinearSolver(sundials_mem, LSA, A);
1621  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinearSolver()");
1622 
1623  // Set the linear system evaluation function
1624  flag = ARKStepSetLinSysFn(sundials_mem, ARKStepSolver::LinSysSetup);
1625  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinSysFn()");
1626 }
1627 
1629 {
1630  // Free any existing matrix and linear solver
1631  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1632  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
1633 
1634  // Create linear solver
1635  LSA = SUNLinSol_SPGMR(*Y, PREC_NONE, 0, Sundials::GetContext());
1636  MFEM_VERIFY(LSA, "error in SUNLinSol_SPGMR()");
1637 
1638  // Attach linear solver
1639  flag = ARKStepSetLinearSolver(sundials_mem, LSA, NULL);
1640  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinearSolver()");
1641 }
1642 
1644 {
1645  // Free any existing matrix and linear solver
1646  if (M != NULL) { SUNMatDestroy(M); M = NULL; }
1647  if (LSM != NULL) { SUNLinSolFree(LSM); LSM = NULL; }
1648 
1649  // Wrap linear solver as SUNLinearSolver and SUNMatrix
1651  MFEM_VERIFY(LSM, "error in SUNLinSolNewEmpty()");
1652 
1653  LSM->content = this;
1654  LSM->ops->gettype = LSGetType;
1655  LSM->ops->solve = ARKStepSolver::MassSysSolve;
1656  LSA->ops->free = LSFree;
1657 
1659  MFEM_VERIFY(M, "error in SUNMatNewEmpty()");
1660 
1661  M->content = this;
1662  M->ops->getid = SUNMatGetID;
1663  M->ops->matvec = ARKStepSolver::MassMult1;
1664  M->ops->destroy = MatDestroy;
1665 
1666  // Attach the linear solver and matrix
1667  flag = ARKStepSetMassLinearSolver(sundials_mem, LSM, M, tdep);
1668  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinearSolver()");
1669 
1670  // Set the linear system function
1671  flag = ARKStepSetMassFn(sundials_mem, ARKStepSolver::MassSysSetup);
1672  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMassFn()");
1673 }
1674 
1676 {
1677  // Free any existing matrix and linear solver
1678  if (M != NULL) { SUNMatDestroy(A); M = NULL; }
1679  if (LSM != NULL) { SUNLinSolFree(LSM); LSM = NULL; }
1680 
1681  // Create linear solver
1682  LSM = SUNLinSol_SPGMR(*Y, PREC_NONE, 0, Sundials::GetContext());
1683  MFEM_VERIFY(LSM, "error in SUNLinSol_SPGMR()");
1684 
1685  // Attach linear solver
1686  flag = ARKStepSetMassLinearSolver(sundials_mem, LSM, NULL, tdep);
1687  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMassLinearSolver()");
1688 
1689  // Attach matrix multiplication function
1690  flag = ARKStepSetMassTimes(sundials_mem, NULL, ARKStepSolver::MassMult2,
1691  this);
1692  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMassTimes()");
1693 }
1694 
1696 {
1697  step_mode = itask;
1698 }
1699 
1700 void ARKStepSolver::SetSStolerances(double reltol, double abstol)
1701 {
1702  flag = ARKStepSStolerances(sundials_mem, reltol, abstol);
1703  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSStolerances()");
1704 }
1705 
1706 void ARKStepSolver::SetMaxStep(double dt_max)
1707 {
1708  flag = ARKStepSetMaxStep(sundials_mem, dt_max);
1709  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMaxStep()");
1710 }
1711 
1713 {
1714  flag = ARKStepSetOrder(sundials_mem, order);
1715  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetOrder()");
1716 }
1717 
1719 {
1720  flag = ARKStepSetTableNum(sundials_mem, ARKODE_DIRK_NONE, table_id);
1721  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetTableNum()");
1722 }
1723 
1725 {
1726  flag = ARKStepSetTableNum(sundials_mem, table_id, ARKODE_ERK_NONE);
1727  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetTableNum()");
1728 }
1729 
1731  ARKODE_DIRKTableID itable_id)
1732 {
1733  flag = ARKStepSetTableNum(sundials_mem, itable_id, etable_id);
1734  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetTableNum()");
1735 }
1736 
1738 {
1739  flag = ARKStepSetFixedStep(sundials_mem, dt);
1740  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetFixedStep()");
1741 }
1742 
1744 {
1745  long int nsteps, expsteps, accsteps, step_attempts;
1746  long int nfe_evals, nfi_evals;
1747  long int nlinsetups, netfails;
1748  double hinused, hlast, hcur, tcur;
1749  long int nniters, nncfails;
1750 
1751  // Get integrator stats
1752  flag = ARKStepGetTimestepperStats(sundials_mem,
1753  &expsteps,
1754  &accsteps,
1755  &step_attempts,
1756  &nfe_evals,
1757  &nfi_evals,
1758  &nlinsetups,
1759  &netfails);
1760  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepGetTimestepperStats()");
1761 
1762  flag = ARKStepGetStepStats(sundials_mem,
1763  &nsteps,
1764  &hinused,
1765  &hlast,
1766  &hcur,
1767  &tcur);
1768 
1769  // Get nonlinear solver stats
1770  flag = ARKStepGetNonlinSolvStats(sundials_mem,
1771  &nniters,
1772  &nncfails);
1773  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepGetNonlinSolvStats()");
1774 
1775  mfem::out <<
1776  "ARKStep:\n"
1777  "num steps: " << nsteps << "\n"
1778  "num exp rhs evals: " << nfe_evals << "\n"
1779  "num imp rhs evals: " << nfi_evals << "\n"
1780  "num lin setups: " << nlinsetups << "\n"
1781  "num nonlin sol iters: " << nniters << "\n"
1782  "num nonlin conv fail: " << nncfails << "\n"
1783  "num steps attempted: " << step_attempts << "\n"
1784  "num acc limited steps: " << accsteps << "\n"
1785  "num exp limited stepfails: " << expsteps << "\n"
1786  "num error test fails: " << netfails << "\n"
1787  "initial dt: " << hinused << "\n"
1788  "last dt: " << hlast << "\n"
1789  "current dt: " << hcur << "\n"
1790  "current t: " << tcur << "\n" << endl;
1791 
1792  return;
1793 }
1794 
1796 {
1797  delete Y;
1798  SUNMatDestroy(A);
1799  SUNLinSolFree(LSA);
1800  SUNNonlinSolFree(NLS);
1801  ARKStepFree(&sundials_mem);
1802 }
1803 
1804 // ---------------------------------------------------------------------------
1805 // KINSOL interface
1806 // ---------------------------------------------------------------------------
1807 
1808 // Wrapper for evaluating the nonlinear residual F(u) = 0
1809 int KINSolver::Mult(const N_Vector u, N_Vector fu, void *user_data)
1810 {
1811  const SundialsNVector mfem_u(u);
1812  SundialsNVector mfem_fu(fu);
1813  KINSolver *self = static_cast<KINSolver*>(user_data);
1814 
1815  // Compute the non-linear action F(u).
1816  self->oper->Mult(mfem_u, mfem_fu);
1817 
1818  // Return success
1819  return 0;
1820 }
1821 
1822 // Wrapper for computing Jacobian-vector products
1823 int KINSolver::GradientMult(N_Vector v, N_Vector Jv, N_Vector u,
1824  booleantype *new_u, void *user_data)
1825 {
1826  const SundialsNVector mfem_v(v);
1827  SundialsNVector mfem_Jv(Jv);
1828  KINSolver *self = static_cast<KINSolver*>(user_data);
1829 
1830  // Update Jacobian information if needed
1831  if (*new_u)
1832  {
1833  const SundialsNVector mfem_u(u);
1834  self->jacobian = &self->oper->GetGradient(mfem_u);
1835  *new_u = SUNFALSE;
1836  }
1837 
1838  // Compute the Jacobian-vector product
1839  self->jacobian->Mult(mfem_v, mfem_Jv);
1840 
1841  // Return success
1842  return 0;
1843 }
1844 
1845 // Wrapper for evaluating linear systems J u = b
1846 int KINSolver::LinSysSetup(N_Vector u, N_Vector, SUNMatrix J,
1847  void *, N_Vector, N_Vector )
1848 {
1849  const SundialsNVector mfem_u(u);
1850  KINSolver *self = static_cast<KINSolver*>(GET_CONTENT(J));
1851 
1852  // Update the Jacobian
1853  self->jacobian = &self->oper->GetGradient(mfem_u);
1854 
1855  // Set the Jacobian solve operator
1856  self->prec->SetOperator(*self->jacobian);
1857 
1858  // Return success
1859  return (0);
1860 }
1861 
1862 // Wrapper for solving linear systems J u = b
1863 int KINSolver::LinSysSolve(SUNLinearSolver LS, SUNMatrix, N_Vector u,
1864  N_Vector b, realtype)
1865 {
1866  SundialsNVector mfem_u(u), mfem_b(b);
1867  KINSolver *self = static_cast<KINSolver*>(GET_CONTENT(LS));
1868 
1869  // Solve for u = [J(u)]^{-1} b, maybe approximately.
1870  self->prec->Mult(mfem_b, mfem_u);
1871 
1872  // Return success
1873  return (0);
1874 }
1875 
1876 int KINSolver::PrecSetup(N_Vector uu,
1877  N_Vector uscale,
1878  N_Vector fval,
1879  N_Vector fscale,
1880  void *user_data)
1881 {
1882  SundialsNVector mfem_u(uu);
1883  KINSolver *self = static_cast<KINSolver *>(user_data);
1884 
1885  // Update the Jacobian
1886  self->jacobian = &self->oper->GetGradient(mfem_u);
1887 
1888  // Set the Jacobian solve operator
1889  self->prec->SetOperator(*self->jacobian);
1890 
1891  return 0;
1892 }
1893 
1894 int KINSolver::PrecSolve(N_Vector uu,
1895  N_Vector uscale,
1896  N_Vector fval,
1897  N_Vector fscale,
1898  N_Vector vv,
1899  void *user_data)
1900 {
1901  KINSolver *self = static_cast<KINSolver *>(user_data);
1902  SundialsNVector mfem_v(vv);
1903 
1904  self->wrk = 0.0;
1905 
1906  // Solve for u = P^{-1} v
1907  self->prec->Mult(mfem_v, self->wrk);
1908 
1909  mfem_v = self->wrk;
1910 
1911  return 0;
1912 }
1913 
1914 KINSolver::KINSolver(int strategy, bool oper_grad)
1915  : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
1916  f_scale(NULL), jacobian(NULL), maa(0)
1917 {
1918  Y = new SundialsNVector();
1919  y_scale = new SundialsNVector();
1920  f_scale = new SundialsNVector();
1921 
1922  // Default abs_tol and print_level
1923  abs_tol = pow(UNIT_ROUNDOFF, 1.0/3.0);
1924  print_level = 0;
1925 }
1926 
1927 #ifdef MFEM_USE_MPI
1928 KINSolver::KINSolver(MPI_Comm comm, int strategy, bool oper_grad)
1929  : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
1930  f_scale(NULL), jacobian(NULL), maa(0)
1931 {
1932  Y = new SundialsNVector(comm);
1933  y_scale = new SundialsNVector(comm);
1934  f_scale = new SundialsNVector(comm);
1935 
1936  // Default abs_tol and print_level
1937  abs_tol = pow(UNIT_ROUNDOFF, 1.0/3.0);
1938  print_level = 0;
1939 }
1940 #endif
1941 
1942 
1944 {
1945  // Initialize the base class
1947  jacobian = NULL;
1948 
1949  // Get the vector length
1950  long local_size = height;
1951 #ifdef MFEM_USE_MPI
1952  long global_size;
1953 #endif
1954 
1955  if (Parallel())
1956  {
1957 #ifdef MFEM_USE_MPI
1958  MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
1959  Y->GetComm());
1960 #endif
1961  }
1962 
1963  if (sundials_mem)
1964  {
1965  // Check if the problem size has changed since the last SetOperator call
1966  int resize = 0;
1967  if (!Parallel())
1968  {
1969  resize = (Y->Size() != local_size);
1970  }
1971  else
1972  {
1973 #ifdef MFEM_USE_MPI
1974  int l_resize = (Y->Size() != local_size) ||
1975  (saved_global_size != global_size);
1976  MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR,
1977  Y->GetComm());
1978 #endif
1979  }
1980 
1981  // Free existing solver memory and re-create with new vector size
1982  if (resize)
1983  {
1984  KINFree(&sundials_mem);
1985  sundials_mem = NULL;
1986  }
1987  }
1988 
1989  if (!sundials_mem)
1990  {
1991  if (!Parallel())
1992  {
1993  Y->SetSize(local_size);
1994  }
1995 #ifdef MFEM_USE_MPI
1996  else
1997  {
1998  Y->SetSize(local_size, global_size);
1999  y_scale->SetSize(local_size, global_size);
2000  f_scale->SetSize(local_size, global_size);
2001  saved_global_size = global_size;
2002  }
2003 #endif
2004 
2005  // Create the solver memory
2007  MFEM_VERIFY(sundials_mem, "Error in KINCreate().");
2008 
2009  // Set number of acceleration vectors
2010  if (maa > 0)
2011  {
2012  flag = KINSetMAA(sundials_mem, maa);
2013  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetMAA()");
2014  }
2015 
2016  // Initialize KINSOL
2017  flag = KINInit(sundials_mem, KINSolver::Mult, *Y);
2018  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINInit()");
2019 
2020  // Attach the KINSolver as user-defined data
2021  flag = KINSetUserData(sundials_mem, this);
2022  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetUserData()");
2023 
2024  // Set the linear solver
2025  if (prec || jfnk)
2026  {
2028  }
2029  else
2030  {
2031  // Free any existing linear solver
2032  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
2033  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
2034 
2035  LSA = SUNLinSol_SPGMR(*Y, PREC_NONE, 0, Sundials::GetContext());
2036  MFEM_VERIFY(LSA, "error in SUNLinSol_SPGMR()");
2037 
2038  flag = KINSetLinearSolver(sundials_mem, LSA, NULL);
2039  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetLinearSolver()");
2040 
2041  // Set Jacobian-vector product function
2042  if (use_oper_grad)
2043  {
2044  flag = KINSetJacTimesVecFn(sundials_mem, KINSolver::GradientMult);
2045  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetJacTimesVecFn()");
2046  }
2047  }
2048  }
2049 }
2050 
2052 {
2053  if (jfnk)
2054  {
2055  SetJFNKSolver(solver);
2056  }
2057  else
2058  {
2059  // Store the solver
2060  prec = &solver;
2061 
2062  // Free any existing linear solver
2063  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
2064  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
2065 
2066  // Wrap KINSolver as SUNLinearSolver and SUNMatrix
2068  MFEM_VERIFY(LSA, "error in SUNLinSolNewEmpty()");
2069 
2070  LSA->content = this;
2071  LSA->ops->gettype = LSGetType;
2072  LSA->ops->solve = KINSolver::LinSysSolve;
2073  LSA->ops->free = LSFree;
2074 
2076  MFEM_VERIFY(A, "error in SUNMatNewEmpty()");
2077 
2078  A->content = this;
2079  A->ops->getid = MatGetID;
2080  A->ops->destroy = MatDestroy;
2081 
2082  // Attach the linear solver and matrix
2083  flag = KINSetLinearSolver(sundials_mem, LSA, A);
2084  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINSetLinearSolver()");
2085 
2086  // Set the Jacobian evaluation function
2087  flag = KINSetJacFn(sundials_mem, KINSolver::LinSysSetup);
2088  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINSetJacFn()");
2089  }
2090 }
2091 
2093 {
2094  // Store the solver
2095  prec = &solver;
2096 
2097  wrk.SetSize(height);
2098 
2099  // Free any existing linear solver
2100  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
2101  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
2102 
2103  // Setup FGMRES
2104  LSA = SUNLinSol_SPFGMR(*Y, prec ? PREC_RIGHT : PREC_NONE, maxli,
2106  MFEM_VERIFY(LSA, "error in SUNLinSol_SPFGMR()");
2107 
2108  flag = SUNLinSol_SPFGMRSetMaxRestarts(LSA, maxlrs);
2109  MFEM_VERIFY(flag == SUNLS_SUCCESS, "error in SUNLinSol_SPFGMR()");
2110 
2111  flag = KINSetLinearSolver(sundials_mem, LSA, NULL);
2112  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINSetLinearSolver()");
2113 
2114  if (prec)
2115  {
2116  flag = KINSetPreconditioner(sundials_mem,
2119  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINSetPreconditioner()");
2120  }
2121 }
2122 
2124 {
2125  flag = KINSetScaledStepTol(sundials_mem, sstol);
2126  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetScaledStepTol()");
2127 }
2128 
2129 void KINSolver::SetMaxSetupCalls(int max_calls)
2130 {
2131  flag = KINSetMaxSetupCalls(sundials_mem, max_calls);
2132  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetMaxSetupCalls()");
2133 }
2134 
2135 void KINSolver::SetMAA(int m_aa)
2136 {
2137  // Store internally as maa must be set before calling KINInit() to
2138  // set the maximum acceleration space size.
2139  maa = m_aa;
2140  if (sundials_mem)
2141  {
2142  flag = KINSetMAA(sundials_mem, maa);
2143  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetMAA()");
2144  }
2145 }
2146 
2148 {
2149  MFEM_ABORT("this method is not supported! Use SetPrintLevel(int) instead.");
2150 }
2151 
2152 // Compute the scaling vectors and solve nonlinear system
2153 void KINSolver::Mult(const Vector&, Vector &x) const
2154 {
2155  // residual norm tolerance
2156  double tol;
2157 
2158  // Uses c = 1, corresponding to x_scale.
2159  c = 1.0;
2160 
2161  if (!iterative_mode) { x = 0.0; }
2162 
2163  // For relative tolerance, r = 1 / |residual(x)|, corresponding to fx_scale.
2164  if (rel_tol > 0.0)
2165  {
2166 
2167  oper->Mult(x, r);
2168 
2169  // Note that KINSOL uses infinity norms.
2170  double norm = r.Normlinf();
2171 #ifdef MFEM_USE_MPI
2172  if (Parallel())
2173  {
2174  double lnorm = norm;
2175  MPI_Allreduce(&lnorm, &norm, 1, MPI_DOUBLE, MPI_MAX,
2176  Y->GetComm());
2177  }
2178 #endif
2179  if (abs_tol > rel_tol * norm)
2180  {
2181  r = 1.0;
2182  tol = abs_tol;
2183  }
2184  else
2185  {
2186  r = 1.0 / norm;
2187  tol = rel_tol;
2188  }
2189  }
2190  else
2191  {
2192  r = 1.0;
2193  tol = abs_tol;
2194  }
2195 
2196  // Set the residual norm tolerance
2197  flag = KINSetFuncNormTol(sundials_mem, tol);
2198  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetFuncNormTol()");
2199 
2200  // Solve the nonlinear system by calling the other Mult method
2201  KINSolver::Mult(x, c, r);
2202 }
2203 
2204 // Solve the nonlinear system using the provided scaling vectors
2206  const Vector &x_scale, const Vector &fx_scale) const
2207 {
2208  flag = KINSetNumMaxIters(sundials_mem, max_iter);
2209  MFEM_ASSERT(flag == KIN_SUCCESS, "KINSetNumMaxIters() failed!");
2210 
2211  Y->MakeRef(x, 0, x.Size());
2212  y_scale->MakeRef(const_cast<Vector&>(x_scale), 0, x_scale.Size());
2213  f_scale->MakeRef(const_cast<Vector&>(fx_scale), 0, fx_scale.Size());
2214 
2215  int rank = -1;
2216  if (!Parallel())
2217  {
2218  rank = 0;
2219  }
2220  else
2221  {
2222 #ifdef MFEM_USE_MPI
2223  MPI_Comm_rank(Y->GetComm(), &rank);
2224 #endif
2225  }
2226 
2227  if (rank == 0)
2228  {
2229  flag = KINSetPrintLevel(sundials_mem, print_level);
2230  MFEM_VERIFY(flag == KIN_SUCCESS, "KINSetPrintLevel() failed!");
2231 
2232 #ifdef SUNDIALS_BUILD_WITH_MONITORING
2233  if (jfnk && print_level)
2234  {
2235  flag = SUNLinSolSetInfoFile_SPFGMR(LSA, stdout);
2236  MFEM_VERIFY(flag == SUNLS_SUCCESS,
2237  "error in SUNLinSolSetInfoFile_SPFGMR()");
2238 
2239  flag = SUNLinSolSetPrintLevel_SPFGMR(LSA, 1);
2240  MFEM_VERIFY(flag == SUNLS_SUCCESS,
2241  "error in SUNLinSolSetPrintLevel_SPFGMR()");
2242  }
2243 #endif
2244  }
2245 
2246  if (!iterative_mode) { x = 0.0; }
2247 
2248  // Solve the nonlinear system
2249  flag = KINSol(sundials_mem, *Y, global_strategy, *y_scale, *f_scale);
2250  converged = (flag >= 0);
2251 
2252  // Make sure host is up to date
2253  Y->HostRead();
2254 
2255  // Get number of nonlinear iterations
2256  long int tmp_nni;
2257  flag = KINGetNumNonlinSolvIters(sundials_mem, &tmp_nni);
2258  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINGetNumNonlinSolvIters()");
2259  final_iter = (int) tmp_nni;
2260 
2261  // Get the residual norm
2262  flag = KINGetFuncNorm(sundials_mem, &final_norm);
2263  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINGetFuncNorm()");
2264 }
2265 
2267 {
2268  delete Y;
2269  delete y_scale;
2270  delete f_scale;
2271  SUNMatDestroy(A);
2272  SUNLinSolFree(LSA);
2273  KINFree(&sundials_mem);
2274 }
2275 
2276 } // namespace mfem
2277 
2278 #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.
Definition: sundials.cpp:1206
void SetSize(int s, long glob_size=0)
Resize the vector to size s.
Definition: sundials.cpp:507
virtual double GetTime() const
Read the currently set time.
Definition: operator.hpp:330
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition: sundials.cpp:696
static bool IsAvailable()
Return true if an actual device (e.g. GPU) has been configured.
Definition: device.hpp:243
MFEM_DEPRECATED void * CVodeCreate(int lmm, SUNContext)
Definition: sundials.cpp:81
void SetJFNKSolver(Solver &solver)
Definition: sundials.cpp:2092
MFEM_DEPRECATED SUNLinearSolver SUNLinSolNewEmpty(SUNContext)
Definition: sundials.cpp:58
SundialsNVector * q
Quadrature vector.
Definition: sundials.hpp:532
void GetForwardSolution(double tB, mfem::Vector &yy)
Get Interpolated Forward solution y at backward integration time tB.
Definition: sundials.cpp:1012
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
Definition: sundials.cpp:1288
virtual int SUNImplicitSolveB(Vector &x, const Vector &b, double tol)
Solve the ODE linear system as setup by the method SUNImplicitSetup().
Definition: operator.hpp:586
const Operator * oper
Definition: solvers.hpp:120
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:862
std::function< int(Vector y, Vector w, CVODESolver *)> EWTFunction
Typedef declaration for error weight functions.
Definition: sundials.hpp:405
static int LinSysSolveB(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system A x = b.
Definition: sundials.cpp:1167
MFEM_DEPRECATED void * KINCreate(SUNContext)
Definition: sundials.cpp:96
MFEM_DEPRECATED void * ARKStepCreate(ARKRhsFn fe, ARKRhsFn fi, realtype t0, N_Vector y0, SUNContext)
Definition: sundials.cpp:88
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. ...
Definition: operator.hpp:547
~SundialsNVector()
Calls SUNDIALS N_VDestroy function if the N_Vector is owned by &#39;this&#39;.
Definition: sundials.cpp:493
void SetScaledStepTol(double sstol)
Set KINSOL&#39;s scaled step tolerance.
Definition: sundials.cpp:2123
RootFunction root_func
A class member to facilitate pointing to a user-specified root function.
Definition: sundials.hpp:402
void SetSize(int s)
Resize the vector to size s.
Definition: vector.hpp:513
SundialsNVector * Y
State vector.
Definition: sundials.hpp:326
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Definition: sundials.cpp:1695
void Delete()
Delete the owned pointers and reset the Memory object.
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
Definition: sundials.cpp:1718
void SetIRKTableNum(ARKODE_DIRKTableID table_id)
Choose a specific Butcher table for a diagonally implicit RK method.
Definition: sundials.cpp:1724
Base abstract class for first order time dependent operators.
Definition: operator.hpp:289
virtual Operator & GetGradient(const Vector &x) const
Evaluate the gradient operator at the point x. The default behavior in class Operator is to generate ...
Definition: operator.hpp:99
virtual void Step(Vector &x, double &t, double &dt)
Integrate the ODE with CVODE using the specified step mode.
Definition: sundials.cpp:784
static constexpr double default_abs_tolB
Default scalar backward absolute tolerance.
Definition: sundials.hpp:540
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.
void SetData(double *d)
Definition: sundials.cpp:513
SundialsNVector * f_scale
scaling vectors
Definition: sundials.hpp:840
void SetSVtolerances(double reltol, Vector abstol)
Set the scalar relative and vector of absolute tolerances.
Definition: sundials.cpp:868
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.
Definition: sundials.cpp:1233
void MakeRef(Vector &base, int offset, int s)
Reset the Vector to be a reference to a sub-vector of base.
Definition: sundials.hpp:250
int indexB
backward problem index
Definition: sundials.hpp:514
virtual ~CVODESolver()
Destroy the associated CVODE memory and SUNDIALS objects.
Definition: sundials.cpp:951
bool use_oper_grad
use the Jv prod function
Definition: sundials.hpp:839
int Size() const
Returns the size of the vector.
Definition: vector.hpp:200
void SetDevicePtrOwner(bool own) const
Set/clear the ownership flag for the device pointer. Ownership indicates whether the pointer will be ...
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
Definition: sundials.cpp:992
MFEM_DEPRECATED SUNLinearSolver SUNLinSol_SPGMR(N_Vector y, int pretype, int maxl, SUNContext)
Definition: sundials.cpp:65
virtual void SetPrintLevel(int print_lvl)
Set the print level for the KINSetPrintLevel function.
Definition: sundials.hpp:959
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 .
Definition: sundials.cpp:1368
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:1706
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
Definition: sundials.cpp:1712
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.
Definition: operator.hpp:333
void * sundials_mem
SUNDIALS mem structure.
Definition: sundials.hpp:321
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:830
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 .
Definition: sundials.cpp:1894
long saved_global_size
Global vector length on last initialization.
Definition: sundials.hpp:324
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:659
virtual void Step(Vector &x, double &t, double &dt)
Definition: sundials.cpp:1260
SundialsNVector()
Creates an empty SundialsNVector.
Definition: sundials.cpp:434
void UseSundialsLinearSolverB()
Use built in SUNDIALS Newton solver.
Definition: sundials.cpp:1133
static int MassSysSetup(realtype t, SUNMatrix M, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
Definition: sundials.cpp:1402
MFEM_DEPRECATED SUNMatrix SUNMatNewEmpty(SUNContext)
Definition: sundials.cpp:51
void SetSStolerancesB(double reltol, double abstol)
Tolerance specification functions for the adjoint problem.
Definition: sundials.cpp:1180
Explicit RK method.
Definition: sundials.hpp:667
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
Definition: solvers.hpp:130
CVODESSolver(int lmm)
Definition: sundials.cpp:964
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
Definition: sundials.cpp:2129
void _SetDataAndSize_()
Set data and length from the internal N_Vector x.
Definition: sundials.cpp:386
ARKStepSolver(Type type=EXPLICIT)
Construct a serial wrapper to SUNDIALS&#39; ARKode integrator.
Definition: sundials.cpp:1445
void SetIMEXTableNum(ARKODE_ERKTableID etable_id, ARKODE_DIRKTableID itable_id)
Choose a specific Butcher table for an IMEX RK method.
Definition: sundials.cpp:1730
Implicit RK method.
Definition: sundials.hpp:668
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
Definition: sundials.cpp:1059
MFEM_DEPRECATED N_Vector N_VNewEmpty_Parallel(MPI_Comm comm, sunindextype local_length, sunindextype global_length, SUNContext)
Definition: sundials.cpp:105
MFEM_DEPRECATED SUNMemoryHelper SUNMemoryHelper_NewEmpty(SUNContext)
Definition: sundials.cpp:129
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:661
int global_strategy
KINSOL solution strategy.
Definition: sundials.hpp:838
SUNLinearSolver LSA
Linear solver for A.
Definition: sundials.hpp:330
static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
Compute the matrix-vector product .
Definition: sundials.cpp:1423
MFEM_DEPRECATED N_Vector N_VNewEmpty_Serial(sunindextype vec_length, SUNContext)
Definition: sundials.cpp:44
virtual void SetEvalMode(const EvalMode new_eval_mode)
Set the evaluation mode of the time-dependent operator.
Definition: operator.hpp:356
SUNMatrix M
Mass matrix M.
Definition: sundials.hpp:329
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: sundials.cpp:2051
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:67
bool Parallel() const
Definition: sundials.hpp:335
bool jfnk
enable JFNK
Definition: sundials.hpp:843
void SetSVtolerancesB(double reltol, Vector abstol)
Tolerance specification functions for the adjoint problem.
Definition: sundials.cpp:1186
MFEM_DEPRECATED SUNLinearSolver SUNLinSol_SPFGMR(N_Vector y, int pretype, int maxl, SUNContext)
Definition: sundials.cpp:73
int lmm_type
Linear multistep method type.
Definition: sundials.hpp:377
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:374
void InitQuadIntegration(mfem::Vector &q0, double reltolQ=1e-3, double abstolQ=1e-8)
Definition: sundials.cpp:1071
void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
Initialize Quadrature Integration (Adjoint)
Definition: sundials.cpp:1086
double b
Definition: lissajous.cpp:42
Interface to the KINSOL library – nonlinear solver methods.
Definition: sundials.hpp:835
static constexpr double default_rel_tolB
Default scalar backward relative tolerance.
Definition: sundials.hpp:538
SundialsNVector * yy
State vector.
Definition: sundials.hpp:534
Vector interface for SUNDIALS N_Vectors.
Definition: sundials.hpp:173
int GetAdjointHeight()
Returns the size of the adjoint problem state space.
Definition: operator.hpp:594
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:880
Type
Types of ARKODE solvers.
Definition: sundials.hpp:665
Type rk_type
Runge-Kutta type.
Definition: sundials.hpp:673
int ARKODE_ERKTableID
Definition: sundials.hpp:50
virtual ~ARKStepSolver()
Destroy the associated ARKode memory and SUNDIALS objects.
Definition: sundials.cpp:1795
SundialsNVector * qB
State vector.
Definition: sundials.hpp:535
bool use_implicit
True for implicit or imex integration.
Definition: sundials.hpp:675
SUNNonlinearSolver NLS
Nonlinear solver.
Definition: sundials.hpp:332
void SetData(double *d)
Definition: vector.hpp:150
constexpr ARKODE_DIRKTableID ARKODE_DIRK_NONE
Definition: sundials.hpp:53
CVODESolver(int lmm)
Construct a serial wrapper to SUNDIALS&#39; CVODE integrator.
Definition: sundials.cpp:682
virtual const double * HostRead() const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), false).
Definition: vector.hpp:453
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Wrapper to compute the nonlinear residual .
Definition: sundials.cpp:1809
int step_mode
ARKStep step mode (ARK_NORMAL or ARK_ONE_STEP).
Definition: sundials.hpp:674
SundialsNVector * y_scale
Definition: sundials.hpp:840
Wrapper for hypre&#39;s parallel vector class.
Definition: hypre.hpp:161
static int LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2)
Setup the linear system .
Definition: sundials.cpp:1846
SUNLinearSolver LSM
Linear solver for M.
Definition: sundials.hpp:331
int ncheck
number of checkpoints used so far
Definition: sundials.hpp:513
void * SUNContext
Definition: sundials.hpp:58
static int RHS1(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Definition: sundials.cpp:1331
N_Vector x
The actual SUNDIALS object.
Definition: sundials.hpp:179
MFEM_DEPRECATED N_Vector N_VMake_MPIPlusX(MPI_Comm comm, N_Vector local_vector, SUNContext)
Definition: sundials.cpp:140
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to CVODE.
Definition: sundials.cpp:811
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
Definition: sundials.cpp:1943
virtual ~KINSolver()
Destroy the associated KINSOL memory.
Definition: sundials.cpp:2266
void SetHostPtrOwner(bool own) const
Set/clear the ownership flag for the host pointer. Ownership indicates whether the pointer will be de...
int max_iter
Limit for the number of iterations the solver is allowed to do.
Definition: solvers.hpp:151
void SetRootFinder(int components, RootFunction func)
Initialize Root Finder.
Definition: sundials.cpp:650
Vector wrk
Work vector needed for the JFNK PC.
Definition: sundials.hpp:844
EWTFunction ewt_func
A class member to facilitate pointing to a user-specified error weight function.
Definition: sundials.hpp:408
void SetFixedStep(double dt)
Use a fixed time step size (disable temporal adaptivity).
Definition: sundials.cpp:1737
static int root(realtype t, N_Vector y, realtype *gout, void *user_data)
Prototype to define root finding for CVODE.
Definition: sundials.cpp:638
bool IsKnown(const void *h_ptr)
Return true if the pointer is known by the memory manager.
int ARKODE_DIRKTableID
Definition: sundials.hpp:51
bool reinit
Flag to signal memory reinitialization is need.
Definition: sundials.hpp:323
int maa
number of acceleration vectors
Definition: sundials.hpp:842
static N_Vector MakeNVector(bool use_device)
Create a N_Vector.
Definition: sundials.cpp:525
void UseMFEMLinearSolverB()
Set Linear Solver for the backward problem.
Definition: sundials.cpp:1101
constexpr ARKODE_ERKTableID ARKODE_ERK_NONE
Definition: sundials.hpp:52
void SetDataAndSize(double *d, int s)
Set the Vector data and size.
Definition: vector.hpp:157
std::function< int(realtype t, Vector y, Vector gout, CVODESolver *)> RootFunction
Typedef for root finding functions.
Definition: sundials.hpp:399
void _SetNvecDataAndSize_(long glob_size=0)
Set data and length of internal N_Vector x from &#39;this&#39;.
Definition: sundials.cpp:313
Vector & Set(const double a, const Vector &x)
(*this) = a * x
Definition: vector.cpp:267
int height
Dimension of the output / number of rows in the matrix.
Definition: operator.hpp:27
void SetMAA(int maa)
Set the number of acceleration vectors to use with KIN_FP or KIN_PICARD.
Definition: sundials.cpp:2135
virtual void QuadratureIntegration(const Vector &y, Vector &qdot) const
Provide the operator integration of a quadrature equation.
Definition: operator.hpp:526
int maxli
Maximum linear iterations.
Definition: sundials.hpp:845
MemoryManager mm
The (single) global memory manager object.
long GetNumSteps()
Get the number of internal steps taken so far.
Definition: sundials.cpp:892
double rel_tol
Relative tolerance.
Definition: solvers.hpp:154
T * ReadWrite(Memory< T > &mem, int size, bool on_dev=true)
Get a pointer for read+write access to mem with the mfem::Device&#39;s DeviceMemoryClass, if on_dev = true, or the mfem::Device&#39;s HostMemoryClass, otherwise.
Definition: device.hpp:353
MPI_Comm GetComm() const
Returns the MPI communicator for the internal N_Vector x.
Definition: sundials.hpp:230
void SetDataAndSize(double *d, int s, long glob_size=0)
Set the vector data and size.
Definition: sundials.cpp:519
MFEM_DEPRECATED N_Vector N_VNewWithMemHelp_Cuda(sunindextype length, booleantype use_managed_mem, SUNMemoryHelper helper, SUNContext)
Definition: sundials.cpp:119
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 .
Definition: sundials.cpp:658
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:1743
static int RHS2(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Definition: sundials.cpp:1351
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 .
Definition: sundials.cpp:1823
const Operator * jacobian
stores oper-&gt;GetGradient()
Definition: sundials.hpp:841
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
Definition: sundials.cpp:1026
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.
Definition: sundials.cpp:1148
SUNLinearSolver LSB
Linear solver for A.
Definition: sundials.hpp:531
Singleton class for SUNContext and SundialsMemHelper objects.
Definition: sundials.hpp:137
SundialsNVector * yB
State vector.
Definition: sundials.hpp:533
static constexpr double default_rel_tol
Default scalar relative tolerance.
Definition: sundials.hpp:342
void SetWFTolerances(EWTFunction func)
Set multiplicative error weights.
Definition: sundials.cpp:1198
SUNMatrix AB
Linear system A = I - gamma J, M - gamma J, or J.
Definition: sundials.hpp:530
static bool UseManagedMemory()
Definition: sundials.hpp:304
int step_mode
CVODE step mode (CV_NORMAL or CV_ONE_STEP).
Definition: sundials.hpp:378
static int MassMult2(N_Vector x, N_Vector v, realtype t, void *mtimes_data)
Compute the matrix-vector product at time t.
Definition: sundials.cpp:1433
KINSolver(int strategy, bool oper_grad=true)
Construct a serial wrapper to SUNDIALS&#39; KINSOL nonlinear solver.
Definition: sundials.cpp:1914
Implicit-explicit ARK method.
Definition: sundials.hpp:669
RefCoord t[3]
static constexpr double default_abs_tol
Default scalar absolute tolerance.
Definition: sundials.hpp:344
int flag
Last flag returned from a call to SUNDIALS.
Definition: sundials.hpp:322
static SundialsMemHelper & GetMemHelper()
Provides access to the SundialsMemHelper object.
Definition: sundials.cpp:170
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to ARKode.
Definition: sundials.cpp:1597
Vector data type.
Definition: vector.hpp:60
virtual bool UseDevice() const
Return the device flag of the Memory object used by the Vector.
Definition: vector.hpp:121
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:1387
static SUNContext & GetContext()
Provides access to the SUNContext object.
Definition: sundials.cpp:165
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:672
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:906
void Init(TimeDependentAdjointOperator &f_)
Definition: sundials.cpp:1021
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition: ode.hpp:26
static int PrecSetup(N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, void *user_data)
Setup the preconditioner.
Definition: sundials.cpp:1876
RefCoord s[3]
void SetMaxOrder(int max_order)
Set the maximum method order.
Definition: sundials.cpp:900
virtual ~CVODESSolver()
Destroy the associated CVODES memory and SUNDIALS objects.
Definition: sundials.cpp:1316
static int ewt(N_Vector y, N_Vector w, void *user_data)
Error control function.
Definition: sundials.cpp:1249
double u(const Vector &xvec)
Definition: lor_mms.hpp:24
void UseMFEMMassLinearSolver(int tdep)
Attach mass matrix linear system setup, solve, and matrix-vector product methods from the TimeDepende...
Definition: sundials.cpp:1643
Base class for solvers.
Definition: operator.hpp:655
void SetMaxNStepsB(int mxstepsB)
Set the maximum number of backward steps.
Definition: sundials.cpp:1065
OutStream out(std::cout)
Global stream used by the library for standard output. Initially it uses the same std::streambuf as s...
Definition: globals.hpp:66
Abstract operator.
Definition: operator.hpp:24
void Init(TimeDependentOperator &f_)
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults...
Definition: sundials.cpp:1461
T * HostReadWrite(Memory< T > &mem, int size)
Shortcut to ReadWrite(Memory&lt;T&gt; &amp;mem, int size, false)
Definition: device.hpp:360
static int RHS(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Number of components in gout.
Definition: sundials.cpp:621
int maxlrs
Maximum linear solver restarts.
Definition: sundials.hpp:846
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:1700
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
Definition: sundials.cpp:1628
virtual void Init(TimeDependentOperator &f_)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:1863
static int MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:1412
bool MPIPlusX() const
Definition: sundials.hpp:286
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 .
Definition: operator.hpp:566
Settings for the output behavior of the IterativeSolver.
Definition: solvers.hpp:78
virtual void Step(Vector &x, double &t, double &dt)
Integrate the ODE with ARKode using the specified step mode.
Definition: sundials.cpp:1557
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.
Definition: sundials.cpp:1219
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
Definition: sundials.cpp:1002
double abs_tol
Absolute tolerance.
Definition: solvers.hpp:157
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition: sundials.cpp:842
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1803
void SetMaxNSteps(int steps)
Set the maximum number of time steps.
Definition: sundials.cpp:886
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
Definition: sundials.cpp:857
void UseSundialsMassLinearSolver(int tdep)
Attach the SUNDIALS GMRES linear solver and the mass matrix matrix-vector product method from the Tim...
Definition: sundials.cpp:1675