Reaktoro
A unified framework for modeling chemically reactive systems
ODE.hpp
1 // Reaktoro is a unified framework for modeling chemically reactive systems.
2 //
3 // Copyright (C) 2014-2018 Allan Leal
4 //
5 // This library is free software; you can redistribute it and/or
6 // modify it under the terms of the GNU Lesser General Public
7 // License as published by the Free Software Foundation; either
8 // version 2.1 of the License, or (at your option) any later version.
9 //
10 // This library is distributed in the hope that it will be useful,
11 // but WITHOUT ANY WARRANTY; without even the implied warranty of
12 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
13 // Lesser General Public License for more details.
14 //
15 // You should have received a copy of the GNU Lesser General Public License
16 // along with this library. If not, see <http://www.gnu.org/licenses/>.
17 
18 #pragma once
19 
20 // C++ includes
21 #include <functional>
22 #include <memory>
23 
24 // Reaktoro includes
25 #include <Reaktoro/Math/Matrix.hpp>
26 
27 namespace Reaktoro {
28 
30 using ODEFunction = std::function<int(double, VectorConstRef, VectorRef)>;
31 
33 using ODEJacobian = std::function<int(double, VectorConstRef, MatrixRef)>;
34 
36 enum class ODEStepMode { Adams, BDF };
37 
39 enum class ODEIterationMode { Functional, Newton };
40 
43 struct ODEOptions
44 {
46  ODEStepMode step = ODEStepMode::BDF;
47 
49  ODEIterationMode iteration = ODEIterationMode::Newton;
50 
57 
60  double initial_step = 0.0;
61 
64  double stop_time = 0.0;
65 
67  double min_step = 0.0;
68 
70  double max_step = 0.0;
71 
73  double reltol = 1e-4;
74 
76  double abstol = 1.0e-6;
77 
80  unsigned max_order_bdf = 5;
81 
84  unsigned max_order_adams = 5;
85 
87  unsigned max_num_steps = 500;
88 
90  unsigned max_hnil_warnings = 10;
91 
94 
97 
100 
103 
105  Vector abstols;
106 };
107 
111 {
112 public:
114  ODEProblem();
115 
117  ODEProblem(const ODEProblem& other);
118 
120  virtual ~ODEProblem();
121 
123  auto operator=(ODEProblem other) -> ODEProblem&;
124 
126  auto setNumEquations(unsigned num) -> void;
127 
129  auto setFunction(const ODEFunction& f) -> void;
130 
132  auto setJacobian(const ODEJacobian& J) -> void;
133 
135  auto initialized() const -> bool;
136 
138  auto numEquations() const -> unsigned;
139 
141  auto function() const -> const ODEFunction&;
142 
144  auto jacobian() const -> const ODEJacobian&;
145 
151  auto function(double t, VectorConstRef y, VectorRef f) const -> int;
152 
158  auto jacobian(double t, VectorConstRef y, MatrixRef J) const -> int;
159 
160 private:
161  struct Impl;
162 
163  std::unique_ptr<Impl> pimpl;
164 };
165 
169 {
170 public:
172  ODESolver();
173 
175  ODESolver(const ODESolver& other) = delete;
176 
178  virtual ~ODESolver();
179 
181  auto operator=(ODESolver other) -> ODESolver&;
182 
185  auto setOptions(const ODEOptions& options) -> void;
186 
189  auto setProblem(const ODEProblem& problem) -> void;
190 
195  auto initialize(double tstart, VectorConstRef y) -> void;
196 
200  auto integrate(double& t, VectorRef y) -> void;
201 
206  auto integrate(double& t, VectorRef y, double tfinal) -> void;
207 
212  auto solve(double& t, double dt, VectorRef y) -> void;
213 
214 private:
215  struct Impl;
216 
217  std::unique_ptr<Impl> pimpl;
218 };
219 
220 } // namespace Reaktoro
Eigen::Ref< Eigen::MatrixXd > MatrixRef
Alias to Eigen type Ref<MatrixXd>.
Definition: Matrix.hpp:43
auto initialize(double tstart, VectorConstRef y) -> void
Initializes the ODE solver.
Definition: ODE.cpp:404
std::function< int(double, VectorConstRef, VectorRef)> ODEFunction
The function signature of the right-hand side function of a system of ordinary differential equations...
Definition: ODE.hpp:30
ODEStepMode
The linear multistep method to be used in ODESolver.
Definition: ODE.hpp:36
double abstol
The scalar absolute error tolerance.
Definition: ODE.hpp:76
unsigned max_num_steps
The maximum allowed number of steps before reaching the final time.
Definition: ODE.hpp:87
unsigned max_hnil_warnings
The maximum number of warnings for t + h = t, with h being too small compared to t.
Definition: ODE.hpp:90
auto setNumEquations(unsigned num) -> void
Set the number of ordinary differential equations.
Definition: ODE.cpp:336
unsigned max_num_nonlinear_iterations
The maximum number of nonlinear iterations.
Definition: ODE.hpp:96
unsigned max_order_bdf
The maximum order for the BDF integration scheme.
Definition: ODE.hpp:80
ODEProblem()
Construct a default ODEProblem instance.
Definition: ODE.cpp:319
A struct that defines the options for the ODESolver.
Definition: ODE.hpp:44
A class that defines a system of ordinary differential equations (ODE) problem.
Definition: ODE.hpp:111
auto operator=(ODEProblem other) -> ODEProblem &
Assign a ODEProblem instance to this instance.
Definition: ODE.cpp:330
double reltol
The scalar relative error tolerance.
Definition: ODE.hpp:73
auto setFunction(const ODEFunction &f) -> void
Set the right-hand side function of the system of ordinary differential equations.
Definition: ODE.cpp:341
auto setProblem(const ODEProblem &problem) -> void
Set the ODE problem.
Definition: ODE.cpp:399
A wrapper class for CVODE, a library for solving ordinary differential equations.
Definition: ODE.hpp:169
The namespace containing all components of the Reaktoro library.
Definition: ChemicalScalar.hpp:24
auto solve(double &t, double dt, VectorRef y) -> void
Solve the ODE equations from a given start time to a final one.
Definition: ODE.cpp:419
auto integrate(double &t, VectorRef y) -> void
Integrate the ODE performing a single step.
Definition: ODE.cpp:409
std::function< int(double, VectorConstRef, MatrixRef)> ODEJacobian
The function signature of the Jacobian of the right-hand side function of a system of ordinary differ...
Definition: ODE.hpp:33
double stop_time
The value of the independent variable t past which the solution is not to proceed.
Definition: ODE.hpp:64
auto jacobian() const -> const ODEJacobian &
Return the Jacobian of the right-hand side function of the system of ordinary differential equations.
Definition: ODE.cpp:366
auto initialized() const -> bool
Return true if the problem has bee initialized.
Definition: ODE.cpp:351
unsigned max_order_adams
The maximum order for the Adams integration scheme.
Definition: ODE.hpp:84
double min_step
The lower bound on the magnitude of the step size.
Definition: ODE.hpp:67
virtual ~ODESolver()
Destroy this ODESolver instance.
Definition: ODE.cpp:385
ODESolver()
Construct a default ODESolver instance.
Definition: ODE.cpp:381
ODEIterationMode iteration
The type of nonlinear solver iteration used in the integration.
Definition: ODE.hpp:49
auto numEquations() const -> unsigned
Return the number of ordinary differential equations.
Definition: ODE.cpp:356
unsigned max_num_convergence_failures
The maximum number of convergence failures.
Definition: ODE.hpp:99
Vector abstols
The vector of absolute error tolerances for each component.
Definition: ODE.hpp:105
auto operator=(ODESolver other) -> ODESolver &
Assign a ODESolver instance to this instance.
Definition: ODE.cpp:388
bool stability_limit_detection
The flag that enables the STAbility Limit Detection (STALD) algorithm.
Definition: ODE.hpp:56
ODEStepMode step
The linear multistep method used in the integration.
Definition: ODE.hpp:46
double initial_step
The initial step size to be used in the integration.
Definition: ODE.hpp:60
virtual ~ODEProblem()
Destroy this ODEProblem instance.
Definition: ODE.cpp:327
unsigned max_num_error_test_failures
The maximum number of error test failures.
Definition: ODE.hpp:93
ODESolver(const ODESolver &other)=delete
Construct a copy of an ODESolver instance.
ODEIterationMode
The type of nonlinear solver iteration to be used in ODESolver.
Definition: ODE.hpp:39
double nonlinear_convergence_coefficient
The coefficient in the nonlinear convergence test.
Definition: ODE.hpp:102
Eigen::Ref< Eigen::VectorXd > VectorRef
< Alias to Eigen type VectorXd.
Definition: Matrix.hpp:29
auto setOptions(const ODEOptions &options) -> void
Set the options for the ODE solver.
Definition: ODE.cpp:394
auto setJacobian(const ODEJacobian &J) -> void
Set the Jacobian of the right-hand side function of the system of ordinary differential equations.
Definition: ODE.cpp:346
Eigen::Ref< const Eigen::VectorXd > VectorConstRef
< Alias to Eigen type Ref<VectorXd>.
Definition: Matrix.hpp:31
double max_step
The upper bound on the magnitude of the step size.
Definition: ODE.hpp:70