Cantera 2.6.0
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
OneDim Class Reference

Container class for multiple-domain 1D problems. More...

#include <OneDim.h>

Inheritance diagram for OneDim:
[legend]
Collaboration diagram for OneDim:
[legend]

Public Member Functions

 OneDim (std::vector< Domain1D * > domains)
 Construct a OneDim container for the domains in the list domains. More...
 
 OneDim (const OneDim &)=delete
 
OneDimoperator= (const OneDim &)=delete
 
void addDomain (Domain1D *d)
 Add a domain. Domains are added left-to-right. More...
 
MultiJacjacobian ()
 Return a reference to the Jacobian evaluator. More...
 
MultiNewtonnewton ()
 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...
 
Domain1Ddomain (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...
 
Domain1Dleft ()
 Pointer to left-most domain (first added). More...
 
Domain1Dright ()
 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...
 
Domain1DpointDomain (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_inttransientMask ()
 
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)
 
AnyMap serialize (const double *soln) const
 
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_fpjacobianTimeStats ()
 Return CPU time spent evaluating Jacobians in each call to solve() More...
 
const vector_fpevalTimeStats ()
 Return CPU time spent on non-Jacobian function evaluations in each call to solve() More...
 
const vector_intjacobianCountStats ()
 Return number of Jacobian evaluations made in each call to solve() More...
 
const vector_intevalCountStats ()
 Return number of non-Jacobian function evaluations made in each call to solve() More...
 
const vector_inttimeStepStats ()
 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 (doublereal *x, doublereal *xnew)
 

Protected Attributes

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< MultiJacm_jac
 Jacobian evaluator. More...
 
std::unique_ptr< MultiNewtonm_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
 
int m_ss_jac_age
 
int m_ts_jac_age
 
Func1m_interrupt
 Function called at the start of every call to eval. More...
 
Func1m_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 Attributes

int m_nevals
 
doublereal m_evaltime
 
std::vector< size_t > m_gridpts
 
vector_int m_jacEvals
 
vector_fp m_jacElapsed
 
vector_int m_funcEvals
 
vector_fp m_funcElapsed
 
vector_int m_timeSteps
 Number of time steps taken in each call to solve() (for example, for each successive grid refinement) More...
 

Detailed Description

Container class for multiple-domain 1D problems.

Each domain is represented by an instance of Domain1D.

Definition at line 26 of file OneDim.h.

Constructor & Destructor Documentation

◆ OneDim() [1/2]

OneDim ( )

Definition at line 20 of file OneDim.cpp.

◆ OneDim() [2/2]

OneDim ( std::vector< Domain1D * >  domains)

Construct a OneDim container for the domains in the list domains.

Definition at line 33 of file OneDim.cpp.

References OneDim::addDomain(), OneDim::init(), OneDim::m_newt, and OneDim::resize().

◆ ~OneDim()

~OneDim ( )
virtual

Definition at line 52 of file OneDim.cpp.

Member Function Documentation

◆ addDomain()

void addDomain ( Domain1D d)

Add a domain. Domains are added left-to-right.

Definition at line 80 of file OneDim.cpp.

References OneDim::resize(), and Domain1D::setContainer().

Referenced by OneDim::OneDim().

◆ jacobian()

MultiJac & jacobian ( )

Return a reference to the Jacobian evaluator.

Definition at line 102 of file OneDim.cpp.

References OneDim::m_jac.

Referenced by Domain1D::needJacUpdate().

◆ newton()

MultiNewton & newton ( )

Return a reference to the Newton iterator.

Definition at line 106 of file OneDim.cpp.

References OneDim::m_newt.

Referenced by OneDim::timeStep().

◆ solve()

int solve ( doublereal *  x0,
doublereal *  x1,
int  loglevel 
)

Solve F(x) = 0, where F(x) is the multi-domain residual function.

Parameters
x0Starting estimate of solution.
x1Final solution satisfying F(x1) = 0.
loglevelControls amount of diagnostic output.

Definition at line 225 of file OneDim.cpp.

References OneDim::eval(), OneDim::m_jac, OneDim::m_jac_ok, OneDim::m_newt, OneDim::m_rdt, and Cantera::npos.

Referenced by Sim1D::newtonSolve(), and OneDim::timeStep().

◆ nDomains()

size_t nDomains ( ) const
inline

◆ domain()

Domain1D & domain ( size_t  i) const
inline

◆ domainIndex()

size_t domainIndex ( const std::string &  name)

Definition at line 56 of file OneDim.cpp.

◆ checkDomainIndex()

void checkDomainIndex ( size_t  n) const
inline

Check that the specified domain index is in range.

Throws an exception if n is greater than nDomains()-1

Definition at line 68 of file OneDim.h.

◆ checkDomainArraySize()

void checkDomainArraySize ( size_t  nn) const
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.

Definition at line 78 of file OneDim.h.

◆ start()

size_t start ( size_t  i) const
inline

The index of the start of domain i in the solution vector.

Definition at line 86 of file OneDim.h.

Referenced by MultiNewton::boundStep(), OneDim::component(), Sim1D::finalize(), MultiNewton::norm2(), Sim1D::refine(), Sim1D::restore(), and Sim1D::Sim1D().

◆ size()

size_t size ( ) const
inline

Total solution vector length;.

Definition at line 97 of file OneDim.h.

References OneDim::m_size.

Referenced by MultiNewton::norm2(), OneDim::resize(), Sim1D::resize(), and MultiNewton::step().

◆ left()

Domain1D * left ( )
inline

Pointer to left-most domain (first added).

Definition at line 102 of file OneDim.h.

Referenced by OneDim::init(), OneDim::initTimeInteg(), and OneDim::setSteadyMode().

◆ right()

Domain1D * right ( )
inline

Pointer to right-most domain (last added).

Definition at line 107 of file OneDim.h.

Referenced by OneDim::pointDomain().

◆ nVars()

size_t nVars ( size_t  jg)
inline

Number of solution components at global point jg.

Definition at line 112 of file OneDim.h.

◆ loc()

size_t loc ( size_t  jg)
inline

Location in the solution vector of the first component of global point jg.

Definition at line 118 of file OneDim.h.

◆ component()

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.

Definition at line 66 of file OneDim.cpp.

References Domain1D::componentName(), OneDim::domain(), Domain1D::nComponents(), OneDim::nDomains(), Cantera::npos, and OneDim::start().

Referenced by Sim1D::setInitialGuess().

◆ bandwidth()

size_t bandwidth ( ) const
inline

Jacobian bandwidth.

Definition at line 127 of file OneDim.h.

References OneDim::m_bw.

◆ init()

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 337 of file OneDim.cpp.

References Domain1D::init(), OneDim::left(), and Domain1D::right().

Referenced by OneDim::OneDim().

◆ points()

size_t points ( )
inline

Total number of points.

Definition at line 139 of file OneDim.h.

◆ ssnorm()

doublereal ssnorm ( doublereal *  x,
doublereal *  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 291 of file OneDim.cpp.

References OneDim::eval(), OneDim::m_size, and Cantera::npos.

Referenced by MultiNewton::solve(), and OneDim::timeStep().

◆ rdt()

doublereal rdt ( ) const
inline

Reciprocal of the time step.

Definition at line 151 of file OneDim.h.

References OneDim::m_rdt.

Referenced by OneDim::eval(), and MultiNewton::solve().

◆ initTimeInteg()

void initTimeInteg ( doublereal  dt,
doublereal *  x 
)

Prepare for time stepping beginning with solution x and timestep dt.

Definition at line 301 of file OneDim.cpp.

References Domain1D::initTimeInteg(), OneDim::left(), OneDim::m_jac, OneDim::m_rdt, Domain1D::right(), and Cantera::Tiny.

Referenced by OneDim::timeStep().

◆ transient()

bool transient ( ) const
inline

True if transient mode.

Definition at line 159 of file OneDim.h.

References OneDim::m_rdt.

◆ steady()

bool steady ( ) const
inline

True if steady mode.

Definition at line 164 of file OneDim.h.

References OneDim::m_rdt.

◆ setSteadyMode()

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 320 of file OneDim.cpp.

References OneDim::left(), OneDim::m_jac, OneDim::m_rdt, Domain1D::right(), and Domain1D::setSteadyMode().

◆ eval()

void eval ( size_t  j,
double *  x,
double *  r,
doublereal  rdt = -1.0,
int  count = 1 
)

Evaluate the multi-domain residual function.

Parameters
jif j != npos, only evaluate residual for points j-1, j, and j + 1; otherwise, evaluate at all grid points.
xsolution vector
ron return, contains the residual vector
rdtReciprocal of the time step. if omitted, then the default value is used.
countSet to zero to omit this call from the statistics

Definition at line 259 of file OneDim.cpp.

References Func1::eval(), OneDim::m_interrupt, OneDim::m_rdt, OneDim::m_size, Cantera::npos, and OneDim::rdt().

Referenced by OneDim::solve(), MultiNewton::solve(), OneDim::ssnorm(), and MultiNewton::step().

◆ pointDomain()

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 247 of file OneDim.cpp.

References Domain1D::left(), Domain1D::loc(), and OneDim::right().

◆ resize()

void resize ( )
virtual

Call after one or more grids has changed size, for example after being refined.

Reimplemented in Sim1D.

Definition at line 169 of file OneDim.cpp.

References Domain1D::bandwidth(), Domain1D::loc(), OneDim::m_bw, OneDim::m_jac, OneDim::m_jac_ok, OneDim::m_newt, OneDim::m_size, Domain1D::nComponents(), OneDim::nDomains(), Domain1D::nPoints(), Cantera::npos, OneDim::saveStats(), and OneDim::size().

Referenced by OneDim::addDomain(), OneDim::OneDim(), and Sim1D::resize().

◆ transientMask()

vector_int & transientMask ( )
inline

Definition at line 200 of file OneDim.h.

◆ timeStep()

doublereal timeStep ( int  nsteps,
double  dt,
double *  x,
double *  r,
int  loglevel 
)

Take time steps using Backward Euler.

Parameters
nstepsnumber of steps
dtinitial step size
xcurrent solution vector
rsolution vector after time stepping
loglevelcontrols amount of printed diagnostics
Returns
size of last timestep taken

Definition at line 349 of file OneDim.cpp.

References Cantera::debuglog(), Func1::eval(), OneDim::initTimeInteg(), OneDim::m_nsteps, OneDim::m_nsteps_max, OneDim::m_size, OneDim::m_tfactor, OneDim::m_time_step_callback, OneDim::m_tmax, OneDim::m_tmin, OneDim::newton(), MultiNewton::setOptions(), OneDim::solve(), OneDim::ssnorm(), and Cantera::writelog().

◆ resetBadValues()

void resetBadValues ( double *  x)

Definition at line 417 of file OneDim.cpp.

◆ writeStats()

void writeStats ( int  printTime = 1)

Write statistics about the number of iterations and Jacobians at each grid level.

Parameters
printTimeBoolean 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 121 of file OneDim.cpp.

References OneDim::m_timeSteps, OneDim::saveStats(), and Cantera::writelog().

◆ save()

void save ( const std::string &  fname,
std::string  id,
const std::string &  desc,
doublereal *  sol,
int  loglevel 
)

Definition at line 424 of file OneDim.cpp.

◆ serialize()

AnyMap serialize ( const double *  soln) const

Definition at line 464 of file OneDim.cpp.

◆ setMinTimeStep()

void setMinTimeStep ( doublereal  tmin)
inline

Definition at line 234 of file OneDim.h.

◆ setMaxTimeStep()

void setMaxTimeStep ( doublereal  tmax)
inline

Definition at line 237 of file OneDim.h.

◆ setTimeStepFactor()

void setTimeStepFactor ( doublereal  tfactor)
inline

Definition at line 240 of file OneDim.h.

◆ setMaxTimeStepCount()

void setMaxTimeStepCount ( int  nmax)
inline

Set the maximum number of timeteps allowed before successful steady-state solve.

Definition at line 246 of file OneDim.h.

References OneDim::m_nsteps_max.

◆ maxTimeStepCount()

int maxTimeStepCount ( ) const
inline

Return the maximum number of timeteps allowed before successful steady-state solve.

Definition at line 252 of file OneDim.h.

References OneDim::m_nsteps_max.

◆ setJacAge()

void setJacAge ( int  ss_age,
int  ts_age = -1 
)

Definition at line 111 of file OneDim.cpp.

◆ saveStats()

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:

  • number of grid points
  • number of Jacobian evaluations
  • CPU time spent evaluating Jacobians
  • number of non-Jacobian function evaluations
  • CPU time spent evaluating functions
  • number of time steps

Definition at line 138 of file OneDim.cpp.

References OneDim::m_jac, OneDim::m_nsteps, and OneDim::m_timeSteps.

Referenced by OneDim::evalCountStats(), OneDim::evalTimeStats(), OneDim::gridSizeStats(), OneDim::jacobianCountStats(), OneDim::jacobianTimeStats(), Domain1D::needJacUpdate(), OneDim::resize(), OneDim::timeStepStats(), and OneDim::writeStats().

◆ clearStats()

void clearStats ( )

Clear saved statistics.

Definition at line 156 of file OneDim.cpp.

References OneDim::m_nsteps, and OneDim::m_timeSteps.

◆ gridSizeStats()

const std::vector< size_t > & gridSizeStats ( )
inline

Return total grid size in each call to solve()

Definition at line 276 of file OneDim.h.

References OneDim::saveStats().

◆ jacobianTimeStats()

const vector_fp & jacobianTimeStats ( )
inline

Return CPU time spent evaluating Jacobians in each call to solve()

Definition at line 282 of file OneDim.h.

References OneDim::saveStats().

◆ evalTimeStats()

const vector_fp & evalTimeStats ( )
inline

Return CPU time spent on non-Jacobian function evaluations in each call to solve()

Definition at line 289 of file OneDim.h.

References OneDim::saveStats().

◆ jacobianCountStats()

const vector_int & jacobianCountStats ( )
inline

Return number of Jacobian evaluations made in each call to solve()

Definition at line 295 of file OneDim.h.

References OneDim::saveStats().

◆ evalCountStats()

const vector_int & evalCountStats ( )
inline

Return number of non-Jacobian function evaluations made in each call to solve()

Definition at line 302 of file OneDim.h.

References OneDim::saveStats().

◆ timeStepStats()

const vector_int & timeStepStats ( )
inline

Return number of time steps taken in each call to solve()

Definition at line 308 of file OneDim.h.

References OneDim::m_timeSteps, and OneDim::saveStats().

◆ setInterrupt()

void setInterrupt ( Func1 interrupt)
inline

Set a function that will be called every time eval is called.

Can be used to provide keyboard interrupt support in the high-level language interfaces.

Definition at line 316 of file OneDim.h.

References OneDim::m_interrupt.

◆ setTimeStepCallback()

void setTimeStepCallback ( Func1 callback)
inline

Set a function that will be called after each successful timestep.

The function will be called with the size of the timestep as the argument. Intended to be used for observing solver progress for debugging purposes.

Definition at line 324 of file OneDim.h.

References OneDim::m_time_step_callback.

◆ evalSSJacobian()

void evalSSJacobian ( doublereal *  x,
doublereal *  xnew 
)
protected

Definition at line 237 of file OneDim.cpp.

Member Data Documentation

◆ m_tmin

doublereal m_tmin
protected

minimum timestep size

Definition at line 331 of file OneDim.h.

Referenced by OneDim::timeStep().

◆ m_tmax

doublereal m_tmax
protected

maximum timestep size

Definition at line 332 of file OneDim.h.

Referenced by OneDim::timeStep().

◆ m_tfactor

doublereal m_tfactor
protected

factor time step is multiplied by if time stepping fails ( < 1 )

Definition at line 335 of file OneDim.h.

Referenced by OneDim::timeStep().

◆ m_jac

std::unique_ptr<MultiJac> m_jac
protected

◆ m_newt

std::unique_ptr<MultiNewton> m_newt
protected

Newton iterator.

Definition at line 338 of file OneDim.h.

Referenced by OneDim::newton(), OneDim::OneDim(), OneDim::resize(), and OneDim::solve().

◆ m_rdt

doublereal m_rdt
protected

◆ m_jac_ok

bool m_jac_ok
protected

if true, Jacobian is current

Definition at line 340 of file OneDim.h.

Referenced by OneDim::resize(), and OneDim::solve().

◆ m_bw

size_t m_bw
protected

Jacobian bandwidth.

Definition at line 342 of file OneDim.h.

Referenced by OneDim::bandwidth(), and OneDim::resize().

◆ m_size

size_t m_size
protected

solution vector size

Definition at line 343 of file OneDim.h.

Referenced by OneDim::eval(), OneDim::resize(), OneDim::size(), OneDim::ssnorm(), and OneDim::timeStep().

◆ m_dom

std::vector<Domain1D*> m_dom
protected

Definition at line 345 of file OneDim.h.

◆ m_connect

std::vector<Domain1D*> m_connect
protected

Definition at line 345 of file OneDim.h.

◆ m_bulk

std::vector<Domain1D*> m_bulk
protected

Definition at line 345 of file OneDim.h.

◆ m_init

bool m_init
protected

Definition at line 347 of file OneDim.h.

◆ m_nvars

std::vector<size_t> m_nvars
protected

Definition at line 348 of file OneDim.h.

◆ m_loc

std::vector<size_t> m_loc
protected

Definition at line 349 of file OneDim.h.

◆ m_mask

vector_int m_mask
protected

Definition at line 350 of file OneDim.h.

◆ m_pts

size_t m_pts
protected

Definition at line 351 of file OneDim.h.

◆ m_ss_jac_age

int m_ss_jac_age
protected

Definition at line 354 of file OneDim.h.

◆ m_ts_jac_age

int m_ts_jac_age
protected

Definition at line 354 of file OneDim.h.

◆ m_interrupt

Func1* m_interrupt
protected

Function called at the start of every call to eval.

Definition at line 357 of file OneDim.h.

Referenced by OneDim::eval(), and OneDim::setInterrupt().

◆ m_time_step_callback

Func1* m_time_step_callback
protected

User-supplied function called after each successful timestep.

Definition at line 360 of file OneDim.h.

Referenced by OneDim::setTimeStepCallback(), and OneDim::timeStep().

◆ m_nsteps

int m_nsteps
protected

Number of time steps taken in the current call to solve()

Definition at line 363 of file OneDim.h.

Referenced by OneDim::clearStats(), OneDim::saveStats(), and OneDim::timeStep().

◆ m_nsteps_max

int m_nsteps_max
protected

Maximum number of timesteps allowed per call to solve()

Definition at line 366 of file OneDim.h.

Referenced by OneDim::maxTimeStepCount(), OneDim::setMaxTimeStepCount(), and OneDim::timeStep().

◆ m_nevals

int m_nevals
private

Definition at line 370 of file OneDim.h.

◆ m_evaltime

doublereal m_evaltime
private

Definition at line 371 of file OneDim.h.

◆ m_gridpts

std::vector<size_t> m_gridpts
private

Definition at line 372 of file OneDim.h.

◆ m_jacEvals

vector_int m_jacEvals
private

Definition at line 373 of file OneDim.h.

◆ m_jacElapsed

vector_fp m_jacElapsed
private

Definition at line 374 of file OneDim.h.

◆ m_funcEvals

vector_int m_funcEvals
private

Definition at line 375 of file OneDim.h.

◆ m_funcElapsed

vector_fp m_funcElapsed
private

Definition at line 376 of file OneDim.h.

◆ m_timeSteps

vector_int m_timeSteps
private

Number of time steps taken in each call to solve() (for example, for each successive grid refinement)

Definition at line 380 of file OneDim.h.

Referenced by OneDim::clearStats(), OneDim::saveStats(), OneDim::timeStepStats(), and OneDim::writeStats().


The documentation for this class was generated from the following files: