33 OneDim(vector<shared_ptr<Domain1D>>& domains);
61 int solve(
double* x0,
double* x1,
int loglevel);
79 if (n >=
m_dom.size()) {
80 throw IndexError(
"OneDim::checkDomainIndex",
"domains", n,
89 if (
m_dom.size() > nn) {
97 if (
m_dom[i]->nComponents()) {
98 return m_dom[i]->loc();
113 return m_dom[0].get();
118 return m_dom.back().get();
133 std::tuple<string, size_t, string>
component(
size_t i);
157 double ssnorm(
double* x,
double* r);
169 return (
m_rdt != 0.0);
174 return (
m_rdt == 0.0);
196 void eval(
size_t j,
double* x,
double* r,
double rdt=-1.0,
int count = 1);
224 double timeStep(
int nsteps,
double dt,
double* x,
double* r,
int loglevel);
280 void setJacAge(
int ss_age,
int ts_age=-1);
Base class for one-dimensional domains.
Base class for 'functor' classes that evaluate a function of one variable.
An array index is out of range.
Class MultiJac evaluates the Jacobian of a system of equations defined by a residual function supplie...
Newton iterator for multi-domain, one-dimensional problems.
Container class for multiple-domain 1D problems.
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.
void checkDomainIndex(size_t n) const
Check that the specified domain index is in range.
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;.
size_t loc(size_t jg)
Location in the solution vector of the first component of global point jg.
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.
vector< int > & transientMask()
Access the vector indicating which equations contain a transient term.
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...
size_t bandwidth() const
Jacobian bandwidth.
double m_rdt
reciprocal of time step
void setMaxTimeStep(double tmax)
Set the maximum time step allowed during time stepping.
shared_ptr< vector< double > > m_state
Solution vector.
void setMinTimeStep(double tmin)
Set the minimum time step allowed during time stepping.
vector< double > m_funcElapsed
Time [s] spent on residual function evaluations on this grid (not counting evaluations used to constr...
const vector< int > & evalCountStats()
Return number of non-Jacobian function evaluations made in each call to solve()
vector< shared_ptr< Domain1D > > m_connect
All connector and boundary domains.
vector< int > m_mask
Transient mask. See transientMask().
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.
const vector< int > & jacobianCountStats()
Return number of Jacobian evaluations made in each call to solve()
bool transient() const
True if transient mode.
void checkDomainArraySize(size_t nn) const
Check that an array size is at least nDomains().
const vector< size_t > & gridSizeStats()
Return total grid size in each call to solve()
void setMaxTimeStepCount(int nmax)
Set the maximum number of timeteps allowed before successful steady-state solve.
const vector< int > & timeStepStats()
Return number of time steps taken in each call to solve()
void setTimeStepFactor(double tfactor)
Sets a factor by which the time step is reduced if the time stepping fails.
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 nVars(size_t jg)
Number of solution components at global point jg.
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 maxTimeStepCount() const
Return 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.
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
size_t points()
Total number of points.
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.
const vector< double > & jacobianTimeStats()
Return CPU time spent evaluating Jacobians in each call to solve()
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.
const vector< double > & evalTimeStats()
Return CPU time spent on non-Jacobian function evaluations in each call to solve()
Domain1D * left()
Pointer to left-most domain (first added).
bool steady() const
True if steady mode.
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.
MultiJac & jacobian()
Return a reference to the Jacobian evaluator of an OneDim object.
Namespace for the Cantera kernel.