MFEM  v4.1.0
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-2020, 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_MPI
24 #include <nvector/nvector_parallel.h>
25 #endif
26 
27 // SUNDIALS linear solvers
28 #include <sunlinsol/sunlinsol_spgmr.h>
29 
30 // Access SUNDIALS object's content pointer
31 #define GET_CONTENT(X) ( X->content )
32 
33 using namespace std;
34 
35 namespace mfem
36 {
37 
38 // ---------------------------------------------------------------------------
39 // SUNMatrix interface functions
40 // ---------------------------------------------------------------------------
41 
42 // Return the matrix ID
43 static SUNMatrix_ID MatGetID(SUNMatrix A)
44 {
45  return (SUNMATRIX_CUSTOM);
46 }
47 
48 static void MatDestroy(SUNMatrix A)
49 {
50  if (A->content) { A->content = NULL; }
51  if (A->ops) { free(A->ops); A->ops = NULL; }
52  free(A); A = NULL;
53  return;
54 }
55 
56 // ---------------------------------------------------------------------------
57 // SUNLinearSolver interface functions
58 // ---------------------------------------------------------------------------
59 
60 // Return the linear solver type
61 static SUNLinearSolver_Type LSGetType(SUNLinearSolver LS)
62 {
63  return (SUNLINEARSOLVER_MATRIX_ITERATIVE);
64 }
65 
66 static int LSFree(SUNLinearSolver LS)
67 {
68  if (LS->content) { LS->content = NULL; }
69  if (LS->ops) { free(LS->ops); LS->ops = NULL; }
70  free(LS); LS = NULL;
71  return (0);
72 }
73 
74 // ---------------------------------------------------------------------------
75 // CVODE interface
76 // ---------------------------------------------------------------------------
77 
78 int CVODESolver::RHS(realtype t, const N_Vector y, N_Vector ydot,
79  void *user_data)
80 {
81  // Get data from N_Vectors
82  const Vector mfem_y(y);
83  Vector mfem_ydot(ydot);
84  CVODESolver *self = static_cast<CVODESolver*>(user_data);
85 
86  // Compute y' = f(t, y)
87  self->f->SetTime(t);
88  self->f->Mult(mfem_y, mfem_ydot);
89 
90  // Return success
91  return (0);
92 }
93 
94 int CVODESolver::LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix A,
95  booleantype jok, booleantype *jcur,
96  realtype gamma, void *user_data, N_Vector tmp1,
97  N_Vector tmp2, N_Vector tmp3)
98 {
99  // Get data from N_Vectors
100  const Vector mfem_y(y);
101  const Vector mfem_fy(fy);
102  CVODESolver *self = static_cast<CVODESolver*>(GET_CONTENT(A));
103 
104  // Compute the linear system
105  self->f->SetTime(t);
106  return (self->f->SUNImplicitSetup(mfem_y, mfem_fy, jok, jcur, gamma));
107 }
108 
109 int CVODESolver::LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x,
110  N_Vector b, realtype tol)
111 {
112  Vector mfem_x(x);
113  const Vector mfem_b(b);
114  CVODESolver *self = static_cast<CVODESolver*>(GET_CONTENT(LS));
115 
116  // Solve the linear system
117  return (self->f->SUNImplicitSolve(mfem_b, mfem_x, tol));
118 }
119 
120 CVODESolver::CVODESolver(int lmm)
121  : lmm_type(lmm), step_mode(CV_NORMAL)
122 {
123  // Allocate an empty serial N_Vector
124  y = N_VNewEmpty_Serial(0);
125  MFEM_VERIFY(y, "error in N_VNewEmpty_Serial()");
126 }
127 
128 #ifdef MFEM_USE_MPI
129 CVODESolver::CVODESolver(MPI_Comm comm, int lmm)
130  : lmm_type(lmm), step_mode(CV_NORMAL)
131 {
132  if (comm == MPI_COMM_NULL)
133  {
134 
135  // Allocate an empty serial N_Vector
136  y = N_VNewEmpty_Serial(0);
137  MFEM_VERIFY(y, "error in N_VNewEmpty_Serial()");
138 
139  }
140  else
141  {
142 
143  // Allocate an empty parallel N_Vector
144  y = N_VNewEmpty_Parallel(comm, 0, 0); // calls MPI_Allreduce()
145  MFEM_VERIFY(y, "error in N_VNewEmpty_Parallel()");
146 
147  }
148 }
149 #endif
150 
152 {
153  // Initialize the base class
154  ODESolver::Init(f_);
155 
156  // Get the vector length
157  long local_size = f_.Height();
158 #ifdef MFEM_USE_MPI
159  long global_size;
160 #endif
161 
162  if (Parallel())
163  {
164 #ifdef MFEM_USE_MPI
165  MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
166  NV_COMM_P(y));
167 #endif
168  }
169 
170  // Get current time
171  double t = f_.GetTime();
172 
173  if (sundials_mem)
174  {
175  // Check if the problem size has changed since the last Init() call
176  int resize = 0;
177  if (!Parallel())
178  {
179  resize = (NV_LENGTH_S(y) != local_size);
180  }
181  else
182  {
183 #ifdef MFEM_USE_MPI
184  int l_resize = (NV_LOCLENGTH_P(y) != local_size) ||
185  (saved_global_size != global_size);
186  MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR, NV_COMM_P(y));
187 #endif
188  }
189 
190  // Free existing solver memory and re-create with new vector size
191  if (resize)
192  {
193  CVodeFree(&sundials_mem);
194  sundials_mem = NULL;
195  }
196  }
197 
198  if (!sundials_mem)
199  {
200  // Temporarly set N_Vector wrapper data to create CVODE. The correct
201  // initial condition will be set using CVodeReInit() when Step() is
202  // called.
203  if (!Parallel())
204  {
205  NV_LENGTH_S(y) = local_size;
206  NV_DATA_S(y) = new double[local_size](); // value-initialize
207  }
208  else
209  {
210 #ifdef MFEM_USE_MPI
211  NV_LOCLENGTH_P(y) = local_size;
212  NV_GLOBLENGTH_P(y) = global_size;
213  saved_global_size = global_size;
214  NV_DATA_P(y) = new double[local_size](); // value-initialize
215 #endif
216  }
217 
218  // Create CVODE
219  sundials_mem = CVodeCreate(lmm_type);
220  MFEM_VERIFY(sundials_mem, "error in CVodeCreate()");
221 
222  // Initialize CVODE
223  flag = CVodeInit(sundials_mem, CVODESolver::RHS, t, y);
224  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeInit()");
225 
226  // Attach the CVODESolver as user-defined data
227  flag = CVodeSetUserData(sundials_mem, this);
228  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetUserData()");
229 
230  // Set default tolerances
231  flag = CVodeSStolerances(sundials_mem, default_rel_tol, default_abs_tol);
232  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetSStolerances()");
233 
234  // Attach MFEM linear solver by default
236 
237  // Delete the allocated data in y.
238  if (!Parallel())
239  {
240  delete [] NV_DATA_S(y);
241  NV_DATA_S(y) = NULL;
242  }
243  else
244  {
245 #ifdef MFEM_USE_MPI
246  delete [] NV_DATA_P(y);
247  NV_DATA_P(y) = NULL;
248 #endif
249  }
250  }
251 
252  // Set the reinit flag to call CVodeReInit() in the next Step() call.
253  reinit = true;
254 }
255 
256 void CVODESolver::Step(Vector &x, double &t, double &dt)
257 {
258  if (!Parallel())
259  {
260  NV_DATA_S(y) = x.GetData();
261  MFEM_VERIFY(NV_LENGTH_S(y) == x.Size(), "");
262  }
263  else
264  {
265 #ifdef MFEM_USE_MPI
266  NV_DATA_P(y) = x.GetData();
267  MFEM_VERIFY(NV_LOCLENGTH_P(y) == x.Size(), "");
268 #endif
269  }
270 
271  // Reinitialize CVODE memory if needed
272  if (reinit)
273  {
274  flag = CVodeReInit(sundials_mem, t, y);
275  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeReInit()");
276 
277  // reset flag
278  reinit = false;
279  }
280 
281  // Integrate the system
282  double tout = t + dt;
283  flag = CVode(sundials_mem, tout, y, &t, step_mode);
284  MFEM_VERIFY(flag >= 0, "error in CVode()");
285 
286  // Return the last incremental step size
287  flag = CVodeGetLastStep(sundials_mem, &dt);
288  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetLastStep()");
289 }
290 
292 {
293  // Free any existing matrix and linear solver
294  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
295  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
296 
297  // Wrap linear solver as SUNLinearSolver and SUNMatrix
298  LSA = SUNLinSolNewEmpty();
299  MFEM_VERIFY(LSA, "error in SUNLinSolNewEmpty()");
300 
301  LSA->content = this;
302  LSA->ops->gettype = LSGetType;
303  LSA->ops->solve = CVODESolver::LinSysSolve;
304  LSA->ops->free = LSFree;
305 
306  A = SUNMatNewEmpty();
307  MFEM_VERIFY(A, "error in SUNMatNewEmpty()");
308 
309  A->content = this;
310  A->ops->getid = MatGetID;
311  A->ops->destroy = MatDestroy;
312 
313  // Attach the linear solver and matrix
314  flag = CVodeSetLinearSolver(sundials_mem, LSA, A);
315  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinearSolver()");
316 
317  // Set the linear system evaluation function
318  flag = CVodeSetLinSysFn(sundials_mem, CVODESolver::LinSysSetup);
319  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinSysFn()");
320 }
321 
323 {
324  // Free any existing matrix and linear solver
325  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
326  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
327 
328  // Create linear solver
329  LSA = SUNLinSol_SPGMR(y, PREC_NONE, 0);
330  MFEM_VERIFY(LSA, "error in SUNLinSol_SPGMR()");
331 
332  // Attach linear solver
333  flag = CVodeSetLinearSolver(sundials_mem, LSA, NULL);
334  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetLinearSolver()");
335 }
336 
338 {
339  step_mode = itask;
340 }
341 
342 void CVODESolver::SetSStolerances(double reltol, double abstol)
343 {
344  flag = CVodeSStolerances(sundials_mem, reltol, abstol);
345  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSStolerances()");
346 }
347 
348 void CVODESolver::SetMaxStep(double dt_max)
349 {
350  flag = CVodeSetMaxStep(sundials_mem, dt_max);
351  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetMaxStep()");
352 }
353 
354 void CVODESolver::SetMaxOrder(int max_order)
355 {
356  flag = CVodeSetMaxOrd(sundials_mem, max_order);
357  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeSetMaxOrd()");
358 }
359 
361 {
362  long int nsteps, nfevals, nlinsetups, netfails;
363  int qlast, qcur;
364  double hinused, hlast, hcur, tcur;
365  long int nniters, nncfails;
366 
367  // Get integrator stats
368  flag = CVodeGetIntegratorStats(sundials_mem,
369  &nsteps,
370  &nfevals,
371  &nlinsetups,
372  &netfails,
373  &qlast,
374  &qcur,
375  &hinused,
376  &hlast,
377  &hcur,
378  &tcur);
379  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetIntegratorStats()");
380 
381  // Get nonlinear solver stats
382  flag = CVodeGetNonlinSolvStats(sundials_mem,
383  &nniters,
384  &nncfails);
385  MFEM_VERIFY(flag == CV_SUCCESS, "error in CVodeGetNonlinSolvStats()");
386 
387  mfem::out <<
388  "CVODE:\n"
389  "num steps: " << nsteps << "\n"
390  "num rhs evals: " << nfevals << "\n"
391  "num lin setups: " << nlinsetups << "\n"
392  "num nonlin sol iters: " << nniters << "\n"
393  "num nonlin conv fail: " << nncfails << "\n"
394  "num error test fails: " << netfails << "\n"
395  "last order: " << qlast << "\n"
396  "current order: " << qcur << "\n"
397  "initial dt: " << hinused << "\n"
398  "last dt: " << hlast << "\n"
399  "current dt: " << hcur << "\n"
400  "current t: " << tcur << "\n" << endl;
401 
402  return;
403 }
404 
406 {
407  N_VDestroy(y);
408  SUNMatDestroy(A);
409  SUNLinSolFree(LSA);
410  SUNNonlinSolFree(NLS);
411  CVodeFree(&sundials_mem);
412 }
413 
414 // ---------------------------------------------------------------------------
415 // ARKStep interface
416 // ---------------------------------------------------------------------------
417 
418 int ARKStepSolver::RHS1(realtype t, const N_Vector y, N_Vector ydot,
419  void *user_data)
420 {
421  // Get data from N_Vectors
422  const Vector mfem_y(y);
423  Vector mfem_ydot(ydot);
424  ARKStepSolver *self = static_cast<ARKStepSolver*>(user_data);
425 
426  // Compute f(t, y) in y' = f(t, y) or fe(t, y) in y' = fe(t, y) + fi(t, y)
427  self->f->SetTime(t);
428  if (self->rk_type == IMEX)
429  {
430  self->f->SetEvalMode(TimeDependentOperator::ADDITIVE_TERM_1);
431  }
432  self->f->Mult(mfem_y, mfem_ydot);
433 
434  // Return success
435  return (0);
436 }
437 
438 int ARKStepSolver::RHS2(realtype t, const N_Vector y, N_Vector ydot,
439  void *user_data)
440 {
441  // Get data from N_Vectors
442  const Vector mfem_y(y);
443  Vector mfem_ydot(ydot);
444  ARKStepSolver *self = static_cast<ARKStepSolver*>(user_data);
445 
446  // Compute fi(t, y) in y' = fe(t, y) + fi(t, y)
447  self->f->SetTime(t);
448  self->f->SetEvalMode(TimeDependentOperator::ADDITIVE_TERM_2);
449  self->f->Mult(mfem_y, mfem_ydot);
450 
451  // Return success
452  return (0);
453 }
454 
455 int ARKStepSolver::LinSysSetup(realtype t, N_Vector y, N_Vector fy, SUNMatrix A,
456  SUNMatrix M, booleantype jok, booleantype *jcur,
457  realtype gamma, void *user_data, N_Vector tmp1,
458  N_Vector tmp2, N_Vector tmp3)
459 {
460  // Get data from N_Vectors
461  const Vector mfem_y(y);
462  const Vector mfem_fy(fy);
463  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(A));
464 
465  // Compute the linear system
466  self->f->SetTime(t);
467  if (self->rk_type == IMEX)
468  {
469  self->f->SetEvalMode(TimeDependentOperator::ADDITIVE_TERM_2);
470  }
471  return (self->f->SUNImplicitSetup(mfem_y, mfem_fy, jok, jcur, gamma));
472 }
473 
474 int ARKStepSolver::LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x,
475  N_Vector b, realtype tol)
476 {
477  Vector mfem_x(x);
478  const Vector mfem_b(b);
479  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(LS));
480 
481  // Solve the linear system
482  if (self->rk_type == IMEX)
483  {
485  }
486  return (self->f->SUNImplicitSolve(mfem_b, mfem_x, tol));
487 }
488 
489 int ARKStepSolver::MassSysSetup(realtype t, SUNMatrix M, void *user_data,
490  N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
491 {
492  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(M));
493 
494  // Compute the mass matrix system
495  self->f->SetTime(t);
496  return (self->f->SUNMassSetup());
497 }
498 
499 int ARKStepSolver::MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x,
500  N_Vector b, realtype tol)
501 {
502  Vector mfem_x(x);
503  const Vector mfem_b(b);
504  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(LS));
505 
506  // Solve the mass matrix system
507  return (self->f->SUNMassSolve(mfem_b, mfem_x, tol));
508 }
509 
510 int ARKStepSolver::MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
511 {
512  const Vector mfem_x(x);
513  Vector mfem_v(v);
514  ARKStepSolver *self = static_cast<ARKStepSolver*>(GET_CONTENT(M));
515 
516  // Compute the mass matrix-vector product
517  return (self->f->SUNMassMult(mfem_x, mfem_v));
518 }
519 
520 int ARKStepSolver::MassMult2(N_Vector x, N_Vector v, realtype t,
521  void* mtimes_data)
522 {
523  const Vector mfem_x(x);
524  Vector mfem_v(v);
525  ARKStepSolver *self = static_cast<ARKStepSolver*>(mtimes_data);
526 
527  // Compute the mass matrix-vector product
528  self->f->SetTime(t);
529  return (self->f->SUNMassMult(mfem_x, mfem_v));
530 }
531 
533  : rk_type(type), step_mode(ARK_NORMAL),
534  use_implicit(type == IMPLICIT || type == IMEX)
535 {
536  // Allocate an empty serial N_Vector
537  y = N_VNewEmpty_Serial(0);
538  MFEM_VERIFY(y, "error in N_VNewEmpty_Serial()");
539 }
540 
541 #ifdef MFEM_USE_MPI
542 ARKStepSolver::ARKStepSolver(MPI_Comm comm, Type type)
543  : rk_type(type), step_mode(ARK_NORMAL),
544  use_implicit(type == IMPLICIT || type == IMEX)
545 {
546  if (comm == MPI_COMM_NULL)
547  {
548 
549  // Allocate an empty serial N_Vector
550  y = N_VNewEmpty_Serial(0);
551  MFEM_VERIFY(y, "error in N_VNewEmpty_Serial()");
552 
553  }
554  else
555  {
556 
557  // Allocate an empty parallel N_Vector
558  y = N_VNewEmpty_Parallel(comm, 0, 0); // calls MPI_Allreduce()
559  MFEM_VERIFY(y, "error in N_VNewEmpty_Parallel()");
560 
561  }
562 }
563 #endif
564 
566 {
567  // Initialize the base class
568  ODESolver::Init(f_);
569 
570  // Get the vector length
571  long local_size = f_.Height();
572 #ifdef MFEM_USE_MPI
573  long global_size;
574 #endif
575 
576  if (Parallel())
577  {
578 #ifdef MFEM_USE_MPI
579  MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
580  NV_COMM_P(y));
581 #endif
582  }
583 
584  // Get current time
585  double t = f_.GetTime();
586 
587  if (sundials_mem)
588  {
589  // Check if the problem size has changed since the last Init() call
590  int resize = 0;
591  if (!Parallel())
592  {
593  resize = (NV_LENGTH_S(y) != local_size);
594  }
595  else
596  {
597 #ifdef MFEM_USE_MPI
598  int l_resize = (NV_LOCLENGTH_P(y) != local_size) ||
599  (saved_global_size != global_size);
600  MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR, NV_COMM_P(y));
601 #endif
602  }
603 
604  // Free existing solver memory and re-create with new vector size
605  if (resize)
606  {
607  ARKStepFree(&sundials_mem);
608  sundials_mem = NULL;
609  }
610  }
611 
612  if (!sundials_mem)
613  {
614 
615  // Temporarly set N_Vector wrapper data to create ARKStep. The correct
616  // initial condition will be set using ARKStepReInit() when Step() is
617  // called.
618  if (!Parallel())
619  {
620  NV_LENGTH_S(y) = local_size;
621  NV_DATA_S(y) = new double[local_size](); // value-initialize
622  }
623  else
624  {
625 #ifdef MFEM_USE_MPI
626  NV_LOCLENGTH_P(y) = local_size;
627  NV_GLOBLENGTH_P(y) = global_size;
628  saved_global_size = global_size;
629  NV_DATA_P(y) = new double[local_size](); // value-initialize
630 #endif
631  }
632 
633  // Create ARKStep memory
634  if (rk_type == IMPLICIT)
635  {
636  sundials_mem = ARKStepCreate(NULL, ARKStepSolver::RHS1, t, y);
637  }
638  else if (rk_type == EXPLICIT)
639  {
640  sundials_mem = ARKStepCreate(ARKStepSolver::RHS1, NULL, t, y);
641  }
642  else
643  {
645  t, y);
646  }
647  MFEM_VERIFY(sundials_mem, "error in ARKStepCreate()");
648 
649  // Attach the ARKStepSolver as user-defined data
650  flag = ARKStepSetUserData(sundials_mem, this);
651  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetUserData()");
652 
653  // Set default tolerances
654  flag = ARKStepSStolerances(sundials_mem, default_rel_tol, default_abs_tol);
655  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetSStolerances()");
656 
657  // If implicit, attach MFEM linear solver by default
659 
660  // Delete the allocated data in y.
661  if (!Parallel())
662  {
663  delete [] NV_DATA_S(y);
664  NV_DATA_S(y) = NULL;
665  }
666  else
667  {
668 #ifdef MFEM_USE_MPI
669  delete [] NV_DATA_P(y);
670  NV_DATA_P(y) = NULL;
671 #endif
672  }
673  }
674 
675  // Set the reinit flag to call ARKStepReInit() in the next Step() call.
676  reinit = true;
677 }
678 
679 void ARKStepSolver::Step(Vector &x, double &t, double &dt)
680 {
681  if (!Parallel())
682  {
683  NV_DATA_S(y) = x.GetData();
684  MFEM_VERIFY(NV_LENGTH_S(y) == x.Size(), "");
685  }
686  else
687  {
688 #ifdef MFEM_USE_MPI
689  NV_DATA_P(y) = x.GetData();
690  MFEM_VERIFY(NV_LOCLENGTH_P(y) == x.Size(), "");
691 #endif
692  }
693 
694  // Reinitialize ARKStep memory if needed
695  if (reinit)
696  {
697  if (rk_type == IMPLICIT)
698  {
699  flag = ARKStepReInit(sundials_mem, NULL, ARKStepSolver::RHS1, t, y);
700  }
701  else if (rk_type == EXPLICIT)
702  {
703  flag = ARKStepReInit(sundials_mem, ARKStepSolver::RHS1, NULL, t, y);
704  }
705  else
706  {
707  flag = ARKStepReInit(sundials_mem,
709  }
710  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepReInit()");
711 
712  // reset flag
713  reinit = false;
714  }
715 
716  // Integrate the system
717  double tout = t + dt;
718  flag = ARKStepEvolve(sundials_mem, tout, y, &t, step_mode);
719  MFEM_VERIFY(flag >= 0, "error in ARKStepEvolve()");
720 
721  // Return the last incremental step size
722  flag = ARKStepGetLastStep(sundials_mem, &dt);
723  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepGetLastStep()");
724 }
725 
727 {
728  // Free any existing matrix and linear solver
729  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
730  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
731 
732  // Wrap linear solver as SUNLinearSolver and SUNMatrix
733  LSA = SUNLinSolNewEmpty();
734  MFEM_VERIFY(LSA, "error in SUNLinSolNewEmpty()");
735 
736  LSA->content = this;
737  LSA->ops->gettype = LSGetType;
738  LSA->ops->solve = ARKStepSolver::LinSysSolve;
739  LSA->ops->free = LSFree;
740 
741  A = SUNMatNewEmpty();
742  MFEM_VERIFY(A, "error in SUNMatNewEmpty()");
743 
744  A->content = this;
745  A->ops->getid = MatGetID;
746  A->ops->destroy = MatDestroy;
747 
748  // Attach the linear solver and matrix
749  flag = ARKStepSetLinearSolver(sundials_mem, LSA, A);
750  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinearSolver()");
751 
752  // Set the linear system evaluation function
753  flag = ARKStepSetLinSysFn(sundials_mem, ARKStepSolver::LinSysSetup);
754  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinSysFn()");
755 }
756 
758 {
759  // Free any existing matrix and linear solver
760  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
761  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
762 
763  // Create linear solver
764  LSA = SUNLinSol_SPGMR(y, PREC_NONE, 0);
765  MFEM_VERIFY(LSA, "error in SUNLinSol_SPGMR()");
766 
767  // Attach linear solver
768  flag = ARKStepSetLinearSolver(sundials_mem, LSA, NULL);
769  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinearSolver()");
770 }
771 
773 {
774  // Free any existing matrix and linear solver
775  if (M != NULL) { SUNMatDestroy(M); M = NULL; }
776  if (LSM != NULL) { SUNLinSolFree(LSM); LSM = NULL; }
777 
778  // Wrap linear solver as SUNLinearSolver and SUNMatrix
779  LSM = SUNLinSolNewEmpty();
780  MFEM_VERIFY(LSM, "error in SUNLinSolNewEmpty()");
781 
782  LSM->content = this;
783  LSM->ops->gettype = LSGetType;
784  LSM->ops->solve = ARKStepSolver::MassSysSolve;
785  LSA->ops->free = LSFree;
786 
787  M = SUNMatNewEmpty();
788  MFEM_VERIFY(M, "error in SUNMatNewEmpty()");
789 
790  M->content = this;
791  M->ops->getid = SUNMatGetID;
792  M->ops->matvec = ARKStepSolver::MassMult1;
793  M->ops->destroy = MatDestroy;
794 
795  // Attach the linear solver and matrix
796  flag = ARKStepSetMassLinearSolver(sundials_mem, LSM, M, tdep);
797  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetLinearSolver()");
798 
799  // Set the linear system function
800  flag = ARKStepSetMassFn(sundials_mem, ARKStepSolver::MassSysSetup);
801  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMassFn()");
802 }
803 
805 {
806  // Free any existing matrix and linear solver
807  if (M != NULL) { SUNMatDestroy(A); M = NULL; }
808  if (LSM != NULL) { SUNLinSolFree(LSM); LSM = NULL; }
809 
810  // Create linear solver
811  LSM = SUNLinSol_SPGMR(y, PREC_NONE, 0);
812  MFEM_VERIFY(LSM, "error in SUNLinSol_SPGMR()");
813 
814  // Attach linear solver
815  flag = ARKStepSetMassLinearSolver(sundials_mem, LSM, NULL, tdep);
816  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMassLinearSolver()");
817 
818  // Attach matrix multiplication function
819  flag = ARKStepSetMassTimes(sundials_mem, NULL, ARKStepSolver::MassMult2,
820  this);
821  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMassTimes()");
822 }
823 
825 {
826  step_mode = itask;
827 }
828 
829 void ARKStepSolver::SetSStolerances(double reltol, double abstol)
830 {
831  flag = ARKStepSStolerances(sundials_mem, reltol, abstol);
832  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSStolerances()");
833 }
834 
835 void ARKStepSolver::SetMaxStep(double dt_max)
836 {
837  flag = ARKStepSetMaxStep(sundials_mem, dt_max);
838  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetMaxStep()");
839 }
840 
841 void ARKStepSolver::SetOrder(int order)
842 {
843  flag = ARKStepSetOrder(sundials_mem, order);
844  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetOrder()");
845 }
846 
847 void ARKStepSolver::SetERKTableNum(int table_num)
848 {
849  flag = ARKStepSetTableNum(sundials_mem, -1, table_num);
850  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetTableNum()");
851 }
852 
853 void ARKStepSolver::SetIRKTableNum(int table_num)
854 {
855  flag = ARKStepSetTableNum(sundials_mem, table_num, -1);
856  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetTableNum()");
857 }
858 
859 void ARKStepSolver::SetIMEXTableNum(int etable_num, int itable_num)
860 {
861  flag = ARKStepSetTableNum(sundials_mem, itable_num, etable_num);
862  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetTableNum()");
863 }
864 
866 {
867  flag = ARKStepSetFixedStep(sundials_mem, dt);
868  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepSetFixedStep()");
869 }
870 
872 {
873  long int nsteps, expsteps, accsteps, step_attempts;
874  long int nfe_evals, nfi_evals;
875  long int nlinsetups, netfails;
876  double hinused, hlast, hcur, tcur;
877  long int nniters, nncfails;
878 
879  // Get integrator stats
880  flag = ARKStepGetTimestepperStats(sundials_mem,
881  &expsteps,
882  &accsteps,
883  &step_attempts,
884  &nfe_evals,
885  &nfi_evals,
886  &nlinsetups,
887  &netfails);
888  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepGetTimestepperStats()");
889 
890  flag = ARKStepGetStepStats(sundials_mem,
891  &nsteps,
892  &hinused,
893  &hlast,
894  &hcur,
895  &tcur);
896 
897  // Get nonlinear solver stats
898  flag = ARKStepGetNonlinSolvStats(sundials_mem,
899  &nniters,
900  &nncfails);
901  MFEM_VERIFY(flag == ARK_SUCCESS, "error in ARKStepGetNonlinSolvStats()");
902 
903  mfem::out <<
904  "ARKStep:\n"
905  "num steps: " << nsteps << "\n"
906  "num exp rhs evals: " << nfe_evals << "\n"
907  "num imp rhs evals: " << nfi_evals << "\n"
908  "num lin setups: " << nlinsetups << "\n"
909  "num nonlin sol iters: " << nniters << "\n"
910  "num nonlin conv fail: " << nncfails << "\n"
911  "num steps attempted: " << step_attempts << "\n"
912  "num acc limited steps: " << accsteps << "\n"
913  "num exp limited stepfails: " << expsteps << "\n"
914  "num error test fails: " << netfails << "\n"
915  "initial dt: " << hinused << "\n"
916  "last dt: " << hlast << "\n"
917  "current dt: " << hcur << "\n"
918  "current t: " << tcur << "\n" << endl;
919 
920  return;
921 }
922 
924 {
925  N_VDestroy(y);
926  SUNMatDestroy(A);
927  SUNLinSolFree(LSA);
928  SUNNonlinSolFree(NLS);
929  ARKStepFree(&sundials_mem);
930 }
931 
932 // ---------------------------------------------------------------------------
933 // KINSOL interface
934 // ---------------------------------------------------------------------------
935 
936 // Wrapper for evaluating the nonlinear residual F(u) = 0
937 int KINSolver::Mult(const N_Vector u, N_Vector fu, void *user_data)
938 {
939  const Vector mfem_u(u);
940  Vector mfem_fu(fu);
941  KINSolver *self = static_cast<KINSolver*>(user_data);
942 
943  // Compute the non-linear action F(u).
944  self->oper->Mult(mfem_u, mfem_fu);
945 
946  // Return success
947  return 0;
948 }
949 
950 // Wrapper for computing Jacobian-vector products
951 int KINSolver::GradientMult(N_Vector v, N_Vector Jv, N_Vector u,
952  booleantype *new_u, void *user_data)
953 {
954  const Vector mfem_v(v);
955  Vector mfem_Jv(Jv);
956  KINSolver *self = static_cast<KINSolver*>(user_data);
957 
958  // Update Jacobian information if needed
959  if (*new_u)
960  {
961  const Vector mfem_u(u);
962  self->jacobian = &self->oper->GetGradient(mfem_u);
963  *new_u = SUNFALSE;
964  }
965 
966  // Compute the Jacobian-vector product
967  self->jacobian->Mult(mfem_v, mfem_Jv);
968 
969  // Return success
970  return 0;
971 }
972 
973 // Wrapper for evaluating linear systems J u = b
974 int KINSolver::LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J,
975  void *user_data, N_Vector tmp1, N_Vector tmp2)
976 {
977  const Vector mfem_u(u);
978  KINSolver *self = static_cast<KINSolver*>(GET_CONTENT(J));
979 
980  // Update the Jacobian
981  self->jacobian = &self->oper->GetGradient(mfem_u);
982 
983  // Set the Jacobian solve operator
984  self->prec->SetOperator(*self->jacobian);
985 
986  // Return success
987  return (0);
988 }
989 
990 // Wrapper for solving linear systems J u = b
991 int KINSolver::LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u,
992  N_Vector b, realtype tol)
993 {
994  Vector mfem_u(u), mfem_b(b);
995  KINSolver *self = static_cast<KINSolver*>(GET_CONTENT(LS));
996 
997  // Solve for u = [J(u)]^{-1} b, maybe approximately.
998  self->prec->Mult(mfem_b, mfem_u);
999 
1000  // Return success
1001  return (0);
1002 }
1003 
1004 KINSolver::KINSolver(int strategy, bool oper_grad)
1005  : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
1006  f_scale(NULL), jacobian(NULL), maa(0)
1007 {
1008  // Allocate empty serial N_Vectors
1009  y = N_VNewEmpty_Serial(0);
1010  y_scale = N_VNewEmpty_Serial(0);
1011  f_scale = N_VNewEmpty_Serial(0);
1012  MFEM_VERIFY(y && y_scale && f_scale, "Error in N_VNewEmpty_Serial().");
1013 
1014  // Default abs_tol and print_level
1015  abs_tol = pow(UNIT_ROUNDOFF, 1.0/3.0);
1016  print_level = 0;
1017 }
1018 
1019 #ifdef MFEM_USE_MPI
1020 KINSolver::KINSolver(MPI_Comm comm, int strategy, bool oper_grad)
1021  : global_strategy(strategy), use_oper_grad(oper_grad), y_scale(NULL),
1022  f_scale(NULL), jacobian(NULL), maa(0)
1023 {
1024  if (comm == MPI_COMM_NULL)
1025  {
1026 
1027  // Allocate empty serial N_Vectors
1028  y = N_VNewEmpty_Serial(0);
1029  y_scale = N_VNewEmpty_Serial(0);
1030  f_scale = N_VNewEmpty_Serial(0);
1031  MFEM_VERIFY(y && y_scale && f_scale, "error in N_VNewEmpty_Serial()");
1032 
1033  }
1034  else
1035  {
1036 
1037  // Allocate empty parallel N_Vectors
1038  y = N_VNewEmpty_Parallel(comm, 0, 0);
1039  y_scale = N_VNewEmpty_Parallel(comm, 0, 0);
1040  f_scale = N_VNewEmpty_Parallel(comm, 0, 0);
1041  MFEM_VERIFY(y && y_scale && f_scale, "error in N_VNewEmpty_Parallel()");
1042 
1043  }
1044 
1045  // Default abs_tol and print_level
1046  abs_tol = pow(UNIT_ROUNDOFF, 1.0/3.0);
1047  print_level = 0;
1048 }
1049 #endif
1050 
1051 
1053 {
1054  // Initialize the base class
1056  jacobian = NULL;
1057 
1058  // Get the vector length
1059  long local_size = height;
1060 #ifdef MFEM_USE_MPI
1061  long global_size;
1062 #endif
1063 
1064  if (Parallel())
1065  {
1066 #ifdef MFEM_USE_MPI
1067  MPI_Allreduce(&local_size, &global_size, 1, MPI_LONG, MPI_SUM,
1068  NV_COMM_P(y));
1069 #endif
1070  }
1071 
1072  if (sundials_mem)
1073  {
1074  // Check if the problem size has changed since the last SetOperator call
1075  int resize = 0;
1076  if (!Parallel())
1077  {
1078  resize = (NV_LENGTH_S(y) != local_size);
1079  }
1080  else
1081  {
1082 #ifdef MFEM_USE_MPI
1083  int l_resize = (NV_LOCLENGTH_P(y) != local_size) ||
1084  (saved_global_size != global_size);
1085  MPI_Allreduce(&l_resize, &resize, 1, MPI_INT, MPI_LOR, NV_COMM_P(y));
1086 #endif
1087  }
1088 
1089  // Free existing solver memory and re-create with new vector size
1090  if (resize)
1091  {
1092  KINFree(&sundials_mem);
1093  sundials_mem = NULL;
1094  }
1095  }
1096 
1097  if (!sundials_mem)
1098  {
1099  // Set actual size and data in the N_Vector y.
1100  if (!Parallel())
1101  {
1102  NV_LENGTH_S(y) = local_size;
1103  NV_DATA_S(y) = new double[local_size](); // value-initialize
1104  NV_LENGTH_S(y_scale) = local_size;
1105  NV_DATA_S(y_scale) = NULL;
1106  NV_LENGTH_S(f_scale) = local_size;
1107  NV_DATA_S(f_scale) = NULL;
1108  }
1109  else
1110  {
1111 #ifdef MFEM_USE_MPI
1112  NV_LOCLENGTH_P(y) = local_size;
1113  NV_GLOBLENGTH_P(y) = global_size;
1114  NV_DATA_P(y) = new double[local_size](); // value-initialize
1115  NV_LOCLENGTH_P(y_scale) = local_size;
1116  NV_GLOBLENGTH_P(y_scale) = global_size;
1117  NV_DATA_P(y_scale) = NULL;
1118  NV_LOCLENGTH_P(f_scale) = local_size;
1119  NV_GLOBLENGTH_P(f_scale) = global_size;
1120  NV_DATA_P(f_scale) = NULL;
1121  saved_global_size = global_size;
1122 #endif
1123  }
1124 
1125  // Create the solver memory
1126  sundials_mem = KINCreate();
1127  MFEM_VERIFY(sundials_mem, "Error in KINCreate().");
1128 
1129  // Set number of acceleration vectors
1130  if (maa > 0)
1131  {
1132  flag = KINSetMAA(sundials_mem, maa);
1133  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetMAA()");
1134  }
1135 
1136  // Initialize KINSOL
1137  flag = KINInit(sundials_mem, KINSolver::Mult, y);
1138  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINInit()");
1139 
1140  // Attach the KINSolver as user-defined data
1141  flag = KINSetUserData(sundials_mem, this);
1142  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetUserData()");
1143 
1144  // Set the linear solver
1145  if (prec)
1146  {
1148  }
1149  else
1150  {
1151  // Free any existing linear solver
1152  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1153  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
1154 
1155  LSA = SUNLinSol_SPGMR(y, PREC_NONE, 0);
1156  MFEM_VERIFY(LSA, "error in SUNLinSol_SPGMR()");
1157 
1158  flag = KINSetLinearSolver(sundials_mem, LSA, NULL);
1159  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetLinearSolver()");
1160 
1161  // Set Jacobian-vector product function
1162  if (use_oper_grad)
1163  {
1164  flag = KINSetJacTimesVecFn(sundials_mem, KINSolver::GradientMult);
1165  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetJacTimesVecFn()");
1166  }
1167  }
1168 
1169  // Delete the allocated data in y.
1170  if (!Parallel())
1171  {
1172  delete [] NV_DATA_S(y);
1173  NV_DATA_S(y) = NULL;
1174  }
1175  else
1176  {
1177 #ifdef MFEM_USE_MPI
1178  delete [] NV_DATA_P(y);
1179  NV_DATA_P(y) = NULL;
1180 #endif
1181  }
1182  }
1183 }
1184 
1186 {
1187  // Store the solver
1188  prec = &solver;
1189 
1190  // Free any existing linear solver
1191  if (A != NULL) { SUNMatDestroy(A); A = NULL; }
1192  if (LSA != NULL) { SUNLinSolFree(LSA); LSA = NULL; }
1193 
1194  // Wrap KINSolver as SUNLinearSolver and SUNMatrix
1195  LSA = SUNLinSolNewEmpty();
1196  MFEM_VERIFY(LSA, "error in SUNLinSolNewEmpty()");
1197 
1198  LSA->content = this;
1199  LSA->ops->gettype = LSGetType;
1200  LSA->ops->solve = KINSolver::LinSysSolve;
1201  LSA->ops->free = LSFree;
1202 
1203  A = SUNMatNewEmpty();
1204  MFEM_VERIFY(A, "error in SUNMatNewEmpty()");
1205 
1206  A->content = this;
1207  A->ops->getid = MatGetID;
1208  A->ops->destroy = MatDestroy;
1209 
1210  // Attach the linear solver and matrix
1211  flag = KINSetLinearSolver(sundials_mem, LSA, A);
1212  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINSetLinearSolver()");
1213 
1214  // Set the Jacobian evaluation function
1215  flag = KINSetJacFn(sundials_mem, KINSolver::LinSysSetup);
1216  MFEM_VERIFY(flag == KIN_SUCCESS, "error in KINSetJacFn()");
1217 }
1218 
1220 {
1221  flag = KINSetScaledStepTol(sundials_mem, sstol);
1222  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetScaledStepTol()");
1223 }
1224 
1225 void KINSolver::SetMaxSetupCalls(int max_calls)
1226 {
1227  flag = KINSetMaxSetupCalls(sundials_mem, max_calls);
1228  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetMaxSetupCalls()");
1229 }
1230 
1231 void KINSolver::SetMAA(int m_aa)
1232 {
1233  // Store internally as maa must be set before calling KINInit() to
1234  // set the maximum acceleration space size.
1235  maa = m_aa;
1236  if (sundials_mem)
1237  {
1238  flag = KINSetMAA(sundials_mem, maa);
1239  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetMAA()");
1240  }
1241 }
1242 
1243 // Compute the scaling vectors and solve nonlinear system
1244 void KINSolver::Mult(const Vector &b, Vector &x) const
1245 {
1246  // residual norm tolerance
1247  double tol;
1248 
1249  // Uses c = 1, corresponding to x_scale.
1250  c = 1.0;
1251 
1252  if (!iterative_mode) { x = 0.0; }
1253 
1254  // For relative tolerance, r = 1 / |residual(x)|, corresponding to fx_scale.
1255  if (rel_tol > 0.0)
1256  {
1257 
1258  oper->Mult(x, r);
1259 
1260  // Note that KINSOL uses infinity norms.
1261  double norm = r.Normlinf();
1262 #ifdef MFEM_USE_MPI
1263  if (Parallel())
1264  {
1265  double lnorm = norm;
1266  MPI_Allreduce(&lnorm, &norm, 1, MPI_DOUBLE, MPI_MAX, NV_COMM_P(y));
1267  }
1268 #endif
1269  if (abs_tol > rel_tol * norm)
1270  {
1271  r = 1.0;
1272  tol = abs_tol;
1273  }
1274  else
1275  {
1276  r = 1.0 / norm;
1277  tol = rel_tol;
1278  }
1279  }
1280  else
1281  {
1282  r = 1.0;
1283  tol = abs_tol;
1284  }
1285 
1286  // Set the residual norm tolerance
1287  flag = KINSetFuncNormTol(sundials_mem, tol);
1288  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINSetFuncNormTol()");
1289 
1290  // Solve the nonlinear system by calling the other Mult method
1291  KINSolver::Mult(x, c, r);
1292 }
1293 
1294 // Solve the nonlinear system using the provided scaling vectors
1296  const Vector &x_scale, const Vector &fx_scale) const
1297 {
1298  flag = KINSetNumMaxIters(sundials_mem, max_iter);
1299  MFEM_ASSERT(flag == KIN_SUCCESS, "KINSetNumMaxIters() failed!");
1300 
1301  if (!Parallel())
1302  {
1303 
1304  NV_DATA_S(y) = x.GetData();
1305  MFEM_VERIFY(NV_LENGTH_S(y) == x.Size(), "");
1306  NV_DATA_S(y_scale) = x_scale.GetData();
1307  NV_DATA_S(f_scale) = fx_scale.GetData();
1308 
1309  flag = KINSetPrintLevel(sundials_mem, print_level);
1310  MFEM_VERIFY(flag == KIN_SUCCESS, "KINSetPrintLevel() failed!");
1311  }
1312  else
1313  {
1314 
1315 #ifdef MFEM_USE_MPI
1316  NV_DATA_P(y) = x.GetData();
1317  MFEM_VERIFY(NV_LOCLENGTH_P(y) == x.Size(), "");
1318  NV_DATA_P(y_scale) = x_scale.GetData();
1319  NV_DATA_P(f_scale) = fx_scale.GetData();
1320 
1321  int rank;
1322  MPI_Comm_rank(NV_COMM_P(y), &rank);
1323  if (rank == 0)
1324  {
1325  flag = KINSetPrintLevel(sundials_mem, print_level);
1326  MFEM_VERIFY(flag == KIN_SUCCESS, "KINSetPrintLevel() failed!");
1327  }
1328 #endif
1329 
1330  }
1331 
1332  if (!iterative_mode) { x = 0.0; }
1333 
1334  // Solve the nonlinear system
1336  converged = (flag >= 0);
1337 
1338  // Get number of nonlinear iterations
1339  long int tmp_nni;
1340  flag = KINGetNumNonlinSolvIters(sundials_mem, &tmp_nni);
1341  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINGetNumNonlinSolvIters()");
1342  final_iter = (int) tmp_nni;
1343 
1344  // Get the residual norm
1345  flag = KINGetFuncNorm(sundials_mem, &final_norm);
1346  MFEM_ASSERT(flag == KIN_SUCCESS, "error in KINGetFuncNorm()");
1347 }
1348 
1350 {
1351  N_VDestroy(y);
1352  N_VDestroy(y_scale);
1353  N_VDestroy(f_scale);
1354  SUNMatDestroy(A);
1355  SUNLinSolFree(LSA);
1356  KINFree(&sundials_mem);
1357 }
1358 
1359 } // namespace mfem
1360 
1361 #endif
virtual double GetTime() const
Read the currently set time.
Definition: operator.hpp:300
void Init(TimeDependentOperator &f_)
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition: sundials.cpp:151
const Operator * oper
Definition: solvers.hpp:43
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:342
void SetScaledStepTol(double sstol)
Set KINSOL&#39;s scaled step tolerance.
Definition: sundials.cpp:1219
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Definition: sundials.cpp:824
Base abstract class for first order time dependent operators.
Definition: operator.hpp:259
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:89
virtual void Step(Vector &x, double &t, double &dt)
Integrate the ODE with CVODE using the specified step mode.
Definition: sundials.cpp:256
virtual void Init(TimeDependentOperator &f)
Associate a TimeDependentOperator with the ODE solver.
Definition: ode.cpp:18
virtual void SetTime(const double _t)
Set the current time.
Definition: operator.hpp:303
virtual ~CVODESolver()
Destroy the associated CVODE memory and SUNDIALS objects.
Definition: sundials.cpp:405
bool use_oper_grad
use the Jv prod function
Definition: sundials.hpp:374
int Size() const
Returns the size of the vector.
Definition: vector.hpp:157
void SetIRKTableNum(int table_num)
Choose a specific Butcher table for a diagonally implicit RK method.
Definition: sundials.cpp:853
void SetIMEXTableNum(int etable_num, int itable_num)
Choose a specific Butcher table for an IMEX RK method.
Definition: sundials.cpp:859
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:455
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:835
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
Definition: sundials.cpp:841
virtual void Mult(const Vector &x, Vector &y) const =0
Operator application: y=A(x).
void * sundials_mem
SUNDIALS mem structure.
Definition: sundials.hpp:48
double Normlinf() const
Returns the l_infinity norm of the vector.
Definition: vector.cpp:729
long saved_global_size
Global vector length on last initialization.
Definition: sundials.hpp:51
bool iterative_mode
If true, use the second argument of Mult() as an initial guess.
Definition: operator.hpp:504
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:489
double * GetData() const
Return a pointer to the beginning of the Vector data.
Definition: vector.hpp:166
Explicit RK method.
Definition: sundials.hpp:203
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
Definition: sundials.cpp:1225
ARKStepSolver(Type type=EXPLICIT)
Construct a serial wrapper to SUNDIALS&#39; ARKode integrator.
Definition: sundials.cpp:532
Implicit RK method.
Definition: sundials.hpp:204
Interface to ARKode&#39;s ARKStep module – additive Runge-Kutta methods.
Definition: sundials.hpp:197
int global_strategy
KINSOL solution strategy.
Definition: sundials.hpp:373
SUNLinearSolver LSA
Linear solver for A.
Definition: sundials.hpp:56
static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
Compute the matrix-vector product .
Definition: sundials.cpp:510
virtual void SetEvalMode(const EvalMode new_eval_mode)
Set the evaluation mode of the time-dependent operator.
Definition: operator.hpp:326
SUNMatrix M
Mass matrix M.
Definition: sundials.hpp:55
virtual void SetSolver(Solver &solver)
Set the linear solver for inverting the Jacobian.
Definition: sundials.cpp:1185
int Height() const
Get the height (size of output) of the Operator. Synonym with NumRows().
Definition: operator.hpp:57
bool Parallel() const
Definition: sundials.hpp:61
int lmm_type
Linear multistep method type.
Definition: sundials.hpp:95
Interface to the CVODE library – linear multi-step methods.
Definition: sundials.hpp:92
double b
Definition: lissajous.cpp:42
Interface to the KINSOL library – nonlinear solver methods.
Definition: sundials.hpp:370
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition: sundials.cpp:348
Type
Types of ARKODE solvers.
Definition: sundials.hpp:201
Type rk_type
Runge-Kutta type.
Definition: sundials.hpp:209
virtual ~ARKStepSolver()
Destroy the associated ARKode memory and SUNDIALS objects.
Definition: sundials.cpp:923
bool use_implicit
True for implicit or imex integration.
Definition: sundials.hpp:211
SUNNonlinearSolver NLS
Nonlinear solver.
Definition: sundials.hpp:58
CVODESolver(int lmm)
Construct a serial wrapper to SUNDIALS&#39; CVODE integrator.
Definition: sundials.cpp:120
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Wrapper to compute the nonlinear residual .
Definition: sundials.cpp:937
int step_mode
ARKStep step mode (ARK_NORMAL or ARK_ONE_STEP).
Definition: sundials.hpp:210
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:974
SUNLinearSolver LSM
Linear solver for M.
Definition: sundials.hpp:57
static int RHS1(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Definition: sundials.cpp:418
N_Vector f_scale
scaling vectors
Definition: sundials.hpp:375
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to CVODE.
Definition: sundials.cpp:291
virtual void SetOperator(const Operator &op)
Set the nonlinear Operator of the system and initialize KINSOL.
Definition: sundials.cpp:1052
virtual ~KINSolver()
Destroy the associated KINSOL memory.
Definition: sundials.cpp:1349
void SetFixedStep(double dt)
Use a fixed time step size (disable temporal adaptivity).
Definition: sundials.cpp:865
bool reinit
Flag to signal memory reinitialization is need.
Definition: sundials.hpp:50
int maa
number of acceleration vectors
Definition: sundials.hpp:377
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:1231
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:94
void PrintInfo() const
Print various ARKStep statistics.
Definition: sundials.cpp:871
static int RHS2(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Definition: sundials.cpp:438
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:951
const Operator * jacobian
stores oper-&gt;GetGradient()
Definition: sundials.hpp:376
static constexpr double default_rel_tol
Default scalar relative tolerance.
Definition: sundials.hpp:68
N_Vector y_scale
Definition: sundials.hpp:375
int step_mode
CVODE step mode (CV_NORMAL or CV_ONE_STEP).
Definition: sundials.hpp:96
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:520
KINSolver(int strategy, bool oper_grad=true)
Construct a serial wrapper to SUNDIALS&#39; KINSOL nonlinear solver.
Definition: sundials.cpp:1004
Implicit-explicit ARK method.
Definition: sundials.hpp:205
static constexpr double default_abs_tol
Default scalar absolute tolerance.
Definition: sundials.hpp:70
int flag
Last flag returned from a call to SUNDIALS.
Definition: sundials.hpp:49
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to ARKode.
Definition: sundials.cpp:726
Vector data type.
Definition: vector.hpp:48
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:474
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:109
void PrintInfo() const
Print various CVODE statistics.
Definition: sundials.cpp:360
TimeDependentOperator * f
Pointer to the associated TimeDependentOperator.
Definition: ode.hpp:26
void SetMaxOrder(int max_order)
Set the maximum method order.
Definition: sundials.cpp:354
void UseMFEMMassLinearSolver(int tdep)
Attach mass matrix linear system setup, solve, and matrix-vector product methods from the TimeDepende...
Definition: sundials.cpp:772
Base class for solvers.
Definition: operator.hpp:500
N_Vector y
State vector.
Definition: sundials.hpp:53
void SetERKTableNum(int table_num)
Choose a specific Butcher table for an explicit RK method.
Definition: sundials.cpp:847
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:565
static int RHS(realtype t, const N_Vector y, N_Vector ydot, void *user_data)
Wrapper to compute the ODE rhs function.
Definition: sundials.cpp:78
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition: sundials.cpp:829
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
Definition: sundials.cpp:757
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:991
static int MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x, N_Vector b, realtype tol)
Solve the linear system .
Definition: sundials.cpp:499
virtual void Step(Vector &x, double &t, double &dt)
Integrate the ODE with ARKode using the specified step mode.
Definition: sundials.cpp:679
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition: sundials.cpp:322
virtual void SetOperator(const Operator &op)
Also calls SetOperator for the preconditioner if there is one.
Definition: solvers.cpp:1328
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
Definition: sundials.cpp:337
void UseSundialsMassLinearSolver(int tdep)
Attach the SUNDIALS GMRES linear solver and the mass matrix matrix-vector product method from the Tim...
Definition: sundials.cpp:804