Cantera
3.1.0a1
|
Container class for multiple-domain 1D problems. More...
#include <OneDim.h>
Container class for multiple-domain 1D problems.
Each domain is represented by an instance of Domain1D.
Public Member Functions | |
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... | |
virtual void | resize () |
Call after one or more grids has changed size, for example after being refined. 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 Member Functions | |
void | evalSSJacobian (double *x, double *xnew) |
Protected Attributes | |
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 Attributes | |
int | m_nevals = 0 |
double | m_evaltime = 0 |
vector< size_t > | m_gridpts |
vector< int > | m_jacEvals |
vector< double > | m_jacElapsed |
vector< int > | m_funcEvals |
vector< double > | m_funcElapsed |
vector< int > | m_timeSteps |
Number of time steps taken in each call to solve() (for example, for each successive grid refinement) More... | |
Construct a OneDim container for the domains in the list domains.
Definition at line 24 of file OneDim.cpp.
void addDomain | ( | shared_ptr< Domain1D > | d | ) |
Add a domain. Domains are added left-to-right.
Definition at line 64 of file OneDim.cpp.
MultiNewton & newton | ( | ) |
Return a reference to the Newton iterator.
Definition at line 91 of file OneDim.cpp.
int solve | ( | double * | x0, |
double * | x1, | ||
int | loglevel | ||
) |
Solve F(x) = 0, where F(x) is the multi-domain residual function.
x0 | Starting estimate of solution. |
x1 | Final solution satisfying F(x1) = 0. |
loglevel | Controls amount of diagnostic output. |
Definition at line 212 of file OneDim.cpp.
|
inline |
|
inline |
Check that the specified domain index is in range.
Throws an exception if n is greater than nDomains()-1
|
inline |
Check that an array size is at least nDomains().
Throws an exception if nn is less than nDomains(). Used before calls which take an array pointer.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
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.
Definition at line 50 of file OneDim.cpp.
void init | ( | ) |
Initialize all domains.
On the first call, this methods calls the init method of each domain, proceeding from left to right. Subsequent calls do nothing.
Definition at line 324 of file OneDim.cpp.
double ssnorm | ( | double * | x, |
double * | r | ||
) |
Steady-state max norm (infinity norm) of the residual evaluated using solution x.
On return, array r contains the steady-state residual values. Used only for diagnostic output.
Definition at line 278 of file OneDim.cpp.
void initTimeInteg | ( | double | dt, |
double * | x | ||
) |
Prepare for time stepping beginning with solution x and timestep dt.
Definition at line 288 of file OneDim.cpp.
void setSteadyMode | ( | ) |
Prepare to solve the steady-state problem.
After invoking this method, subsequent calls to solve() will solve the steady-state problem. Sets the reciprocal of the time step to zero, and, if it was previously non- zero, signals that a new Jacobian will be needed.
Definition at line 307 of file OneDim.cpp.
void eval | ( | size_t | j, |
double * | x, | ||
double * | r, | ||
double | rdt = -1.0 , |
||
int | count = 1 |
||
) |
Evaluate the multi-domain residual function.
j | if j != npos, only evaluate residual for points j-1, j, and j + 1; otherwise, evaluate at all grid points. |
x | solution vector |
r | on return, contains the residual vector |
rdt | Reciprocal of the time step. if omitted, then the default value is used. |
count | Set to zero to omit this call from the statistics |
Definition at line 246 of file OneDim.cpp.
Domain1D * pointDomain | ( | size_t | i | ) |
Return a pointer to the domain global point i belongs to.
The domains are scanned right-to-left, and the first one with starting location less or equal to i is returned.
Definition at line 234 of file OneDim.cpp.
|
virtual |
Call after one or more grids has changed size, for example after being refined.
Reimplemented in Sim1D.
Definition at line 154 of file OneDim.cpp.
double timeStep | ( | int | nsteps, |
double | dt, | ||
double * | x, | ||
double * | r, | ||
int | loglevel | ||
) |
Take time steps using Backward Euler.
nsteps | number of steps |
dt | initial step size |
x | current solution vector |
r | solution vector after time stepping |
loglevel | controls amount of printed diagnostics |
Definition at line 336 of file OneDim.cpp.
void writeStats | ( | int | printTime = 1 | ) |
Write statistics about the number of iterations and Jacobians at each grid level.
printTime | Boolean that indicates whether time should be printed out The default is true. It's turned off for test problems where we don't want to print any times |
Definition at line 106 of file OneDim.cpp.
|
inline |
|
inline |
void saveStats | ( | ) |
Save statistics on function and Jacobian evaluation, and reset the counters.
Statistics are saved only if the number of Jacobian evaluations is greater than zero. The statistics saved are:
Definition at line 123 of file OneDim.cpp.
void clearStats | ( | ) |
Clear saved statistics.
Definition at line 141 of file OneDim.cpp.
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
inline |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
|
private |