Cantera  3.2.0a5
Loading...
Searching...
No Matches
OneDim.cpp
Go to the documentation of this file.
1//! @file OneDim.cpp
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
11
12#include <fstream>
13#include <ctime>
14
15using namespace std;
16
17namespace Cantera
18{
19
20OneDim::OneDim(vector<shared_ptr<Domain1D>>& domains)
21{
22 for (auto& dom : domains) {
23 addDomain(dom);
24 }
25 init();
26 resize();
27}
28
29size_t OneDim::domainIndex(const string& name) const
30{
31 for (size_t n = 0; n < m_dom.size(); n++) {
32 if (domain(n).id() == name) {
33 return n;
34 }
35 }
36 throw CanteraError("OneDim::domainIndex", "Domain '{}' not found", name);
37}
38
39std::tuple<string, size_t, string> OneDim::component(size_t i) const {
40 if (i >= m_size) {
41 throw IndexError("OneDim::component", "components", i, m_size);
42 }
43 const auto& [n, j, k] = m_componentInfo[i];
44 Domain1D& dom = domain(n);
45 return make_tuple(dom.id(), j, dom.componentName(k));
46}
47
48string OneDim::componentName(size_t i) const {
49 const auto& [dom, pt, comp] = component(i);
50 return fmt::format("domain {}, component {} at point {}", dom, comp, pt);
51}
52
53pair<string, string> OneDim::componentTableHeader() const
54{
55 return {"", "Domain Pt. Component"};
56}
57
58string OneDim::componentTableLabel(size_t i) const
59{
60 const auto& [dom, pt, comp] = component(i);
61 return fmt::format("{:8s} {:3d} {:<12s}", dom, pt, comp);
62}
63
64double OneDim::upperBound(size_t i) const
65{
66 const auto& [n, j, k] = m_componentInfo[i];
67 Domain1D& dom = domain(n);
68 return dom.upperBound(k);
69}
70
71double OneDim::lowerBound(size_t i) const
72{
73 const auto& [n, j, k] = m_componentInfo[i];
74 Domain1D& dom = domain(n);
75 return dom.lowerBound(k);
76}
77
78void OneDim::addDomain(shared_ptr<Domain1D> d)
79{
80 // if 'd' is not the first domain, link it to the last domain
81 // added (the rightmost one)
82 size_t n = m_dom.size();
83 if (n > 0) {
84 m_dom.back()->append(d.get());
85 }
86
87 // every other domain is a connector
88 if (n % 2 == 0) {
89 m_connect.push_back(d);
90 } else {
91 m_bulk.push_back(d);
92 }
93
94 // add it also to the global domain list, and set its container and position
95 m_dom.push_back(d);
96 d->setData(m_state);
97 d->setContainer(this, m_dom.size()-1);
98 resize();
99}
100
101double OneDim::weightedNorm(const double* step) const
102{
103 double sum = 0.0;
104 const double* x = m_state->data();
105 size_t nd = nDomains();
106 for (size_t n = 0; n < nd; n++) {
107 Domain1D& dom = domain(n);
108 double d_sum = 0.0;
109 size_t nv = dom.nComponents();
110 size_t np = dom.nPoints();
111 size_t dstart = start(n);
112
113 for (size_t i = 0; i < nv; i++) {
114 double esum = 0.0;
115 for (size_t j = 0; j < np; j++) {
116 esum += fabs(x[dstart + nv*j + i]);
117 }
118 double ewt = dom.rtol(i)*esum/np + dom.atol(i);
119 for (size_t j = 0; j < np; j++) {
120 double f = step[dstart + nv*j + i]/ewt;
121 d_sum += f*f;
122 }
123 }
124 sum += d_sum;
125 }
126 return sqrt(sum / size());
127}
128
130{
131 warn_deprecated("OneDim::jacobian",
132 "Replaced by linearSolver(). To be removed after Cantera 3.2.");
133 auto multijac = dynamic_pointer_cast<MultiJac>(m_jac);
134 if (multijac) {
135 return *multijac;
136 } else {
137 throw CanteraError("OneDim::jacobian", "Active Jacobian is not a MultiJac");
138 }
139}
140
141void OneDim::writeStats(int printTime)
142{
143 saveStats();
144 writelog("\nStatistics:\n\n Grid Timesteps Functions Time Jacobians Time\n");
145 size_t n = m_gridpts.size();
146 for (size_t i = 0; i < n; i++) {
147 if (printTime) {
148 writelog("{:5d} {:5d} {:6d} {:9.4f} {:5d} {:9.4f}\n",
150 m_jacEvals[i], m_jacElapsed[i]);
151 } else {
152 writelog("{:5d} {:5d} {:6d} NA {:5d} NA\n",
154 }
155 }
156}
157
159{
160 if (m_jac) {
161 int nev = m_jac->nEvals();
162 if (nev > 0 && m_nevals > 0) {
163 m_gridpts.push_back(m_pts);
164 m_jacEvals.push_back(m_jac->nEvals());
165 m_jacElapsed.push_back(m_jac->elapsedTime());
166 m_funcEvals.push_back(m_nevals);
167 m_nevals = 0;
168 m_funcElapsed.push_back(m_evaltime);
169 m_evaltime = 0.0;
171 m_nsteps = 0;
172 }
173 }
174}
175
177{
178 m_gridpts.clear();
179 m_jacEvals.clear();
180 m_jacElapsed.clear();
181 m_funcEvals.clear();
182 m_funcElapsed.clear();
183 m_timeSteps.clear();
184 m_nevals = 0;
185 m_evaltime = 0.0;
186 m_nsteps = 0;
187}
188
190{
191 m_bw = 0;
192 m_nvars.clear();
193 m_loc.clear();
194 m_componentInfo.clear();
195 size_t lc = 0;
196
197 // save the statistics for the last grid
198 saveStats();
199 m_pts = 0;
200 for (size_t i = 0; i < nDomains(); i++) {
201 const auto& d = m_dom[i];
202
203 size_t np = d->nPoints();
204 size_t nv = d->nComponents();
205 for (size_t n = 0; n < np; n++) {
206 m_nvars.push_back(nv);
207 m_loc.push_back(lc);
208 lc += nv;
209 m_pts++;
210 for (size_t k = 0; k < nv; k++) {
211 m_componentInfo.emplace_back(i, n, k);
212 }
213 }
214
215 // update the Jacobian bandwidth
216
217 // bandwidth of the local block
218 size_t bw1 = d->bandwidth();
219 if (bw1 == npos) {
220 bw1 = std::max<size_t>(2*d->nComponents(), 1) - 1;
221 }
222 m_bw = std::max(m_bw, bw1);
223
224 // bandwidth of the block coupling the first point of this
225 // domain to the last point of the previous domain
226 if (i > 0) {
227 size_t bw2 = m_dom[i-1]->bandwidth();
228 if (bw2 == npos) {
229 bw2 = m_dom[i-1]->nComponents();
230 }
231 bw2 += d->nComponents() - 1;
232 m_bw = std::max(m_bw, bw2);
233 }
234 m_size = d->loc() + d->size();
235 }
236 if (!m_jac) {
237 m_jac = newSystemJacobian("banded-direct");
238 }
240}
241
243{
244 Domain1D* d = right();
245 while (d) {
246 if (d->loc() <= i) {
247 return d;
248 }
249 d = d->left();
250 }
251 return 0;
252}
253
254void OneDim::eval(size_t j, double* x, double* r, double rdt, int count)
255{
256 clock_t t0 = clock();
257 if (m_interrupt) {
259 }
260 fill(r, r + m_size, 0.0);
261 if (j == npos) {
262 fill(m_mask.begin(), m_mask.end(), 0);
263 }
264 if (rdt < 0.0) {
265 rdt = m_rdt;
266 }
267
268 // iterate over the bulk domains first
269 for (const auto& d : m_bulk) {
270 d->eval(j, x, r, m_mask.data(), rdt);
271 }
272
273 // then over the connector domains
274 for (const auto& d : m_connect) {
275 d->eval(j, x, r, m_mask.data(), rdt);
276 }
277
278 // increment counter and time
279 if (count) {
280 clock_t t1 = clock();
281 m_evaltime += double(t1 - t0)/CLOCKS_PER_SEC;
282 m_nevals++;
283 }
284}
285
286void OneDim::evalJacobian(double* x0)
287{
288 m_jac->reset();
289 clock_t t0 = clock();
290 m_work1.resize(size());
291 m_work2.resize(size());
292 eval(npos, x0, m_work1.data(), 0.0, 0);
293 size_t ipt = 0;
294 for (size_t j = 0; j < points(); j++) {
295 size_t nv = nVars(j);
296 for (size_t n = 0; n < nv; n++) {
297 // perturb x(n); preserve sign(x(n))
298 double xsave = x0[ipt];
299 double dx = fabs(xsave) * m_jacobianRelPerturb + m_jacobianAbsPerturb;
300 if (xsave < 0) {
301 dx = -dx;
302 }
303 x0[ipt] = xsave + dx;
304 double rdx = 1.0 / (x0[ipt] - xsave);
305
306 // calculate perturbed residual
307 eval(j, x0, m_work2.data(), 0.0, 0);
308
309 // compute nth column of Jacobian
310 for (size_t i = j - 1; i != j+2; i++) {
311 if (i != npos && i < points()) {
312 size_t mv = nVars(i);
313 size_t iloc = loc(i);
314 for (size_t m = 0; m < mv; m++) {
315 double delta = m_work2[m+iloc] - m_work1[m+iloc];
316 if (std::abs(delta) > m_jacobianThreshold || m+iloc == ipt) {
317 m_jac->setValue(m + iloc, ipt, delta * rdx);
318 }
319 }
320 }
321 }
322 x0[ipt] = xsave;
323 ipt++;
324 }
325 }
326
327 m_jac->updateElapsed(double(clock() - t0) / CLOCKS_PER_SEC);
328 m_jac->incrementEvals();
329 m_jac->setAge(0);
330}
331
332void OneDim::initTimeInteg(double dt, double* x)
333{
335 // iterate over all domains, preparing each one to begin time stepping
336 Domain1D* d = left();
337 while (d) {
338 d->initTimeInteg(dt, x);
339 d = d->right();
340 }
341}
342
344{
345 if (m_rdt == 0) {
346 return;
347 }
349 // iterate over all domains, preparing them for steady-state solution
350 Domain1D* d = left();
351 while (d) {
352 d->setSteadyMode();
353 d = d->right();
354 }
355}
356
358{
359 if (!m_init) {
360 Domain1D* d = left();
361 while (d) {
362 d->init();
363 d = d->right();
364 }
365 }
366 m_init = true;
367}
368
370{
371 for (auto dom : m_dom) {
372 dom->resetBadValues(x);
373 }
374}
375
376}
Base class for exceptions thrown by Cantera classes.
Base class for one-dimensional domains.
Definition Domain1D.h:29
size_t nComponents() const
Number of components at each grid point.
Definition Domain1D.h:174
double rtol(size_t n)
Relative tolerance of the nth component.
Definition Domain1D.h:294
Domain1D * left() const
Return a pointer to the left neighbor.
Definition Domain1D.h:698
string id() const
Returns the identifying tag for this domain.
Definition Domain1D.h:718
size_t nPoints() const
Number of grid points in this domain.
Definition Domain1D.h:202
double lowerBound(size_t n) const
Lower bound on the nth component.
Definition Domain1D.h:329
double upperBound(size_t n) const
Upper bound on the nth component.
Definition Domain1D.h:324
Domain1D * right() const
Return a pointer to the right neighbor.
Definition Domain1D.h:703
virtual string componentName(size_t n) const
Name of component n. May be overloaded.
Definition Domain1D.cpp:75
virtual void init()
Initialize.
Definition Domain1D.h:152
double atol(size_t n)
Absolute tolerance of the nth component.
Definition Domain1D.h:299
void initTimeInteg(double dt, const double *x0)
Performs the setup required before starting a time-stepping solution.
Definition Domain1D.h:355
void setSteadyMode()
Set the internally-stored reciprocal of the time step to 0.0, which is used to indicate that the prob...
Definition Domain1D.h:364
virtual size_t loc(size_t j=0) const
Location of the start of the local solution vector in the global solution vector.
Definition Domain1D.h:657
virtual double eval(double t) const
Evaluate the function.
Definition Func1.cpp:28
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:25
size_t start(size_t i) const
The index of the start of domain i in the solution vector.
Definition OneDim.h:104
void init()
Initialize all domains.
Definition OneDim.cpp:357
double weightedNorm(const double *step) const override
Compute the weighted norm of a step vector.
Definition OneDim.cpp:101
void resize() override
Call to set the size of internal data structures after first defining the system or if the problem si...
Definition OneDim.cpp:189
string componentName(size_t i) const override
Get the name of the i-th component of the state vector.
Definition OneDim.cpp:48
void saveStats()
Save statistics on function and Jacobian evaluation, and reset the counters.
Definition OneDim.cpp:158
pair< string, string > componentTableHeader() const override
Get header lines describing the column names included in a component label.
Definition OneDim.cpp:53
size_t loc(size_t jg)
Location in the solution vector of the first component of global point jg.
Definition OneDim.h:130
double upperBound(size_t i) const override
Get the upper bound for global component i in the state vector.
Definition OneDim.cpp:64
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:254
void addDomain(shared_ptr< Domain1D > d)
Add a domain. Domains are added left-to-right.
Definition OneDim.cpp:78
string componentTableLabel(size_t i) const override
Get elements of the component name, aligned with the column headings given by componentTableHeader().
Definition OneDim.cpp:58
size_t nDomains() const
Number of domains.
Definition OneDim.h:73
vector< double > m_jacElapsed
Time [s] spent evaluating Jacobians on this grid.
Definition OneDim.h:294
Domain1D * right()
Pointer to right-most domain (last added).
Definition OneDim.h:120
size_t domainIndex(const string &name) const
Get the index of the domain named name.
Definition OneDim.cpp:29
OneDim()=default
Default constructor.
void setSteadyMode() override
Prepare to solve the steady-state problem.
Definition OneDim.cpp:343
vector< double > m_funcElapsed
Time [s] spent on residual function evaluations on this grid (not counting evaluations used to constr...
Definition OneDim.h:302
vector< shared_ptr< Domain1D > > m_connect
All connector and boundary domains.
Definition OneDim.h:262
vector< std::tuple< size_t, size_t, size_t > > m_componentInfo
Domain, grid point, and component indices for each element of the global state vector.
Definition OneDim.h:280
vector< shared_ptr< Domain1D > > m_bulk
All bulk/flow domains.
Definition OneDim.h:265
vector< int > m_funcEvals
Number of residual function evaluations on this grid (not counting evaluations used to construct Jaco...
Definition OneDim.h:298
vector< size_t > m_loc
Location in the state vector of the first component of each point, across all domains.
Definition OneDim.h:276
void evalJacobian(double *x0) override
Evaluates the Jacobian at x0 using finite differences.
Definition OneDim.cpp:286
double m_evaltime
Total time [s] spent in eval()
Definition OneDim.h:290
size_t nVars(size_t jg)
Number of solution components at global point jg.
Definition OneDim.h:125
std::tuple< string, size_t, string > component(size_t i) const
Return the domain, local point index, and component name for the i-th component of the global solutio...
Definition OneDim.cpp:39
Domain1D * pointDomain(size_t i)
Return a pointer to the domain global point i belongs to.
Definition OneDim.cpp:242
void resetBadValues(double *x) override
Reset values such as negative species concentrations.
Definition OneDim.cpp:369
vector< size_t > m_gridpts
Number of grid points in this grid.
Definition OneDim.h:292
size_t points()
Total number of points.
Definition OneDim.h:152
vector< size_t > m_nvars
Number of variables at each point, across all domains.
Definition OneDim.h:272
int m_nevals
Number of calls to eval()
Definition OneDim.h:289
bool m_init
Indicates whether one-time initialization for each domain has been completed.
Definition OneDim.h:268
void writeStats(int printTime=1)
Write statistics about the number of iterations and Jacobians at each grid level.
Definition OneDim.cpp:141
void clearStats()
Clear saved statistics.
Definition OneDim.cpp:176
size_t m_pts
Total number of points.
Definition OneDim.h:283
vector< int > m_jacEvals
Number of Jacobian evaluations on this grid.
Definition OneDim.h:293
void initTimeInteg(double dt, double *x) override
Prepare for time stepping beginning with solution x and timestep dt.
Definition OneDim.cpp:332
vector< int > m_timeSteps
Number of time steps taken in each call to solve() (for example, for each successive grid refinement)
Definition OneDim.h:306
Domain1D & domain(size_t i) const
Return a reference to domain i.
Definition OneDim.h:78
vector< shared_ptr< Domain1D > > m_dom
All domains comprising the system.
Definition OneDim.h:259
Domain1D * left()
Pointer to left-most domain (first added).
Definition OneDim.h:115
double lowerBound(size_t i) const override
Get the lower bound for global component i in the state vector.
Definition OneDim.cpp:71
int m_nsteps
Number of time steps taken in the current call to solve()
size_t m_size
Solution vector size
virtual void resize()
Call to set the size of internal data structures after first defining the system or if the problem si...
double m_jacobianAbsPerturb
Absolute perturbation of each component in finite difference Jacobian.
size_t size() const
Total solution vector length;.
double rdt() const
Reciprocal of the time step.
virtual void initTimeInteg(double dt, double *x)
Prepare for time stepping beginning with solution x and timestep dt.
double m_rdt
Reciprocal of time step.
double m_jacobianThreshold
Threshold for ignoring small elements in Jacobian.
shared_ptr< SystemJacobian > m_jac
Jacobian evaluator.
shared_ptr< vector< double > > m_state
Solution vector.
vector< int > m_mask
Transient mask.
Func1 * m_interrupt
Function called at the start of every call to eval.
size_t m_bw
Jacobian bandwidth.
virtual void setSteadyMode()
Prepare to solve the steady-state problem.
double m_jacobianRelPerturb
Relative perturbation of each component in finite difference Jacobian.
vector< double > m_work1
Work arrays used during Jacobian evaluation.
MultiJac & jacobian()
Return a reference to the Jacobian evaluator of an OneDim object.
Definition OneDim.cpp:129
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Definition global.h:176
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
shared_ptr< SystemJacobian > newSystemJacobian(const string &type)
Create a SystemJacobian object of the specified type.
void warn_deprecated(const string &source, const AnyBase &node, const string &message)
A deprecation warning for syntax in an input file.
Definition AnyMap.cpp:1997