Cantera
2.4.0
|
One-dimensional simulations. More...
#include <Sim1D.h>
Public Member Functions | |
Sim1D () | |
Default constructor. More... | |
Sim1D (std::vector< Domain1D *> &domains) | |
Standard constructor. More... | |
void | save (const std::string &fname, const std::string &id, const std::string &desc, int loglevel=1) |
void | saveResidual (const std::string &fname, const std::string &id, const std::string &desc, int loglevel=1) |
void | showSolution (std::ostream &s) |
Print to stream s the current solution for all domains. More... | |
void | showSolution () |
const doublereal * | solution () |
void | setTimeStep (double stepsize, size_t n, const int *tsteps) |
void | solve (int loglevel=0, bool refine_grid=true) |
void | eval (doublereal rdt=-1.0, int count=1) |
void | getResidual (double rdt, double *resid) |
int | refine (int loglevel=0) |
Refine the grid in all domains. More... | |
int | setFixedTemperature (doublereal t) |
Add node for fixed temperature point of freely propagating flame. More... | |
void | setRefineCriteria (int dom=-1, doublereal ratio=10.0, doublereal slope=0.8, doublereal curve=0.8, doublereal prune=-0.1) |
Set grid refinement criteria. More... | |
vector_fp | getRefineCriteria (int dom) |
Get the grid refinement criteria. More... | |
void | setMaxGridPoints (int dom, int npoints) |
Set the maximum number of grid points in the domain. More... | |
size_t | maxGridPoints (size_t dom) |
Get the maximum number of grid points in this domain. More... | |
void | setGridMin (int dom, double gridmin) |
Set the minimum grid spacing in the specified domain(s). More... | |
void | restore (const std::string &fname, const std::string &id, int loglevel=2) |
Initialize the solution with a previously-saved solution. More... | |
void | restoreTimeSteppingSolution () |
Set the current solution vector to the last successful time-stepping solution. More... | |
void | restoreSteadySolution () |
Set the current solution vector and grid to the last successful steady- state solution. More... | |
void | getInitialSoln () |
void | setSolution (const doublereal *soln) |
const doublereal * | solution () const |
doublereal | jacobian (int i, int j) |
void | evalSSJacobian () |
void | solveAdjoint (const double *b, double *lambda) |
Solve the equation \( J^T \lambda = b \). More... | |
virtual void | resize () |
Call after one or more grids has changed size, e.g. after being refined. More... | |
void | setSteadyCallback (Func1 *callback) |
Set a function that will be called after each successful steady-state solve, before regridding. More... | |
Setting initial values | |
These methods are used to set the initial values of solution components. | |
void | setInitialGuess (const std::string &component, vector_fp &locs, vector_fp &vals) |
Set initial guess for one component for all domains. More... | |
void | setValue (size_t dom, size_t comp, size_t localPoint, doublereal value) |
Set a single value in the solution vector. More... | |
doublereal | value (size_t dom, size_t comp, size_t localPoint) const |
Get one entry in the solution vector. More... | |
doublereal | workValue (size_t dom, size_t comp, size_t localPoint) const |
void | setProfile (size_t dom, size_t comp, const vector_fp &pos, const vector_fp &values) |
Specify a profile for one component of one domain. More... | |
void | setFlatProfile (size_t dom, size_t comp, doublereal v) |
Set component 'comp' of domain 'dom' to value 'v' at all points. More... | |
Public Member Functions inherited from OneDim | |
OneDim (std::vector< Domain1D *> domains) | |
Construct a OneDim container for the domains in the list domains. More... | |
OneDim (const OneDim &)=delete | |
OneDim & | operator= (const OneDim &)=delete |
void | addDomain (Domain1D *d) |
Add a domain. Domains are added left-to-right. More... | |
MultiJac & | jacobian () |
Return a reference to the Jacobian evaluator. More... | |
MultiNewton & | newton () |
Return a reference to the Newton iterator. More... | |
int | solve (doublereal *x0, doublereal *x1, int loglevel) |
Solve F(x) = 0, where F(x) is the multi-domain residual function. More... | |
size_t | nDomains () const |
Number of domains. More... | |
Domain1D & | domain (size_t i) const |
Return a reference to domain i. More... | |
size_t | domainIndex (const std::string &name) |
void | checkDomainIndex (size_t n) const |
Check that the specified domain index is in range. More... | |
void | checkDomainArraySize (size_t nn) const |
Check that an array size is at least nDomains(). More... | |
size_t | start (size_t i) const |
The index of the start of domain i in the solution vector. More... | |
size_t | size () const |
Total solution vector length;. More... | |
Domain1D * | left () |
Pointer to left-most domain (first added). More... | |
Domain1D * | right () |
Pointer to right-most domain (last added). More... | |
size_t | nVars (size_t jg) |
Number of solution components at global point jg. More... | |
size_t | loc (size_t jg) |
Location in the solution vector of the first component of global point jg. More... | |
std::tuple< std::string, size_t, std::string > | component (size_t i) |
Return the domain, local point index, and component name for the i-th component of the global solution vector. More... | |
size_t | bandwidth () const |
Jacobian bandwidth. More... | |
void | init () |
size_t | points () |
Total number of points. More... | |
doublereal | ssnorm (doublereal *x, doublereal *r) |
Steady-state max norm (infinity norm) of the residual evaluated using solution x. More... | |
doublereal | rdt () const |
Reciprocal of the time step. More... | |
void | initTimeInteg (doublereal dt, doublereal *x) |
Prepare for time stepping beginning with solution x and timestep dt. More... | |
bool | transient () const |
True if transient mode. More... | |
bool | steady () const |
True if steady mode. More... | |
void | setSteadyMode () |
void | eval (size_t j, double *x, double *r, doublereal rdt=-1.0, int count=1) |
Evaluate the multi-domain residual function. More... | |
Domain1D * | pointDomain (size_t i) |
Return a pointer to the domain global point i belongs to. More... | |
vector_int & | transientMask () |
double | timeStep (int nsteps, double dt, double *x, double *r, int loglevel) |
void | resetBadValues (double *x) |
void | writeStats (int printTime=1) |
Write statistics about the number of iterations and Jacobians at each grid level. More... | |
void | save (const std::string &fname, std::string id, const std::string &desc, doublereal *sol, int loglevel) |
void | setMinTimeStep (doublereal tmin) |
void | setMaxTimeStep (doublereal tmax) |
void | setTimeStepFactor (doublereal tfactor) |
void | setMaxTimeStepCount (int nmax) |
Set the maximum number of timeteps allowed before successful steady-state solve. More... | |
int | maxTimeStepCount () const |
Return the maximum number of timeteps allowed before successful steady-state solve. More... | |
void | setJacAge (int ss_age, int ts_age=-1) |
void | saveStats () |
Save statistics on function and Jacobian evaluation, and reset the counters. More... | |
void | clearStats () |
Clear saved statistics. More... | |
const std::vector< size_t > & | gridSizeStats () |
Return total grid size in each call to solve() More... | |
const vector_fp & | jacobianTimeStats () |
Return CPU time spent evaluating Jacobians in each call to solve() More... | |
const vector_fp & | evalTimeStats () |
Return CPU time spent on non-Jacobian function evaluations in each call to solve() More... | |
const vector_int & | jacobianCountStats () |
Return number of Jacobian evaluations made in each call to solve() More... | |
const vector_int & | evalCountStats () |
Return number of non-Jacobian function evaluations made in each call to solve() More... | |
const vector_int & | timeStepStats () |
Return number of time steps taken in each call to solve() More... | |
void | setInterrupt (Func1 *interrupt) |
Set a function that will be called every time eval is called. More... | |
void | setTimeStepCallback (Func1 *callback) |
Set a function that will be called after each successful timestep. More... | |
Protected Attributes | |
vector_fp | m_x |
the solution vector More... | |
vector_fp | m_xlast_ts |
the solution vector after the last successful timestepping More... | |
vector_fp | m_xlast_ss |
the solution vector after the last successful steady-state solve (stored before grid refinement) More... | |
std::vector< vector_fp > | m_grid_last_ss |
the grids for each domain after the last successful steady-state solve (stored before grid refinement) More... | |
vector_fp | m_xnew |
a work array used to hold the residual or the new solution More... | |
doublereal | m_tstep |
timestep More... | |
vector_int | m_steps |
array of number of steps to take before re-attempting the steady-state solution More... | |
Func1 * | m_steady_callback |
User-supplied function called after a successful steady-state solve. More... | |
Protected Attributes inherited from OneDim | |
doublereal | m_tmin |
minimum timestep size More... | |
doublereal | m_tmax |
maximum timestep size More... | |
doublereal | m_tfactor |
factor time step is multiplied by if time stepping fails ( < 1 ) More... | |
std::unique_ptr< MultiJac > | m_jac |
Jacobian evaluator. More... | |
std::unique_ptr< MultiNewton > | m_newt |
Newton iterator. More... | |
doublereal | m_rdt |
reciprocal of time step More... | |
bool | m_jac_ok |
if true, Jacobian is current More... | |
size_t | m_bw |
Jacobian bandwidth. More... | |
size_t | m_size |
solution vector size More... | |
std::vector< Domain1D * > | m_dom |
std::vector< Domain1D * > | m_connect |
std::vector< Domain1D * > | m_bulk |
bool | m_init |
std::vector< size_t > | m_nvars |
std::vector< size_t > | m_loc |
vector_int | m_mask |
size_t | m_pts |
doublereal | m_solve_time |
int | m_ss_jac_age |
int | m_ts_jac_age |
Func1 * | m_interrupt |
Function called at the start of every call to eval. More... | |
Func1 * | m_time_step_callback |
User-supplied function called after each successful timestep. More... | |
int | m_nsteps |
Number of time steps taken in the current call to solve() More... | |
int | m_nsteps_max |
Maximum number of timesteps allowed per call to solve() More... | |
Private Member Functions | |
void | finalize () |
Calls method _finalize in each domain. More... | |
int | newtonSolve (int loglevel) |
Wrapper around the Newton solver. More... | |
Additional Inherited Members | |
Protected Member Functions inherited from OneDim | |
void | evalSSJacobian (doublereal *x, doublereal *xnew) |
One-dimensional simulations.
Class Sim1D extends class OneDim by storing the solution vector, and by adding a hybrid Newton/time-stepping solver.
|
inline |
Standard constructor.
domains | A vector of pointers to the domains to be linked together. The domain pointers must be entered in left-to-right order — i.e., the pointer to the leftmost domain is domain[0], the pointer to the domain to its right is domain[1], etc. |
Definition at line 21 of file Sim1D.cpp.
References Domain1D::_getInitialSoln(), OneDim::domain(), Sim1D::m_steps, Sim1D::m_tstep, Sim1D::m_x, OneDim::nDomains(), Sim1D::resize(), and OneDim::start().
Set initial guess for one component for all domains.
component | component name |
locs | A vector of relative positions, beginning with 0.0 at the left of the domain, and ending with 1.0 at the right of the domain. |
vals | A vector of values corresponding to the relative position locations. |
Definition at line 37 of file Sim1D.cpp.
References OneDim::component(), Domain1D::componentName(), OneDim::domain(), Domain1D::nComponents(), OneDim::nDomains(), and Sim1D::setProfile().
void setValue | ( | size_t | dom, |
size_t | comp, | ||
size_t | localPoint, | ||
doublereal | value | ||
) |
Set a single value in the solution vector.
dom | domain number, beginning with 0 for the leftmost domain. |
comp | component number |
localPoint | grid point within the domain, beginning with 0 for the leftmost grid point in the domain. |
value | the value. |
Definition at line 50 of file Sim1D.cpp.
References AssertThrowMsg, OneDim::domain(), Domain1D::loc(), Sim1D::m_x, and Sim1D::value().
Referenced by Sim1D::setFlatProfile(), and Sim1D::setProfile().
doublereal value | ( | size_t | dom, |
size_t | comp, | ||
size_t | localPoint | ||
) | const |
Get one entry in the solution vector.
dom | domain number, beginning with 0 for the leftmost domain. |
comp | component number |
localPoint | grid point within the domain, beginning with 0 for the leftmost grid point in the domain. |
Definition at line 58 of file Sim1D.cpp.
References AssertThrowMsg, OneDim::domain(), Domain1D::loc(), and Sim1D::m_x.
Referenced by Sim1D::refine(), Sim1D::setFixedTemperature(), and Sim1D::setValue().
Specify a profile for one component of one domain.
dom | domain number, beginning with 0 for the leftmost domain. |
comp | component number |
pos | A vector of relative positions, beginning with 0.0 at the left of the domain, and ending with 1.0 at the right of the domain. |
values | A vector of values corresponding to the relative position locations. |
Note that the vector pos and values can have lengths different than the number of grid points, but their lengths must be equal. The values at the grid points will be linearly interpolated based on the (pos, values) specification.
Definition at line 74 of file Sim1D.cpp.
References OneDim::domain(), Cantera::linearInterp(), Domain1D::nPoints(), and Sim1D::setValue().
Referenced by Sim1D::setInitialGuess().
void setFlatProfile | ( | size_t | dom, |
size_t | comp, | ||
doublereal | v | ||
) |
Set component 'comp' of domain 'dom' to value 'v' at all points.
Definition at line 140 of file Sim1D.cpp.
References OneDim::domain(), Domain1D::nPoints(), and Sim1D::setValue().
void showSolution | ( | std::ostream & | s | ) |
Print to stream s the current solution for all domains.
Definition at line 148 of file Sim1D.cpp.
References OneDim::domain(), Domain1D::domainType(), Sim1D::m_x, OneDim::nDomains(), and OneDim::start().
int refine | ( | int | loglevel = 0 | ) |
Refine the grid in all domains.
Definition at line 341 of file Sim1D.cpp.
References OneDim::domain(), Sim1D::finalize(), Sim1D::m_grid_last_ss, Sim1D::m_x, Sim1D::m_xlast_ss, Domain1D::nComponents(), OneDim::nDomains(), Domain1D::nPoints(), Domain1D::refiner(), Sim1D::resize(), Domain1D::setupGrid(), OneDim::start(), Sim1D::value(), and Cantera::writelog().
int setFixedTemperature | ( | doublereal | t | ) |
Add node for fixed temperature point of freely propagating flame.
Definition at line 427 of file Sim1D.cpp.
References OneDim::domain(), Domain1D::domainType(), Sim1D::finalize(), StFlow::m_tfixed, Sim1D::m_x, StFlow::m_zfixed, Domain1D::nComponents(), OneDim::nDomains(), Domain1D::nPoints(), Sim1D::resize(), Domain1D::setupGrid(), and Sim1D::value().
void setRefineCriteria | ( | int | dom = -1 , |
doublereal | ratio = 10.0 , |
||
doublereal | slope = 0.8 , |
||
doublereal | curve = 0.8 , |
||
doublereal | prune = -0.1 |
||
) |
Set grid refinement criteria.
If dom >= 0, then the settings apply only to the specified domain. If dom < 0, the settings are applied to each domain.
Definition at line 514 of file Sim1D.cpp.
References OneDim::domain(), OneDim::nDomains(), Domain1D::refiner(), and Refiner::setCriteria().
vector_fp getRefineCriteria | ( | int | dom | ) |
Get the grid refinement criteria.
dom must be greater than or equal to zero (i.e., the domain must be specified).
Definition at line 528 of file Sim1D.cpp.
References OneDim::domain(), Refiner::getCriteria(), and Domain1D::refiner().
void setMaxGridPoints | ( | int | dom, |
int | npoints | ||
) |
Set the maximum number of grid points in the domain.
If dom >= 0, then the settings apply only to the specified domain. If dom < 0, the settings are applied to each domain.
Definition at line 552 of file Sim1D.cpp.
References OneDim::domain(), OneDim::nDomains(), Domain1D::refiner(), and Refiner::setMaxPoints().
size_t maxGridPoints | ( | size_t | dom | ) |
Get the maximum number of grid points in this domain.
dom | domain number, beginning with 0 for the leftmost domain. |
Definition at line 565 of file Sim1D.cpp.
References OneDim::domain(), Refiner::maxPoints(), and Domain1D::refiner().
void setGridMin | ( | int | dom, |
double | gridmin | ||
) |
Set the minimum grid spacing in the specified domain(s).
dom | Domain index. If dom == -1, the specified spacing is applied to all domains. |
gridmin | The minimum allowable grid spacing [m] |
Definition at line 539 of file Sim1D.cpp.
References OneDim::domain(), OneDim::nDomains(), Domain1D::refiner(), and Refiner::setGridMin().
void restore | ( | const std::string & | fname, |
const std::string & | id, | ||
int | loglevel = 2 |
||
) |
Initialize the solution with a previously-saved solution.
Definition at line 107 of file Sim1D.cpp.
References XML_Node::build(), OneDim::domain(), Sim1D::finalize(), XML_Node::findID(), XML_Node::getChildren(), Cantera::intValue(), OneDim::loc(), Sim1D::m_x, Sim1D::m_xlast_ts, OneDim::nDomains(), Domain1D::resize(), Sim1D::resize(), Domain1D::restore(), and Cantera::writelog().
void restoreTimeSteppingSolution | ( | ) |
Set the current solution vector to the last successful time-stepping solution.
This can be used to examine the solver progress after a failed integration.
Definition at line 168 of file Sim1D.cpp.
References Sim1D::m_x, and Sim1D::m_xlast_ts.
void restoreSteadySolution | ( | ) |
Set the current solution vector and grid to the last successful steady- state solution.
This can be used to examine the solver progress after a failure during grid refinement.
Definition at line 177 of file Sim1D.cpp.
References OneDim::domain(), Sim1D::m_grid_last_ss, Sim1D::m_x, Sim1D::m_xlast_ss, OneDim::nDomains(), and Domain1D::setupGrid().
void solveAdjoint | ( | const double * | b, |
double * | lambda | ||
) |
Solve the equation \( J^T \lambda = b \).
Here, \( J = \partial f/\partial x \) is the Jacobian matrix of the system of equations \( f(x,p)=0 \). This can be used to efficiently solve for the sensitivities of a scalar objective function \( g(x,p) \) to a vector of parameters \( p \) by solving:
\[ J^T \lambda = \left( \frac{\partial g}{\partial x} \right)^T \]
for \( \lambda \) and then computing:
\[ \left.\frac{dg}{dp}\right|_{f=0} = \frac{\partial g}{\partial p} - \lambda^T \frac{\partial f}{\partial p} \]
|
virtual |
Call after one or more grids has changed size, e.g. after being refined.
Reimplemented from OneDim.
Definition at line 605 of file Sim1D.cpp.
References Sim1D::m_x, Sim1D::m_xnew, OneDim::resize(), and OneDim::size().
Referenced by Sim1D::refine(), Sim1D::restore(), Sim1D::setFixedTemperature(), and Sim1D::Sim1D().
|
inline |
Set a function that will be called after each successful steady-state solve, before regridding.
Intended to be used for observing solver progress for debugging purposes.
Definition at line 218 of file Sim1D.h.
References Sim1D::m_steady_callback.
|
private |
Calls method _finalize in each domain.
Definition at line 197 of file Sim1D.cpp.
References Domain1D::_finalize(), OneDim::domain(), Sim1D::m_x, OneDim::nDomains(), and OneDim::start().
Referenced by Sim1D::refine(), Sim1D::restore(), and Sim1D::setFixedTemperature().
|
private |
Wrapper around the Newton solver.
Definition at line 213 of file Sim1D.cpp.
References Sim1D::m_x, Sim1D::m_xnew, and OneDim::solve().
|
protected |
the solution vector
Definition at line 224 of file Sim1D.h.
Referenced by Sim1D::finalize(), Sim1D::newtonSolve(), Sim1D::refine(), Sim1D::resize(), Sim1D::restore(), Sim1D::restoreSteadySolution(), Sim1D::restoreTimeSteppingSolution(), Sim1D::setFixedTemperature(), Sim1D::setValue(), Sim1D::showSolution(), Sim1D::Sim1D(), and Sim1D::value().
|
protected |
the solution vector after the last successful timestepping
Definition at line 227 of file Sim1D.h.
Referenced by Sim1D::restore(), and Sim1D::restoreTimeSteppingSolution().
|
protected |
the solution vector after the last successful steady-state solve (stored before grid refinement)
Definition at line 231 of file Sim1D.h.
Referenced by Sim1D::refine(), and Sim1D::restoreSteadySolution().
|
protected |
the grids for each domain after the last successful steady-state solve (stored before grid refinement)
Definition at line 235 of file Sim1D.h.
Referenced by Sim1D::refine(), and Sim1D::restoreSteadySolution().
|
protected |
a work array used to hold the residual or the new solution
Definition at line 238 of file Sim1D.h.
Referenced by Sim1D::newtonSolve(), and Sim1D::resize().
|
protected |
|
protected |
array of number of steps to take before re-attempting the steady-state solution
Definition at line 245 of file Sim1D.h.
Referenced by Sim1D::Sim1D().
|
protected |
User-supplied function called after a successful steady-state solve.
Definition at line 248 of file Sim1D.h.
Referenced by Sim1D::setSteadyCallback().