Cantera  3.2.0
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 n Number of equations
58 * @param abstol array of N absolute tolerance values
59 */
60 virtual void setTolerances(double reltol, size_t n,
61 double* abstol) {
62 warn("setTolerances");
63 }
64
65 //! Set error tolerances.
66 /*!
67 * @param reltol scalar relative tolerance
68 * @param abstol scalar absolute tolerance
69 */
70 virtual void setTolerances(double reltol, double abstol) {
71 warn("setTolerances");
72 }
73
74 //! Set the sensitivity error tolerances
75 /*!
76 * @param reltol scalar relative tolerance
77 * @param abstol scalar absolute tolerance
78 */
79 virtual void setSensitivityTolerances(double reltol, double abstol)
80 { }
81
82 //! Set the linear solver type.
83 /*!
84 * @param linSolverType Type of the linear solver
85 */
86 virtual void setLinearSolverType(const string& linSolverType) {
87 warn("setLinearSolverType");
88 }
89
90 //! Configure how many event/root functions the integrator should monitor
91 /**
92 * Callers toggle this to enable/disable root finding on the fly (when
93 * ReactorNet enforces advance limits).
94 * @param nroots Number of root functions
95 */
96 virtual void setRootFunctionCount([[maybe_unused]] size_t nroots) {}
97
98 //! Current value of the independent variable tracked by the integrator
99 virtual double currentTime() const = 0;
100
101 //! Set preconditioner used by the linear solver
102 /*!
103 * @param preconditioner preconditioner object used for the linear solver
104 */
105 virtual void setPreconditioner(shared_ptr<SystemJacobian> preconditioner) {
107 if (preconditioner->preconditionerSide() == "none") {
108 m_prec_side = PreconditionerSide::NO_PRECONDITION;
109 } else if (preconditioner->preconditionerSide() == "left") {
111 } else if (preconditioner->preconditionerSide() == "right") {
113 } else if (preconditioner->preconditionerSide() == "both") {
115 } else {
116 throw CanteraError("Integrator::setPreconditioner",
117 "Invalid option '{}'", preconditioner->preconditionerSide());
118 }
119 }
120
121 //! Solve a linear system Ax=b where A is the preconditioner
122 /*!
123 * @param[in] stateSize length of the rhs and output vectors
124 * @param[in] rhs right hand side vector used in linear system
125 * @param[out] output output vector for solution
126 */
127 virtual void preconditionerSolve(size_t stateSize, double* rhs, double* output) {
128 m_preconditioner->solve(stateSize, rhs, output);
129 }
130
131 //! Return the side of the system on which the preconditioner is applied
133 return m_prec_side;
134 }
135
136 //! Return preconditioner reference to object
137 virtual shared_ptr<SystemJacobian> preconditioner() {
138 return m_preconditioner;
139 }
140
141 //! Return the integrator problem type
142 virtual string linearSolverType() const {
143 warn("linearSolverType");
144 return "";
145 }
146
147 /**
148 * Initialize the integrator for a new problem. Call after all options have
149 * been set.
150 * @param t0 initial time
151 * @param func RHS evaluator object for system of equations.
152 */
153 virtual void initialize(double t0, FuncEval& func) {
154 warn("initialize");
155 }
156
157 virtual void reinitialize(double t0, FuncEval& func) {
158 warn("reinitialize");
159 }
160
161 //! Integrate the system of equations.
162 /*!
163 * @param tout Integrate to this time. Note that this is the
164 * absolute time value, not a time interval.
165 */
166 virtual void integrate(double tout) {
167 warn("integrate");
168 }
169
170 /**
171 * Integrate the system of equations.
172 * @param tout integrate to this time. Note that this is the
173 * absolute time value, not a time interval.
174 */
175 virtual double step(double tout) {
176 warn("step");
177 return 0.0;
178 }
179
180 //! The current value of the solution of equation k.
181 virtual double& solution(size_t k) {
182 warn("solution");
183 return m_dummy;
184 }
185
186 //! The current value of the solution of the system of equations.
187 virtual double* solution() {
188 warn("solution");
189 return 0;
190 }
191
192 //! n-th derivative of the output function at time tout.
193 virtual double* derivative(double tout, int n) {
194 warn("derivative");
195 return 0;
196 }
197
198 //! Order used during the last solution step
199 virtual int lastOrder() const {
200 warn("lastOrder");
201 return 0;
202 }
203
204 //! The number of equations.
205 virtual int nEquations() const {
206 warn("nEquations");
207 return 0;
208 }
209
210 //! The number of function evaluations.
211 virtual int nEvals() const {
212 warn("nEvals");
213 return 0;
214 }
215
216 virtual int maxOrder() const {
217 warn("maxOrder");
218 return 0;
219 }
220
221 //! Set the maximum integration order that will be used.
222 virtual void setMaxOrder(int n) {
223 warn("setMaxorder");
224 }
225
226 //! Set the solution method
227 virtual void setMethod(MethodType t) {
228 warn("setMethodType");
229 }
230
231 //! Set the maximum step size
232 virtual void setMaxStepSize(double hmax) {
233 warn("setMaxStepSize");
234 }
235
236 //! Set the minimum step size
237 virtual void setMinStepSize(double hmin) {
238 warn("setMinStepSize");
239 }
240
241 //! Set the maximum permissible number of error test failures
242 virtual void setMaxErrTestFails(int n) {
243 warn("setMaxErrTestFails");
244 }
245
246 //! Set the maximum number of time-steps the integrator can take
247 //! before reaching the next output time
248 //! @param nmax The maximum number of steps, setting this value
249 //! to zero disables this option.
250 virtual void setMaxSteps(int nmax) {
251 warn("setMaxStep");
252 }
253
254 //! Returns the maximum number of time-steps the integrator can take
255 //! before reaching the next output time
256 virtual int maxSteps() {
257 warn("maxSteps");
258 return 0;
259 }
260
261 virtual void setBandwidth(int N_Upper, int N_Lower) {
262 warn("setBandwidth");
263 }
264
265 virtual int nSensParams() {
266 warn("nSensParams()");
267 return 0;
268 }
269
270 virtual double sensitivity(size_t k, size_t p) {
271 warn("sensitivity");
272 return 0.0;
273 }
274
275 //! Get solver stats from integrator
276 virtual AnyMap solverStats() const {
277 AnyMap stats;
278 warn("solverStats");
279 return stats;
280 }
281
282 virtual int maxNonlinIterations() const {
283 warn("maxNonlinIterations");
284 return 0;
285 }
286 virtual void setMaxNonlinIterations(int n) {
287 warn("setMaxNonlinIterations");
288 }
289
290 virtual int maxNonlinConvFailures() const {
291 warn("maxNonlinConvFailures");
292 return 0;
293 }
294 virtual void setMaxNonlinConvFailures(int n) {
295 warn("setMaxNonlinConvFailures");
296 }
297
298 virtual bool algebraicInErrorTest() const {
299 warn("algebraicInErrorTest");
300 return true;
301 }
302 virtual void includeAlgebraicInErrorTest(bool yesno) {
303 warn("includeAlgebraicInErrorTest");
304 }
305
306protected:
307 //! Pointer to preconditioner object used in integration which is
308 //! set by setPreconditioner and initialized inside of
309 //! ReactorNet::initialize()
310 shared_ptr<SystemJacobian> m_preconditioner;
311 //! Type of preconditioning used in applyOptions
312 PreconditionerSide m_prec_side = PreconditionerSide::NO_PRECONDITION;
313 // methods for DAE solvers
314
315private:
316 double m_dummy;
317 void warn(const string& msg) const {
318 writelog(">>>> Warning: method "+msg+" of base class "
319 +"Integrator called. Nothing done.\n");
320 }
321};
322
323// defined in Integrators.cpp
324
325//! Create new Integrator object
326//! @param itype Integration mode; either @c CVODE or @c IDA
327//! @ingroup odeGroup
328Integrator* newIntegrator(const string& itype);
329
330} // namespace
331
332#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:232
virtual void setTolerances(double reltol, size_t n, double *abstol)
Set error tolerances.
Definition Integrator.h:60
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:250
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
Definition Integrator.h:256
shared_ptr< SystemJacobian > m_preconditioner
Pointer to preconditioner object used in integration which is set by setPreconditioner and initialize...
Definition Integrator.h:310
Integrator()
Default Constructor.
Definition Integrator.h:47
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
Definition Integrator.h:242
virtual void setRootFunctionCount(size_t nroots)
Configure how many event/root functions the integrator should monitor.
Definition Integrator.h:96
virtual double * derivative(double tout, int n)
n-th derivative of the output function at time tout.
Definition Integrator.h:193
virtual AnyMap solverStats() const
Get solver stats from integrator.
Definition Integrator.h:276
virtual void initialize(double t0, FuncEval &func)
Initialize the integrator for a new problem.
Definition Integrator.h:153
virtual void setLinearSolverType(const string &linSolverType)
Set the linear solver type.
Definition Integrator.h:86
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:132
virtual void setTolerances(double reltol, double abstol)
Set error tolerances.
Definition Integrator.h:70
virtual void setMinStepSize(double hmin)
Set the minimum step size.
Definition Integrator.h:237
virtual double step(double tout)
Integrate the system of equations.
Definition Integrator.h:175
virtual void setSensitivityTolerances(double reltol, double abstol)
Set the sensitivity error tolerances.
Definition Integrator.h:79
virtual string linearSolverType() const
Return the integrator problem type.
Definition Integrator.h:142
virtual void integrate(double tout)
Integrate the system of equations.
Definition Integrator.h:166
PreconditionerSide m_prec_side
Type of preconditioning used in applyOptions.
Definition Integrator.h:312
virtual int lastOrder() const
Order used during the last solution step.
Definition Integrator.h:199
virtual int nEvals() const
The number of function evaluations.
Definition Integrator.h:211
virtual int nEquations() const
The number of equations.
Definition Integrator.h:205
virtual void setPreconditioner(shared_ptr< SystemJacobian > preconditioner)
Set preconditioner used by the linear solver.
Definition Integrator.h:105
virtual double * solution()
The current value of the solution of the system of equations.
Definition Integrator.h:187
virtual ~Integrator()
Destructor.
Definition Integrator.h:51
virtual shared_ptr< SystemJacobian > preconditioner()
Return preconditioner reference to object.
Definition Integrator.h:137
virtual void setMethod(MethodType t)
Set the solution method.
Definition Integrator.h:227
virtual void preconditionerSolve(size_t stateSize, double *rhs, double *output)
Solve a linear system Ax=b where A is the preconditioner.
Definition Integrator.h:127
virtual double & solution(size_t k)
The current value of the solution of equation k.
Definition Integrator.h:181
virtual void setMaxOrder(int n)
Set the maximum integration order that will be used.
Definition Integrator.h:222
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:176
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