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