Cantera  4.0.0a1
Loading...
Searching...
No Matches
Integrator.h
Go to the documentation of this file.
1/**
2 * @file Integrator.h
3 */
4
5// This file is part of Cantera. See License.txt in the top-level directory or
6// at https://cantera.org/license.txt for license and copyright information.
7
8#ifndef CT_INTEGRATOR_H
9#define CT_INTEGRATOR_H
10#include "FuncEval.h"
11#include "SystemJacobian.h"
12
13#include "cantera/base/global.h"
14#include "cantera/base/AnyMap.h"
15
16namespace Cantera
17{
18
19/**
20 * Specifies the method used to integrate the system of equations.
21 * Not all methods are supported by all integrators.
22 */
24 BDF_Method, //!< Backward Differentiation
25 Adams_Method //! Adams
26};
27
28//! Specifies the method used for iteration.
29/*!
30 * Not all methods are supported by all integrators.
31 */
33 //! Newton Iteration
35 //! Functional Iteration
37};
38
39//! Abstract base class for ODE system integrators.
40/*!
41 * @ingroup odeGroup
42 */
44{
45public:
46 //! Default Constructor
48 }
49
50 //! Destructor
51 virtual ~Integrator() {
52 }
53
54 //! Set error tolerances.
55 /*!
56 * @param reltol scalar relative tolerance
57 * @param abstol vector of absolute tolerance values, length nEquations().
58 */
59 virtual void setTolerances(double reltol, span<const double> abstol) {
60 warn("setTolerances");
61 }
62
63 //! Set error tolerances.
64 /*!
65 * @param reltol scalar relative tolerance
66 * @param abstol scalar absolute tolerance
67 */
68 virtual void setTolerances(double reltol, double abstol) {
69 warn("setTolerances");
70 }
71
72 //! Set the sensitivity error tolerances
73 /*!
74 * @param reltol scalar relative tolerance
75 * @param abstol scalar absolute tolerance
76 */
77 virtual void setSensitivityTolerances(double reltol, double abstol)
78 { }
79
80 //! Set the linear solver type.
81 /*!
82 * @param linSolverType Type of the linear solver
83 */
84 virtual void setLinearSolverType(const string& linSolverType) {
85 warn("setLinearSolverType");
86 }
87
88 //! Configure how many event/root functions the integrator should monitor
89 /**
90 * Callers toggle this to enable/disable root finding on the fly (when
91 * ReactorNet enforces advance limits).
92 * @param nroots Number of root functions
93 */
94 virtual void setRootFunctionCount([[maybe_unused]] size_t nroots) {}
95
96 //! Current value of the independent variable tracked by the integrator
97 virtual double currentTime() const = 0;
98
99 //! Set preconditioner used by the linear solver
100 /*!
101 * @param preconditioner preconditioner object used for the linear solver
102 */
103 virtual void setPreconditioner(shared_ptr<SystemJacobian> preconditioner) {
105 if (preconditioner->preconditionerSide() == "none") {
106 m_prec_side = PreconditionerSide::NO_PRECONDITION;
107 } else if (preconditioner->preconditionerSide() == "left") {
109 } else if (preconditioner->preconditionerSide() == "right") {
111 } else if (preconditioner->preconditionerSide() == "both") {
113 } else {
114 throw CanteraError("Integrator::setPreconditioner",
115 "Invalid option '{}'", preconditioner->preconditionerSide());
116 }
117 }
118
119 //! Solve a linear system Ax=b where A is the preconditioner
120 /*!
121 * @param[in] rhs right hand side vector used in linear system
122 * @param[out] output output vector for solution
123 */
124 virtual void preconditionerSolve(span<const double> rhs, span<double> output) {
125 m_preconditioner->solve(rhs, output);
126 }
127
128 //! Return the side of the system on which the preconditioner is applied
130 return m_prec_side;
131 }
132
133 //! Return preconditioner reference to object
134 virtual shared_ptr<SystemJacobian> preconditioner() {
135 return m_preconditioner;
136 }
137
138 //! Return the integrator problem type
139 virtual string linearSolverType() const {
140 warn("linearSolverType");
141 return "";
142 }
143
144 /**
145 * Initialize the integrator for a new problem. Call after all options have
146 * been set.
147 * @param t0 initial time
148 * @param func RHS evaluator object for system of equations.
149 */
150 virtual void initialize(double t0, FuncEval& func) {
151 warn("initialize");
152 }
153
154 virtual void reinitialize(double t0, FuncEval& func) {
155 warn("reinitialize");
156 }
157
158 //! Integrate the system of equations.
159 /*!
160 * @param tout Integrate to this time. Note that this is the
161 * absolute time value, not a time interval.
162 */
163 virtual void integrate(double tout) {
164 warn("integrate");
165 }
166
167 /**
168 * Integrate the system of equations.
169 * @param tout integrate to this time. Note that this is the
170 * absolute time value, not a time interval.
171 */
172 virtual double step(double tout) {
173 warn("step");
174 return 0.0;
175 }
176
177 //! The current value of the solution of equation k.
178 virtual double& solution(size_t k) {
179 warn("solution");
180 return m_dummy;
181 }
182
183 //! The current value of the solution of the system of equations.
184 virtual span<double> solution() {
185 warn("solution");
186 return {};
187 }
188
189 //! n-th derivative of the output function at time tout.
190 virtual span<double> derivative(double tout, int n) {
191 warn("derivative");
192 return {};
193 }
194
195 //! Order used during the last solution step
196 virtual int lastOrder() const {
197 warn("lastOrder");
198 return 0;
199 }
200
201 //! The number of equations.
202 virtual int nEquations() const {
203 warn("nEquations");
204 return 0;
205 }
206
207 //! The number of function evaluations.
208 virtual int nEvals() const {
209 warn("nEvals");
210 return 0;
211 }
212
213 virtual int maxOrder() const {
214 warn("maxOrder");
215 return 0;
216 }
217
218 //! Set the maximum integration order that will be used.
219 virtual void setMaxOrder(int n) {
220 warn("setMaxorder");
221 }
222
223 //! Set the solution method
224 virtual void setMethod(MethodType t) {
225 warn("setMethodType");
226 }
227
228 //! Set the maximum step size
229 virtual void setMaxStepSize(double hmax) {
230 warn("setMaxStepSize");
231 }
232
233 //! Set the minimum step size
234 virtual void setMinStepSize(double hmin) {
235 warn("setMinStepSize");
236 }
237
238 //! Set the maximum permissible number of error test failures
239 virtual void setMaxErrTestFails(int n) {
240 warn("setMaxErrTestFails");
241 }
242
243 //! Set the maximum number of time-steps the integrator can take
244 //! before reaching the next output time
245 //! @param nmax The maximum number of steps, setting this value
246 //! to zero disables this option.
247 virtual void setMaxSteps(int nmax) {
248 warn("setMaxStep");
249 }
250
251 //! Returns the maximum number of time-steps the integrator can take
252 //! before reaching the next output time
253 virtual int maxSteps() {
254 warn("maxSteps");
255 return 0;
256 }
257
258 virtual void setBandwidth(int N_Upper, int N_Lower) {
259 warn("setBandwidth");
260 }
261
262 virtual int nSensParams() {
263 warn("nSensParams()");
264 return 0;
265 }
266
267 virtual double sensitivity(size_t k, size_t p) {
268 warn("sensitivity");
269 return 0.0;
270 }
271
272 //! Get solver stats from integrator
273 virtual AnyMap solverStats() const {
274 AnyMap stats;
275 warn("solverStats");
276 return stats;
277 }
278
279 virtual int maxNonlinIterations() const {
280 warn("maxNonlinIterations");
281 return 0;
282 }
283 virtual void setMaxNonlinIterations(int n) {
284 warn("setMaxNonlinIterations");
285 }
286
287 virtual int maxNonlinConvFailures() const {
288 warn("maxNonlinConvFailures");
289 return 0;
290 }
291 virtual void setMaxNonlinConvFailures(int n) {
292 warn("setMaxNonlinConvFailures");
293 }
294
295 virtual bool algebraicInErrorTest() const {
296 warn("algebraicInErrorTest");
297 return true;
298 }
299 virtual void includeAlgebraicInErrorTest(bool yesno) {
300 warn("includeAlgebraicInErrorTest");
301 }
302
303protected:
304 //! Pointer to preconditioner object used in integration which is
305 //! set by setPreconditioner and initialized inside of
306 //! ReactorNet::initialize()
307 shared_ptr<SystemJacobian> m_preconditioner;
308 //! Type of preconditioning used in applyOptions
309 PreconditionerSide m_prec_side = PreconditionerSide::NO_PRECONDITION;
310 // methods for DAE solvers
311
312private:
313 double m_dummy;
314 void warn(const string& msg) const {
315 writelog(">>>> Warning: method "+msg+" of base class "
316 +"Integrator called. Nothing done.\n");
317 }
318};
319
320// defined in Integrators.cpp
321
322//! Create new Integrator object
323//! @param itype Integration mode; either @c CVODE or @c IDA
324//! @ingroup odeGroup
325Integrator* newIntegrator(const string& itype);
326
327} // namespace
328
329#endif
Declarations for class SystemJacobian.
A map of string keys to values whose type can vary at runtime.
Definition AnyMap.h:431
Base class for exceptions thrown by Cantera classes.
Virtual base class for ODE/DAE right-hand-side function evaluators.
Definition FuncEval.h:32
Abstract base class for ODE system integrators.
Definition Integrator.h:44
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
Definition Integrator.h:229
virtual void setMaxSteps(int nmax)
Set the maximum number of time-steps the integrator can take before reaching the next output time.
Definition Integrator.h:247
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
Definition Integrator.h:253
shared_ptr< SystemJacobian > m_preconditioner
Pointer to preconditioner object used in integration which is set by setPreconditioner and initialize...
Definition Integrator.h:307
virtual span< double > derivative(double tout, int n)
n-th derivative of the output function at time tout.
Definition Integrator.h:190
Integrator()
Default Constructor.
Definition Integrator.h:47
virtual void preconditionerSolve(span< const double > rhs, span< double > output)
Solve a linear system Ax=b where A is the preconditioner.
Definition Integrator.h:124
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
Definition Integrator.h:239
virtual void setRootFunctionCount(size_t nroots)
Configure how many event/root functions the integrator should monitor.
Definition Integrator.h:94
virtual AnyMap solverStats() const
Get solver stats from integrator.
Definition Integrator.h:273
virtual void initialize(double t0, FuncEval &func)
Initialize the integrator for a new problem.
Definition Integrator.h:150
virtual void setTolerances(double reltol, span< const double > abstol)
Set error tolerances.
Definition Integrator.h:59
virtual void setLinearSolverType(const string &linSolverType)
Set the linear solver type.
Definition Integrator.h:84
virtual double currentTime() const =0
Current value of the independent variable tracked by the integrator.
virtual PreconditionerSide preconditionerSide()
Return the side of the system on which the preconditioner is applied.
Definition Integrator.h:129
virtual void setTolerances(double reltol, double abstol)
Set error tolerances.
Definition Integrator.h:68
virtual void setMinStepSize(double hmin)
Set the minimum step size.
Definition Integrator.h:234
virtual double step(double tout)
Integrate the system of equations.
Definition Integrator.h:172
virtual span< double > solution()
The current value of the solution of the system of equations.
Definition Integrator.h:184
virtual void setSensitivityTolerances(double reltol, double abstol)
Set the sensitivity error tolerances.
Definition Integrator.h:77
virtual string linearSolverType() const
Return the integrator problem type.
Definition Integrator.h:139
virtual void integrate(double tout)
Integrate the system of equations.
Definition Integrator.h:163
PreconditionerSide m_prec_side
Type of preconditioning used in applyOptions.
Definition Integrator.h:309
virtual int lastOrder() const
Order used during the last solution step.
Definition Integrator.h:196
virtual int nEvals() const
The number of function evaluations.
Definition Integrator.h:208
virtual int nEquations() const
The number of equations.
Definition Integrator.h:202
virtual void setPreconditioner(shared_ptr< SystemJacobian > preconditioner)
Set preconditioner used by the linear solver.
Definition Integrator.h:103
virtual ~Integrator()
Destructor.
Definition Integrator.h:51
virtual shared_ptr< SystemJacobian > preconditioner()
Return preconditioner reference to object.
Definition Integrator.h:134
virtual void setMethod(MethodType t)
Set the solution method.
Definition Integrator.h:224
virtual double & solution(size_t k)
The current value of the solution of equation k.
Definition Integrator.h:178
virtual void setMaxOrder(int n)
Set the maximum integration order that will be used.
Definition Integrator.h:219
This file contains definitions for utility functions and text for modules, inputfiles and logging,...
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition global.h:171
Integrator * newIntegrator(const string &itype)
Create new Integrator object.
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
MethodType
Specifies the method used to integrate the system of equations.
Definition Integrator.h:23
@ BDF_Method
Backward Differentiation.
Definition Integrator.h:24
@ Adams_Method
Adams.
Definition Integrator.h:25
PreconditionerSide
Specifies the side of the system on which the preconditioner is applied.
@ LEFT_PRECONDITION
No preconditioning.
@ BOTH_PRECONDITION
Right side preconditioning.
@ RIGHT_PRECONDITION
Left side preconditioning.
IterType
Specifies the method used for iteration.
Definition Integrator.h:32
@ Functional_Iter
Functional Iteration.
Definition Integrator.h:36
@ Newton_Iter
Newton Iteration.
Definition Integrator.h:34