Reference documentation for deal.II version 9.3.2
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
ida.h
Go to the documentation of this file.
1 //-----------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2020 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_ida_h
17 #define dealii_sundials_ida_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/mpi.h>
22 #ifdef DEAL_II_WITH_SUNDIALS
23 # if DEAL_II_SUNDIALS_VERSION_LT(4, 0, 0)
24 
26 # include <deal.II/base/exceptions.h>
27 # include <deal.II/base/logstream.h>
29 # ifdef DEAL_II_WITH_PETSC
32 # endif
33 # include <deal.II/lac/vector.h>
35 
36 # ifdef DEAL_II_SUNDIALS_WITH_IDAS
37 # include <idas/idas.h>
38 # else
39 # include <ida/ida.h>
40 # endif
41 
42 # include <sundials/sundials_config.h>
43 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
44 # include <ida/ida_spbcgs.h>
45 # include <ida/ida_spgmr.h>
46 # include <ida/ida_sptfqmr.h>
47 # endif
48 # include <boost/signals2.hpp>
49 
50 # include <nvector/nvector_serial.h>
51 # include <sundials/sundials_math.h>
52 # include <sundials/sundials_types.h>
53 
54 # include <memory>
55 
56 
58 
59 // Shorthand notation for IDA error codes.
60 # define AssertIDA(code) Assert(code >= 0, ExcIDAError(code))
61 
62 namespace SUNDIALS
63 {
233  template <typename VectorType = Vector<double>>
234  class IDA
235  {
236  public:
240  class AdditionalData
241  {
242  public:
252  enum InitialConditionCorrection
253  {
257  none = 0,
258 
265  use_y_diff = 1,
266 
270  use_y_dot = 2
271  };
272 
304  AdditionalData( // Initial parameters
305  const double initial_time = 0.0,
306  const double final_time = 1.0,
307  const double initial_step_size = 1e-2,
308  const double output_period = 1e-1,
309  // Running parameters
310  const double minimum_step_size = 1e-6,
311  const unsigned int maximum_order = 5,
312  const unsigned int maximum_non_linear_iterations = 10,
313  // Error parameters
314  const double absolute_tolerance = 1e-6,
315  const double relative_tolerance = 1e-5,
316  const bool ignore_algebraic_terms_for_errors = true,
317  // Initial conditions parameters
318  const InitialConditionCorrection &ic_type = use_y_diff,
319  const InitialConditionCorrection &reset_type = use_y_diff,
320  const unsigned int maximum_non_linear_iterations_ic = 5)
321  : initial_time(initial_time)
322  , final_time(final_time)
323  , initial_step_size(initial_step_size)
324  , minimum_step_size(minimum_step_size)
325  , absolute_tolerance(absolute_tolerance)
326  , relative_tolerance(relative_tolerance)
327  , maximum_order(maximum_order)
328  , output_period(output_period)
329  , ignore_algebraic_terms_for_errors(ignore_algebraic_terms_for_errors)
330  , ic_type(ic_type)
331  , reset_type(reset_type)
332  , maximum_non_linear_iterations_ic(maximum_non_linear_iterations_ic)
333  , maximum_non_linear_iterations(maximum_non_linear_iterations)
334  {}
335 
378  void
379  add_parameters(ParameterHandler &prm)
380  {
381  prm.add_parameter("Initial time", initial_time);
382  prm.add_parameter("Final time", final_time);
383  prm.add_parameter("Time interval between each output", output_period);
384 
385  prm.enter_subsection("Running parameters");
386  prm.add_parameter("Initial step size", initial_step_size);
387  prm.add_parameter("Minimum step size", minimum_step_size);
388  prm.add_parameter("Maximum order of BDF", maximum_order);
389  prm.add_parameter("Maximum number of nonlinear iterations",
390  maximum_non_linear_iterations);
391  prm.leave_subsection();
392 
393  prm.enter_subsection("Error control");
394  prm.add_parameter("Absolute error tolerance", absolute_tolerance);
395  prm.add_parameter("Relative error tolerance", relative_tolerance);
396  prm.add_parameter(
397  "Ignore algebraic terms for error computations",
398  ignore_algebraic_terms_for_errors,
399  "Indicate whether or not to suppress algebraic variables "
400  "in the local error test.");
401  prm.leave_subsection();
402 
403  prm.enter_subsection("Initial condition correction parameters");
404  static std::string ic_type_str = "use_y_diff";
405  prm.add_parameter(
406  "Correction type at initial time",
407  ic_type_str,
408  "This is one of the following three options for the "
409  "initial condition calculation. \n"
410  " none: do not try to make initial conditions consistent. \n"
411  " use_y_diff: compute the algebraic components of y and differential\n"
412  " components of y_dot, given the differential components of y. \n"
413  " This option requires that the user specifies differential and \n"
414  " algebraic components in the function get_differential_components.\n"
415  " use_y_dot: compute all components of y, given y_dot.",
416  Patterns::Selection("none|use_y_diff|use_y_dot"));
417  prm.add_action("Correction type at initial time",
418  [&](const std::string &value) {
419  if (value == "use_y_diff")
420  ic_type = use_y_diff;
421  else if (value == "use_y_dot")
422  ic_type = use_y_dot;
423  else if (value == "none")
424  ic_type = none;
425  else
426  AssertThrow(false, ExcInternalError());
427  });
428 
429  static std::string reset_type_str = "use_y_diff";
430  prm.add_parameter(
431  "Correction type after restart",
432  reset_type_str,
433  "This is one of the following three options for the "
434  "initial condition calculation. \n"
435  " none: do not try to make initial conditions consistent. \n"
436  " use_y_diff: compute the algebraic components of y and differential\n"
437  " components of y_dot, given the differential components of y. \n"
438  " This option requires that the user specifies differential and \n"
439  " algebraic components in the function get_differential_components.\n"
440  " use_y_dot: compute all components of y, given y_dot.",
441  Patterns::Selection("none|use_y_diff|use_y_dot"));
442  prm.add_action("Correction type after restart",
443  [&](const std::string &value) {
444  if (value == "use_y_diff")
445  reset_type = use_y_diff;
446  else if (value == "use_y_dot")
447  reset_type = use_y_dot;
448  else if (value == "none")
449  reset_type = none;
450  else
451  AssertThrow(false, ExcInternalError());
452  });
453  prm.add_parameter("Maximum number of nonlinear iterations",
454  maximum_non_linear_iterations_ic);
455  prm.leave_subsection();
456  }
457 
461  double initial_time;
462 
466  double final_time;
467 
471  double initial_step_size;
472 
476  double minimum_step_size;
477 
481  double absolute_tolerance;
482 
486  double relative_tolerance;
487 
491  unsigned int maximum_order;
492 
496  double output_period;
497 
501  bool ignore_algebraic_terms_for_errors;
502 
517  InitialConditionCorrection ic_type;
518 
529  InitialConditionCorrection reset_type;
530 
534  unsigned maximum_non_linear_iterations_ic;
535 
539  unsigned int maximum_non_linear_iterations;
540  };
541 
583  IDA(const AdditionalData &data = AdditionalData(),
584  const MPI_Comm & mpi_comm = MPI_COMM_WORLD);
585 
589  ~IDA();
590 
595  unsigned int
596  solve_dae(VectorType &solution, VectorType &solution_dot);
597 
619  void
620  reset(const double t, const double h, VectorType &y, VectorType &yp);
621 
625  std::function<void(VectorType &)> reinit_vector;
626 
637  std::function<int(const double t,
638  const VectorType &y,
639  const VectorType &y_dot,
640  VectorType & res)>
641  residual;
642 
673  std::function<int(const double t,
674  const VectorType &y,
675  const VectorType &y_dot,
676  const double alpha)>
677  setup_jacobian;
678 
707  std::function<int(const VectorType &rhs, VectorType &dst)>
708  solve_jacobian_system;
709 
724  std::function<void(const double t,
725  const VectorType & sol,
726  const VectorType & sol_dot,
727  const unsigned int step_number)>
728  output_step;
729 
746  std::function<bool(const double t, VectorType &sol, VectorType &sol_dot)>
747  solver_should_restart;
748 
763  std::function<IndexSet()> differential_components;
764 
771  std::function<VectorType &()> get_local_tolerances;
772 
776  DeclException1(ExcIDAError,
777  int,
778  << "One of the SUNDIALS IDA internal functions "
779  << " returned a negative error code: " << arg1
780  << ". Please consult SUNDIALS manual.");
781 
782 
783  private:
788  DeclException1(ExcFunctionNotProvided,
789  std::string,
790  << "Please provide an implementation for the function \""
791  << arg1 << "\"");
792 
798  void
799  set_functions_to_trigger_an_assert();
800 
804  AdditionalData data;
805 
809  void *ida_mem;
810 
814  N_Vector yy;
815 
819  N_Vector yp;
820 
824  N_Vector abs_tolls;
825 
829  N_Vector diff_id;
830 
836  MPI_Comm communicator;
837 
842 
843 # ifdef DEAL_II_WITH_PETSC
844 # ifdef PETSC_USE_COMPLEX
845  static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
846  "Sundials does not support complex scalar types, "
847  "but PETSc is configured to use a complex scalar type!");
848 
849  static_assert(
850  !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
851  "Sundials does not support complex scalar types, "
852  "but PETSc is configured to use a complex scalar type!");
853 # endif // PETSC_USE_COMPLEX
854 # endif // DEAL_II_WITH_PETSC
855  };
856 
857 } // namespace SUNDIALS
858 
860 
861 # endif
862 #endif // DEAL_II_WITH_SUNDIALS
863 
864 #endif
void add_parameter(const std::string &entry, ParameterType &parameter, const std::string &documentation="", const Patterns::PatternBase &pattern= *Patterns::Tools::Convert< ParameterType >::to_pattern(), const bool has_to_be_set=false)
void add_action(const std::string &entry, const std::function< void(const std::string &value)> &action)
void enter_subsection(const std::string &subsection)
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:396
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:397
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1575
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
DeclException1(ExcARKodeError, int,<< "One of the SUNDIALS ARKode internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")