Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
arkode.h
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2019 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 //---------------------------------------------------------------
15 
16 #ifndef dealii_sundials_arkode_h
17 #define dealii_sundials_arkode_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/mpi.h>
22 
23 #ifdef DEAL_II_WITH_SUNDIALS
24 
25 # include <deal.II/base/conditional_ostream.h>
26 # include <deal.II/base/exceptions.h>
27 # include <deal.II/base/logstream.h>
28 # include <deal.II/base/mpi.h>
29 # include <deal.II/base/parameter_handler.h>
30 # ifdef DEAL_II_WITH_PETSC
31 # include <deal.II/lac/petsc_block_vector.h>
32 # include <deal.II/lac/petsc_vector.h>
33 # endif
34 # include <deal.II/lac/vector.h>
35 # include <deal.II/lac/vector_memory.h>
36 
37 # include <arkode/arkode.h>
38 # include <arkode/arkode_impl.h>
39 # include <nvector/nvector_serial.h>
40 # ifdef DEAL_II_WITH_MPI
41 # include <nvector/nvector_parallel.h>
42 # endif
43 # include <boost/signals2.hpp>
44 
45 # include <sundials/sundials_math.h>
46 # include <sundials/sundials_types.h>
47 
48 # include <memory>
49 
50 
51 DEAL_II_NAMESPACE_OPEN
52 
53 // Shorthand notation for ARKODE error codes.
54 # define AssertARKode(code) Assert(code >= 0, ExcARKodeError(code))
55 
59 namespace SUNDIALS
60 {
313  template <typename VectorType = Vector<double>>
314  class ARKode
315  {
316  public:
321  {
322  public:
350  // Initial parameters
351  const double initial_time = 0.0,
352  const double final_time = 1.0,
353  const double initial_step_size = 1e-2,
354  const double output_period = 1e-1,
355  // Running parameters
356  const double minimum_step_size = 1e-6,
357  const unsigned int maximum_order = 5,
358  const unsigned int maximum_non_linear_iterations = 10,
359  const bool implicit_function_is_linear = false,
360  const bool implicit_function_is_time_independent = false,
361  // Error parameters
362  const double absolute_tolerance = 1e-6,
363  const double relative_tolerance = 1e-5)
376  {}
377 
420  void
422  {
423  prm.add_parameter("Initial time", initial_time);
424  prm.add_parameter("Final time", final_time);
425  prm.add_parameter("Time interval between each output", output_period);
426 
427  prm.enter_subsection("Running parameters");
428  prm.add_parameter("Initial step size", initial_step_size);
429  prm.add_parameter("Minimum step size", minimum_step_size);
430  prm.add_parameter("Maximum order of ARK", maximum_order);
431  prm.add_parameter("Maximum number of nonlinear iterations",
433  prm.add_parameter("Implicit function is linear",
435  prm.add_parameter("Implicit function is time independent",
437  prm.leave_subsection();
438 
439  prm.enter_subsection("Error control");
440  prm.add_parameter("Absolute error tolerance", absolute_tolerance);
441  prm.add_parameter("Relative error tolerance", relative_tolerance);
442  prm.leave_subsection();
443  }
444 
448  double initial_time;
449 
453  double final_time;
454 
459 
464 
469 
474 
478  unsigned int maximum_order;
479 
484 
490 
495 
501  };
502 
515  const MPI_Comm mpi_comm = MPI_COMM_WORLD);
516 
520  ~ARKode();
521 
526  unsigned int
527  solve_ode(VectorType &solution);
528 
543  void
544  reset(const double t, const double h, const VectorType &y);
545 
550  std::function<void(VectorType &)> reinit_vector;
551 
568  std::function<
569  int(const double t, const VectorType &y, VectorType &explicit_f)>
571 
588  std::function<int(const double t, const VectorType &y, VectorType &res)>
590 
675  std::function<int(const int convfail,
676  const double t,
677  const double gamma,
678  const VectorType &ypred,
679  const VectorType &fpred,
680  bool & j_is_current)>
682 
727  std::function<int(const double t,
728  const double gamma,
729  const VectorType &ycur,
730  const VectorType &fcur,
731  const VectorType &rhs,
732  VectorType & dst)>
734 
735 
769  std::function<int(const double t)> setup_mass;
770 
789  std::function<int(const VectorType &rhs, VectorType &dst)>
791 
806  std::function<void(const double t,
807  const VectorType & sol,
808  const unsigned int step_number)>
810 
828  std::function<bool(const double t, VectorType &sol)> solver_should_restart;
829 
836  std::function<VectorType &()> get_local_tolerances;
837 
841  DeclException1(ExcARKodeError,
842  int,
843  << "One of the SUNDIALS ARKode internal functions "
844  << " returned a negative error code: " << arg1
845  << ". Please consult SUNDIALS manual.");
846 
847 
848  private:
853  DeclException1(ExcFunctionNotProvided,
854  std::string,
855  << "Please provide an implementation for the function \""
856  << arg1 << "\"");
857 
863  void
865 
870 
874  void *arkode_mem;
875 
879  N_Vector yy;
880 
884  N_Vector abs_tolls;
885 
891  MPI_Comm communicator;
892 
897 
898 # ifdef DEAL_II_WITH_PETSC
899 # ifdef PETSC_USE_COMPLEX
900  static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
901  "Sundials does not support complex scalar types, "
902  "but PETSc is configured to use a complex scalar type!");
903 
904  static_assert(
905  !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
906  "Sundials does not support complex scalar types, "
907  "but PETSc is configured to use a complex scalar type!");
908 # endif // PETSC_USE_COMPLEX
909 # endif // DEAL_II_WITH_PETSC
910  };
911 
912 } // namespace SUNDIALS
913 
914 DEAL_II_NAMESPACE_CLOSE
915 #endif
916 
917 
918 #endif
void reset(const double t, const double h, const VectorType &y)
Definition: arkode.cc:354
void * arkode_mem
Definition: arkode.h:874
GrowingVectorMemory< VectorType > mem
Definition: arkode.h:896
std::function< int(const VectorType &rhs, VectorType &dst)> solve_mass_system
Definition: arkode.h:790
MPI_Comm communicator
Definition: arkode.h:891
N_Vector yy
Definition: arkode.h:879
std::function< void(const double t, const VectorType &sol, const unsigned int step_number)> output_step
Definition: arkode.h:809
N_Vector abs_tolls
Definition: arkode.h:884
std::function< VectorType &()> get_local_tolerances
Definition: arkode.h:836
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern())
unsigned int maximum_non_linear_iterations
Definition: arkode.h:489
std::function< int(const double t)> setup_mass
Definition: arkode.h:769
std::function< bool(const double t, VectorType &sol)> solver_should_restart
Definition: arkode.h:828
void enter_subsection(const std::string &subsection)
std::function< int(const double t, const double gamma, const VectorType &ycur, const VectorType &fcur, const VectorType &rhs, VectorType &dst)> solve_jacobian_system
Definition: arkode.h:733
AdditionalData data
Definition: arkode.h:869
std::function< int(const double t, const VectorType &y, VectorType &res)> implicit_function
Definition: arkode.h:589
DeclException1(ExcARKodeError, int,<< "One of the SUNDIALS ARKode internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")
void add_parameters(ParameterHandler &prm)
Definition: arkode.h:421
std::function< int(const double t, const VectorType &y, VectorType &explicit_f)> explicit_function
Definition: arkode.h:570
std::function< void(VectorType &)> reinit_vector
Definition: arkode.h:550
unsigned int solve_ode(VectorType &solution)
Definition: arkode.cc:272
std::function< int(const int convfail, const double t, const double gamma, const VectorType &ypred, const VectorType &fpred, bool &j_is_current)> setup_jacobian
Definition: arkode.h:681
void set_functions_to_trigger_an_assert()
Definition: arkode.cc:495
AdditionalData(const double initial_time=0.0, const double final_time=1.0, const double initial_step_size=1e-2, const double output_period=1e-1, const double minimum_step_size=1e-6, const unsigned int maximum_order=5, const unsigned int maximum_non_linear_iterations=10, const bool implicit_function_is_linear=false, const bool implicit_function_is_time_independent=false, const double absolute_tolerance=1e-6, const double relative_tolerance=1e-5)
Definition: arkode.h:349
ARKode(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
Definition: arkode.cc:240