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