Reference documentation for deal.II version 9.4.1
\(\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\}}\)
Loading...
Searching...
No Matches
ida.h
Go to the documentation of this file.
1// ---------------------------------------------------------------------
2//
3// Copyright (C) 2017 - 2022 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
17#ifndef dealii_sundials_ida_h
18#define dealii_sundials_ida_h
19
20#include <deal.II/base/config.h>
21
22#include <deal.II/base/mpi.h>
23#ifdef DEAL_II_WITH_SUNDIALS
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
42
43# include <boost/signals2.hpp>
44
45# include <nvector/nvector_serial.h>
46# include <sundials/sundials_config.h>
47# include <sundials/sundials_math.h>
48# include <sundials/sundials_types.h>
49
50# include <memory>
51
52
54
55// Shorthand notation for IDA error codes.
56# define AssertIDA(code) Assert(code >= 0, ExcIDAError(code))
57
58namespace SUNDIALS
59{
233 template <typename VectorType = Vector<double>>
234 class IDA
235 {
236 public:
241 {
242 public:
253 {
257 none = 0,
258
266
270 use_y_dot = 2
271 };
272
307 AdditionalData( // Initial parameters
308 const double initial_time = 0.0,
309 const double final_time = 1.0,
310 const double initial_step_size = 1e-2,
311 const double output_period = 1e-1,
312 // Running parameters
313 const double minimum_step_size = 1e-6,
314 const unsigned int maximum_order = 5,
315 const unsigned int maximum_non_linear_iterations = 10,
316 const double ls_norm_factor = 0,
317 // Error parameters
318 const double absolute_tolerance = 1e-6,
319 const double relative_tolerance = 1e-5,
320 const bool ignore_algebraic_terms_for_errors = true,
321 // Initial conditions parameters
324 const unsigned int maximum_non_linear_iterations_ic = 5)
339 {}
340
383 void
385 {
386 prm.add_parameter("Initial time", initial_time);
387 prm.add_parameter("Final time", final_time);
388 prm.add_parameter("Time interval between each output", output_period);
389
390 prm.enter_subsection("Running parameters");
391 prm.add_parameter("Initial step size", initial_step_size);
392 prm.add_parameter("Minimum step size", minimum_step_size);
393 prm.add_parameter("Maximum order of BDF", maximum_order);
394 prm.add_parameter("Maximum number of nonlinear iterations",
396 prm.leave_subsection();
397
398 prm.enter_subsection("Error control");
399 prm.add_parameter("Absolute error tolerance", absolute_tolerance);
400 prm.add_parameter("Relative error tolerance", relative_tolerance);
401 prm.add_parameter(
402 "Ignore algebraic terms for error computations",
404 "Indicate whether or not to suppress algebraic variables "
405 "in the local error test.");
406 prm.leave_subsection();
407
408 prm.enter_subsection("Initial condition correction parameters");
409 static std::string ic_type_str = "use_y_diff";
410 prm.add_parameter(
411 "Correction type at initial time",
412 ic_type_str,
413 "This is one of the following three options for the "
414 "initial condition calculation. \n"
415 " none: do not try to make initial conditions consistent. \n"
416 " use_y_diff: compute the algebraic components of y and differential\n"
417 " components of y_dot, given the differential components of y. \n"
418 " This option requires that the user specifies differential and \n"
419 " algebraic components in the function get_differential_components.\n"
420 " use_y_dot: compute all components of y, given y_dot.",
421 Patterns::Selection("none|use_y_diff|use_y_dot"));
422 prm.add_action("Correction type at initial time",
423 [&](const std::string &value) {
424 if (value == "use_y_diff")
426 else if (value == "use_y_dot")
428 else if (value == "none")
429 ic_type = none;
430 else
432 });
433
434 static std::string reset_type_str = "use_y_diff";
435 prm.add_parameter(
436 "Correction type after restart",
437 reset_type_str,
438 "This is one of the following three options for the "
439 "initial condition calculation. \n"
440 " none: do not try to make initial conditions consistent. \n"
441 " use_y_diff: compute the algebraic components of y and differential\n"
442 " components of y_dot, given the differential components of y. \n"
443 " This option requires that the user specifies differential and \n"
444 " algebraic components in the function get_differential_components.\n"
445 " use_y_dot: compute all components of y, given y_dot.",
446 Patterns::Selection("none|use_y_diff|use_y_dot"));
447 prm.add_action("Correction type after restart",
448 [&](const std::string &value) {
449 if (value == "use_y_diff")
451 else if (value == "use_y_dot")
453 else if (value == "none")
455 else
457 });
458 prm.add_parameter("Maximum number of nonlinear iterations",
460 prm.add_parameter(
461 "Factor to use when converting from the integrator tolerance to the linear solver tolerance",
463 prm.leave_subsection();
464 }
465
470
475
480
485
490
495
499 unsigned int maximum_order;
500
505
510
526
538
543
548
554 };
555
599
607 IDA(const AdditionalData &data, const MPI_Comm &mpi_comm);
608
612 ~IDA();
613
618 unsigned int
619 solve_dae(VectorType &solution, VectorType &solution_dot);
620
642 void
643 reset(const double t, const double h, VectorType &y, VectorType &yp);
644
648 std::function<void(VectorType &)> reinit_vector;
649
660 std::function<int(const double t,
661 const VectorType &y,
662 const VectorType &y_dot,
663 VectorType & res)>
665
699 std::function<int(const double t,
700 const VectorType &y,
701 const VectorType &y_dot,
702 const double alpha)>
704
737 std::function<int(const VectorType &rhs, VectorType &dst)>
739
781 std::function<
782 int(const VectorType &rhs, VectorType &dst, const double tolerance)>
784
799 std::function<void(const double t,
800 const VectorType & sol,
801 const VectorType & sol_dot,
802 const unsigned int step_number)>
804
821 std::function<bool(const double t, VectorType &sol, VectorType &sol_dot)>
823
839
846 std::function<VectorType &()> get_local_tolerances;
847
851 DeclException1(ExcIDAError,
852 int,
853 << "One of the SUNDIALS IDA internal functions "
854 << " returned a negative error code: " << arg1
855 << ". Please consult SUNDIALS manual.");
856
857
858 private:
863 DeclException1(ExcFunctionNotProvided,
864 std::string,
865 << "Please provide an implementation for the function \""
866 << arg1 << "\"");
867
873 void
875
880
884 void *ida_mem;
885
886# if DEAL_II_SUNDIALS_VERSION_GTE(6, 0, 0)
890 SUNContext ida_ctx;
891# endif
892
898
903
904# ifdef DEAL_II_WITH_PETSC
905# ifdef PETSC_USE_COMPLEX
906 static_assert(!std::is_same<VectorType, PETScWrappers::MPI::Vector>::value,
907 "Sundials does not support complex scalar types, "
908 "but PETSc is configured to use a complex scalar type!");
909
910 static_assert(
911 !std::is_same<VectorType, PETScWrappers::MPI::BlockVector>::value,
912 "Sundials does not support complex scalar types, "
913 "but PETSc is configured to use a complex scalar type!");
914# endif // PETSC_USE_COMPLEX
915# endif // DEAL_II_WITH_PETSC
916 };
917} // namespace SUNDIALS
918
920
921#endif // DEAL_II_WITH_SUNDIALS
922
923#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)
InitialConditionCorrection ic_type
Definition: ida.h:525
bool ignore_algebraic_terms_for_errors
Definition: ida.h:509
void add_parameters(ParameterHandler &prm)
Definition: ida.h:384
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 ls_norm_factor=0, 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:307
unsigned maximum_non_linear_iterations_ic
Definition: ida.h:542
InitialConditionCorrection reset_type
Definition: ida.h:537
unsigned int maximum_order
Definition: ida.h:499
unsigned int maximum_non_linear_iterations
Definition: ida.h:547
DEAL_II_DEPRECATED std::function< int(const VectorType &rhs, VectorType &dst)> solve_jacobian_system
Definition: ida.h:738
unsigned int solve_dae(VectorType &solution, VectorType &solution_dot)
Definition: ida.cc:172
MPI_Comm mpi_communicator
Definition: ida.h:897
void set_functions_to_trigger_an_assert()
Definition: ida.cc:466
const AdditionalData data
Definition: ida.h:879
void reset(const double t, const double h, VectorType &y, VectorType &yp)
Definition: ida.cc:229
std::function< IndexSet()> differential_components
Definition: ida.h:838
DeclException1(ExcIDAError, int,<< "One of the SUNDIALS IDA internal functions "<< " returned a negative error code: "<< arg1<< ". Please consult SUNDIALS manual.")
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, VectorType &res)> residual
Definition: ida.h:664
std::function< VectorType &()> get_local_tolerances
Definition: ida.h:846
void * ida_mem
Definition: ida.h:884
std::function< void(const double t, const VectorType &sol, const VectorType &sol_dot, const unsigned int step_number)> output_step
Definition: ida.h:803
std::function< bool(const double t, VectorType &sol, VectorType &sol_dot)> solver_should_restart
Definition: ida.h:822
std::function< int(const double t, const VectorType &y, const VectorType &y_dot, const double alpha)> setup_jacobian
Definition: ida.h:703
DeclException1(ExcFunctionNotProvided, std::string,<< "Please provide an implementation for the function \""<< arg1<< "\"")
GrowingVectorMemory< VectorType > mem
Definition: ida.h:902
std::function< void(VectorType &)> reinit_vector
Definition: ida.h:648
std::function< int(const VectorType &rhs, VectorType &dst, const double tolerance)> solve_with_jacobian
Definition: ida.h:783
SUNContext ida_ctx
Definition: ida.h:890
#define DEAL_II_DEPRECATED
Definition: config.h:164
#define DEAL_II_NAMESPACE_OPEN
Definition: config.h:442
#define DEAL_II_NAMESPACE_CLOSE
Definition: config.h:443
static ::ExceptionBase & ExcInternalError()
#define AssertThrow(cond, exc)
Definition: exceptions.h:1583