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