Cantera 2.6.0
Integrator.h
Go to the documentation of this file.
1/**
2 * @file Integrator.h
3 */
4
5/**
6 * @defgroup odeGroup ODE Integrators
7 */
8
9// This file is part of Cantera. See License.txt in the top-level directory or
10// at https://cantera.org/license.txt for license and copyright information.
11
12#ifndef CT_INTEGRATOR_H
13#define CT_INTEGRATOR_H
14#include "FuncEval.h"
15
16#include "cantera/base/global.h"
17
18namespace Cantera
19{
20
21const int DIAG = 1;
22const int DENSE = 2;
23const int NOJAC = 4;
24const int JAC = 8;
25const int GMRES = 16;
26const int BAND = 32;
27
28/**
29 * Specifies the method used to integrate the system of equations.
30 * Not all methods are supported by all integrators.
31 */
33 BDF_Method, //!< Backward Differentiation
34 Adams_Method //! Adams
35};
36
37//! Specifies the method used for iteration.
38/*!
39 * Not all methods are supported by all integrators.
40 */
42 //! Newton Iteration
44 //! Functional Iteration
46};
47
48//! Abstract base class for ODE system integrators.
49/*!
50 * @ingroup odeGroup
51 */
53{
54public:
55 //! Default Constructor
57 }
58
59 //! Destructor
60 virtual ~Integrator() {
61 }
62
63 //! Set error tolerances.
64 /*!
65 * @param reltol scalar relative tolerance
66 * @param n Number of equations
67 * @param abstol array of N absolute tolerance values
68 */
69 virtual void setTolerances(doublereal reltol, size_t n,
70 doublereal* abstol) {
71 warn("setTolerances");
72 }
73
74 //! Set error tolerances.
75 /*!
76 * @param reltol scalar relative tolerance
77 * @param abstol scalar absolute tolerance
78 */
79 virtual void setTolerances(doublereal reltol, doublereal abstol) {
80 warn("setTolerances");
81 }
82
83 //! Set the sensitivity error tolerances
84 /*!
85 * @param reltol scalar relative tolerance
86 * @param abstol scalar absolute tolerance
87 */
88 virtual void setSensitivityTolerances(doublereal reltol, doublereal abstol)
89 { }
90
91 //! Set the problem type.
92 /*!
93 * @param probtype Type of the problem
94 */
95 virtual void setProblemType(int probtype) {
96 warn("setProblemType");
97 }
98
99 /**
100 * Initialize the integrator for a new problem. Call after all options have
101 * been set.
102 * @param t0 initial time
103 * @param func RHS evaluator object for system of equations.
104 */
105 virtual void initialize(doublereal t0, FuncEval& func) {
106 warn("initialize");
107 }
108
109 virtual void reinitialize(doublereal t0, FuncEval& func) {
110 warn("reinitialize");
111 }
112
113 //! Integrate the system of equations.
114 /*!
115 * @param tout Integrate to this time. Note that this is the
116 * absolute time value, not a time interval.
117 */
118 virtual void integrate(doublereal tout) {
119 warn("integrate");
120 }
121
122 /**
123 * Integrate the system of equations.
124 * @param tout integrate to this time. Note that this is the
125 * absolute time value, not a time interval.
126 */
127 virtual doublereal step(doublereal tout) {
128 warn("step");
129 return 0.0;
130 }
131
132 //! The current value of the solution of equation k.
133 virtual doublereal& solution(size_t k) {
134 warn("solution");
135 return m_dummy;
136 }
137
138 //! The current value of the solution of the system of equations.
139 virtual doublereal* solution() {
140 warn("solution");
141 return 0;
142 }
143
144 //! n-th derivative of the output function at time tout.
145 virtual double* derivative(double tout, int n) {
146 warn("derivative");
147 return 0;
148 }
149
150 //! Order used during the last solution step
151 virtual int lastOrder() const {
152 warn("lastOrder");
153 return 0;
154 }
155
156 //! The number of equations.
157 virtual int nEquations() const {
158 warn("nEquations");
159 return 0;
160 }
161
162 //! The number of function evaluations.
163 virtual int nEvals() const {
164 warn("nEvals");
165 return 0;
166 }
167
168 //! Set the maximum integration order that will be used.
169 virtual void setMaxOrder(int n) {
170 warn("setMaxorder");
171 }
172
173 //! Set the solution method
174 virtual void setMethod(MethodType t) {
175 warn("setMethodType");
176 }
177
178 //! Set the maximum step size
179 virtual void setMaxStepSize(double hmax) {
180 warn("setMaxStepSize");
181 }
182
183 //! Set the minimum step size
184 virtual void setMinStepSize(double hmin) {
185 warn("setMinStepSize");
186 }
187
188 //! Set the maximum permissible number of error test failures
189 virtual void setMaxErrTestFails(int n) {
190 warn("setMaxErrTestFails");
191 }
192
193 //! Set the maximum number of time-steps the integrator can take
194 //! before reaching the next output time
195 //! @param nmax The maximum number of steps, setting this value
196 //! to zero disables this option.
197 virtual void setMaxSteps(int nmax) {
198 warn("setMaxStep");
199 }
200
201 //! Returns the maximum number of time-steps the integrator can take
202 //! before reaching the next output time
203 virtual int maxSteps() {
204 warn("maxSteps");
205 return 0;
206 }
207
208 virtual void setBandwidth(int N_Upper, int N_Lower) {
209 warn("setBandwidth");
210 }
211
212 virtual int nSensParams() {
213 warn("nSensParams()");
214 return 0;
215 }
216
217 virtual double sensitivity(size_t k, size_t p) {
218 warn("sensitivity");
219 return 0.0;
220 }
221
222private:
223 doublereal m_dummy;
224 void warn(const std::string& msg) const {
225 writelog(">>>> Warning: method "+msg+" of base class "
226 +"Integrator called. Nothing done.\n");
227 }
228};
229
230// defined in ODE_integrators.cpp
231Integrator* newIntegrator(const std::string& itype);
232
233} // namespace
234
235#endif
Virtual base class for ODE right-hand-side function evaluators.
Definition: FuncEval.h:27
Abstract base class for ODE system integrators.
Definition: Integrator.h:53
virtual doublereal step(doublereal tout)
Integrate the system of equations.
Definition: Integrator.h:127
virtual void setTolerances(doublereal reltol, doublereal abstol)
Set error tolerances.
Definition: Integrator.h:79
virtual void setMaxStepSize(double hmax)
Set the maximum step size.
Definition: Integrator.h:179
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:197
virtual void setProblemType(int probtype)
Set the problem type.
Definition: Integrator.h:95
virtual int maxSteps()
Returns the maximum number of time-steps the integrator can take before reaching the next output time...
Definition: Integrator.h:203
Integrator()
Default Constructor.
Definition: Integrator.h:56
virtual void setMaxErrTestFails(int n)
Set the maximum permissible number of error test failures.
Definition: Integrator.h:189
virtual void setTolerances(doublereal reltol, size_t n, doublereal *abstol)
Set error tolerances.
Definition: Integrator.h:69
virtual double * derivative(double tout, int n)
n-th derivative of the output function at time tout.
Definition: Integrator.h:145
virtual void initialize(doublereal t0, FuncEval &func)
Initialize the integrator for a new problem.
Definition: Integrator.h:105
virtual void setMinStepSize(double hmin)
Set the minimum step size.
Definition: Integrator.h:184
virtual void integrate(doublereal tout)
Integrate the system of equations.
Definition: Integrator.h:118
virtual doublereal * solution()
The current value of the solution of the system of equations.
Definition: Integrator.h:139
virtual int lastOrder() const
Order used during the last solution step.
Definition: Integrator.h:151
virtual int nEvals() const
The number of function evaluations.
Definition: Integrator.h:163
virtual int nEquations() const
The number of equations.
Definition: Integrator.h:157
virtual void setSensitivityTolerances(doublereal reltol, doublereal abstol)
Set the sensitivity error tolerances.
Definition: Integrator.h:88
virtual ~Integrator()
Destructor.
Definition: Integrator.h:60
virtual void setMethod(MethodType t)
Set the solution method.
Definition: Integrator.h:174
virtual doublereal & solution(size_t k)
The current value of the solution of equation k.
Definition: Integrator.h:133
virtual void setMaxOrder(int n)
Set the maximum integration order that will be used.
Definition: Integrator.h:169
This file contains definitions for utility functions and text for modules, inputfiles,...
void writelog(const std::string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition: global.h:164
Namespace for the Cantera kernel.
Definition: AnyMap.h:29
MethodType
Specifies the method used to integrate the system of equations.
Definition: Integrator.h:32
@ BDF_Method
Backward Differentiation.
Definition: Integrator.h:33
@ Adams_Method
Adams.
Definition: Integrator.h:34
IterType
Specifies the method used for iteration.
Definition: Integrator.h:41
@ Functional_Iter
Functional Iteration.
Definition: Integrator.h:45
@ Newton_Iter
Newton Iteration.
Definition: Integrator.h:43