Cantera  3.2.0a2
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 26 of file OneDim.h.

Public Member Functions

 OneDim ()=default
 Default constructor.
 
 OneDim (vector< shared_ptr< Domain1D > > &domains)
 Construct a OneDim container for the domains in the list domains.
 
void addDomain (shared_ptr< Domain1D > d)
 Add a domain. Domains are added left-to-right.
 
shared_ptr< SystemJacobiangetJacobian ()
 
double weightedNorm (const double *step) const override
 Compute the weighted norm of a step vector.
 
MultiJacjacobian ()
 Return a reference to the Jacobian evaluator of an OneDim object.
 
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.
 
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) const
 Return the domain, local point index, and component name for the i-th component of the global solution vector.
 
string componentName (size_t i) const override
 Get the name of the i-th component of the state vector.
 
pair< string, string > componentTableHeader () const override
 Get header lines describing the column names included in a component label.
 
string componentTableLabel (size_t i) const override
 Get elements of the component name, aligned with the column headings given by componentTableHeader().
 
double upperBound (size_t i) const override
 Get the upper bound for global component i in the state vector.
 
double lowerBound (size_t i) const override
 Get the lower bound for global component i in the state vector.
 
void init ()
 Initialize all domains.
 
size_t points ()
 Total number of points.
 
void initTimeInteg (double dt, double *x) override
 Prepare for time stepping beginning with solution x and timestep dt.
 
void setSteadyMode () override
 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 eval (double *x, double *r, double rdt=-1.0, int count=1) override
 Evaluate the residual function.
 
void evalJacobian (double *x0) override
 Evaluates the Jacobian at x0 using finite differences.
 
Domain1DpointDomain (size_t i)
 Return a pointer to the domain global point i belongs to.
 
void resize () override
 Call to set the size of internal data structures after first defining the system or if the problem size changes, for example after grid refinement.
 
void resetBadValues (double *x) override
 Reset values such as negative species concentrations.
 
void writeStats (int printTime=1)
 Write statistics about the number of iterations and Jacobians at each grid level.
 
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()
 
- Public Member Functions inherited from SteadyStateSystem
 SteadyStateSystem (const SteadyStateSystem &)=delete
 
SteadyStateSystemoperator= (const SteadyStateSystem &)=delete
 
virtual void eval (double *x, double *r, double rdt=-1.0, int count=1)=0
 Evaluate the residual function.
 
virtual void evalJacobian (double *x0)=0
 Evaluates the Jacobian at x0 using finite differences.
 
virtual double weightedNorm (const double *step) const =0
 Compute the weighted norm of step.
 
int solve (double *x0, double *x1, int loglevel)
 Solve \( F(x) = 0 \), where \( F(x) \) is the residual function.
 
void setInitialGuess (const double *x)
 Set the initial guess. Should be called before solve().
 
void solve (int loglevel=0)
 Solve the steady-state problem, taking internal timesteps as necessary until the Newton solver can converge for the steady problem.
 
void getState (double *x) const
 Get the converged steady-state solution after calling solve().
 
double ssnorm (double *x, double *r)
 Steady-state max norm (infinity norm) of the residual evaluated using solution x.
 
size_t size () const
 Total solution vector length;.
 
virtual void resize ()
 Call to set the size of internal data structures after first defining the system or if the problem size changes, for example after grid refinement.
 
size_t bandwidth () const
 Jacobian bandwidth.
 
virtual string componentName (size_t i) const
 Get the name of the i-th component of the state vector.
 
virtual pair< string, string > componentTableHeader () const
 Get header lines describing the column names included in a component label.
 
virtual string componentTableLabel (size_t i) const
 Get elements of the component name, aligned with the column headings given by componentTableHeader().
 
virtual double upperBound (size_t i) const =0
 Get the upper bound for global component i in the state vector.
 
virtual double lowerBound (size_t i) const =0
 Get the lower bound for global component i in the state vector.
 
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 the linear solver being used to hold the Jacobian matrix and solve linear systems as part of each Newton iteration.
 
double rdt () const
 Reciprocal of the time step.
 
bool transient () const
 True if transient mode.
 
bool steady () const
 True if steady mode.
 
virtual void initTimeInteg (double dt, double *x)
 Prepare for time stepping beginning with solution x and timestep dt.
 
virtual void setSteadyMode ()
 Prepare to solve the steady-state problem.
 
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.
 
virtual void resetBadValues (double *x)
 Reset values such as negative species concentrations.
 
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 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.
 
virtual void writeDebugInfo (const string &header_suffix, const string &message, int loglevel, int attempt_counter)
 Write solver debugging based on the specified log level.
 
virtual void clearDebugFile ()
 Delete debug output file that may be created by writeDebugInfo() when solving with high loglevel.
 
void setTimeStep (double stepsize, size_t n, const int *tsteps)
 Set the number of time steps to try when the steady Newton solver is unsuccessful.
 
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
 Get the maximum number of timeteps allowed before successful steady-state solve.
 

Protected Attributes

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< std::tuple< size_t, size_t, size_t > > m_componentInfo
 Domain, grid point, and component indices for each element of the global state vector.
 
size_t m_pts = 0
 Total number of points.
 
- Protected Attributes inherited from SteadyStateSystem
vector< int > m_steps = { 10 }
 Array of number of steps to take after each unsuccessful steady-state solve before re-attempting the steady-state solution.
 
double m_tstep = 1.0e-5
 Initial timestep.
 
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.
 
vector< double > m_xnew
 Work array used to hold the residual or the new solution.
 
vector< double > m_xlast_ts
 State vector after the last successful set of time steps.
 
unique_ptr< MultiNewtonm_newt
 Newton iterator.
 
double m_rdt = 0.0
 Reciprocal of time step.
 
shared_ptr< SystemJacobianm_jac
 Jacobian evaluator.
 
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< int > m_mask
 Transient mask.
 
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.
 
int m_attempt_counter = 0
 Counter used to manage the number of states stored in the debug log file generated by writeDebugInfo()
 
int m_max_history = 10
 Constant that determines the maximum number of states stored in the debug log file generated by writeDebugInfo()
 
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)
 

Additional Inherited Members

- Protected Member Functions inherited from SteadyStateSystem
void evalSSJacobian (double *x, double *rsd)
 Evaluate the steady-state Jacobian, accessible via linearSolver()
 

Constructor & Destructor Documentation

◆ OneDim() [1/2]

OneDim ( )
default

Default constructor.

◆ OneDim() [2/2]

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

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

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

◆ getJacobian()

shared_ptr< SystemJacobian > getJacobian ( )
inline
Deprecated:
To be removed after Cantera 3.2. Use linearSolver() instead.

Definition at line 39 of file OneDim.h.

◆ weightedNorm()

double weightedNorm ( const double *  step) const
overridevirtual

Compute the weighted norm of a step vector.

The weighted norm of a step vector \( \Delta x \) is defined as

\[ ||\Delta x|| = \sqrt{ \frac{1}{N} \sum_{d,n,j} \left( \frac{\Delta x_{d,n,j}}{w_{d,n}} \right)^2 } \]

where the error weight for solution component \( n \) in domain \( d \) is given by

\[ w_{d,n} = \frac{\epsilon_{{\rm r};d,n}}{J_d} \sum_j |x_{d,n,j}| + \epsilon_{{\rm a};d,n} \]

Here, \( \epsilon_{{\rm r};d,n} \) is the relative error tolerance for component \( n \) in domain \( d \), and multiplies the average magnitude of solution component \( n \) in the domain. The second term, \( \epsilon_{{\rm a};d,n} \), is the absolute error tolerance for component \( n \) in domain \( d \). \( N \) is the total number of state variables across all domains and \( J_d \) is the number of grid points in domain \( d \).

Implements SteadyStateSystem.

Definition at line 98 of file OneDim.cpp.

◆ nDomains()

size_t nDomains ( ) const
inline

Number of domains.

Definition at line 73 of file OneDim.h.

◆ domain()

Domain1D & domain ( size_t  i) const
inline

Return a reference to domain i.

Definition at line 78 of file OneDim.h.

◆ domainIndex()

size_t domainIndex ( const string &  name)

Get the index of the domain named name.

Definition at line 29 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 87 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 96 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 104 of file OneDim.h.

◆ left()

Domain1D * left ( )
inline

Pointer to left-most domain (first added).

Definition at line 115 of file OneDim.h.

◆ right()

Domain1D * right ( )
inline

Pointer to right-most domain (last added).

Definition at line 120 of file OneDim.h.

◆ nVars()

size_t nVars ( size_t  jg)
inline

Number of solution components at global point jg.

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

◆ component()

std::tuple< string, size_t, string > component ( size_t  i) const

Return the domain, local point index, and component name for the i-th component of the global solution vector.

Definition at line 39 of file OneDim.cpp.

◆ componentName()

string componentName ( size_t  i) const
overridevirtual

Get the name of the i-th component of the state vector.

Reimplemented from SteadyStateSystem.

Definition at line 45 of file OneDim.cpp.

◆ componentTableHeader()

pair< string, string > componentTableHeader ( ) const
overridevirtual

Get header lines describing the column names included in a component label.

Headings should be aligned with items formatted in the output of componentTableLabel(). Two lines are allowed. If only one is needed, the first line should be blank.

Reimplemented from SteadyStateSystem.

Definition at line 50 of file OneDim.cpp.

◆ componentTableLabel()

string componentTableLabel ( size_t  i) const
overridevirtual

Get elements of the component name, aligned with the column headings given by componentTableHeader().

Reimplemented from SteadyStateSystem.

Definition at line 55 of file OneDim.cpp.

◆ upperBound()

double upperBound ( size_t  i) const
overridevirtual

Get the upper bound for global component i in the state vector.

Implements SteadyStateSystem.

Definition at line 61 of file OneDim.cpp.

◆ lowerBound()

double lowerBound ( size_t  i) const
overridevirtual

Get the lower bound for global component i in the state vector.

Implements SteadyStateSystem.

Definition at line 68 of file OneDim.cpp.

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

◆ points()

size_t points ( )
inline

Total number of points.

Definition at line 152 of file OneDim.h.

◆ initTimeInteg()

void initTimeInteg ( double  dt,
double *  x 
)
overridevirtual

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

Reimplemented from SteadyStateSystem.

Definition at line 329 of file OneDim.cpp.

◆ setSteadyMode()

void setSteadyMode ( )
overridevirtual

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.

Reimplemented from SteadyStateSystem.

Definition at line 340 of file OneDim.cpp.

◆ eval() [1/2]

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

◆ eval() [2/2]

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

Evaluate the residual function.

Parameters
[in]xState vector
[out]rOn return, contains the residual vector
rdtReciprocal of the time step. if omitted, then the internally stored value (accessible using the rdt() method) is used.
countSet to zero to omit this call from the statistics

Implements SteadyStateSystem.

Reimplemented in Sim1D.

Definition at line 172 of file OneDim.h.

◆ evalJacobian()

void evalJacobian ( double *  x0)
overridevirtual

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

Implements SteadyStateSystem.

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

◆ resize()

void resize ( )
overridevirtual

Call to set the size of internal data structures after first defining the system or if the problem size changes, for example after grid refinement.

The m_size variable should be updated before calling this method

Reimplemented from SteadyStateSystem.

Reimplemented in Sim1D.

Definition at line 186 of file OneDim.cpp.

◆ resetBadValues()

void resetBadValues ( double *  x)
overridevirtual

Reset values such as negative species concentrations.

Reimplemented from SteadyStateSystem.

Definition at line 366 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 138 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 155 of file OneDim.cpp.

◆ clearStats()

void clearStats ( )

Clear saved statistics.

Definition at line 173 of file OneDim.cpp.

◆ gridSizeStats()

const vector< size_t > & gridSizeStats ( )
inline

Return total grid size in each call to solve()

Definition at line 215 of file OneDim.h.

◆ jacobianTimeStats()

const vector< double > & jacobianTimeStats ( )
inline

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

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

◆ jacobianCountStats()

const vector< int > & jacobianCountStats ( )
inline

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

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

◆ timeStepStats()

const vector< int > & timeStepStats ( )
inline

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

Definition at line 247 of file OneDim.h.

Member Data Documentation

◆ m_dom

vector<shared_ptr<Domain1D> > m_dom
protected

All domains comprising the system.

Definition at line 254 of file OneDim.h.

◆ m_connect

vector<shared_ptr<Domain1D> > m_connect
protected

All connector and boundary domains.

Definition at line 257 of file OneDim.h.

◆ m_bulk

vector<shared_ptr<Domain1D> > m_bulk
protected

All bulk/flow domains.

Definition at line 260 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 263 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 267 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 271 of file OneDim.h.

◆ m_componentInfo

vector<std::tuple<size_t, size_t, size_t> > m_componentInfo
protected

Domain, grid point, and component indices for each element of the global state vector.

Length size()

Definition at line 275 of file OneDim.h.

◆ m_pts

size_t m_pts = 0
protected

Total number of points.

Definition at line 278 of file OneDim.h.

◆ m_nevals

int m_nevals = 0
private

Number of calls to eval()

Definition at line 284 of file OneDim.h.

◆ m_evaltime

double m_evaltime = 0
private

Total time [s] spent in eval()

Definition at line 285 of file OneDim.h.

◆ m_gridpts

vector<size_t> m_gridpts
private

Number of grid points in this grid.

Definition at line 287 of file OneDim.h.

◆ m_jacEvals

vector<int> m_jacEvals
private

Number of Jacobian evaluations on this grid.

Definition at line 288 of file OneDim.h.

◆ m_jacElapsed

vector<double> m_jacElapsed
private

Time [s] spent evaluating Jacobians on this grid.

Definition at line 289 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 293 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 297 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 301 of file OneDim.h.


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