Cantera  3.2.0a2
Loading...
Searching...
No Matches
SteadyReactorSolver Class Reference

Adapter class to enable using the SteadyStateSystem solver with ReactorNet. More...

#include <ReactorNet.h>

Inheritance diagram for SteadyReactorSolver:
[legend]

Detailed Description

Adapter class to enable using the SteadyStateSystem solver with ReactorNet.

See also
ReactorNet::solveSteady
Since
New in Cantera 3.2.

Definition at line 440 of file ReactorNet.h.

Public Member Functions

 SteadyReactorSolver (ReactorNet *net, double *x0)
 
void eval (double *x, double *r, double rdt=-1.0, int count=1) override
 Evaluate the residual function.
 
void initTimeInteg (double dt, double *x) override
 Prepare for time stepping beginning with solution x and timestep dt.
 
void evalJacobian (double *x0) override
 Evaluates the Jacobian at x0 using finite differences.
 
double weightedNorm (const double *step) const override
 Compute the weighted norm of step.
 
string componentName (size_t i) const override
 Get the name of the i-th component of the state vector.
 
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 resetBadValues (double *x) override
 Reset values such as negative species concentrations.
 
void writeDebugInfo (const string &header_suffix, const string &message, int loglevel, int attempt_counter) override
 Write solver debugging based on the specified log level.
 
- 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.
 

Private Attributes

ReactorNetm_net = nullptr
 
vector< double > m_initialState
 Initial value of each state variable.
 
vector< size_t > m_algebraic
 Indices of variables that are held constant in the time-stepping mode of the steady-state solver.
 

Additional Inherited Members

- Protected Member Functions inherited from SteadyStateSystem
void evalSSJacobian (double *x, double *rsd)
 Evaluate the steady-state Jacobian, accessible via linearSolver()
 
- 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.
 

Constructor & Destructor Documentation

◆ SteadyReactorSolver()

SteadyReactorSolver ( ReactorNet net,
double *  x0 
)

Definition at line 744 of file ReactorNet.cpp.

Member Function Documentation

◆ eval()

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

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.

Definition at line 766 of file ReactorNet.cpp.

◆ 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 782 of file ReactorNet.cpp.

◆ 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 788 of file ReactorNet.cpp.

◆ weightedNorm()

double weightedNorm ( const double *  step) const
overridevirtual

Compute the weighted norm of step.

The purpose is to measure the "size" of the step vector \( \Delta x \) in a way that takes into account the relative importance or scale of different solution components. Each component of the step vector is normalized by a weight that depends on both the current magnitude of the solution vector and specified tolerances. This makes the norm dimensionless and scaled appropriately, avoiding issues where some components dominate due to differences in their scales. See OneDim::weightedNorm() for a representative implementation.

Implements SteadyStateSystem.

Definition at line 823 of file ReactorNet.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 835 of file ReactorNet.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 840 of file ReactorNet.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 845 of file ReactorNet.cpp.

◆ resetBadValues()

void resetBadValues ( double *  x)
overridevirtual

Reset values such as negative species concentrations.

Reimplemented from SteadyStateSystem.

Definition at line 850 of file ReactorNet.cpp.

◆ writeDebugInfo()

void writeDebugInfo ( const string &  header_suffix,
const string &  message,
int  loglevel,
int  attempt_counter 
)
overridevirtual

Write solver debugging based on the specified log level.

See also
Sim1D::writeDebugInfo for a specific implementation of this capability.

Reimplemented from SteadyStateSystem.

Definition at line 855 of file ReactorNet.cpp.

Member Data Documentation

◆ m_net

ReactorNet* m_net = nullptr
private

Definition at line 456 of file ReactorNet.h.

◆ m_initialState

vector<double> m_initialState
private

Initial value of each state variable.

Definition at line 459 of file ReactorNet.h.

◆ m_algebraic

vector<size_t> m_algebraic
private

Indices of variables that are held constant in the time-stepping mode of the steady-state solver.

Definition at line 463 of file ReactorNet.h.


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