Cantera  4.0.0a1
Loading...
Searching...
No Matches
SteadyStateSystem.h
Go to the documentation of this file.
1//! @file SteadyStateSystem.h
2
3// This file is part of Cantera. See License.txt in the top-level directory or
4// at https://cantera.org/license.txt for license and copyright information.
5
6#ifndef CT_STEADYSTATESYSTEM_H
7#define CT_STEADYSTATESYSTEM_H
8
10#include "SystemJacobian.h"
11
12namespace Cantera
13{
14
15class Func1;
16class MultiNewton;
17
18//! Base class for representing a system of differential-algebraic equations and solving
19//! for its steady-state response.
20//!
21//! @since New in %Cantera 3.2.
22//! @ingroup numerics
24{
25public:
27 virtual ~SteadyStateSystem();
28 SteadyStateSystem(const SteadyStateSystem&) = delete;
29 SteadyStateSystem& operator=(const SteadyStateSystem&) = delete;
30
31 //! Evaluate the residual function
32 //!
33 //! @param[in] x State vector
34 //! @param[out] r On return, contains the residual vector
35 //! @param rdt Reciprocal of the time step. if omitted, then the internally stored
36 //! value (accessible using the rdt() method) is used.
37 //! @param count Set to zero to omit this call from the statistics
38 virtual void eval(span<const double> x, span<double> r,
39 double rdt=-1.0, int count=1) = 0;
40
41 //! Evaluates the Jacobian at x0 using finite differences.
42 //!
43 //! The Jacobian is computed by perturbing each component of `x0`, evaluating the
44 //! residual function, and then estimating the partial derivatives numerically using
45 //! finite differences to determine the corresponding column of the Jacobian.
46 //!
47 //! @param x0 State vector at which to evaluate the Jacobian
48 virtual void evalJacobian(span<const double> x0) = 0;
49
50 //! Compute the weighted norm of `step`.
51 //!
52 //! The purpose is to measure the "size" of the step vector @f$ \Delta x @f$ in a
53 //! way that takes into account the relative importance or scale of different
54 //! solution components. Each component of the step vector is normalized by a weight
55 //! that depends on both the current magnitude of the solution vector and specified
56 //! tolerances. This makes the norm dimensionless and scaled appropriately, avoiding
57 //! issues where some components dominate due to differences in their scales. See
58 //! OneDim::weightedNorm() for a representative implementation.
59 virtual double weightedNorm(span<const double> step) const = 0;
60
61 //! Set the initial guess. Should be called before solve().
62 void setInitialGuess(span<const double> x);
63
64 //! Solve the steady-state problem, taking internal timesteps as necessary until
65 //! the Newton solver can converge for the steady problem.
66 //! @param loglevel Controls amount of diagnostic output.
67 void solve(int loglevel=0);
68
69 //! Get the converged steady-state solution after calling solve().
70 void getState(span<double> x) const;
71
72 //! Steady-state max norm (infinity norm) of the residual evaluated using solution
73 //! x. On return, array r contains the steady-state residual values.
74 double ssnorm(span<const double> x, span<double> r);
75
76 //! Total solution vector length;
77 size_t size() const {
78 return m_size;
79 }
80
81 //! Call to set the size of internal data structures after first defining the system
82 //! or if the problem size changes, for example after grid refinement.
83 //!
84 //! The #m_size variable should be updated before calling this method
85 virtual void resize();
86
87 //! Jacobian bandwidth.
88 size_t bandwidth() const {
89 return m_bw;
90 }
91
92 //! Get the name of the i-th component of the state vector
93 virtual string componentName(size_t i) const { return fmt::format("{}", i); }
94
95 //! Get header lines describing the column names included in a component label.
96 //! Headings should be aligned with items formatted in the output of
97 //! componentTableLabel(). Two lines are allowed. If only one is needed, the first
98 //! line should be blank.
99 virtual pair<string, string> componentTableHeader() const {
100 return {"", "component name"};
101 }
102
103 //! Get elements of the component name, aligned with the column headings given
104 //! by componentTableHeader().
105 virtual string componentTableLabel(size_t i) const { return componentName(i); }
106
107 //! Get the upper bound for global component *i* in the state vector
108 virtual double upperBound(size_t i) const = 0;
109
110 //! Get the lower bound for global component *i* in the state vector
111 virtual double lowerBound(size_t i) const = 0;
112
113 //! Return a reference to the Newton iterator.
115
116 //! Set the linear solver used to hold the Jacobian matrix and solve linear systems
117 //! as part of each Newton iteration.
118 void setLinearSolver(shared_ptr<SystemJacobian> solver);
119
120 //! Get the the linear solver being used to hold the Jacobian matrix and solve
121 //! linear systems as part of each Newton iteration.
122 shared_ptr<SystemJacobian> linearSolver() const { return m_jac; }
123
124 //! Reciprocal of the time step.
125 double rdt() const {
126 return m_rdt;
127 }
128
129 //! True if transient mode.
130 bool transient() const {
131 return (m_rdt != 0.0);
132 }
133
134 //! True if steady mode.
135 bool steady() const {
136 return (m_rdt == 0.0);
137 }
138
139 //! Prepare for time stepping beginning with solution *x* and timestep *dt*.
140 virtual void initTimeInteg(double dt, span<const double> x);
141
142 //! Prepare to solve the steady-state problem. After invoking this method,
143 //! subsequent calls to solve() will solve the steady-state problem. Sets the
144 //! reciprocal of the time step to zero, and, if it was previously non-zero, signals
145 //! that a new Jacobian will be needed.
146 virtual void setSteadyMode();
147
148 //! Access the vector indicating which equations contain a transient term.
149 //! Elements are 1 for equations with a transient terms and 0 otherwise.
150 vector<int>& transientMask() {
151 return m_mask;
152 }
153
154 //! Take time steps using Backward Euler.
155 //!
156 //! @param nsteps number of steps
157 //! @param dt initial step size
158 //! @param x current solution vector
159 //! @param r solution vector after time stepping
160 //! @param loglevel controls amount of printed diagnostics
161 //! @returns size of last timestep taken
162 double timeStep(int nsteps, double dt, span<double> x, span<double> r,
163 int loglevel);
164
165 //! Reset values such as negative species concentrations
166 virtual void resetBadValues(span<double> x) {}
167
168 //! @name Options
169 //! @{
170
171 //! Set the number of time steps to try when the steady Newton solver is
172 //! unsuccessful.
173 //! @param stepsize Initial time step size [s]
174 //! @param tsteps A sequence of time step counts to take after subsequent failures
175 //! of the steady-state solver. The last value in `tsteps` will be used again
176 //! after further unsuccessful solution attempts.
177 void setTimeStep(double stepsize, span<const int> tsteps);
178
179 //! Set the minimum time step allowed during time stepping
180 void setMinTimeStep(double tmin) {
181 m_tmin = tmin;
182 }
183
184 //! Set the maximum time step allowed during time stepping
185 void setMaxTimeStep(double tmax) {
186 m_tmax = tmax;
187 }
188
189 //! Sets a factor by which the time step is reduced if the time stepping fails.
190 //! The default value is 0.5.
191 //!
192 //! @param tfactor factor time step is multiplied by if time stepping fails
193 void setTimeStepFactor(double tfactor) {
194 m_tfactor = tfactor;
195 }
196
197 //! Set the maximum number of timeteps allowed before successful steady-state solve
198 void setMaxTimeStepCount(int nmax) {
199 m_nsteps_max = nmax;
200 }
201
202 //! Get the maximum number of timeteps allowed before successful steady-state solve
203 int maxTimeStepCount() const {
204 return m_nsteps_max;
205 }
206 //! @}
207
208 //! Set the maximum number of steps that can be taken using the same Jacobian
209 //! before it must be re-evaluated.
210 //! @param ss_age Age limit during steady-state mode
211 //! @param ts_age Age limit during time stepping mode. If not specified, the
212 //! steady-state age is also used during time stepping.
213 void setJacAge(int ss_age, int ts_age=-1);
214
215 //! Set a function that will be called every time #eval is called.
216 //! Can be used to provide keyboard interrupt support in the high-level
217 //! language interfaces.
218 void setInterrupt(Func1* interrupt) {
219 m_interrupt = interrupt;
220 }
221
222 //! Set a function that will be called after each successful timestep. The
223 //! function will be called with the size of the timestep as the argument.
224 //! Intended to be used for observing solver progress for debugging
225 //! purposes.
226 void setTimeStepCallback(Func1* callback) {
227 m_time_step_callback = callback;
228 }
229
230 //! Configure perturbations used to evaluate finite difference Jacobian
231 //! @param relative Relative perturbation (multiplied by the absolute value of
232 //! each component). Default `1.0e-5`.
233 //! @param absolute Absolute perturbation (independent of component value).
234 //! Default `1.0e-10`.
235 //! @param threshold Threshold below which to exclude elements from the Jacobian
236 //! Default `0.0`.
237 void setJacobianPerturbation(double relative, double absolute, double threshold) {
238 m_jacobianRelPerturb = relative;
239 m_jacobianAbsPerturb = absolute;
240 m_jacobianThreshold = threshold;
241 }
242
243 //! Write solver debugging based on the specified log level.
244 //!
245 //! @see Sim1D::writeDebugInfo for a specific implementation of this capability.
246 virtual void writeDebugInfo(const string& header_suffix, const string& message,
247 int loglevel, int attempt_counter) {}
248
249 //! Delete debug output file that may be created by writeDebugInfo() when solving
250 //! with high `loglevel`.
251 virtual void clearDebugFile() {}
252
253protected:
254 //! Evaluate the steady-state Jacobian, accessible via linearSolver()
255 //! @param[in] x Current state vector, length size()
256 void evalSSJacobian(span<const double> x);
257
258 //! Array of number of steps to take after each unsuccessful steady-state solve
259 //! before re-attempting the steady-state solution. For subsequent time stepping
260 //! calls, the final value is reused. See setTimeStep().
261 vector<int> m_steps = { 10 };
262
263 double m_tstep = 1.0e-5; //!< Initial timestep
264 double m_tmin = 1e-16; //!< Minimum timestep size
265 double m_tmax = 1e+08; //!< Maximum timestep size
266
267 //! Factor time step is multiplied by if time stepping fails ( < 1 )
268 double m_tfactor = 0.5;
269
270 shared_ptr<vector<double>> m_state; //!< Solution vector
271
272 //! Work array used to hold the residual or the new solution
273 vector<double> m_xnew;
274
275 //! State vector after the last successful set of time steps
276 vector<double> m_xlast_ts;
277
278 unique_ptr<MultiNewton> m_newt; //!< Newton iterator
279 double m_rdt = 0.0; //!< Reciprocal of time step
280
281 shared_ptr<SystemJacobian> m_jac; //!< Jacobian evaluator
282 bool m_jac_ok = false; //!< If `true`, Jacobian is current
283
284 size_t m_bw = 0; //!< Jacobian bandwidth
285 size_t m_size = 0; //!< %Solution vector size
286
287 //! Work arrays used during Jacobian evaluation
288 vector<double> m_work1, m_work2;
289
290 vector<int> m_mask; //!< Transient mask. See transientMask()
291 int m_ss_jac_age = 20; //!< Maximum age of the Jacobian in steady-state mode.
292 int m_ts_jac_age = 20; //!< Maximum age of the Jacobian in time-stepping mode.
293
294 //! Counter used to manage the number of states stored in the debug log file
295 //! generated by writeDebugInfo()
297
298 //! Constant that determines the maximum number of states stored in the debug log
299 //! file generated by writeDebugInfo()
301
302 //! Function called at the start of every call to #eval.
303 Func1* m_interrupt = nullptr;
304
305 //! User-supplied function called after each successful timestep.
307
308 int m_nsteps = 0; //!< Number of time steps taken in the current call to solve()
309 int m_nsteps_max = 500; //!< Maximum number of timesteps allowed per call to solve()
310
311 //! Threshold for ignoring small elements in Jacobian
313 //! Relative perturbation of each component in finite difference Jacobian
314 double m_jacobianRelPerturb = 1e-5;
315 //! Absolute perturbation of each component in finite difference Jacobian
316 double m_jacobianAbsPerturb = 1e-10;
317};
318
319}
320
321#endif
Declarations for class SystemJacobian.
Base class for 'functor' classes that evaluate a function of one variable.
Definition Func1.h:75
Newton iterator for multi-domain, one-dimensional problems.
Definition MultiNewton.h:24
Base class for representing a system of differential-algebraic equations and solving for its steady-s...
int m_nsteps
Number of time steps taken in the current call to solve()
size_t m_size
Solution vector size
int m_nsteps_max
Maximum number of timesteps allowed per call to solve()
virtual void resize()
Call to set the size of internal data structures after first defining the system or if the problem si...
shared_ptr< SystemJacobian > linearSolver() const
Get the the linear solver being used to hold the Jacobian matrix and solve linear systems as part of ...
double timeStep(int nsteps, double dt, span< double > x, span< double > r, int loglevel)
Take time steps using Backward Euler.
virtual void resetBadValues(span< double > x)
Reset values such as negative species concentrations.
vector< double > m_xnew
Work array used to hold the residual or the new solution.
unique_ptr< MultiNewton > m_newt
Newton iterator.
void getState(span< double > x) const
Get the converged steady-state solution after calling solve().
double m_jacobianAbsPerturb
Absolute perturbation of each component in finite difference Jacobian.
size_t size() const
Total solution vector length;.
vector< int > & transientMask()
Access the vector indicating which equations contain a transient term.
virtual double upperBound(size_t i) const =0
Get the upper bound for global component i in the state vector.
double rdt() const
Reciprocal of the time step.
virtual string componentName(size_t i) const
Get the name of the i-th component of the state vector.
virtual string componentTableLabel(size_t i) const
Get elements of the component name, aligned with the column headings given by componentTableHeader().
virtual void initTimeInteg(double dt, span< const double > x)
Prepare for time stepping beginning with solution x and timestep dt.
void evalSSJacobian(span< const double > x)
Evaluate the steady-state Jacobian, accessible via linearSolver()
size_t bandwidth() const
Jacobian bandwidth.
virtual double weightedNorm(span< const double > step) const =0
Compute the weighted norm of step.
double m_rdt
Reciprocal of time step.
double m_jacobianThreshold
Threshold for ignoring small elements in Jacobian.
void setMaxTimeStep(double tmax)
Set the maximum time step allowed during time stepping.
shared_ptr< SystemJacobian > m_jac
Jacobian evaluator.
void setInitialGuess(span< const double > x)
Set the initial guess. Should be called before solve().
shared_ptr< vector< double > > m_state
Solution vector.
void setMinTimeStep(double tmin)
Set the minimum time step allowed during time stepping.
virtual void eval(span< const double > x, span< double > r, double rdt=-1.0, int count=1)=0
Evaluate the residual function.
vector< int > m_mask
Transient mask.
void setTimeStepCallback(Func1 *callback)
Set a function that will be called after each successful timestep.
Func1 * m_interrupt
Function called at the start of every call to eval.
virtual double lowerBound(size_t i) const =0
Get the lower bound for global component i in the state vector.
bool transient() const
True if transient mode.
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
void solve(int loglevel=0)
Solve the steady-state problem, taking internal timesteps as necessary until the Newton solver can co...
double ssnorm(span< const double > x, span< double > r)
Steady-state max norm (infinity norm) of the residual evaluated using solution x.
void setTimeStepFactor(double tfactor)
Sets a factor by which the time step is reduced if the time stepping fails.
void setTimeStep(double stepsize, span< const int > tsteps)
Set the number of time steps to try when the steady Newton solver is unsuccessful.
vector< int > m_steps
Array of number of steps to take after each unsuccessful steady-state solve before re-attempting the ...
size_t m_bw
Jacobian bandwidth.
double m_tfactor
Factor time step is multiplied by if time stepping fails ( < 1 )
bool m_jac_ok
If true, Jacobian is current.
double m_tstep
Initial timestep.
int maxTimeStepCount() const
Get the maximum number of timeteps allowed before successful steady-state solve.
int m_ts_jac_age
Maximum age of the Jacobian in time-stepping mode.
void setJacAge(int ss_age, int ts_age=-1)
Set the maximum number of steps that can be taken using the same Jacobian before it must be re-evalua...
virtual void writeDebugInfo(const string &header_suffix, const string &message, int loglevel, int attempt_counter)
Write solver debugging based on the specified log level.
double m_tmin
Minimum timestep size.
void setJacobianPerturbation(double relative, double absolute, double threshold)
Configure perturbations used to evaluate finite difference Jacobian.
double m_tmax
Maximum timestep size.
int m_ss_jac_age
Maximum age of the Jacobian in steady-state mode.
virtual void evalJacobian(span< const double > x0)=0
Evaluates the Jacobian at x0 using finite differences.
virtual void setSteadyMode()
Prepare to solve the steady-state problem.
virtual void clearDebugFile()
Delete debug output file that may be created by writeDebugInfo() when solving with high loglevel.
int m_attempt_counter
Counter used to manage the number of states stored in the debug log file generated by writeDebugInfo(...
void setLinearSolver(shared_ptr< SystemJacobian > solver)
Set the linear solver used to hold the Jacobian matrix and solve linear systems as part of each Newto...
int m_max_history
Constant that determines the maximum number of states stored in the debug log file generated by write...
Func1 * m_time_step_callback
User-supplied function called after each successful timestep.
virtual pair< string, string > componentTableHeader() const
Get header lines describing the column names included in a component label.
double m_jacobianRelPerturb
Relative perturbation of each component in finite difference Jacobian.
vector< double > m_xlast_ts
State vector after the last successful set of time steps.
bool steady() const
True if steady mode.
vector< double > m_work1
Work arrays used during Jacobian evaluation.
MultiNewton & newton()
Return a reference to the Newton iterator.
void setInterrupt(Func1 *interrupt)
Set a function that will be called every time eval is called.
This file contains definitions of constants, types and terms that are used in internal routines and a...
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595