Cantera
3.1.0a1
|
One-dimensional simulations. More...
#include <Sim1D.h>
One-dimensional simulations.
Class Sim1D extends class OneDim by storing the solution vector, and by adding a hybrid Newton/time-stepping solver.
Public Member Functions | |
Sim1D () | |
Default constructor. More... | |
Sim1D (vector< shared_ptr< Domain1D >> &domains) | |
Standard constructor. More... | |
void | setTimeStep (double stepsize, size_t n, const int *tsteps) |
void | solve (int loglevel=0, bool refine_grid=true) |
void | eval (double 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 (double t) |
Add node for fixed temperature point of freely propagating flame. More... | |
double | fixedTemperature () |
Return temperature at the point used to fix the flame location. More... | |
double | fixedTemperatureLocation () |
Return location of the point where temperature is fixed. More... | |
void | setRefineCriteria (int dom=-1, double ratio=10.0, double slope=0.8, double curve=0.8, double prune=-0.1) |
Set grid refinement criteria. More... | |
vector< double > | 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 | 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 () |
double | jacobian (int i, int j) |
void | evalSSJacobian () |
void | solveAdjoint (const double *b, double *lambda) |
Solve the equation \( J^T \lambda = b \). More... | |
void | resize () override |
Call after one or more grids has changed size, for example 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 string &component, vector< double > &locs, vector< double > &vals) |
Set initial guess for one component for all domains. More... | |
void | setValue (size_t dom, size_t comp, size_t localPoint, double value) |
Set a single value in the solution vector. More... | |
double | value (size_t dom, size_t comp, size_t localPoint) const |
Get one entry in the solution vector. More... | |
double | workValue (size_t dom, size_t comp, size_t localPoint) const |
void | setProfile (size_t dom, size_t comp, const vector< double > &pos, const vector< double > &values) |
Specify a profile for one component of one domain. More... | |
void | setFlatProfile (size_t dom, size_t comp, double v) |
Set component 'comp' of domain 'dom' to value 'v' at all points. More... | |
Logging, saving and restoring of solutions | |
void | show (std::ostream &s) |
Output information on current solution for all domains to stream. More... | |
void | show () |
Show logging information on current solution for all domains. More... | |
void | save (const string &fname, const string &name, const string &desc, bool overwrite=false, int compression=0, const string &basis="") |
Save current simulation data to a container file or CSV format. More... | |
void | saveResidual (const string &fname, const string &name, const string &desc, bool overwrite=false, int compression=0) |
Save the residual of the current solution to a container file. More... | |
AnyMap | restore (const string &fname, const string &name) |
Retrieve data and settings from a previously saved simulation. More... | |
Public Member Functions inherited from OneDim | |
OneDim (vector< shared_ptr< 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 (shared_ptr< Domain1D > d) |
Add a domain. Domains are added left-to-right. More... | |
MultiJac & | jacobian () |
Return a reference to the Jacobian evaluator of an OneDim object. More... | |
MultiNewton & | newton () |
Return a reference to the Newton iterator. More... | |
int | solve (double *x0, double *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 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< 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 solution vector. More... | |
size_t | bandwidth () const |
Jacobian bandwidth. More... | |
void | init () |
Initialize all domains. More... | |
size_t | points () |
Total number of points. More... | |
double | ssnorm (double *x, double *r) |
Steady-state max norm (infinity norm) of the residual evaluated using solution x. More... | |
double | rdt () const |
Reciprocal of the time step. More... | |
void | initTimeInteg (double dt, double *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 () |
Prepare to solve the steady-state problem. More... | |
void | eval (size_t j, double *x, double *r, double 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) |
Take time steps using Backward Euler. More... | |
void | resetBadValues (double *x) |
void | writeStats (int printTime=1) |
Write statistics about the number of iterations and Jacobians at each grid level. More... | |
void | setMinTimeStep (double tmin) |
void | setMaxTimeStep (double tmax) |
void | setTimeStepFactor (double 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 vector< size_t > & | gridSizeStats () |
Return total grid size in each call to solve() More... | |
const vector< double > & | jacobianTimeStats () |
Return CPU time spent evaluating Jacobians in each call to solve() More... | |
const vector< double > & | 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< double > | m_xlast_ts |
the solution vector after the last successful timestepping More... | |
vector< double > | m_xlast_ss |
the solution vector after the last successful steady-state solve (stored before grid refinement) More... | |
vector< vector< double > > | m_grid_last_ss |
the grids for each domain after the last successful steady-state solve (stored before grid refinement) More... | |
vector< double > | m_xnew |
a work array used to hold the residual or the new solution More... | |
double | 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 | |
double | m_tmin = 1e-16 |
minimum timestep size More... | |
double | m_tmax = 1e+08 |
maximum timestep size More... | |
double | m_tfactor = 0.5 |
factor time step is multiplied by if time stepping fails ( < 1 ) More... | |
shared_ptr< vector< double > > | m_state |
Solution vector. More... | |
unique_ptr< MultiJac > | m_jac |
Jacobian evaluator. More... | |
unique_ptr< MultiNewton > | m_newt |
Newton iterator. More... | |
double | m_rdt = 0.0 |
reciprocal of time step More... | |
bool | m_jac_ok = false |
if true, Jacobian is current More... | |
size_t | m_bw = 0 |
Jacobian bandwidth. More... | |
size_t | m_size = 0 |
solution vector size More... | |
vector< shared_ptr< Domain1D > > | m_dom |
vector< shared_ptr< Domain1D > > | m_connect |
vector< shared_ptr< Domain1D > > | m_bulk |
bool | m_init = false |
vector< size_t > | m_nvars |
vector< size_t > | m_loc |
vector< int > | m_mask |
size_t | m_pts = 0 |
int | m_ss_jac_age = 20 |
int | m_ts_jac_age = 20 |
Func1 * | m_interrupt = nullptr |
Function called at the start of every call to eval. More... | |
Func1 * | m_time_step_callback = nullptr |
User-supplied function called after each successful timestep. More... | |
int | m_nsteps = 0 |
Number of time steps taken in the current call to solve() More... | |
int | m_nsteps_max = 500 |
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 (double *x, double *xnew) |
|
inline |
Standard constructor.
domains | A vector of shared pointers to the domains to be linked together. The domain pointers must be entered in left-to-right order — that is, the pointer to the leftmost domain is domain[0], the pointer to the domain to its right is domain[1], etc. |
void setInitialGuess | ( | const string & | component, |
vector< double > & | locs, | ||
vector< double > & | vals | ||
) |
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. |
void setValue | ( | size_t | dom, |
size_t | comp, | ||
size_t | localPoint, | ||
double | value | ||
) |
double value | ( | size_t | dom, |
size_t | comp, | ||
size_t | localPoint | ||
) | const |
void setProfile | ( | size_t | dom, |
size_t | comp, | ||
const vector< double > & | pos, | ||
const vector< double > & | values | ||
) |
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.
void setFlatProfile | ( | size_t | dom, |
size_t | comp, | ||
double | v | ||
) |
void show | ( | std::ostream & | s | ) |
void show | ( | ) |
void save | ( | const string & | fname, |
const string & | name, | ||
const string & | desc, | ||
bool | overwrite = false , |
||
int | compression = 0 , |
||
const string & | basis = "" |
||
) |
Save current simulation data to a container file or CSV format.
In order to save the content of a Sim1D object, individual domains are converted to SolutionArray objects and saved using the SolutionArray::save() method. For HDF and YAML output, all domains are written to a single container file with shared header information. Simulation settings of individual domains are preserved as meta data of the corresponding SolutionArray objects. For CSV files, only state and auxiliary data of the main 1D domain are saved.
The complete state of the current object can be restored from HDF and YAML container files using the restore() method, while individual domains can be loaded using SolutionArray::restore() for further analysis. While CSV do not contain complete information, they can still be used for setting initial states of individual simulation objects for some Cantera API's.
fname | Name of output file (CSV, YAML or HDF) |
name | Identifier of storage location within the container file; this node/group contains header information and multiple subgroups holding domain-specific SolutionArray data (YAML/HDF only) |
desc | Custom comment describing the dataset to be stored (YAML/HDF only) |
overwrite | Force overwrite if file/name exists; optional (default=false) |
compression | Compression level (0-9); optional (default=0; HDF only) |
basis | Output mass ("Y"/"mass") or mole ("X"/"mole") fractions; if not specified (default=""), the native basis of the underlying ThermoPhase manager is used - |
void saveResidual | ( | const string & | fname, |
const string & | name, | ||
const string & | desc, | ||
bool | overwrite = false , |
||
int | compression = 0 |
||
) |
Save the residual of the current solution to a container file.
fname | Name of output container file |
name | Identifier of solution within the container file |
desc | Description of the solution |
overwrite | Force overwrite if name exists; optional (default=false) |
compression | Compression level (optional; HDF only) |
AnyMap restore | ( | const string & | fname, |
const string & | name | ||
) |
Retrieve data and settings from a previously saved simulation.
This method restores a simulation object from YAML or HDF data previously saved using the save() method.
fname | Name of container file (YAML or HDF) |
name | Identifier of location within the container file; this node/group contains header information and subgroups with domain-specific SolutionArray data |
int refine | ( | int | loglevel = 0 | ) |
int setFixedTemperature | ( | double | t | ) |
double fixedTemperature | ( | ) |
double fixedTemperatureLocation | ( | ) |
void setRefineCriteria | ( | int | dom = -1 , |
double | ratio = 10.0 , |
||
double | slope = 0.8 , |
||
double | curve = 0.8 , |
||
double | 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.
vector< double > getRefineCriteria | ( | int | dom | ) |
Get the grid refinement criteria.
dom must be greater than or equal to zero (that is, the domain must be specified).
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.
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. |
void setGridMin | ( | int | dom, |
double | gridmin | ||
) |
void restoreTimeSteppingSolution | ( | ) |
void restoreSteadySolution | ( | ) |
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} \]
|
overridevirtual |
|
inline |
|
private |
|
private |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |