MFEM v4.9.0
Finite element discretization library
Loading...
Searching...
No Matches
sundials.hpp
Go to the documentation of this file.
1// Copyright (c) 2010-2025, 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#ifndef MFEM_SUNDIALS
13#define MFEM_SUNDIALS
14
15#include "../config/config.hpp"
16
17#ifdef MFEM_USE_SUNDIALS
18
19#ifdef MFEM_USE_MPI
20#include <mpi.h>
21#include "hypre.hpp"
22#endif
23
24#include "ode.hpp"
25#include "solvers.hpp"
26
27#include <sundials/sundials_config.h>
28// Check for appropriate SUNDIALS version
29#if !defined(SUNDIALS_VERSION_MAJOR) || (SUNDIALS_VERSION_MAJOR < 5)
30#error MFEM requires SUNDIALS version 5.0.0 or newer!
31#endif
32#if defined(MFEM_USE_CUDA) && ((SUNDIALS_VERSION_MAJOR == 5) && (SUNDIALS_VERSION_MINOR < 4))
33#error MFEM requires SUNDIALS version 5.4.0 or newer when MFEM_USE_CUDA=TRUE!
34#endif
35#if defined(MFEM_USE_HIP) && ((SUNDIALS_VERSION_MAJOR == 5) && (SUNDIALS_VERSION_MINOR < 7))
36#error MFEM requires SUNDIALS version 5.7.0 or newer when MFEM_USE_HIP=TRUE!
37#endif
38#if defined(MFEM_USE_CUDA) && !defined(SUNDIALS_NVECTOR_CUDA)
39#error MFEM_USE_CUDA=TRUE requires SUNDIALS to be built with CUDA support
40#endif
41#if defined(MFEM_USE_HIP) && !defined(SUNDIALS_NVECTOR_HIP)
42#error MFEM_USE_HIP=TRUE requires SUNDIALS to be built with HIP support
43#endif
44#include <sundials/sundials_matrix.h>
45#include <sundials/sundials_linearsolver.h>
46#include <arkode/arkode_arkstep.h>
47#include <cvodes/cvodes.h>
48#include <kinsol/kinsol.h>
49#if defined(MFEM_USE_CUDA)
50#include <sunmemory/sunmemory_cuda.h>
51#elif defined(MFEM_USE_HIP)
52#include <sunmemory/sunmemory_hip.h>
53#endif
54
55#include <functional>
56
57#define MFEM_SUNDIALS_VERSION \
58 (SUNDIALS_VERSION_MAJOR*10000 + SUNDIALS_VERSION_MINOR*100 + \
59 SUNDIALS_VERSION_PATCH)
60
61#if (SUNDIALS_VERSION_MAJOR < 6)
62
63/// (DEPRECATED) Map SUNDIALS version >= 6 datatypes and constants to
64/// version < 6 for backwards compatibility
69constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8 = FEHLBERG_13_7_8;
70
71/// (DEPRECATED) There is no SUNContext in SUNDIALS version < 6 so set it to
72/// arbitrary type for more compact backwards compatibility
73using SUNContext = void*;
74
75/// 'sunrealtype' was first introduced in v6.0.0
76typedef realtype sunrealtype;
77/// 'sunbooleantype' was first introduced in v6.0.0
78typedef booleantype sunbooleantype;
79
80/// New constant names introduced in v6.0.0
82
83// KIN_ORTH_MGS was introduced in SUNDIALS v6; here, we define it just so that
84// it can be used as the default option in the second parameter of
85// KINSolver::EnableAndersonAcc -- the actual value of the parameter will be
86// ignored when using SUNDIALS < v6.
87#define KIN_ORTH_MGS 0
88
89#endif // #if SUNDIALS_VERSION_MAJOR < 6
90
91#if (SUNDIALS_VERSION_MAJOR < 7)
92
93/** @brief The enum constant SUN_SUCCESS was added in v7 as a replacement of
94 various *_SUCCESS macros that were removed in v7. */
95enum { SUN_SUCCESS = 0 };
96
97#endif // #if SUNDIALS_VERSION_MAJOR < 7
98
99
100namespace mfem
101{
102
103#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
104
105// ---------------------------------------------------------------------------
106// SUNMemory interface class (used when CUDA or HIP is enabled)
107// ---------------------------------------------------------------------------
108class SundialsMemHelper
109{
110 /// The actual SUNDIALS object
111 SUNMemoryHelper h;
112
113public:
114
115 /// Default constructor -- object must be moved to
116 SundialsMemHelper() = default;
117
118 /// Require a SUNContext as an argument (rather than calling Sundials::GetContext)
119 /// to avoid undefined behavior during the construction of the Sundials singleton.
121
122 /// Implement move assignment
124
125 /// Disable copy construction
126 SundialsMemHelper(const SundialsMemHelper& that_helper) = delete;
127
128 ~SundialsMemHelper() { if (h) { SUNMemoryHelper_Destroy(h); } }
129
130 /// Disable copy assignment
132
133 /// Implement move assignment
135
136 /// Typecasting to SUNDIALS' SUNMemoryHelper type
137 operator SUNMemoryHelper() const { return h; }
138
139 static int SundialsMemHelper_Alloc(SUNMemoryHelper helper, SUNMemory* memptr,
140 size_t memsize, SUNMemoryType mem_type
141#if (SUNDIALS_VERSION_MAJOR >= 6)
142 , void* queue
143#endif
144 );
145
146 static int SundialsMemHelper_Dealloc(SUNMemoryHelper helper, SUNMemory sunmem
147#if (SUNDIALS_VERSION_MAJOR >= 6)
148 , void* queue
149#endif
150 );
151
152};
153
154#else // MFEM_USE_CUDA || MFEM_USE_HIP
155
156// ---------------------------------------------------------------------------
157// Dummy SUNMemory interface class (used when CUDA or HIP is not enabled)
158// ---------------------------------------------------------------------------
160{
161public:
162
163 SundialsMemHelper() = default;
164
166 {
167 // Do nothing
168 }
169};
170
171#endif // MFEM_USE_CUDA || MFEM_USE_HIP
172
173
174/// Singleton class for SUNContext and SundialsMemHelper objects
176{
177public:
178
179 /// Disable copy construction
180 Sundials(Sundials &other) = delete;
181
182 /// Disable copy assignment
183 void operator=(const Sundials &other) = delete;
184
185 /// @brief Initializes SUNContext and SundialsMemHelper objects. Should be called at
186 /// the beginning of the calling program (after Mpi::Init if applicable)
187 ///
188 /// Multiple calls to Sundials::Init() have no effect.
189 static void Init();
190
191 /// @brief Finalize sundials (called automatically at program exit if
192 /// Sundials::Init() has been called).
193 ///
194 /// Multiple calls to Sundials::Finalize() have no effect. This function can be
195 /// called manually to more precisely control when sundials is finalized.
196 static void Finalize();
197
198 /// Provides access to the SUNContext object
199 static SUNContext &GetContext();
200
201 /// Provides access to the SundialsMemHelper object
203
204private:
205 /// Returns a reference to the singleton instance of the class.
206 static Sundials &Instance();
207
208 /// Constructor called by Sundials::Instance (does nothing for version < 6)
209 Sundials();
210
211 /// Destructor called at end of calling program (does nothing for version < 6)
212 ~Sundials();
213
214 bool isInitialized;
215 SUNContext context;
216 SundialsMemHelper memHelper;
217};
218
219
220/// Vector interface for SUNDIALS N_Vectors.
222{
223protected:
225
226 /// The actual SUNDIALS object
227 N_Vector x;
228
229 friend class SundialsSolver;
230
231 /// Set data and length of internal N_Vector x from 'this'.
232 void _SetNvecDataAndSize_(long glob_size = 0);
233
234 /// Set data and length from the internal N_Vector x.
235 void _SetDataAndSize_();
236
237public:
238 /// Creates an empty SundialsNVector.
240
241 /// Creates a SundialsNVector referencing an array of doubles, owned by someone else.
242 /** The pointer @a data_ can be NULL. The data array can be replaced later
243 with SetData(). */
244 SundialsNVector(double *data_, int size_);
245
246 /// Creates a SundialsNVector out of a SUNDIALS N_Vector object.
247 /** The N_Vector @a nv must be destroyed outside. */
248 SundialsNVector(N_Vector nv);
249
250#ifdef MFEM_USE_MPI
251 /// Creates an empty SundialsNVector.
252 SundialsNVector(MPI_Comm comm);
253
254 /// Creates a SundialsNVector with the given local and global sizes.
255 SundialsNVector(MPI_Comm comm, int loc_size, long glob_size);
256
257 /// Creates a SundialsNVector referencing an array of doubles, owned by someone else.
258 /** The pointer @a data_ can be NULL. The data array can be replaced later
259 with SetData(). */
260 SundialsNVector(MPI_Comm comm, double *data_, int loc_size, long glob_size);
261
262 /// Creates a SundialsNVector from a HypreParVector.
263 /** Ownership of the data will not change. */
265#endif
266
267 /// Calls SUNDIALS N_VDestroy function if the N_Vector is owned by 'this'.
269
270 /// Returns the N_Vector_ID for the internal N_Vector.
271 inline N_Vector_ID GetNVectorID() const { return N_VGetVectorID(x); }
272
273 /// Returns the N_Vector_ID for the N_Vector @a x_.
274 inline N_Vector_ID GetNVectorID(N_Vector x_) const { return N_VGetVectorID(x_); }
275
276#ifdef MFEM_USE_MPI
277 /// Returns the MPI communicator for the internal N_Vector x.
278 inline MPI_Comm GetComm() const
279 {
280#if SUNDIALS_VERSION_MAJOR < 7
281 return *static_cast<MPI_Comm*>(N_VGetCommunicator(x));
282#else
283 return N_VGetCommunicator(x);
284#endif
285 }
286
287 /// Returns the MPI global length for the internal N_Vector x.
288 inline long GlobalSize() const { return N_VGetLength(x); }
289#endif
290
291 /// Resize the vector to size @a s.
292 void SetSize(int s, long glob_size = 0);
293
294 /// Set the vector data.
295 /// @warning This method should be called only when OwnsData() is false.
296 void SetData(double *d);
297
298 /// Set the vector data and size.
299 /** The Vector does not assume ownership of the new data. The new size is
300 also used as the new Capacity().
301 @warning This method should be called only when OwnsData() is false. */
302 void SetDataAndSize(double *d, int s, long glob_size = 0);
303
304 /// Reset the Vector to be a reference to a sub-vector of @a base.
305 inline void MakeRef(Vector &base, int offset, int s)
306 {
307 // Ensure that the base is registered/initialized before making an alias
308 base.Read();
309 Vector::MakeRef(base, offset, s);
311 }
312
313 /** @brief Reset the Vector to be a reference to a sub-vector of @a base
314 without changing its current size. */
315 inline void MakeRef(Vector &base, int offset)
316 {
317 // Ensure that the base is registered/initialized before making an alias
318 base.Read();
319 Vector::MakeRef(base, offset);
321 }
322
323 /// Typecasting to SUNDIALS' N_Vector type
324 operator N_Vector() const { return x; }
325
326 /// Changes the ownership of the vector
327 N_Vector StealNVector() { own_NVector = 0; return x; }
328
329 /// Sets ownership of the internal N_Vector
330 void SetOwnership(int own) { own_NVector = own; }
331
332 /// Gets ownership of the internal N_Vector
333 int GetOwnership() const { return own_NVector; }
334
335 /// Copy assignment.
336 /** @note Defining this method overwrites the implicitly defined copy
337 assignment operator. */
338 using Vector::operator=;
339
340#ifdef MFEM_USE_MPI
341 bool MPIPlusX() const
342 { return (GetNVectorID() == SUNDIALS_NVEC_MPIPLUSX); }
343#else
344 bool MPIPlusX() const { return false; }
345#endif
346
347 /// Create a N_Vector.
348 /** @param[in] use_device If true, use the SUNDIALS CUDA or HIP N_Vector. */
349 static N_Vector MakeNVector(bool use_device);
350
351#ifdef MFEM_USE_MPI
352 /// Create a parallel N_Vector.
353 /** @param[in] comm The MPI communicator to use.
354 @param[in] use_device If true, use the SUNDIALS CUDA or HIP N_Vector. */
355 static N_Vector MakeNVector(MPI_Comm comm, bool use_device);
356#endif
357
358#if defined(MFEM_USE_CUDA) || defined(MFEM_USE_HIP)
359 static bool UseManagedMemory()
360 {
362 }
363#else
364 static bool UseManagedMemory()
365 {
366 return false;
367 }
368#endif
369
370};
371
372/// Base class for interfacing with SUNDIALS packages.
374{
375protected:
376 void *sundials_mem; ///< SUNDIALS mem structure.
377 mutable int flag; ///< Last flag returned from a call to SUNDIALS.
378 bool reinit; ///< Flag to signal memory reinitialization is need.
379 long saved_global_size; ///< Global vector length on last initialization.
380
381 SundialsNVector* Y; ///< State vector.
382 SUNMatrix A; /**< Linear system A = I - gamma J,
383 M - gamma J, or J. */
384 SUNMatrix M; ///< Mass matrix M.
385 SUNLinearSolver LSA; ///< Linear solver for A.
386 SUNLinearSolver LSM; ///< Linear solver for M.
387 SUNNonlinearSolver NLS; ///< Nonlinear solver.
388
389#ifdef MFEM_USE_MPI
390 bool Parallel() const
391 { return (Y->MPIPlusX() || Y->GetNVectorID() == SUNDIALS_NVEC_PARALLEL); }
392#else
393 bool Parallel() const { return false; }
394#endif
395
396 /// Default scalar relative tolerance.
397 static constexpr double default_rel_tol = 1e-4;
398 /// Default scalar absolute tolerance.
399 static constexpr double default_abs_tol = 1e-9;
400
401 /** @brief Protected constructor: objects of this type should be constructed
402 only as part of a derived class. */
403 SundialsSolver() : sundials_mem(NULL), flag(0), reinit(false),
404 saved_global_size(0), Y(NULL), A(NULL), M(NULL),
405 LSA(NULL), LSM(NULL), NLS(NULL) { }
406
407 // Helper functions
408 // Serial version
409 void AllocateEmptyNVector(N_Vector &y);
410
411#ifdef MFEM_USE_MPI
412 void AllocateEmptyNVector(N_Vector &y, MPI_Comm comm);
413#endif
414
415public:
416 /// Access the SUNDIALS memory structure.
417 void *GetMem() const { return sundials_mem; }
418
419 /// Returns the last flag returned by a call to a SUNDIALS function.
420 int GetFlag() const { return flag; }
421};
422
423
424// ---------------------------------------------------------------------------
425// Interface to the CVODE library -- linear multi-step methods
426// ---------------------------------------------------------------------------
427
428/// Interface to the CVODE library -- linear multi-step methods.
430{
431protected:
432 int lmm_type; ///< Linear multistep method type.
433 int step_mode; ///< CVODE step mode (CV_NORMAL or CV_ONE_STEP).
434 int root_components; /// Number of components in gout
435
436 /// Wrapper to compute the ODE rhs function.
437 static int RHS(sunrealtype t, const N_Vector y, N_Vector ydot,
438 void *user_data);
439
440 /// Setup the linear system $ A x = b $.
441 static int LinSysSetup(sunrealtype t, N_Vector y, N_Vector fy, SUNMatrix A,
443 sunrealtype gamma, void *user_data, N_Vector tmp1,
444 N_Vector tmp2, N_Vector tmp3);
445
446 /// Solve the linear system $ A x = b $.
447 static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x,
448 N_Vector b, sunrealtype tol);
449
450 /// Prototype to define root finding for CVODE
451 static int root(sunrealtype t, N_Vector y, sunrealtype *gout,
452 void *user_data);
453
454 /// Typedef for root finding functions
455 typedef std::function<int(sunrealtype t, Vector y, Vector gout,
457
458 /// A class member to facilitate pointing to a user-specified root function
460
461 /// Typedef declaration for error weight functions
462 typedef std::function<int(Vector y, Vector w, CVODESolver*)> EWTFunction;
463
464 /** @brief A class member to facilitate pointing to a user-specified error
465 weight function */
467
468public:
469 /// Construct a serial wrapper to SUNDIALS' CVODE integrator.
470 /** @param[in] lmm Specifies the linear multistep method, the options are:
471 - CV_ADAMS - implicit methods for non-stiff systems,
472 - CV_BDF - implicit methods for stiff systems. */
473 CVODESolver(int lmm);
474
475#ifdef MFEM_USE_MPI
476 /// Construct a parallel wrapper to SUNDIALS' CVODE integrator.
477 /** @param[in] comm The MPI communicator used to partition the ODE system
478 @param[in] lmm Specifies the linear multistep method, the options are:
479 - CV_ADAMS - implicit methods for non-stiff systems,
480 - CV_BDF - implicit methods for stiff systems. */
481 CVODESolver(MPI_Comm comm, int lmm);
482#endif
483
484 /** @brief Initialize CVODE: calls CVodeCreate() to create the CVODE
485 memory and set some defaults.
486
487 If the CVODE memory has already been created, it checks if the problem
488 size has changed since the last call to Init(). If the problem is the
489 same then CVodeReInit() will be called in the next call to Step(). If
490 the problem size has changed, the CVODE memory is freed and realloced
491 for the new problem size. */
492 /** @param[in] f_ The TimeDependentOperator that defines the ODE system.
493
494 @note All other methods must be called after Init().
495
496 @note If this method is called a second time with a different problem
497 size, then any non-default user-set options will be lost and will need
498 to be set again. */
499 void Init(TimeDependentOperator &f_) override;
500
501 /// Integrate the ODE with CVODE using the specified step mode.
502 /** @param[in,out] x On output, the solution vector at the requested output
503 time tout = @a t + @a dt.
504 @param[in,out] t On output, the output time reached.
505 @param[in,out] dt On output, the last time step taken.
506
507 @note On input, the values of @a t and @a dt are used to compute desired
508 output time for the integration, tout = @a t + @a dt.
509 */
510 void Step(Vector &x, double &t, double &dt) override;
511
512 /** @brief Attach the linear system setup and solve methods from the
513 TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to
514 CVODE.
515 */
516 void UseMFEMLinearSolver();
517
518 /// Attach SUNDIALS GMRES linear solver to CVODE.
520
521 /// Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
522 /** @param[in] itask The desired step mode. */
523 void SetStepMode(int itask);
524
525 /// Set the scalar relative and scalar absolute tolerances.
526 void SetSStolerances(double reltol, double abstol);
527
528 /// Set the scalar relative and vector of absolute tolerances.
529 void SetSVtolerances(double reltol, Vector abstol);
530
531 /// Initialize Root Finder.
532 void SetRootFinder(int components, RootFunction func);
533
534 /// Set the maximum time step.
535 void SetMaxStep(double dt_max);
536
537 /// Set the maximum number of time steps.
538 void SetMaxNSteps(int steps);
539
540 /// Get the number of internal steps taken so far.
541 long GetNumSteps();
542
543 /** @brief Set the maximum method order.
544
545 CVODE uses adaptive-order integration, based on the local truncation
546 error. The default values for @a max_order are 12 for CV_ADAMS and
547 5 for CV_BDF. Use this if you know a priori that your system is such
548 that higher order integration formulas are unstable.
549
550 @note @a max_order can't be higher than the current maximum order. */
551 void SetMaxOrder(int max_order);
552
553 /// Print various CVODE statistics.
554 void PrintInfo() const;
555
556 /// Destroy the associated CVODE memory and SUNDIALS objects.
557 virtual ~CVODESolver();
558
559};
560
561// ---------------------------------------------------------------------------
562// Interface to the CVODES library -- linear multi-step methods
563// ---------------------------------------------------------------------------
564
566{
567private:
568 using CVODESolver::Init;
569
570protected:
571 int ncheck; ///< number of checkpoints used so far
572 int indexB; ///< backward problem index
573
574 /// Wrapper to compute the ODE RHS Quadrature function.
575 static int RHSQ(sunrealtype t, const N_Vector y, N_Vector qdot,
576 void *user_data);
577
578 /// Wrapper to compute the ODE RHS backward function.
579 static int RHSB(sunrealtype t, N_Vector y,
580 N_Vector yB, N_Vector yBdot, void *user_dataB);
581
582 /// Wrapper to compute the ODE RHS Backwards Quadrature function.
583 static int RHSQB(sunrealtype t, N_Vector y, N_Vector yB,
584 N_Vector qBdot, void *user_dataB);
585
586 /// Error control function
587 static int ewt(N_Vector y, N_Vector w, void *user_data);
588
589 SUNMatrix AB; ///< Linear system A = I - gamma J, M - gamma J, or J.
590 SUNLinearSolver LSB; ///< Linear solver for A.
591 SundialsNVector* q; ///< Quadrature vector.
592 SundialsNVector* yB; ///< State vector.
593 SundialsNVector* yy; ///< State vector.
594 SundialsNVector* qB; ///< State vector.
595
596 /// Default scalar backward relative tolerance
597 static constexpr double default_rel_tolB = 1e-4;
598 /// Default scalar backward absolute tolerance
599 static constexpr double default_abs_tolB = 1e-9;
600 /// Default scalar backward absolute quadrature tolerance
601 static constexpr double default_abs_tolQB = 1e-9;
602
603public:
604 /** Construct a serial wrapper to SUNDIALS' CVODE integrator.
605 @param[in] lmm Specifies the linear multistep method, the options are:
606 CV_ADAMS - implicit methods for non-stiff systems
607 CV_BDF - implicit methods for stiff systems */
608 CVODESSolver(int lmm);
609
610#ifdef MFEM_USE_MPI
611 /** Construct a parallel wrapper to SUNDIALS' CVODE integrator.
612 @param[in] comm The MPI communicator used to partition the ODE system
613 @param[in] lmm Specifies the linear multistep method, the options are:
614 CV_ADAMS - implicit methods for non-stiff systems
615 CV_BDF - implicit methods for stiff systems */
616 CVODESSolver(MPI_Comm comm, int lmm);
617#endif
618
619 /** Initialize CVODE: Calls CVodeInit() and sets some defaults. We define this
620 to force the time dependent operator to be a TimeDependenAdjointOperator.
621 @param[in] f_ the TimeDependentAdjointOperator that defines the ODE system
622
623 @note All other methods must be called after Init(). */
625
626 /// Initialize the adjoint problem
628
629 /** Integrate the ODE with CVODE using the specified step mode.
630
631 @param[out] x Solution vector at the requested output time x=x(t).
632 @param[in,out] t On output, the output time reached.
633 @param[in,out] dt On output, the last time step taken.
634
635 @note On input, the values of t and dt are used to compute desired
636 output time for the integration, tout = t + dt. */
637 void Step(Vector &x, double &t, double &dt) override;
638
639 /// Solve one adjoint time step
640 virtual void StepB(Vector &w, double &t, double &dt);
641
642 /// Set multiplicative error weights
643 void SetWFTolerances(EWTFunction func);
644
645 // Initialize Quadrature Integration
647 double reltolQ = 1e-3,
648 double abstolQ = 1e-8);
649
650 /// Initialize Quadrature Integration (Adjoint)
651 void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB = 1e-3,
652 double abstolQB = 1e-8);
653
654 /// Initialize Adjoint
655 void InitAdjointSolve(int steps, int interpolation);
656
657 /// Set the maximum number of backward steps
658 void SetMaxNStepsB(int mxstepsB);
659
660 /// Get Number of Steps for ForwardSolve
662
663 /// Evaluate Quadrature
664 void EvalQuadIntegration(double t, Vector &q);
665
666 /// Evaluate Quadrature solution
667 void EvalQuadIntegrationB(double t, Vector &dG_dp);
668
669 /// Get Interpolated Forward solution y at backward integration time tB
670 void GetForwardSolution(double tB, mfem::Vector & yy);
671
672 /// Set Linear Solver for the backward problem
674
675 /// Use built in SUNDIALS Newton solver
677
678 /**
679 \brief Tolerance specification functions for the adjoint problem.
680
681 It should be called after InitB() is called.
682
683 \param[in] reltol the scalar relative error tolerance.
684 \param[in] abstol the scalar absolute error tolerance.
685 */
686 void SetSStolerancesB(double reltol, double abstol);
687
688 /**
689 \brief Tolerance specification functions for the adjoint problem.
690
691 It should be called after InitB() is called.
692
693 \param[in] reltol the scalar relative error tolerance
694 \param[in] abstol the vector of absolute error tolerances
695 */
696 void SetSVtolerancesB(double reltol, Vector abstol);
697
698 /// Setup the linear system A x = b
699 static int LinSysSetupB(sunrealtype t, N_Vector y, N_Vector yB, N_Vector fyB,
700 SUNMatrix A,
702 sunrealtype gamma, void *user_data, N_Vector tmp1,
703 N_Vector tmp2, N_Vector tmp3);
704
705 /// Solve the linear system A x = b
706 static int LinSysSolveB(SUNLinearSolver LS, SUNMatrix A, N_Vector x,
707 N_Vector b, sunrealtype tol);
708
709
710 /// Destroy the associated CVODES memory and SUNDIALS objects.
711 virtual ~CVODESSolver();
712};
713
714
715// ---------------------------------------------------------------------------
716// Interface to ARKode's ARKStep module -- Additive Runge-Kutta methods
717// ---------------------------------------------------------------------------
718
719/// Interface to ARKode's ARKStep module -- additive Runge-Kutta methods.
721{
722public:
723 /// Types of ARKODE solvers.
724 enum Type
725 {
726 EXPLICIT, ///< Explicit RK method
727 IMPLICIT, ///< Implicit RK method
728 IMEX ///< Implicit-explicit ARK method
729 };
730
731protected:
732 Type rk_type; ///< Runge-Kutta type.
733 int step_mode; ///< ARKStep step mode (ARK_NORMAL or ARK_ONE_STEP).
734 bool use_implicit; ///< True for implicit or imex integration.
735
736 /** @name Wrappers to compute the ODE RHS functions.
737 RHS1 is explicit RHS and RHS2 the implicit RHS for IMEX integration. When
738 purely implicit or explicit only RHS1 is used. */
739 ///@{
740 static int RHS1(sunrealtype t, const N_Vector y, N_Vector ydot,
741 void *user_data);
742 static int RHS2(sunrealtype t, const N_Vector y, N_Vector ydot,
743 void *user_data);
744 ///@}
745
746 /// Setup the linear system $ A x = b $.
747 static int LinSysSetup(sunrealtype t, N_Vector y, N_Vector fy, SUNMatrix A,
748 SUNMatrix M, sunbooleantype jok, sunbooleantype *jcur,
749 sunrealtype gamma, void *user_data, N_Vector tmp1,
750 N_Vector tmp2, N_Vector tmp3);
751
752 /// Solve the linear system $ A x = b $.
753 static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x,
754 N_Vector b, sunrealtype tol);
755
756 /// Setup the linear system $ M x = b $.
757 static int MassSysSetup(sunrealtype t, SUNMatrix M, void *user_data,
758 N_Vector tmp1, N_Vector tmp2, N_Vector tmp3);
759
760 /// Solve the linear system $ M x = b $.
761 static int MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x,
762 N_Vector b, sunrealtype tol);
763
764 /// Compute the matrix-vector product $ v = M x $.
765 static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v);
766
767 /// Compute the matrix-vector product $v = M_t x $ at time t.
768 static int MassMult2(N_Vector x, N_Vector v, sunrealtype t,
769 void* mtimes_data);
770
771public:
772 /// Construct a serial wrapper to SUNDIALS' ARKode integrator.
773 /** @param[in] type Specifies the RK method type:
774 - EXPLICIT - explicit RK method (default)
775 - IMPLICIT - implicit RK method
776 - IMEX - implicit-explicit ARK method */
778
779#ifdef MFEM_USE_MPI
780 /// Construct a parallel wrapper to SUNDIALS' ARKode integrator.
781 /** @param[in] comm The MPI communicator used to partition the ODE system.
782 @param[in] type Specifies the RK method type:
783 - EXPLICIT - explicit RK method (default)
784 - IMPLICIT - implicit RK method
785 - IMEX - implicit-explicit ARK method */
786 ARKStepSolver(MPI_Comm comm, Type type = EXPLICIT);
787#endif
788
789 /** @brief Initialize ARKode: calls ARKStepCreate() to create the ARKStep
790 memory and set some defaults.
791
792 If the ARKStep has already been created, it checks if the problem size
793 has changed since the last call to Init(). If the problem is the same
794 then ARKStepReInit() will be called in the next call to Step(). If the
795 problem size has changed, the ARKStep memory is freed and realloced
796 for the new problem size. */
797 /** @param[in] f_ The TimeDependentOperator that defines the ODE system
798
799 @note All other methods must be called after Init().
800
801 @note If this method is called a second time with a different problem
802 size, then any non-default user-set options will be lost and will need
803 to be set again. */
804 void Init(TimeDependentOperator &f_) override;
805
806 /// Integrate the ODE with ARKode using the specified step mode.
807 /**
808 @param[in,out] x On output, the solution vector at the requested output
809 time, tout = @a t + @a dt
810 @param[in,out] t On output, the output time reached
811 @param[in,out] dt On output, the last time step taken
812
813 @note On input, the values of @a t and @a dt are used to compute desired
814 output time for the integration, tout = @a t + @a dt.
815 */
816 void Step(Vector &x, real_t &t, real_t &dt) override;
817
818 /** @brief Attach the linear system setup and solve methods from the
819 TimeDependentOperator i.e., SUNImplicitSetup() and SUNImplicitSolve() to
820 ARKode.
821 */
822 void UseMFEMLinearSolver();
823
824 /// Attach a SUNDIALS GMRES linear solver to ARKode.
826
827 /** @brief Attach mass matrix linear system setup, solve, and matrix-vector
828 product methods from the TimeDependentOperator i.e., SUNMassSetup(),
829 SUNMassSolve(), and SUNMassMult() to ARKode.
830
831 @param[in] tdep An integer flag indicating if the mass matrix is time
832 dependent (1) or time independent (0)
833 */
834 void UseMFEMMassLinearSolver(int tdep);
835
836 /** @brief Attach the SUNDIALS GMRES linear solver and the mass matrix
837 matrix-vector product method from the TimeDependentOperator i.e.,
838 SUNMassMult() to ARKode to solve mass matrix systems.
839
840 @param[in] tdep An integer flag indicating if the mass matrix is time
841 dependent (1) or time independent (0)
842 */
843 void UseSundialsMassLinearSolver(int tdep);
844
845 /// Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
846 /** @param[in] itask The desired step mode */
847 void SetStepMode(int itask);
848
849 /// Set the scalar relative and scalar absolute tolerances.
850 void SetSStolerances(double reltol, double abstol);
851
852 /// Set the maximum time step.
853 void SetMaxStep(double dt_max);
854
855 /// Chooses integration order for all explicit / implicit / IMEX methods.
856 /** The default is 4, and the allowed ranges are: [2, 8] for explicit;
857 [2, 5] for implicit; [3, 5] for IMEX. */
858 void SetOrder(int order);
859
860 /// Choose a specific Butcher table for an explicit RK method.
861 /** See ARKODE documentation for all possible options, stability regions, etc.
862 For example, table_id = BOGACKI_SHAMPINE_4_2_3 is 4-stage 3rd order. */
863 void SetERKTableNum(ARKODE_ERKTableID table_id);
864
865 /// Choose a specific Butcher table for a diagonally implicit RK method.
866 /** See ARKODE documentation for all possible options, stability regions, etc.
867 For example, table_id = CASH_5_3_4 is 5-stage 4th order. */
868 void SetIRKTableNum(ARKODE_DIRKTableID table_id);
869
870 /// Choose a specific Butcher table for an IMEX RK method.
871 /** See ARKODE documentation for all possible options, stability regions, etc.
872 For example, etable_id = ARK548L2SA_DIRK_8_4_5 and
873 itable_id = ARK548L2SA_ERK_8_4_5 is 8-stage 5th order. */
874 void SetIMEXTableNum(ARKODE_ERKTableID etable_id, ARKODE_DIRKTableID itable_id);
875
876 /// Use a fixed time step size (disable temporal adaptivity).
877 /** Use of this function is not recommended, since there is no assurance of
878 the validity of the computed solutions. It is primarily provided for
879 code-to-code verification testing purposes. */
880 void SetFixedStep(double dt);
881
882 /// Print various ARKStep statistics.
883 void PrintInfo() const;
884
885 /// Destroy the associated ARKode memory and SUNDIALS objects.
886 virtual ~ARKStepSolver();
887
888};
889
890
891// ---------------------------------------------------------------------------
892// Interface to the KINSOL library -- nonlinear solver methods
893// ---------------------------------------------------------------------------
894
895/// Interface to the KINSOL library -- nonlinear solver methods.
897{
898protected:
899 int global_strategy; ///< KINSOL solution strategy
900 bool use_oper_grad; ///< use the Jv prod function
901 mutable SundialsNVector *y_scale, *f_scale; ///< scaling vectors
902 const Operator *jacobian; ///< stores oper->GetGradient()
903 int aa_n = 0; ///< number of acceleration vectors
904 int aa_delay; ///< Anderson Acceleration delay
905 double aa_damping; ///< Anderson Acceleration damping
906 int aa_orth; ///< Anderson Acceleration orthogonalization routine
907 double fp_damping = 1.0; ///< Fixed Point or Picard damping parameter
908 bool jfnk = false; ///< enable JFNK
909 Vector wrk; ///< Work vector needed for the JFNK PC
910 int maxli = 5; ///< Maximum linear iterations
911 int maxlrs = 0; ///< Maximum linear solver restarts
912
913 /// Wrapper to compute the nonlinear residual $ F(u) = 0 $.
914 static int Mult(const N_Vector u, N_Vector fu, void *user_data);
915
916 /// Wrapper to compute the Jacobian-vector product $ J(u) v = Jv $.
917 static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u,
918 sunbooleantype *new_u, void *user_data);
919
920 /// Setup the linear system $ J u = b $.
921 static int LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J,
922 void *user_data, N_Vector tmp1, N_Vector tmp2);
923
924 /// Solve the linear system $ J u = b $.
925 static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u,
926 N_Vector b, sunrealtype tol);
927
928 /// Setup the preconditioner.
929 static int PrecSetup(N_Vector uu,
930 N_Vector uscale,
931 N_Vector fval,
932 N_Vector fscale,
933 void *user_data);
934
935 /// Solve the preconditioner equation $ Pz = v $.
936 static int PrecSolve(N_Vector uu,
937 N_Vector uscale,
938 N_Vector fval,
939 N_Vector fscale,
940 N_Vector vv,
941 void *user_data);
942
943 void SetJFNKSolver(Solver &solver);
944
945public:
946
947 /// Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.
948 /** @param[in] strategy Specifies the nonlinear solver strategy:
949 KIN_NONE / KIN_LINESEARCH / KIN_PICARD / KIN_FP.
950 @param[in] oper_grad Specifies whether the solver should use its
951 Operator's GetGradient() method to compute the
952 Jacobian of the system. */
953 KINSolver(int strategy, bool oper_grad = true);
954
955#ifdef MFEM_USE_MPI
956 /// Construct a parallel wrapper to SUNDIALS' KINSOL nonlinear solver.
957 /** @param[in] comm The MPI communicator used to partition the system.
958 @param[in] strategy Specifies the nonlinear solver strategy:
959 KIN_NONE / KIN_LINESEARCH / KIN_PICARD / KIN_FP.
960 @param[in] oper_grad Specifies whether the solver should use its
961 Operator's GetGradient() method to compute the
962 Jacobian of the system. */
963 KINSolver(MPI_Comm comm, int strategy, bool oper_grad = true);
964#endif
965
966 /// Destroy the associated KINSOL memory.
967 virtual ~KINSolver();
968
969 /// Set the nonlinear Operator of the system and initialize KINSOL.
970 /** @note If this method is called a second time with a different problem
971 size, then non-default KINSOL-specific options will be lost and will need
972 to be set again. */
973 void SetOperator(const Operator &op) override;
974
975 /// Set the linear solver for inverting the Jacobian.
976 /** @note This function assumes that Operator::GetGradient(const Vector &)
977 is implemented by the Operator specified by
978 SetOperator(const Operator &).
979
980 This method must be called after SetOperator(). */
981 void SetSolver(Solver &solver) override;
982
983 /// Equivalent to SetSolver(solver).
984 void SetPreconditioner(Solver &solver) override { SetSolver(solver); }
985
986 /// Set KINSOL's scaled step tolerance.
987 /** The default tolerance is $ U^\frac{2}{3} $ , where
988 U = machine unit round-off.
989 @note This method must be called after SetOperator(). */
990 void SetScaledStepTol(double sstol);
991
992 /// Set maximum number of nonlinear iterations without a Jacobian update.
993 /** The default is 10.
994 @note This method must be called after SetOperator(). */
995 void SetMaxSetupCalls(int max_calls);
996
997 /// Enable Anderson Acceleration for KIN_FP or KIN_PICARD.
998 /** @note Has to be called once before SetOperator() in order to set up the
999 maximum subspace size. Subsequent calls need @a n less or equal to the
1000 initial subspace size.
1001 @param[in] n Anderson Acceleration subspace size
1002 @param[in] orth Anderson Acceleration orthogonalization routine
1003 @param[in] delay Anderson Acceleration delay
1004 @param[in] damping Anderson Acceleration damping parameter valid from 0 <
1005 d <= 1.0. Default is 1.0 (no damping) */
1006 void EnableAndersonAcc(int n, int orth = KIN_ORTH_MGS, int delay = 0,
1007 double damping = 1.0);
1008
1009 /// Specifies the value of the damping parameter in the fixed point or Picard
1010 /// iteration.
1011 /** param[in] damping fixed point iteration or Picard damping parameter */
1012 void SetDamping(double damping);
1013
1014 /// Set the Jacobian Free Newton Krylov flag. The default is false.
1015 /** This flag indicates to use JFNK as the linear solver for KINSOL. This
1016 means that the Solver object set in SetSolver() or SetPreconditioner() is
1017 used as a preconditioner for an FGMRES algorithm provided by SpFGMR from
1018 SUNDIALS. Furthermore, all Jacobian-vector products in the outer Krylov
1019 method are approximated by a difference quotient and the relative
1020 tolerance for the outer Krylov method is adaptive. See the KINSOL User
1021 Manual for details. */
1022 void SetJFNK(bool use_jfnk) { jfnk = use_jfnk; }
1023
1024 /// Set the maximum number of linear solver iterations
1025 /** @note Only valid in combination with JFNK */
1026 void SetLSMaxIter(int m) { maxli = m; }
1027
1028 /// Set the maximum number of linear solver restarts
1029 /** @note Only valid in combination with JFNK */
1030 void SetLSMaxRestarts(int m) { maxlrs = m; }
1031
1032 /// Set the print level for the KINSetPrintLevel function.
1033 void SetPrintLevel(int print_lvl) override { print_level = print_lvl; }
1034
1035 /// This method is not supported and will throw an error.
1036 void SetPrintLevel(PrintLevel) override;
1037
1038 /// Solve the nonlinear system $ F(x) = 0 $.
1039 /** This method computes the x_scale and fx_scale vectors and calls the
1040 other Mult(Vector&, Vector&, Vector&) const method. The x_scale vector
1041 is a vector of ones and values of fx_scale are determined by comparing
1042 the chosen relative and functional norm (i.e. absolute) tolerances.
1043 @param[in] b Not used, KINSOL always assumes zero RHS
1044 @param[in,out] x On input, initial guess, if @a #iterative_mode = true,
1045 otherwise the initial guess is zero; on output, the
1046 solution */
1047 void Mult(const Vector &b, Vector &x) const override;
1048
1049 /// Solve the nonlinear system $ F(x) = 0 $.
1050 /** Calls KINSol() to solve the nonlinear system. Before calling KINSol(),
1051 this functions uses the data members inherited from class IterativeSolver
1052 to set corresponding KINSOL options.
1053 @param[in,out] x On input, initial guess, if @a #iterative_mode =
1054 true, otherwise the initial guess is zero; on
1055 output, the solution
1056 @param[in] x_scale Elements of a diagonal scaling matrix D, s.t.
1057 D*x has all elements roughly the same when
1058 x is close to a solution
1059 @param[in] fx_scale Elements of a diagonal scaling matrix E, s.t.
1060 D*F(x) has all elements roughly the same when
1061 x is not too close to a solution */
1062 void Mult(Vector &x, const Vector &x_scale, const Vector &fx_scale) const;
1063};
1064
1065} // namespace mfem
1066
1067#endif // MFEM_USE_SUNDIALS
1068
1069#endif // MFEM_SUNDIALS
Interface to ARKode's ARKStep module – additive Runge-Kutta methods.
Definition sundials.hpp:721
static int RHS2(sunrealtype t, const N_Vector y, N_Vector ydot, void *user_data)
void SetMaxStep(double dt_max)
Set the maximum time step.
void PrintInfo() const
Print various ARKStep statistics.
Type rk_type
Runge-Kutta type.
Definition sundials.hpp:732
ARKStepSolver(Type type=EXPLICIT)
Construct a serial wrapper to SUNDIALS' ARKode integrator.
void SetOrder(int order)
Chooses integration order for all explicit / implicit / IMEX methods.
void SetStepMode(int itask)
Select the ARKode step mode: ARK_NORMAL (default) or ARK_ONE_STEP.
Type
Types of ARKODE solvers.
Definition sundials.hpp:725
@ IMPLICIT
Implicit RK method.
Definition sundials.hpp:727
@ IMEX
Implicit-explicit ARK method.
Definition sundials.hpp:728
@ EXPLICIT
Explicit RK method.
Definition sundials.hpp:726
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system .
static int LinSysSetup(sunrealtype t, N_Vector y, N_Vector fy, SUNMatrix A, SUNMatrix M, sunbooleantype jok, sunbooleantype *jcur, sunrealtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
void UseSundialsLinearSolver()
Attach a SUNDIALS GMRES linear solver to ARKode.
void Init(TimeDependentOperator &f_) override
Initialize ARKode: calls ARKStepCreate() to create the ARKStep memory and set some defaults.
static int MassMult2(N_Vector x, N_Vector v, sunrealtype t, void *mtimes_data)
Compute the matrix-vector product at time t.
void SetIRKTableNum(ARKODE_DIRKTableID table_id)
Choose a specific Butcher table for a diagonally implicit RK method.
void SetFixedStep(double dt)
Use a fixed time step size (disable temporal adaptivity).
void SetERKTableNum(ARKODE_ERKTableID table_id)
Choose a specific Butcher table for an explicit RK method.
int step_mode
ARKStep step mode (ARK_NORMAL or ARK_ONE_STEP).
Definition sundials.hpp:733
static int RHS1(sunrealtype t, const N_Vector y, N_Vector ydot, void *user_data)
void UseMFEMMassLinearSolver(int tdep)
Attach mass matrix linear system setup, solve, and matrix-vector product methods from the TimeDepende...
bool use_implicit
True for implicit or imex integration.
Definition sundials.hpp:734
void UseSundialsMassLinearSolver(int tdep)
Attach the SUNDIALS GMRES linear solver and the mass matrix matrix-vector product method from the Tim...
virtual ~ARKStepSolver()
Destroy the associated ARKode memory and SUNDIALS objects.
void SetIMEXTableNum(ARKODE_ERKTableID etable_id, ARKODE_DIRKTableID itable_id)
Choose a specific Butcher table for an IMEX RK method.
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i....
static int MassMult1(SUNMatrix M, N_Vector x, N_Vector v)
Compute the matrix-vector product .
static int MassSysSolve(SUNLinearSolver LS, SUNMatrix M, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system .
static int MassSysSetup(sunrealtype t, SUNMatrix M, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
void Step(Vector &x, real_t &t, real_t &dt) override
Integrate the ODE with ARKode using the specified step mode.
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
static int RHSB(sunrealtype t, N_Vector y, N_Vector yB, N_Vector yBdot, void *user_dataB)
Wrapper to compute the ODE RHS backward function.
void EvalQuadIntegrationB(double t, Vector &dG_dp)
Evaluate Quadrature solution.
void EvalQuadIntegration(double t, Vector &q)
Evaluate Quadrature.
static int LinSysSetupB(sunrealtype t, N_Vector y, N_Vector yB, N_Vector fyB, SUNMatrix A, sunbooleantype jok, sunbooleantype *jcur, sunrealtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system A x = b.
static int LinSysSolveB(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system A x = b.
long GetNumSteps()
Get Number of Steps for ForwardSolve.
static constexpr double default_abs_tolB
Default scalar backward absolute tolerance.
Definition sundials.hpp:599
static constexpr double default_abs_tolQB
Default scalar backward absolute quadrature tolerance.
Definition sundials.hpp:601
void SetMaxNStepsB(int mxstepsB)
Set the maximum number of backward steps.
static int RHSQ(sunrealtype t, const N_Vector y, N_Vector qdot, void *user_data)
Wrapper to compute the ODE RHS Quadrature function.
void Step(Vector &x, double &t, double &dt) override
void InitB(TimeDependentAdjointOperator &f_)
Initialize the adjoint problem.
static int RHSQB(sunrealtype t, N_Vector y, N_Vector yB, N_Vector qBdot, void *user_dataB)
Wrapper to compute the ODE RHS Backwards Quadrature function.
SundialsNVector * q
Quadrature vector.
Definition sundials.hpp:591
int indexB
backward problem index
Definition sundials.hpp:572
void GetForwardSolution(double tB, mfem::Vector &yy)
Get Interpolated Forward solution y at backward integration time tB.
SUNLinearSolver LSB
Linear solver for A.
Definition sundials.hpp:590
void SetSVtolerancesB(double reltol, Vector abstol)
Tolerance specification functions for the adjoint problem.
void UseSundialsLinearSolverB()
Use built in SUNDIALS Newton solver.
void SetWFTolerances(EWTFunction func)
Set multiplicative error weights.
SundialsNVector * yy
State vector.
Definition sundials.hpp:593
void Init(TimeDependentAdjointOperator &f_)
SUNMatrix AB
Linear system A = I - gamma J, M - gamma J, or J.
Definition sundials.hpp:589
int ncheck
number of checkpoints used so far
Definition sundials.hpp:571
SundialsNVector * qB
State vector.
Definition sundials.hpp:594
static int ewt(N_Vector y, N_Vector w, void *user_data)
Error control function.
virtual void StepB(Vector &w, double &t, double &dt)
Solve one adjoint time step.
void InitQuadIntegrationB(mfem::Vector &qB0, double reltolQB=1e-3, double abstolQB=1e-8)
Initialize Quadrature Integration (Adjoint)
static constexpr double default_rel_tolB
Default scalar backward relative tolerance.
Definition sundials.hpp:597
void InitAdjointSolve(int steps, int interpolation)
Initialize Adjoint.
SundialsNVector * yB
State vector.
Definition sundials.hpp:592
void InitQuadIntegration(mfem::Vector &q0, double reltolQ=1e-3, double abstolQ=1e-8)
virtual ~CVODESSolver()
Destroy the associated CVODES memory and SUNDIALS objects.
void SetSStolerancesB(double reltol, double abstol)
Tolerance specification functions for the adjoint problem.
void UseMFEMLinearSolverB()
Set Linear Solver for the backward problem.
Interface to the CVODE library – linear multi-step methods.
Definition sundials.hpp:430
void SetStepMode(int itask)
Select the CVODE step mode: CV_NORMAL (default) or CV_ONE_STEP.
Definition sundials.cpp:913
void SetRootFinder(int components, RootFunction func)
Initialize Root Finder.
Definition sundials.cpp:705
std::function< int(sunrealtype t, Vector y, Vector gout, CVODESolver *)> RootFunction
Typedef for root finding functions.
Definition sundials.hpp:456
void SetSStolerances(double reltol, double abstol)
Set the scalar relative and scalar absolute tolerances.
Definition sundials.cpp:918
void Init(TimeDependentOperator &f_) override
Initialize CVODE: calls CVodeCreate() to create the CVODE memory and set some defaults.
Definition sundials.cpp:752
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix A, N_Vector x, N_Vector b, sunrealtype tol)
Solve the linear system .
Definition sundials.cpp:728
virtual ~CVODESolver()
Destroy the associated CVODE memory and SUNDIALS objects.
static int root(sunrealtype t, N_Vector y, sunrealtype *gout, void *user_data)
Prototype to define root finding for CVODE.
Definition sundials.cpp:692
void SetMaxNSteps(int steps)
Set the maximum number of time steps.
Definition sundials.cpp:942
CVODESolver(int lmm)
Construct a serial wrapper to SUNDIALS' CVODE integrator.
Definition sundials.cpp:738
long GetNumSteps()
Get the number of internal steps taken so far.
Definition sundials.cpp:948
void Step(Vector &x, double &t, double &dt) override
Integrate the ODE with CVODE using the specified step mode.
Definition sundials.cpp:840
EWTFunction ewt_func
A class member to facilitate pointing to a user-specified error weight function.
Definition sundials.hpp:466
static int RHS(sunrealtype t, const N_Vector y, N_Vector ydot, void *user_data)
Number of components in gout.
Definition sundials.cpp:675
void SetMaxStep(double dt_max)
Set the maximum time step.
Definition sundials.cpp:936
int lmm_type
Linear multistep method type.
Definition sundials.hpp:432
void PrintInfo() const
Print various CVODE statistics.
Definition sundials.cpp:962
void UseSundialsLinearSolver()
Attach SUNDIALS GMRES linear solver to CVODE.
Definition sundials.cpp:898
void UseMFEMLinearSolver()
Attach the linear system setup and solve methods from the TimeDependentOperator i....
Definition sundials.cpp:867
RootFunction root_func
A class member to facilitate pointing to a user-specified root function.
Definition sundials.hpp:459
std::function< int(Vector y, Vector w, CVODESolver *)> EWTFunction
Typedef declaration for error weight functions.
Definition sundials.hpp:462
int step_mode
CVODE step mode (CV_NORMAL or CV_ONE_STEP).
Definition sundials.hpp:433
static int LinSysSetup(sunrealtype t, N_Vector y, N_Vector fy, SUNMatrix A, sunbooleantype jok, sunbooleantype *jcur, sunrealtype gamma, void *user_data, N_Vector tmp1, N_Vector tmp2, N_Vector tmp3)
Setup the linear system .
Definition sundials.cpp:713
void SetSVtolerances(double reltol, Vector abstol)
Set the scalar relative and vector of absolute tolerances.
Definition sundials.cpp:924
void SetMaxOrder(int max_order)
Set the maximum method order.
Definition sundials.cpp:956
static MemoryType GetDeviceMemoryType()
Get the current Device MemoryType. This is the MemoryType used by most MFEM classes when allocating m...
Definition device.hpp:277
Wrapper for hypre's parallel vector class.
Definition hypre.hpp:230
int print_level
(DEPRECATED) Legacy print level definition, which is left for compatibility with custom iterative sol...
Definition solvers.hpp:156
Interface to the KINSOL library – nonlinear solver methods.
Definition sundials.hpp:897
SundialsNVector * f_scale
scaling vectors
Definition sundials.hpp:901
KINSolver(int strategy, bool oper_grad=true)
Construct a serial wrapper to SUNDIALS' KINSOL nonlinear solver.
double aa_damping
Anderson Acceleration damping.
Definition sundials.hpp:905
void SetJFNK(bool use_jfnk)
Set the Jacobian Free Newton Krylov flag. The default is false.
void SetJFNKSolver(Solver &solver)
virtual ~KINSolver()
Destroy the associated KINSOL memory.
static int Mult(const N_Vector u, N_Vector fu, void *user_data)
Wrapper to compute the nonlinear residual .
static int LinSysSolve(SUNLinearSolver LS, SUNMatrix J, N_Vector u, N_Vector b, sunrealtype tol)
Solve the linear system .
int aa_delay
Anderson Acceleration delay.
Definition sundials.hpp:904
int global_strategy
KINSOL solution strategy.
Definition sundials.hpp:899
bool jfnk
enable JFNK
Definition sundials.hpp:908
static int PrecSolve(N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, N_Vector vv, void *user_data)
Solve the preconditioner equation .
int maxlrs
Maximum linear solver restarts.
Definition sundials.hpp:911
const Operator * jacobian
stores oper->GetGradient()
Definition sundials.hpp:902
int maxli
Maximum linear iterations.
Definition sundials.hpp:910
Vector wrk
Work vector needed for the JFNK PC.
Definition sundials.hpp:909
void SetScaledStepTol(double sstol)
Set KINSOL's scaled step tolerance.
void SetSolver(Solver &solver) override
Set the linear solver for inverting the Jacobian.
SundialsNVector * y_scale
Definition sundials.hpp:901
void SetOperator(const Operator &op) override
Set the nonlinear Operator of the system and initialize KINSOL.
void SetMaxSetupCalls(int max_calls)
Set maximum number of nonlinear iterations without a Jacobian update.
static int LinSysSetup(N_Vector u, N_Vector fu, SUNMatrix J, void *user_data, N_Vector tmp1, N_Vector tmp2)
Setup the linear system .
int aa_orth
Anderson Acceleration orthogonalization routine.
Definition sundials.hpp:906
static int PrecSetup(N_Vector uu, N_Vector uscale, N_Vector fval, N_Vector fscale, void *user_data)
Setup the preconditioner.
void SetDamping(double damping)
void SetLSMaxRestarts(int m)
Set the maximum number of linear solver restarts.
double fp_damping
Fixed Point or Picard damping parameter.
Definition sundials.hpp:907
void EnableAndersonAcc(int n, int orth=KIN_ORTH_MGS, int delay=0, double damping=1.0)
Enable Anderson Acceleration for KIN_FP or KIN_PICARD.
void SetPrintLevel(int print_lvl) override
Set the print level for the KINSetPrintLevel function.
void SetPreconditioner(Solver &solver) override
Equivalent to SetSolver(solver).
Definition sundials.hpp:984
void SetLSMaxIter(int m)
Set the maximum number of linear solver iterations.
bool use_oper_grad
use the Jv prod function
Definition sundials.hpp:900
int aa_n
number of acceleration vectors
Definition sundials.hpp:903
static int GradientMult(N_Vector v, N_Vector Jv, N_Vector u, sunbooleantype *new_u, void *user_data)
Wrapper to compute the Jacobian-vector product .
Newton's method for solving F(x)=b for a given operator F.
Definition solvers.hpp:781
Abstract class for solving systems of ODEs: dx/dt = f(x,t)
Definition ode.hpp:121
Abstract operator.
Definition operator.hpp:25
Base class for solvers.
Definition operator.hpp:792
SundialsMemHelper & operator=(const SundialsMemHelper &)=delete
Disable copy assignment.
SundialsMemHelper(const SundialsMemHelper &that_helper)=delete
Disable copy construction.
SundialsMemHelper()=default
Default constructor – object must be moved to.
static int SundialsMemHelper_Alloc(SUNMemoryHelper helper, SUNMemory *memptr, size_t memsize, SUNMemoryType mem_type #if(SUNDIALS_VERSION_MAJOR >=6), void *queue #endif)
Definition sundials.cpp:282
static int SundialsMemHelper_Dealloc(SUNMemoryHelper helper, SUNMemory sunmem #if(SUNDIALS_VERSION_MAJOR >=6), void *queue #endif)
Definition sundials.cpp:326
SundialsMemHelper(SUNContext context)
Definition sundials.hpp:165
Vector interface for SUNDIALS N_Vectors.
Definition sundials.hpp:222
N_Vector_ID GetNVectorID(N_Vector x_) const
Returns the N_Vector_ID for the N_Vector x_.
Definition sundials.hpp:274
N_Vector StealNVector()
Changes the ownership of the vector.
Definition sundials.hpp:327
long GlobalSize() const
Returns the MPI global length for the internal N_Vector x.
Definition sundials.hpp:288
MPI_Comm GetComm() const
Returns the MPI communicator for the internal N_Vector x.
Definition sundials.hpp:278
void MakeRef(Vector &base, int offset, int s)
Reset the Vector to be a reference to a sub-vector of base.
Definition sundials.hpp:305
void SetSize(int s, long glob_size=0)
Resize the vector to size s.
Definition sundials.cpp:559
static bool UseManagedMemory()
Definition sundials.hpp:359
static N_Vector MakeNVector(bool use_device)
Create a N_Vector.
Definition sundials.cpp:577
int GetOwnership() const
Gets ownership of the internal N_Vector.
Definition sundials.hpp:333
~SundialsNVector()
Calls SUNDIALS N_VDestroy function if the N_Vector is owned by 'this'.
Definition sundials.cpp:545
void SetDataAndSize(double *d, int s, long glob_size=0)
Set the vector data and size.
Definition sundials.cpp:571
void SetOwnership(int own)
Sets ownership of the internal N_Vector.
Definition sundials.hpp:330
N_Vector_ID GetNVectorID() const
Returns the N_Vector_ID for the internal N_Vector.
Definition sundials.hpp:271
void _SetNvecDataAndSize_(long glob_size=0)
Set data and length of internal N_Vector x from 'this'.
Definition sundials.cpp:364
void _SetDataAndSize_()
Set data and length from the internal N_Vector x.
Definition sundials.cpp:438
N_Vector x
The actual SUNDIALS object.
Definition sundials.hpp:227
void MakeRef(Vector &base, int offset)
Reset the Vector to be a reference to a sub-vector of base without changing its current size.
Definition sundials.hpp:315
void SetData(double *d)
Definition sundials.cpp:565
bool MPIPlusX() const
Definition sundials.hpp:341
SundialsNVector()
Creates an empty SundialsNVector.
Definition sundials.cpp:486
Base class for interfacing with SUNDIALS packages.
Definition sundials.hpp:374
static constexpr double default_abs_tol
Default scalar absolute tolerance.
Definition sundials.hpp:399
SUNMatrix M
Mass matrix M.
Definition sundials.hpp:384
int flag
Last flag returned from a call to SUNDIALS.
Definition sundials.hpp:377
long saved_global_size
Global vector length on last initialization.
Definition sundials.hpp:379
SundialsSolver()
Protected constructor: objects of this type should be constructed only as part of a derived class.
Definition sundials.hpp:403
static constexpr double default_rel_tol
Default scalar relative tolerance.
Definition sundials.hpp:397
void * GetMem() const
Access the SUNDIALS memory structure.
Definition sundials.hpp:417
bool reinit
Flag to signal memory reinitialization is need.
Definition sundials.hpp:378
void AllocateEmptyNVector(N_Vector &y)
SundialsNVector * Y
State vector.
Definition sundials.hpp:381
void * sundials_mem
SUNDIALS mem structure.
Definition sundials.hpp:376
void AllocateEmptyNVector(N_Vector &y, MPI_Comm comm)
SUNLinearSolver LSA
Linear solver for A.
Definition sundials.hpp:385
SUNLinearSolver LSM
Linear solver for M.
Definition sundials.hpp:386
int GetFlag() const
Returns the last flag returned by a call to a SUNDIALS function.
Definition sundials.hpp:420
SUNNonlinearSolver NLS
Nonlinear solver.
Definition sundials.hpp:387
bool Parallel() const
Definition sundials.hpp:390
Singleton class for SUNContext and SundialsMemHelper objects.
Definition sundials.hpp:176
Sundials(Sundials &other)=delete
Disable copy construction.
static SUNContext & GetContext()
Provides access to the SUNContext object.
Definition sundials.cpp:185
void operator=(const Sundials &other)=delete
Disable copy assignment.
static SundialsMemHelper & GetMemHelper()
Provides access to the SundialsMemHelper object.
Definition sundials.cpp:190
static void Finalize()
Finalize sundials (called automatically at program exit if Sundials::Init() has been called).
Definition sundials.cpp:226
static void Init()
Initializes SUNContext and SundialsMemHelper objects. Should be called at the beginning of the callin...
Definition sundials.cpp:174
Base abstract class for first order time dependent operators.
Definition operator.hpp:344
Vector data type.
Definition vector.hpp:82
virtual const real_t * Read(bool on_dev=true) const
Shortcut for mfem::Read(vec.GetMemory(), vec.Size(), on_dev).
Definition vector.hpp:520
void MakeRef(Vector &base, int offset, int size)
Reset the Vector to be a reference to a sub-vector of base.
Definition vector.hpp:660
real_t b
Definition lissajous.cpp:42
real_t u(const Vector &xvec)
Definition lor_mms.hpp:22
float real_t
Definition config.hpp:46
int ARKODE_DIRKTableID
Definition sundials.hpp:66
realtype sunrealtype
'sunrealtype' was first introduced in v6.0.0
Definition sundials.hpp:76
constexpr ARKODE_ERKTableID ARKODE_FEHLBERG_13_7_8
Definition sundials.hpp:69
int ARKODE_ERKTableID
Definition sundials.hpp:65
@ SUN_PREC_NONE
Definition sundials.hpp:81
@ SUN_PREC_LEFT
Definition sundials.hpp:81
@ SUN_PREC_RIGHT
Definition sundials.hpp:81
@ SUN_PREC_BOTH
Definition sundials.hpp:81
booleantype sunbooleantype
'sunbooleantype' was first introduced in v6.0.0
Definition sundials.hpp:78
constexpr ARKODE_ERKTableID ARKODE_ERK_NONE
Definition sundials.hpp:67
constexpr ARKODE_DIRKTableID ARKODE_DIRK_NONE
Definition sundials.hpp:68
@ SUN_SUCCESS
Definition sundials.hpp:95
void * SUNContext
Definition sundials.hpp:73