21 m_newt = make_unique<MultiNewton>(1);
27 m_newt = make_unique<MultiNewton>(1);
28 m_state = make_shared<vector<double>>();
29 for (
auto& dom : domains) {
42 for (
size_t n = 0; n <
m_dom.size(); n++) {
43 if (
domain(n).
id() == name) {
47 throw CanteraError(
"OneDim::domainIndex",
"no domain named >>"+name+
"<<");
68 size_t n =
m_dom.size();
70 m_dom.back()->append(d.get());
83 d->setContainer(
this,
m_dom.size()-1);
109 writelog(
"\nStatistics:\n\n Grid Timesteps Functions Time Jacobians Time\n");
111 for (
size_t i = 0; i < n; i++) {
113 writelog(
"{:5d} {:5d} {:6d} {:9.4f} {:5d} {:9.4f}\n",
117 writelog(
"{:5d} {:5d} {:6d} NA {:5d} NA\n",
126 int nev =
m_jac->nEvals();
164 for (
size_t i = 0; i <
nDomains(); i++) {
165 const auto& d =
m_dom[i];
167 size_t np = d->nPoints();
168 size_t nv = d->nComponents();
169 for (
size_t n = 0; n < np; n++) {
179 size_t bw1 = d->bandwidth();
181 bw1 = std::max<size_t>(2*d->nComponents(), 1) - 1;
188 size_t bw2 =
m_dom[i-1]->bandwidth();
190 bw2 =
m_dom[i-1]->nComponents();
192 bw2 += d->nComponents() - 1;
195 m_size = d->loc() + d->size();
204 m_jac = make_unique<MultiJac>(*
this);
212 m_jac->eval(x, xnew, 0.0);
216 return m_newt->solve(x, xnew, *
this, *
m_jac, loglevel);
221 double rdt_save =
m_rdt;
225 m_jac->eval(x, rsd, 0.0);
241void OneDim::eval(
size_t j,
double* x,
double* r,
double rdt,
int count)
243 clock_t t0 = clock();
256 for (
const auto& d :
m_bulk) {
267 clock_t t1 = clock();
277 for (
size_t i = 0; i <
m_size; i++) {
278 ss = std::max(fabs(r[i]),ss);
285 double rdt_old =
m_rdt;
337 int successiveFailures = 0;
341 writelog(
"\n============================");
342 writelog(
"\n{:<5s} {:<11s} {:<7s}\n",
"step",
"dt (s)",
"log(ss)");
343 writelog(
"============================");
348 writelog(
"\n{:<5d} {:<6.4e} {:>7.4f}", n, dt, log10(ss));
349 }
else if (loglevel > 1) {
351 writelog(
"\nTimestep ({}) dt= {:<11.4e} log(ss)= {:<7.4f}", n, dt, log10(ss));
357 int j0 =
m_jac->nEvals();
360 int status =
solve(x, r, loglevel);
366 writelog(
"\nTimestep ({}) succeeded", n);
368 successiveFailures = 0;
373 if (
m_jac->nEvals() == j0) {
379 dt = std::min(dt,
m_tmax);
382 "Took maximum number of timesteps allowed ({}) without "
388 successiveFailures++;
391 }
else if (loglevel > 1) {
392 writelog(
"\nTimestep ({}) failed", n);
394 if (successiveFailures > 2) {
395 debuglog(
"--> Resetting negative species concentrations", loglevel);
397 successiveFailures = 0;
399 debuglog(
"--> Reducing timestep", loglevel);
402 string err_msg = fmt::format(
403 "Time integration failed. Minimum timestep ({}) reached.",
m_tmin);
413 writelog(
"\n{:<5d} {:<6.4e} {:>7.4f}", n, dt, log10(ss));
414 writelog(
"\n============================");
415 }
else if (loglevel > 1) {
417 writelog(
"\nTimestep ({}) dt= {:<11.4e} log10(ss)= {:<7.4f}\n", n, dt, log10(ss));
427 for (
auto dom :
m_dom) {
428 dom->resetBadValues(x);
Base class for exceptions thrown by Cantera classes.
Base class for one-dimensional domains.
size_t nComponents() const
Number of components at each grid point.
Domain1D * left() const
Return a pointer to the left neighbor.
string id() const
Returns the identifying tag for this domain.
Domain1D * right() const
Return a pointer to the right neighbor.
virtual string componentName(size_t n) const
Name of component n. May be overloaded.
virtual void init()
Initialize.
void initTimeInteg(double dt, const double *x0)
Performs the setup required before starting a time-stepping solution.
void setSteadyMode()
Set the internally-stored reciprocal of the time step to 0.0, which is used to indicate that the prob...
virtual size_t loc(size_t j=0) const
Location of the start of the local solution vector in the global solution vector.
virtual double eval(double t) const
Evaluate the function.
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplie...
Newton iterator for multi-domain, one-dimensional problems.
void setOptions(int maxJacAge=5)
Set options.
int solve(double *x0, double *x1, int loglevel)
Solve F(x) = 0, where F(x) is the multi-domain residual function.
size_t start(size_t i) const
The index of the start of domain i in the solution vector.
int m_nsteps
Number of time steps taken in the current call to solve()
void init()
Initialize all domains.
size_t m_size
solution vector size
int m_nsteps_max
Maximum number of timesteps allowed per call to solve()
virtual void resize()
Call after one or more grids has changed size, for example after being refined.
void resetBadValues(double *x)
Reset values such as negative species concentrations in each domain.
void saveStats()
Save statistics on function and Jacobian evaluation, and reset the counters.
unique_ptr< MultiNewton > m_newt
Newton iterator.
size_t size() const
Total solution vector length;.
void eval(size_t j, double *x, double *r, double rdt=-1.0, int count=1)
Evaluate the multi-domain residual function.
double ssnorm(double *x, double *r)
Steady-state max norm (infinity norm) of the residual evaluated using solution x.
void addDomain(shared_ptr< Domain1D > d)
Add a domain. Domains are added left-to-right.
double rdt() const
Reciprocal of the time step.
void initTimeInteg(double dt, double *x)
Prepare for time stepping beginning with solution x and timestep dt.
size_t nDomains() const
Number of domains.
vector< double > m_jacElapsed
Time [s] spent evaluating Jacobians on this grid.
Domain1D * right()
Pointer to right-most domain (last added).
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...
double m_rdt
reciprocal of time step
shared_ptr< vector< double > > m_state
Solution vector.
vector< double > m_funcElapsed
Time [s] spent on residual function evaluations on this grid (not counting evaluations used to constr...
vector< shared_ptr< Domain1D > > m_connect
All connector and boundary domains.
vector< int > m_mask
Transient mask. See transientMask().
Func1 * m_interrupt
Function called at the start of every call to eval.
vector< shared_ptr< Domain1D > > m_bulk
All bulk/flow domains.
vector< int > m_funcEvals
Number of residual function evaluations on this grid (not counting evaluations used to construct Jaco...
vector< size_t > m_loc
Location in the state vector of the first component of each point, across all domains.
unique_ptr< MultiJac > m_jac
Jacobian evaluator.
double m_evaltime
Total time [s] spent in eval()
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
int m_ts_jac_age
Maximum age of the Jacobian in time-stepping mode.
size_t domainIndex(const string &name)
Get the index of the domain named name.
double timeStep(int nsteps, double dt, double *x, double *r, int loglevel)
Take time steps using Backward Euler.
Domain1D * pointDomain(size_t i)
Return a pointer to the domain global point i belongs to.
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...
vector< size_t > m_gridpts
Number of grid points in this grid.
void evalSSJacobian(double *x, double *rsd)
Evaluate the steady-state Jacobian, accessible via jacobian()
double m_tmin
minimum timestep size
double m_tmax
maximum timestep size
int m_ss_jac_age
Maximum age of the Jacobian in steady-state mode.
vector< size_t > m_nvars
Number of variables at each point, across all domains.
int m_nevals
Number of calls to eval()
bool m_init
Indicates whether one-time initialization for each domain has been completed.
void writeStats(int printTime=1)
Write statistics about the number of iterations and Jacobians at each grid level.
void clearStats()
Clear saved statistics.
void setSteadyMode()
Prepare to solve the steady-state problem.
size_t m_pts
Total number of points.
OneDim()
Default constructor.
vector< int > m_jacEvals
Number of Jacobian evaluations on this grid.
Func1 * m_time_step_callback
User-supplied function called after each successful timestep.
vector< int > m_timeSteps
Number of time steps taken in each call to solve() (for example, for each successive grid refinement)
Domain1D & domain(size_t i) const
Return a reference to domain i.
vector< shared_ptr< Domain1D > > m_dom
All domains comprising the system.
Domain1D * left()
Pointer to left-most domain (first added).
MultiNewton & newton()
Return a reference to the Newton iterator.
MultiJac & jacobian()
Return a reference to the Jacobian evaluator of an OneDim object.
void debuglog(const string &msg, int loglevel)
Write a message to the log only if loglevel > 0.
void writelog(const string &fmt, const Args &... args)
Write a formatted message to the screen.
Namespace for the Cantera kernel.
const size_t npos
index returned by functions to indicate "no position"
const double Tiny
Small number to compare differences of mole fractions against.
offset
Offsets of solution components in the 1D solution array.