Cantera  3.1.0b1
Loading...
Searching...
No Matches
OneDim.h
Go to the documentation of this file.
1/**
2 * @file OneDim.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_ONEDIM_H
9#define CT_ONEDIM_H
10
11#include "Domain1D.h"
12#include "MultiJac.h"
13
14namespace Cantera
15{
16
17class Func1;
18class MultiNewton;
19class AnyMap;
20
21/**
22 * Container class for multiple-domain 1D problems. Each domain is
23 * represented by an instance of Domain1D.
24 * @ingroup onedGroup
25 */
26class OneDim
27{
28public:
29 //! Default constructor
30 OneDim();
31
32 //! Construct a OneDim container for the domains in the list *domains*.
33 OneDim(vector<shared_ptr<Domain1D>>& domains);
34
35 virtual ~OneDim();
36 OneDim(const OneDim&) = delete;
37 OneDim& operator=(const OneDim&) = delete;
38
39 //! Add a domain. Domains are added left-to-right.
40 void addDomain(shared_ptr<Domain1D> d);
41
42 //! Return a reference to the Jacobian evaluator of an OneDim object.
43 //! @ingroup derivGroup
45
46 //! Return a reference to the Newton iterator.
48
49 /**
50 * Solve F(x) = 0, where F(x) is the multi-domain residual function.
51 *
52 * @param x0 Starting estimate of solution.
53 * @param x1 Final solution satisfying F(x1) = 0.
54 * @param loglevel Controls amount of diagnostic output.
55 *
56 * @returns
57 * - 1 for success
58 * - -2 failure (maximum number of damping steps was reached)
59 * - -3 failure (solution was up against the bounds
60 */
61 int solve(double* x0, double* x1, int loglevel);
62
63 //! Number of domains.
64 size_t nDomains() const {
65 return m_dom.size();
66 }
67
68 //! Return a reference to domain i.
69 Domain1D& domain(size_t i) const {
70 return *m_dom[i];
71 }
72
73 //! Get the index of the domain named `name`.
74 size_t domainIndex(const string& name);
75
76 //! Check that the specified domain index is in range.
77 //! Throws an exception if n is greater than nDomains()-1
78 void checkDomainIndex(size_t n) const {
79 if (n >= m_dom.size()) {
80 throw IndexError("OneDim::checkDomainIndex", "domains", n,
81 m_dom.size()-1);
82 }
83 }
84
85 //! Check that an array size is at least nDomains().
86 //! Throws an exception if nn is less than nDomains(). Used before calls
87 //! which take an array pointer.
88 void checkDomainArraySize(size_t nn) const {
89 if (m_dom.size() > nn) {
90 throw ArraySizeError("OneDim::checkDomainArraySize", nn,
91 m_dom.size());
92 }
93 }
94
95 //! The index of the start of domain i in the solution vector.
96 size_t start(size_t i) const {
97 if (m_dom[i]->nComponents()) {
98 return m_dom[i]->loc();
99 } else {
100 // Special case for domains with no solution components to avoid
101 // spurious out-of-bounds memory access
102 return 0;
103 }
104 }
105
106 //! Total solution vector length;
107 size_t size() const {
108 return m_size;
109 }
110
111 //! Pointer to left-most domain (first added).
113 return m_dom[0].get();
114 }
115
116 //! Pointer to right-most domain (last added).
118 return m_dom.back().get();
119 }
120
121 //! Number of solution components at global point `jg`.
122 size_t nVars(size_t jg) {
123 return m_nvars[jg];
124 }
125
126 //! Location in the solution vector of the first component of global point `jg`.
127 size_t loc(size_t jg) {
128 return m_loc[jg];
129 }
130
131 //! Return the domain, local point index, and component name for the i-th
132 //! component of the global solution vector
133 std::tuple<string, size_t, string> component(size_t i);
134
135 //! Jacobian bandwidth.
136 size_t bandwidth() const {
137 return m_bw;
138 }
139
140 /**
141 * Initialize all domains. On the first call, this methods calls the init
142 * method of each domain, proceeding from left to right. Subsequent calls
143 * do nothing.
144 */
145 void init();
146
147 //! Total number of points.
148 size_t points() {
149 return m_pts;
150 }
151
152 /**
153 * Steady-state max norm (infinity norm) of the residual evaluated using
154 * solution x. On return, array r contains the steady-state residual
155 * values. Used only for diagnostic output.
156 */
157 double ssnorm(double* x, double* r);
158
159 //! Reciprocal of the time step.
160 double rdt() const {
161 return m_rdt;
162 }
163
164 //! Prepare for time stepping beginning with solution *x* and timestep *dt*.
165 void initTimeInteg(double dt, double* x);
166
167 //! True if transient mode.
168 bool transient() const {
169 return (m_rdt != 0.0);
170 }
171
172 //! True if steady mode.
173 bool steady() const {
174 return (m_rdt == 0.0);
175 }
176
177 /**
178 * Prepare to solve the steady-state problem. After invoking this method,
179 * subsequent calls to solve() will solve the steady-state problem. Sets
180 * the reciprocal of the time step to zero, and, if it was previously non-
181 * zero, signals that a new Jacobian will be needed.
182 */
183 void setSteadyMode();
184
185 /**
186 * Evaluate the multi-domain residual function
187 *
188 * @param j if j != npos, only evaluate residual for points j-1, j,
189 * and j + 1; otherwise, evaluate at all grid points.
190 * @param x solution vector
191 * @param r on return, contains the residual vector
192 * @param rdt Reciprocal of the time step. if omitted, then
193 * the default value is used.
194 * @param count Set to zero to omit this call from the statistics
195 */
196 void eval(size_t j, double* x, double* r, double rdt=-1.0, int count = 1);
197
198 //! Return a pointer to the domain global point *i* belongs to.
199 /*!
200 * The domains are scanned right-to-left, and the first one with starting
201 * location less or equal to i is returned.
202 */
203 Domain1D* pointDomain(size_t i);
204
205 //! Call after one or more grids has changed size, for example after being refined.
206 virtual void resize();
207
208 //! Access the vector indicating which equations contain a transient term.
209 //! Elements are 1 for equations with a transient terms and 0 otherwise.
210 vector<int>& transientMask() {
211 return m_mask;
212 }
213
214 /**
215 * Take time steps using Backward Euler.
216 *
217 * @param nsteps number of steps
218 * @param dt initial step size
219 * @param x current solution vector
220 * @param r solution vector after time stepping
221 * @param loglevel controls amount of printed diagnostics
222 * @returns size of last timestep taken
223 */
224 double timeStep(int nsteps, double dt, double* x, double* r, int loglevel);
225
226 //! Reset values such as negative species concentrations in each domain.
227 //! @see Domain1D::resetBadValues
228 void resetBadValues(double* x);
229
230 //! Write statistics about the number of iterations and Jacobians at each
231 //! grid level
232 /*!
233 * @param printTime Boolean that indicates whether time should be printed
234 * out The default is true. It's turned off for test
235 * problems where we don't want to print any times
236 */
237 void writeStats(int printTime = 1);
238
239 //! @name Options
240 //! @{
241
242 //! Set the minimum time step allowed during time stepping
243 void setMinTimeStep(double tmin) {
244 m_tmin = tmin;
245 }
246
247 //! Set the maximum time step allowed during time stepping
248 void setMaxTimeStep(double tmax) {
249 m_tmax = tmax;
250 }
251
252 /**
253 * Sets a factor by which the time step is reduced if the time stepping
254 * fails. The default value is 0.5.
255 *
256 * @param tfactor factor time step is multiplied by if time stepping fails
257 */
258 void setTimeStepFactor(double tfactor) {
259 m_tfactor = tfactor;
260 }
261
262 //! Set the maximum number of timeteps allowed before successful
263 //! steady-state solve
264 void setMaxTimeStepCount(int nmax) {
265 m_nsteps_max = nmax;
266 }
267
268 //! Return the maximum number of timeteps allowed before successful
269 //! steady-state solve
270 int maxTimeStepCount() const {
271 return m_nsteps_max;
272 }
273 //! @}
274
275 //! Set the maximum number of steps that can be taken using the same Jacobian
276 //! before it must be re-evaluated.
277 //! @param ss_age Age limit during steady-state mode
278 //! @param ts_age Age limit during time stepping mode. If not specified, the
279 //! steady-state age is also used during time stepping.
280 void setJacAge(int ss_age, int ts_age=-1);
281
282 /**
283 * Save statistics on function and Jacobian evaluation, and reset the
284 * counters. Statistics are saved only if the number of Jacobian
285 * evaluations is greater than zero. The statistics saved are:
286 *
287 * - number of grid points
288 * - number of Jacobian evaluations
289 * - CPU time spent evaluating Jacobians
290 * - number of non-Jacobian function evaluations
291 * - CPU time spent evaluating functions
292 * - number of time steps
293 */
294 void saveStats();
295
296 //! Clear saved statistics
297 void clearStats();
298
299 //! Return total grid size in each call to solve()
300 const vector<size_t>& gridSizeStats() {
301 saveStats();
302 return m_gridpts;
303 }
304
305 //! Return CPU time spent evaluating Jacobians in each call to solve()
306 const vector<double>& jacobianTimeStats() {
307 saveStats();
308 return m_jacElapsed;
309 }
310
311 //! Return CPU time spent on non-Jacobian function evaluations in each call
312 //! to solve()
313 const vector<double>& evalTimeStats() {
314 saveStats();
315 return m_funcElapsed;
316 }
317
318 //! Return number of Jacobian evaluations made in each call to solve()
319 const vector<int>& jacobianCountStats() {
320 saveStats();
321 return m_jacEvals;
322 }
323
324 //! Return number of non-Jacobian function evaluations made in each call to
325 //! solve()
326 const vector<int>& evalCountStats() {
327 saveStats();
328 return m_funcEvals;
329 }
330
331 //! Return number of time steps taken in each call to solve()
332 const vector<int>& timeStepStats() {
333 saveStats();
334 return m_timeSteps;
335 }
336
337 //! Set a function that will be called every time #eval is called.
338 //! Can be used to provide keyboard interrupt support in the high-level
339 //! language interfaces.
340 void setInterrupt(Func1* interrupt) {
341 m_interrupt = interrupt;
342 }
343
344 //! Set a function that will be called after each successful timestep. The
345 //! function will be called with the size of the timestep as the argument.
346 //! Intended to be used for observing solver progress for debugging
347 //! purposes.
348 void setTimeStepCallback(Func1* callback) {
349 m_time_step_callback = callback;
350 }
351
352protected:
353 //! Evaluate the steady-state Jacobian, accessible via jacobian()
354 //! @param[in] x Current state vector, length size()
355 //! @param[out] rsd Storage for the residual, length size()
356 void evalSSJacobian(double* x, double* rsd);
357
358 double m_tmin = 1e-16; //!< minimum timestep size
359 double m_tmax = 1e+08; //!< maximum timestep size
360
361 //! factor time step is multiplied by if time stepping fails ( < 1 )
362 double m_tfactor = 0.5;
363
364 shared_ptr<vector<double>> m_state; //!< Solution vector
365
366 unique_ptr<MultiJac> m_jac; //!< Jacobian evaluator
367 unique_ptr<MultiNewton> m_newt; //!< Newton iterator
368 double m_rdt = 0.0; //!< reciprocal of time step
369 bool m_jac_ok = false; //!< if true, Jacobian is current
370
371 size_t m_bw = 0; //!< Jacobian bandwidth
372 size_t m_size = 0; //!< solution vector size
373
374 //! All domains comprising the system
375 vector<shared_ptr<Domain1D>> m_dom;
376
377 //! All connector and boundary domains
378 vector<shared_ptr<Domain1D>> m_connect;
379
380 //! All bulk/flow domains
381 vector<shared_ptr<Domain1D>> m_bulk;
382
383 //! Indicates whether one-time initialization for each domain has been completed.
384 bool m_init = false;
385
386 //! Number of variables at each point, across all domains. Length points().
387 //! Accessed with nVars().
388 vector<size_t> m_nvars;
389
390 //! Location in the state vector of the first component of each point, across all
391 //! domains. Accessed with loc().
392 vector<size_t> m_loc;
393
394 //! Transient mask. See transientMask().
395 vector<int> m_mask;
396
397 //! Total number of points.
398 size_t m_pts = 0;
399
400 int m_ss_jac_age = 20; //!< Maximum age of the Jacobian in steady-state mode.
401 int m_ts_jac_age = 20; //!< Maximum age of the Jacobian in time-stepping mode.
402
403 //! Function called at the start of every call to #eval.
404 Func1* m_interrupt = nullptr;
405
406 //! User-supplied function called after each successful timestep.
408
409 //! Number of time steps taken in the current call to solve()
410 int m_nsteps = 0;
411
412 //! Maximum number of timesteps allowed per call to solve()
413 int m_nsteps_max = 500;
414
415private:
416 //! @name Statistics
417 //! Solver stats are collected after successfully solving on a particular grid.
418 //! @{
419 int m_nevals = 0; //!< Number of calls to eval()
420 double m_evaltime = 0; //!< Total time [s] spent in eval()
421
422 vector<size_t> m_gridpts; //!< Number of grid points in this grid
423 vector<int> m_jacEvals; //!< Number of Jacobian evaluations on this grid
424 vector<double> m_jacElapsed; //!< Time [s] spent evaluating Jacobians on this grid
425
426 //! Number of residual function evaluations on this grid (not counting evaluations
427 //! used to construct Jacobians).
428 vector<int> m_funcEvals;
429
430 //! Time [s] spent on residual function evaluations on this grid (not counting
431 //! evaluations used to construct Jacobians).
432 vector<double> m_funcElapsed;
433
434 //! Number of time steps taken in each call to solve() (for example, for each
435 //! successive grid refinement)
436 vector<int> m_timeSteps;
437
438 //! @}
439};
440
441}
442
443#endif
Array size error.
Base class for one-dimensional domains.
Definition Domain1D.h:29
Base class for 'functor' classes that evaluate a function of one variable.
Definition Func1.h:75
An array index is out of range.
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplie...
Definition MultiJac.h:24
Newton iterator for multi-domain, one-dimensional problems.
Definition MultiNewton.h:24
Container class for multiple-domain 1D problems.
Definition OneDim.h:27
int solve(double *x0, double *x1, int loglevel)
Solve F(x) = 0, where F(x) is the multi-domain residual function.
Definition OneDim.cpp:208
size_t start(size_t i) const
The index of the start of domain i in the solution vector.
Definition OneDim.h:96
int m_nsteps
Number of time steps taken in the current call to solve()
Definition OneDim.h:410
void init()
Initialize all domains.
Definition OneDim.cpp:319
void checkDomainIndex(size_t n) const
Check that the specified domain index is in range.
Definition OneDim.h:78
size_t m_size
solution vector size
Definition OneDim.h:372
int m_nsteps_max
Maximum number of timesteps allowed per call to solve()
Definition OneDim.h:413
virtual void resize()
Call after one or more grids has changed size, for example after being refined.
Definition OneDim.cpp:154
void resetBadValues(double *x)
Reset values such as negative species concentrations in each domain.
Definition OneDim.cpp:425
void saveStats()
Save statistics on function and Jacobian evaluation, and reset the counters.
Definition OneDim.cpp:123
unique_ptr< MultiNewton > m_newt
Newton iterator.
Definition OneDim.h:367
size_t size() const
Total solution vector length;.
Definition OneDim.h:107
size_t loc(size_t jg)
Location in the solution vector of the first component of global point jg.
Definition OneDim.h:127
void eval(size_t j, double *x, double *r, double rdt=-1.0, int count=1)
Evaluate the multi-domain residual function.
Definition OneDim.cpp:241
double ssnorm(double *x, double *r)
Steady-state max norm (infinity norm) of the residual evaluated using solution x.
Definition OneDim.cpp:273
vector< int > & transientMask()
Access the vector indicating which equations contain a transient term.
Definition OneDim.h:210
void addDomain(shared_ptr< Domain1D > d)
Add a domain. Domains are added left-to-right.
Definition OneDim.cpp:64
double rdt() const
Reciprocal of the time step.
Definition OneDim.h:160
void initTimeInteg(double dt, double *x)
Prepare for time stepping beginning with solution x and timestep dt.
Definition OneDim.cpp:283
size_t nDomains() const
Number of domains.
Definition OneDim.h:64
vector< double > m_jacElapsed
Time [s] spent evaluating Jacobians on this grid.
Definition OneDim.h:424
Domain1D * right()
Pointer to right-most domain (last added).
Definition OneDim.h:117
std::tuple< string, size_t, string > component(size_t i)
Return the domain, local point index, and component name for the i-th component of the global solutio...
Definition OneDim.cpp:50
size_t bandwidth() const
Jacobian bandwidth.
Definition OneDim.h:136
double m_rdt
reciprocal of time step
Definition OneDim.h:368
void setMaxTimeStep(double tmax)
Set the maximum time step allowed during time stepping.
Definition OneDim.h:248
shared_ptr< vector< double > > m_state
Solution vector.
Definition OneDim.h:364
void setMinTimeStep(double tmin)
Set the minimum time step allowed during time stepping.
Definition OneDim.h:243
vector< double > m_funcElapsed
Time [s] spent on residual function evaluations on this grid (not counting evaluations used to constr...
Definition OneDim.h:432
const vector< int > & evalCountStats()
Return number of non-Jacobian function evaluations made in each call to solve()
Definition OneDim.h:326
vector< shared_ptr< Domain1D > > m_connect
All connector and boundary domains.
Definition OneDim.h:378
vector< int > m_mask
Transient mask. See transientMask().
Definition OneDim.h:395
void setTimeStepCallback(Func1 *callback)
Set a function that will be called after each successful timestep.
Definition OneDim.h:348
Func1 * m_interrupt
Function called at the start of every call to eval.
Definition OneDim.h:404
const vector< int > & jacobianCountStats()
Return number of Jacobian evaluations made in each call to solve()
Definition OneDim.h:319
bool transient() const
True if transient mode.
Definition OneDim.h:168
void checkDomainArraySize(size_t nn) const
Check that an array size is at least nDomains().
Definition OneDim.h:88
const vector< size_t > & gridSizeStats()
Return total grid size in each call to solve()
Definition OneDim.h:300
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
Definition OneDim.h:264
const vector< int > & timeStepStats()
Return number of time steps taken in each call to solve()
Definition OneDim.h:332
void setTimeStepFactor(double tfactor)
Sets a factor by which the time step is reduced if the time stepping fails.
Definition OneDim.h:258
vector< shared_ptr< Domain1D > > m_bulk
All bulk/flow domains.
Definition OneDim.h:381
vector< int > m_funcEvals
Number of residual function evaluations on this grid (not counting evaluations used to construct Jaco...
Definition OneDim.h:428
vector< size_t > m_loc
Location in the state vector of the first component of each point, across all domains.
Definition OneDim.h:392
unique_ptr< MultiJac > m_jac
Jacobian evaluator.
Definition OneDim.h:366
double m_evaltime
Total time [s] spent in eval()
Definition OneDim.h:420
size_t nVars(size_t jg)
Number of solution components at global point jg.
Definition OneDim.h:122
size_t m_bw
Jacobian bandwidth.
Definition OneDim.h:371
double m_tfactor
factor time step is multiplied by if time stepping fails ( < 1 )
Definition OneDim.h:362
bool m_jac_ok
if true, Jacobian is current
Definition OneDim.h:369
int maxTimeStepCount() const
Return the maximum number of timeteps allowed before successful steady-state solve.
Definition OneDim.h:270
int m_ts_jac_age
Maximum age of the Jacobian in time-stepping mode.
Definition OneDim.h:401
size_t domainIndex(const string &name)
Get the index of the domain named name.
Definition OneDim.cpp:40
double timeStep(int nsteps, double dt, double *x, double *r, int loglevel)
Take time steps using Backward Euler.
Definition OneDim.cpp:331
Domain1D * pointDomain(size_t i)
Return a pointer to the domain global point i belongs to.
Definition OneDim.cpp:229
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...
Definition OneDim.cpp:96
vector< size_t > m_gridpts
Number of grid points in this grid.
Definition OneDim.h:422
void evalSSJacobian(double *x, double *rsd)
Evaluate the steady-state Jacobian, accessible via jacobian()
Definition OneDim.cpp:219
double m_tmin
minimum timestep size
Definition OneDim.h:358
size_t points()
Total number of points.
Definition OneDim.h:148
double m_tmax
maximum timestep size
Definition OneDim.h:359
int m_ss_jac_age
Maximum age of the Jacobian in steady-state mode.
Definition OneDim.h:400
vector< size_t > m_nvars
Number of variables at each point, across all domains.
Definition OneDim.h:388
int m_nevals
Number of calls to eval()
Definition OneDim.h:419
bool m_init
Indicates whether one-time initialization for each domain has been completed.
Definition OneDim.h:384
void writeStats(int printTime=1)
Write statistics about the number of iterations and Jacobians at each grid level.
Definition OneDim.cpp:106
void clearStats()
Clear saved statistics.
Definition OneDim.cpp:141
void setSteadyMode()
Prepare to solve the steady-state problem.
Definition OneDim.cpp:302
size_t m_pts
Total number of points.
Definition OneDim.h:398
OneDim()
Default constructor.
Definition OneDim.cpp:19
vector< int > m_jacEvals
Number of Jacobian evaluations on this grid.
Definition OneDim.h:423
const vector< double > & jacobianTimeStats()
Return CPU time spent evaluating Jacobians in each call to solve()
Definition OneDim.h:306
Func1 * m_time_step_callback
User-supplied function called after each successful timestep.
Definition OneDim.h:407
vector< int > m_timeSteps
Number of time steps taken in each call to solve() (for example, for each successive grid refinement)
Definition OneDim.h:436
Domain1D & domain(size_t i) const
Return a reference to domain i.
Definition OneDim.h:69
vector< shared_ptr< Domain1D > > m_dom
All domains comprising the system.
Definition OneDim.h:375
const vector< double > & evalTimeStats()
Return CPU time spent on non-Jacobian function evaluations in each call to solve()
Definition OneDim.h:313
Domain1D * left()
Pointer to left-most domain (first added).
Definition OneDim.h:112
bool steady() const
True if steady mode.
Definition OneDim.h:173
MultiNewton & newton()
Return a reference to the Newton iterator.
Definition OneDim.cpp:91
void setInterrupt(Func1 *interrupt)
Set a function that will be called every time eval is called.
Definition OneDim.h:340
MultiJac & jacobian()
Return a reference to the Jacobian evaluator of an OneDim object.
Definition OneDim.cpp:87
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595