Reference documentation for deal.II version 9.2.0
\(\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 
25 # include <deal.II/base/exceptions.h>
26 # include <deal.II/base/logstream.h>
28 # ifdef DEAL_II_WITH_PETSC
31 # endif
32 # include <deal.II/lac/vector.h>
34 
35 # ifdef DEAL_II_SUNDIALS_WITH_IDAS
36 # include <idas/idas.h>
37 # else
38 # include <ida/ida.h>
39 # endif
40 
41 # include <sundials/sundials_config.h>
42 # if DEAL_II_SUNDIALS_VERSION_LT(3, 0, 0)
43 # include <ida/ida_spbcgs.h>
44 # include <ida/ida_spgmr.h>
45 # include <ida/ida_sptfqmr.h>
46 # endif
47 # include <boost/signals2.hpp>
48 
49 # include <nvector/nvector_serial.h>
50 # include <sundials/sundials_math.h>
51 # include <sundials/sundials_types.h>
52 
53 # include <memory>
54 
55 
57 
58 // Shorthand notation for IDA error codes.
59 # define AssertIDA(code) Assert(code >= 0, ExcIDAError(code))
60 
61 namespace SUNDIALS
62 {
234  template <typename VectorType = Vector<double>>
235  class IDA
236  {
237  public:
242  {
243  public:
254  {
258  none = 0,
259 
267 
272  };
273 
305  AdditionalData( // Initial parameters
306  const double initial_time = 0.0,
307  const double final_time = 1.0,
308  const double initial_step_size = 1e-2,
309  const double output_period = 1e-1,
310  // Running parameters
311  const double minimum_step_size = 1e-6,
312  const unsigned int maximum_order = 5,
313  const unsigned int maximum_non_linear_iterations = 10,
314  // Error parameters
315  const double absolute_tolerance = 1e-6,
316  const double relative_tolerance = 1e-5,
317  const bool ignore_algebraic_terms_for_errors = true,
318  // Initial conditions parameters
321  const unsigned int maximum_non_linear_iterations_ic = 5)
331  , ic_type(ic_type)
335  {}
336 
379  void
381  {
382  prm.add_parameter("Initial time", initial_time);
383  prm.add_parameter("Final time", final_time);
384  prm.add_parameter("Time interval between each output", output_period);
385 
386  prm.enter_subsection("Running parameters");
387  prm.add_parameter("Initial step size", initial_step_size);
388  prm.add_parameter("Minimum step size", minimum_step_size);
389  prm.add_parameter("Maximum order of BDF", maximum_order);
390  prm.add_parameter("Maximum number of nonlinear iterations",
392  prm.leave_subsection();
393 
394  prm.enter_subsection("Error control");
395  prm.add_parameter("Absolute error tolerance", absolute_tolerance);
396  prm.add_parameter("Relative error tolerance", relative_tolerance);
397  prm.add_parameter(
398  "Ignore algebraic terms for error computations",
400  "Indicate whether or not to suppress algebraic variables "
401  "in the local error test.");
402  prm.leave_subsection();
403 
404  prm.enter_subsection("Initial condition correction parameters");
405  static std::string ic_type_str = "use_y_diff";
406  prm.add_parameter(
407  "Correction type at initial time",
408  ic_type_str,
409  "This is one of the following three options for the "
410  "initial condition calculation. \n"
411  " none: do not try to make initial conditions consistent. \n"
412  " use_y_diff: compute the algebraic components of y and differential\n"
413  " components of y_dot, given the differential components of y. \n"
414  " This option requires that the user specifies differential and \n"
415  " algebraic components in the function get_differential_components.\n"
416  " use_y_dot: compute all components of y, given y_dot.",
417  Patterns::Selection("none|use_y_diff|use_y_dot"));
418  prm.add_action("Correction type at initial time",
419  [&](const std::string &value) {
420  if (value == "use_y_diff")
422  else if (value == "use_y_dot")
423  ic_type = use_y_dot;
424  else if (value == "none")
425  ic_type = none;
426  else
427  AssertThrow(false, ExcInternalError());
428  });
429 
430  static std::string reset_type_str = "use_y_diff";
431  prm.add_parameter(
432  "Correction type after restart",
433  reset_type_str,
434  "This is one of the following three options for the "
435  "initial condition calculation. \n"
436  " none: do not try to make initial conditions consistent. \n"
437  " use_y_diff: compute the algebraic components of y and differential\n"
438  " components of y_dot, given the differential components of y. \n"
439  " This option requires that the user specifies differential and \n"
440  " algebraic components in the function get_differential_components.\n"
441  " use_y_dot: compute all components of y, given y_dot.",
442  Patterns::Selection("none|use_y_diff|use_y_dot"));
443  prm.add_action("Correction type after restart",
444  [&](const std::string &value) {
445  if (value == "use_y_diff")
447  else if (value == "use_y_dot")
449  else if (value == "none")
450  reset_type = none;
451  else
452  AssertThrow(false, ExcInternalError());
453  });
454  prm.add_parameter("Maximum number of nonlinear iterations",
456  prm.leave_subsection();
457  }
458 
462  double initial_time;
463 
467  double final_time;
468 
473 
478 
483 
488 
492  unsigned int maximum_order;
493 
498 
503 
519 
531 
536 
541  };
542 
585  const MPI_Comm mpi_comm = MPI_COMM_WORLD);
586 
590  ~IDA();
591 
596  unsigned int
597  solve_dae(VectorType &solution, VectorType &solution_dot);
598 
620  void
621  reset(const double t, const double h, VectorType &y, VectorType &yp);
622 
626  std::function<void(VectorType &)> reinit_vector;
627 
638  std::function<int(const double t,
639  const VectorType &y,
640  const VectorType &y_dot,
641  VectorType & res)>
643 
674  std::function<int(const double t,
675  const VectorType &y,
676  const VectorType &y_dot,
677  const double alpha)>
679 
708  std::function<int(const VectorType &rhs, VectorType &dst)>
710 
725  std::function<void(const double t,
726  const VectorType & sol,
727  const VectorType & sol_dot,
728  const unsigned int step_number)>
730 
747  std::function<bool(const double t, VectorType &sol, VectorType &sol_dot)>
749 
759 
766  std::function<VectorType &()> get_local_tolerances;
767 
771  DeclException1(ExcIDAError,
772  int,
773  << "One of the SUNDIALS IDA internal functions "
774  << " returned a negative error code: " << arg1
775  << ". Please consult SUNDIALS manual.");
776 
777 
778  private:
783  DeclException1(ExcFunctionNotProvided,
784  std::string,
785  << "Please provide an implementation for the function \""
786  << arg1 << "\"");
787 
793  void
795 
800 
804  void *ida_mem;
805 
809  N_Vector yy;
810 
814  N_Vector yp;
815 
819  N_Vector abs_tolls;
820 
824  N_Vector diff_id;
825 
832 
837 
838 # ifdef DEAL_II_WITH_PETSC
839 # ifdef PETSC_USE_COMPLEX
841  "Sundials does not support complex scalar types, "
842  "but PETSc is configured to use a complex scalar type!");
843 
844  static_assert(
846  "Sundials does not support complex scalar types, "
847  "but PETSc is configured to use a complex scalar type!");
848 # endif // PETSC_USE_COMPLEX
849 # endif // DEAL_II_WITH_PETSC
850  };
851 
852 } // namespace SUNDIALS
853 
855 
856 #endif // DEAL_II_WITH_SUNDIALS
857 
858 #endif
IndexSet
Definition: index_set.h:74
SUNDIALS::IDA::AdditionalData::maximum_order
unsigned int maximum_order
Definition: ida.h:492
SUNDIALS::IDA::AdditionalData::ic_type
InitialConditionCorrection ic_type
Definition: ida.h:518
SUNDIALS::IDA::communicator
MPI_Comm communicator
Definition: ida.h:831
SUNDIALS::IDA::AdditionalData::reset_type
InitialConditionCorrection reset_type
Definition: ida.h:530
SUNDIALS::IDA::reinit_vector
std::function< void(VectorType &)> reinit_vector
Definition: ida.h:626
SUNDIALS::IDA::AdditionalData::use_y_dot
@ use_y_dot
Definition: ida.h:271
SUNDIALS::IDA::solve_jacobian_system
std::function< int(const VectorType &rhs, VectorType &dst)> solve_jacobian_system
Definition: ida.h:709
SUNDIALS::IDA::AdditionalData
Definition: ida.h:241
SUNDIALS::IDA::AdditionalData::minimum_step_size
double minimum_step_size
Definition: ida.h:477
SUNDIALS::IDA::AdditionalData::relative_tolerance
double relative_tolerance
Definition: ida.h:487
SUNDIALS::IDA::output_step
std::function< void(const double t, const VectorType &sol, const VectorType &sol_dot, const unsigned int step_number)> output_step
Definition: ida.h:729
VectorType
SUNDIALS::IDA::AdditionalData::maximum_non_linear_iterations_ic
unsigned maximum_non_linear_iterations_ic
Definition: ida.h:535
MPI_Comm
Physics::Elasticity::Kinematics::e
SymmetricTensor< 2, dim, Number > e(const Tensor< 2, dim, Number > &F)
SUNDIALS::IDA::solver_should_restart
std::function< bool(const double t, VectorType &sol, VectorType &sol_dot)> solver_should_restart
Definition: ida.h:748
SUNDIALS
Definition: arkode.h:58
SUNDIALS::IDA::mem
GrowingVectorMemory< VectorType > mem
Definition: ida.h:836
SUNDIALS::IDA::~IDA
~IDA()
Definition: ida.cc:171
SUNDIALS::IDA::AdditionalData::AdditionalData
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 double absolute_tolerance=1e-6, const double relative_tolerance=1e-5, const bool ignore_algebraic_terms_for_errors=true, const InitialConditionCorrection &ic_type=use_y_diff, const InitialConditionCorrection &reset_type=use_y_diff, const unsigned int maximum_non_linear_iterations_ic=5)
Definition: ida.h:305
SUNDIALS::IDA::AdditionalData::initial_time
double initial_time
Definition: ida.h:462
GrowingVectorMemory
Definition: vector_memory.h:320
SUNDIALS::IDA::differential_components
std::function< IndexSet()> differential_components
Definition: ida.h:758
bool
conditional_ostream.h
SUNDIALS::IDA::AdditionalData::output_period
double output_period
Definition: ida.h:497
SUNDIALS::IDA::AdditionalData::final_time
double final_time
Definition: ida.h:467
mpi.h
DEAL_II_NAMESPACE_OPEN
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:358
SUNDIALS::IDA::set_functions_to_trigger_an_assert
void set_functions_to_trigger_an_assert()
Definition: ida.cc:443
SUNDIALS::IDA::AdditionalData::initial_step_size
double initial_step_size
Definition: ida.h:472
SUNDIALS::IDA::setup_jacobian
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, const double alpha)> setup_jacobian
Definition: ida.h:678
petsc_block_vector.h
SUNDIALS::IDA::yp
N_Vector yp
Definition: ida.h:814
SUNDIALS::IDA::AdditionalData::ignore_algebraic_terms_for_errors
bool ignore_algebraic_terms_for_errors
Definition: ida.h:502
SUNDIALS::IDA::DeclException1
DeclException1(ExcIDAError, int,<< "One of the SUNDIALS IDA internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")
parameter_handler.h
ParameterHandler::add_action
void add_action(const std::string &entry, const std::function< void(const std::string &value)> &action)
Definition: parameter_handler.cc:830
SUNDIALS::IDA::diff_id
N_Vector diff_id
Definition: ida.h:824
exceptions.h
SUNDIALS::IDA::AdditionalData::InitialConditionCorrection
InitialConditionCorrection
Definition: ida.h:253
SUNDIALS::IDA::AdditionalData::maximum_non_linear_iterations
unsigned int maximum_non_linear_iterations
Definition: ida.h:540
value
static const bool value
Definition: dof_tools_constraints.cc:433
StandardExceptions::ExcInternalError
static ::ExceptionBase & ExcInternalError()
SUNDIALS::IDA::data
AdditionalData data
Definition: ida.h:799
SUNDIALS::IDA::solve_dae
unsigned int solve_dae(VectorType &solution, VectorType &solution_dot)
Definition: ida.cc:189
SUNDIALS::IDA::AdditionalData::none
@ none
Definition: ida.h:258
vector.h
SUNDIALS::IDA::AdditionalData::use_y_diff
@ use_y_diff
Definition: ida.h:266
SUNDIALS::IDA::reset
void reset(const double t, const double h, VectorType &y, VectorType &yp)
Definition: ida.cc:279
ParameterHandler::add_parameter
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)
Definition: parameter_handler.h:2328
SUNDIALS::IDA::IDA
IDA(const AdditionalData &data=AdditionalData(), const MPI_Comm mpi_comm=MPI_COMM_WORLD)
Definition: ida.cc:156
petsc_vector.h
vector_memory.h
SUNDIALS::IDA::residual
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, VectorType &res)> residual
Definition: ida.h:642
ParameterHandler
Definition: parameter_handler.h:845
config.h
SUNDIALS::IDA::abs_tolls
N_Vector abs_tolls
Definition: ida.h:819
ParameterHandler::enter_subsection
void enter_subsection(const std::string &subsection)
Definition: parameter_handler.cc:927
SUNDIALS::IDA::AdditionalData::absolute_tolerance
double absolute_tolerance
Definition: ida.h:482
SUNDIALS::IDA::AdditionalData::add_parameters
void add_parameters(ParameterHandler &prm)
Definition: ida.h:380
SUNDIALS::IDA::get_local_tolerances
std::function< VectorType &()> get_local_tolerances
Definition: ida.h:766
SUNDIALS::IDA::ida_mem
void * ida_mem
Definition: ida.h:804
Patterns::Selection
Definition: patterns.h:383
DEAL_II_NAMESPACE_CLOSE
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:359
ParameterHandler::leave_subsection
void leave_subsection()
Definition: parameter_handler.cc:941
logstream.h
AssertThrow
#define AssertThrow(cond, exc)
Definition: exceptions.h:1531
SUNDIALS::IDA
Definition: ida.h:235
int
SUNDIALS::IDA::yy
N_Vector yy
Definition: ida.h:809