Cantera  3.2.0a1
Loading...
Searching...
No Matches
OneDim Class Reference

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

#include <OneDim.h>

Inheritance diagram for OneDim:
[legend]

Detailed Description

Container class for multiple-domain 1D problems.

Each domain is represented by an instance of Domain1D.

Definition at line 27 of file OneDim.h.

Public Member Functions

 OneDim ()
 Default constructor.
 
 OneDim (vector< shared_ptr< Domain1D > > &domains)
 Construct a OneDim container for the domains in the list domains.
 
 OneDim (const OneDim &)=delete
 
OneDimoperator= (const OneDim &)=delete
 
void addDomain (shared_ptr< Domain1D > d)
 Add a domain. Domains are added left-to-right.
 
MultiJacjacobian ()
 Return a reference to the Jacobian evaluator of an OneDim object.
 
shared_ptr< SystemJacobiangetJacobian ()
 
MultiNewtonnewton ()
 Return a reference to the Newton iterator.
 
void setLinearSolver (shared_ptr< SystemJacobian > solver)
 Set the linear solver used to hold the Jacobian matrix and solve linear systems as part of each Newton iteration.
 
shared_ptr< SystemJacobianlinearSolver () const
 Get the type of the linear solver being used.
 
int solve (double *x0, double *x1, int loglevel)
 Solve F(x) = 0, where F(x) is the multi-domain residual function.
 
size_t nDomains () const
 Number of domains.
 
Domain1Ddomain (size_t i) const
 Return a reference to domain i.
 
size_t domainIndex (const string &name)
 Get the index of the domain named name.
 
void checkDomainIndex (size_t n) const
 Check that the specified domain index is in range.
 
void checkDomainArraySize (size_t nn) const
 Check that an array size is at least nDomains().
 
size_t start (size_t i) const
 The index of the start of domain i in the solution vector.
 
size_t size () const
 Total solution vector length;.
 
Domain1Dleft ()
 Pointer to left-most domain (first added).
 
Domain1Dright ()
 Pointer to right-most domain (last added).
 
size_t nVars (size_t jg)
 Number of solution components at global point jg.
 
size_t loc (size_t jg)
 Location in the solution vector of the first component of global point jg.
 
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.
 
size_t bandwidth () const
 Jacobian bandwidth.
 
void init ()
 Initialize all domains.
 
size_t points ()
 Total number of points.
 
double ssnorm (double *x, double *r)
 Steady-state max norm (infinity norm) of the residual evaluated using solution x.
 
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.
 
bool transient () const
 True if transient mode.
 
bool steady () const
 True if steady mode.
 
void setSteadyMode ()
 Prepare to solve the steady-state problem.
 
void eval (size_t j, double *x, double *r, double rdt=-1.0, int count=1)
 Evaluate the multi-domain residual function.
 
void evalJacobian (double *x0)
 Evaluates the Jacobian at x0 using finite differences.
 
Domain1DpointDomain (size_t i)
 Return a pointer to the domain global point i belongs to.
 
virtual void resize ()
 Call after one or more grids has changed size, for example after being refined.
 
vector< int > & transientMask ()
 Access the vector indicating which equations contain a transient term.
 
double timeStep (int nsteps, double dt, double *x, double *r, int loglevel)
 Take time steps using Backward Euler.
 
void resetBadValues (double *x)
 Reset values such as negative species concentrations in each domain.
 
void writeStats (int printTime=1)
 Write statistics about the number of iterations and Jacobians at each grid level.
 
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-evaluated.
 
void saveStats ()
 Save statistics on function and Jacobian evaluation, and reset the counters.
 
void clearStats ()
 Clear saved statistics.
 
const vector< size_t > & gridSizeStats ()
 Return total grid size in each call to solve()
 
const vector< double > & jacobianTimeStats ()
 Return CPU time spent evaluating Jacobians in each call to solve()
 
const vector< double > & evalTimeStats ()
 Return CPU time spent on non-Jacobian function evaluations in each call to solve()
 
const vector< int > & jacobianCountStats ()
 Return number of Jacobian evaluations made in each call to solve()
 
const vector< int > & evalCountStats ()
 Return number of non-Jacobian function evaluations made in each call to solve()
 
const vector< int > & timeStepStats ()
 Return number of time steps taken in each call to solve()
 
void setInterrupt (Func1 *interrupt)
 Set a function that will be called every time eval is called.
 
void setTimeStepCallback (Func1 *callback)
 Set a function that will be called after each successful timestep.
 
void setJacobianPerturbation (double relative, double absolute, double threshold)
 Configure perturbations used to evaluate finite difference Jacobian.
 
Options
void setMinTimeStep (double tmin)
 Set the minimum time step allowed during time stepping.
 
void setMaxTimeStep (double tmax)
 Set the maximum time step allowed during time stepping.
 
void setTimeStepFactor (double tfactor)
 Sets a factor by which the time step is reduced if the time stepping fails.
 
void setMaxTimeStepCount (int nmax)
 Set the maximum number of timeteps allowed before successful steady-state solve.
 
int maxTimeStepCount () const
 Return the maximum number of timeteps allowed before successful steady-state solve.
 

Protected Member Functions

void evalSSJacobian (double *x, double *rsd)
 Evaluate the steady-state Jacobian, accessible via jacobian()
 

Protected Attributes

double m_tmin = 1e-16
 minimum timestep size
 
double m_tmax = 1e+08
 maximum timestep size
 
double m_tfactor = 0.5
 factor time step is multiplied by if time stepping fails ( < 1 )
 
shared_ptr< vector< double > > m_state
 Solution vector.
 
shared_ptr< SystemJacobianm_jac
 Jacobian evaluator.
 
unique_ptr< MultiNewtonm_newt
 Newton iterator.
 
double m_rdt = 0.0
 reciprocal of time step
 
bool m_jac_ok = false
 if true, Jacobian is current
 
size_t m_bw = 0
 Jacobian bandwidth.
 
size_t m_size = 0
 solution vector size
 
vector< double > m_work1
 Work arrays used during Jacobian evaluation.
 
vector< double > m_work2
 
vector< shared_ptr< Domain1D > > m_dom
 All domains comprising the system.
 
vector< shared_ptr< Domain1D > > m_connect
 All connector and boundary domains.
 
vector< shared_ptr< Domain1D > > m_bulk
 All bulk/flow domains.
 
bool m_init = false
 Indicates whether one-time initialization for each domain has been completed.
 
vector< size_t > m_nvars
 Number of variables at each point, across all domains.
 
vector< size_t > m_loc
 Location in the state vector of the first component of each point, across all domains.
 
vector< int > m_mask
 Transient mask. See transientMask().
 
size_t m_pts = 0
 Total number of points.
 
int m_ss_jac_age = 20
 Maximum age of the Jacobian in steady-state mode.
 
int m_ts_jac_age = 20
 Maximum age of the Jacobian in time-stepping mode.
 
Func1m_interrupt = nullptr
 Function called at the start of every call to eval.
 
Func1m_time_step_callback = nullptr
 User-supplied function called after each successful timestep.
 
int m_nsteps = 0
 Number of time steps taken in the current call to solve()
 
int m_nsteps_max = 500
 Maximum number of timesteps allowed per call to solve()
 
double m_jacobianThreshold = 0.0
 Threshold for ignoring small elements in Jacobian.
 
double m_jacobianRelPerturb = 1e-5
 Relative perturbation of each component in finite difference Jacobian.
 
double m_jacobianAbsPerturb = 1e-10
 Absolute perturbation of each component in finite difference Jacobian.
 

Private Attributes

Statistics

Solver stats are collected after successfully solving on a particular grid.

int m_nevals = 0
 Number of calls to eval()
 
double m_evaltime = 0
 Total time [s] spent in eval()
 
vector< size_t > m_gridpts
 Number of grid points in this grid.
 
vector< int > m_jacEvals
 Number of Jacobian evaluations on this grid.
 
vector< double > m_jacElapsed
 Time [s] spent evaluating Jacobians on this grid.
 
vector< int > m_funcEvals
 Number of residual function evaluations on this grid (not counting evaluations used to construct Jacobians).
 
vector< double > m_funcElapsed
 Time [s] spent on residual function evaluations on this grid (not counting evaluations used to construct Jacobians).
 
vector< int > m_timeSteps
 Number of time steps taken in each call to solve() (for example, for each successive grid refinement)
 

Constructor & Destructor Documentation

◆ OneDim() [1/2]

OneDim ( )

Default constructor.

Definition at line 20 of file OneDim.cpp.

◆ OneDim() [2/2]

OneDim ( vector< shared_ptr< Domain1D > > &  domains)

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

Definition at line 25 of file OneDim.cpp.

◆ ~OneDim()

~OneDim ( )
virtual

Definition at line 37 of file OneDim.cpp.

Member Function Documentation

◆ addDomain()

void addDomain ( shared_ptr< Domain1D d)

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

Definition at line 65 of file OneDim.cpp.

◆ getJacobian()

shared_ptr< SystemJacobian > getJacobian ( )
inline

Definition at line 47 of file OneDim.h.

◆ newton()

MultiNewton & newton ( )

Return a reference to the Newton iterator.

Definition at line 100 of file OneDim.cpp.

◆ setLinearSolver()

void setLinearSolver ( shared_ptr< SystemJacobian solver)

Set the linear solver used to hold the Jacobian matrix and solve linear systems as part of each Newton iteration.

The default is a direct, banded solver.

Definition at line 105 of file OneDim.cpp.

◆ linearSolver()

shared_ptr< SystemJacobian > linearSolver ( ) const
inline

Get the type of the linear solver being used.

Definition at line 59 of file OneDim.h.

◆ solve()

int solve ( double *  x0,
double *  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.
Returns
  • 1 for success
  • -2 failure (maximum number of damping steps was reached)
  • -3 failure (solution was up against the bounds

Definition at line 230 of file OneDim.cpp.

◆ nDomains()

size_t nDomains ( ) const
inline

Number of domains.

Definition at line 76 of file OneDim.h.

◆ domain()

Domain1D & domain ( size_t  i) const
inline

Return a reference to domain i.

Definition at line 81 of file OneDim.h.

◆ domainIndex()

size_t domainIndex ( const string &  name)

Get the index of the domain named name.

Definition at line 41 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 90 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 99 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 107 of file OneDim.h.

◆ size()

size_t size ( ) const
inline

Total solution vector length;.

Definition at line 118 of file OneDim.h.

◆ left()

Domain1D * left ( )
inline

Pointer to left-most domain (first added).

Definition at line 123 of file OneDim.h.

◆ right()

Domain1D * right ( )
inline

Pointer to right-most domain (last added).

Definition at line 128 of file OneDim.h.

◆ nVars()

size_t nVars ( size_t  jg)
inline

Number of solution components at global point jg.

Definition at line 133 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 138 of file OneDim.h.

◆ component()

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

◆ bandwidth()

size_t bandwidth ( ) const
inline

Jacobian bandwidth.

Definition at line 147 of file OneDim.h.

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

◆ points()

size_t points ( )
inline

Total number of points.

Definition at line 159 of file OneDim.h.

◆ ssnorm()

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

◆ rdt()

double rdt ( ) const
inline

Reciprocal of the time step.

Definition at line 171 of file OneDim.h.

◆ initTimeInteg()

void initTimeInteg ( double  dt,
double *  x 
)

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

Definition at line 349 of file OneDim.cpp.

◆ transient()

bool transient ( ) const
inline

True if transient mode.

Definition at line 179 of file OneDim.h.

◆ steady()

bool steady ( ) const
inline

True if steady mode.

Definition at line 184 of file OneDim.h.

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

◆ eval()

void eval ( size_t  j,
double *  x,
double *  r,
double  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 261 of file OneDim.cpp.

◆ evalJacobian()

void evalJacobian ( double *  x0)

Evaluates the Jacobian at x0 using finite differences.

The Jacobian is computed by perturbing each component of x0, evaluating the residual function, and then estimating the partial derivatives numerically using finite differences to determine the corresponding column of the Jacobian.

Parameters
x0State vector at which to evaluate the Jacobian

Definition at line 293 of file OneDim.cpp.

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

◆ resize()

void resize ( )
virtual

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

Reimplemented in Sim1D.

Definition at line 172 of file OneDim.cpp.

◆ transientMask()

vector< int > & transientMask ( )
inline

Access the vector indicating which equations contain a transient term.

Elements are 1 for equations with a transient terms and 0 otherwise.

Definition at line 232 of file OneDim.h.

◆ timeStep()

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

◆ resetBadValues()

void resetBadValues ( double *  x)

Reset values such as negative species concentrations in each domain.

See also
Domain1D::resetBadValues

Definition at line 491 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 124 of file OneDim.cpp.

◆ setMinTimeStep()

void setMinTimeStep ( double  tmin)
inline

Set the minimum time step allowed during time stepping.

Definition at line 265 of file OneDim.h.

◆ setMaxTimeStep()

void setMaxTimeStep ( double  tmax)
inline

Set the maximum time step allowed during time stepping.

Definition at line 270 of file OneDim.h.

◆ setTimeStepFactor()

void setTimeStepFactor ( double  tfactor)
inline

Sets a factor by which the time step is reduced if the time stepping fails.

The default value is 0.5.

Parameters
tfactorfactor time step is multiplied by if time stepping fails

Definition at line 280 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 286 of file OneDim.h.

◆ maxTimeStepCount()

int maxTimeStepCount ( ) const
inline

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

Definition at line 292 of file OneDim.h.

◆ setJacAge()

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-evaluated.

Parameters
ss_ageAge limit during steady-state mode
ts_ageAge limit during time stepping mode. If not specified, the steady-state age is also used during time stepping.

Definition at line 114 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 141 of file OneDim.cpp.

◆ clearStats()

void clearStats ( )

Clear saved statistics.

Definition at line 159 of file OneDim.cpp.

◆ gridSizeStats()

const vector< size_t > & gridSizeStats ( )
inline

Return total grid size in each call to solve()

Definition at line 322 of file OneDim.h.

◆ jacobianTimeStats()

const vector< double > & jacobianTimeStats ( )
inline

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

Definition at line 328 of file OneDim.h.

◆ evalTimeStats()

const vector< double > & evalTimeStats ( )
inline

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

Definition at line 335 of file OneDim.h.

◆ jacobianCountStats()

const vector< int > & jacobianCountStats ( )
inline

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

Definition at line 341 of file OneDim.h.

◆ evalCountStats()

const vector< int > & evalCountStats ( )
inline

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

Definition at line 348 of file OneDim.h.

◆ timeStepStats()

const vector< int > & timeStepStats ( )
inline

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

Definition at line 354 of file OneDim.h.

◆ 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 362 of file OneDim.h.

◆ 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 370 of file OneDim.h.

◆ setJacobianPerturbation()

void setJacobianPerturbation ( double  relative,
double  absolute,
double  threshold 
)
inline

Configure perturbations used to evaluate finite difference Jacobian.

Parameters
relativeRelative perturbation (multiplied by the absolute value of each component). Default 1.0e-5.
absoluteAbsolute perturbation (independent of component value). Default 1.0e-10.
thresholdThreshold below which to exclude elements from the Jacobian Default 0.0.

Definition at line 381 of file OneDim.h.

◆ evalSSJacobian()

void evalSSJacobian ( double *  x,
double *  rsd 
)
protected

Evaluate the steady-state Jacobian, accessible via jacobian()

Parameters
[in]xCurrent state vector, length size()
[out]rsdStorage for the residual, length size()

Definition at line 240 of file OneDim.cpp.

Member Data Documentation

◆ m_tmin

double m_tmin = 1e-16
protected

minimum timestep size

Definition at line 393 of file OneDim.h.

◆ m_tmax

double m_tmax = 1e+08
protected

maximum timestep size

Definition at line 394 of file OneDim.h.

◆ m_tfactor

double m_tfactor = 0.5
protected

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

Definition at line 397 of file OneDim.h.

◆ m_state

shared_ptr<vector<double> > m_state
protected

Solution vector.

Definition at line 399 of file OneDim.h.

◆ m_jac

shared_ptr<SystemJacobian> m_jac
protected

Jacobian evaluator.

Definition at line 401 of file OneDim.h.

◆ m_newt

unique_ptr<MultiNewton> m_newt
protected

Newton iterator.

Definition at line 402 of file OneDim.h.

◆ m_rdt

double m_rdt = 0.0
protected

reciprocal of time step

Definition at line 403 of file OneDim.h.

◆ m_jac_ok

bool m_jac_ok = false
protected

if true, Jacobian is current

Definition at line 404 of file OneDim.h.

◆ m_bw

size_t m_bw = 0
protected

Jacobian bandwidth.

Definition at line 406 of file OneDim.h.

◆ m_size

size_t m_size = 0
protected

solution vector size

Definition at line 407 of file OneDim.h.

◆ m_work1

vector<double> m_work1
protected

Work arrays used during Jacobian evaluation.

Definition at line 410 of file OneDim.h.

◆ m_work2

vector<double> m_work2
protected

Definition at line 410 of file OneDim.h.

◆ m_dom

vector<shared_ptr<Domain1D> > m_dom
protected

All domains comprising the system.

Definition at line 413 of file OneDim.h.

◆ m_connect

vector<shared_ptr<Domain1D> > m_connect
protected

All connector and boundary domains.

Definition at line 416 of file OneDim.h.

◆ m_bulk

vector<shared_ptr<Domain1D> > m_bulk
protected

All bulk/flow domains.

Definition at line 419 of file OneDim.h.

◆ m_init

bool m_init = false
protected

Indicates whether one-time initialization for each domain has been completed.

Definition at line 422 of file OneDim.h.

◆ m_nvars

vector<size_t> m_nvars
protected

Number of variables at each point, across all domains.

Length points(). Accessed with nVars().

Definition at line 426 of file OneDim.h.

◆ m_loc

vector<size_t> m_loc
protected

Location in the state vector of the first component of each point, across all domains.

Accessed with loc().

Definition at line 430 of file OneDim.h.

◆ m_mask

vector<int> m_mask
protected

Transient mask. See transientMask().

Definition at line 433 of file OneDim.h.

◆ m_pts

size_t m_pts = 0
protected

Total number of points.

Definition at line 436 of file OneDim.h.

◆ m_ss_jac_age

int m_ss_jac_age = 20
protected

Maximum age of the Jacobian in steady-state mode.

Definition at line 438 of file OneDim.h.

◆ m_ts_jac_age

int m_ts_jac_age = 20
protected

Maximum age of the Jacobian in time-stepping mode.

Definition at line 439 of file OneDim.h.

◆ m_interrupt

Func1* m_interrupt = nullptr
protected

Function called at the start of every call to eval.

Definition at line 442 of file OneDim.h.

◆ m_time_step_callback

Func1* m_time_step_callback = nullptr
protected

User-supplied function called after each successful timestep.

Definition at line 445 of file OneDim.h.

◆ m_nsteps

int m_nsteps = 0
protected

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

Definition at line 448 of file OneDim.h.

◆ m_nsteps_max

int m_nsteps_max = 500
protected

Maximum number of timesteps allowed per call to solve()

Definition at line 451 of file OneDim.h.

◆ m_jacobianThreshold

double m_jacobianThreshold = 0.0
protected

Threshold for ignoring small elements in Jacobian.

Definition at line 454 of file OneDim.h.

◆ m_jacobianRelPerturb

double m_jacobianRelPerturb = 1e-5
protected

Relative perturbation of each component in finite difference Jacobian.

Definition at line 456 of file OneDim.h.

◆ m_jacobianAbsPerturb

double m_jacobianAbsPerturb = 1e-10
protected

Absolute perturbation of each component in finite difference Jacobian.

Definition at line 458 of file OneDim.h.

◆ m_nevals

int m_nevals = 0
private

Number of calls to eval()

Definition at line 464 of file OneDim.h.

◆ m_evaltime

double m_evaltime = 0
private

Total time [s] spent in eval()

Definition at line 465 of file OneDim.h.

◆ m_gridpts

vector<size_t> m_gridpts
private

Number of grid points in this grid.

Definition at line 467 of file OneDim.h.

◆ m_jacEvals

vector<int> m_jacEvals
private

Number of Jacobian evaluations on this grid.

Definition at line 468 of file OneDim.h.

◆ m_jacElapsed

vector<double> m_jacElapsed
private

Time [s] spent evaluating Jacobians on this grid.

Definition at line 469 of file OneDim.h.

◆ m_funcEvals

vector<int> m_funcEvals
private

Number of residual function evaluations on this grid (not counting evaluations used to construct Jacobians).

Definition at line 473 of file OneDim.h.

◆ m_funcElapsed

vector<double> m_funcElapsed
private

Time [s] spent on residual function evaluations on this grid (not counting evaluations used to construct Jacobians).

Definition at line 477 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 481 of file OneDim.h.


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