Cantera  3.1.0a1
ImplicitSurfChem Class Reference

Advances the surface coverages of the associated set of SurfacePhase objects in time. More...

#include <ImplicitSurfChem.h>

Inheritance diagram for ImplicitSurfChem:
[legend]

Detailed Description

Advances the surface coverages of the associated set of SurfacePhase objects in time.

This function advances a set of SurfPhase objects, each associated with one InterfaceKinetics object, in time. The following equation is used for each surface phase, i.

\[ \dot \theta_k = \dot s_k (\sigma_k / s_0) \]

In this equation,

  • \( \theta_k \) is the site coverage for the kth species.
  • \( \dot s_k \) is the source term for the kth species
  • \( \sigma_k \) is the number of surface sites covered by each species k.
  • \( s_0 \) is the total site density of the interfacial phase.

Additionally, the 0'th equation in the set is discarded. Instead the alternate equation is solved for

\[ \sum_{k=0}^{N-1} \dot \theta_k = 0 \]

This last equation serves to ensure that sum of the \( \theta_k \) values stays constant.

The object uses the CVODE software to advance the surface equations.

The solution vector used by this object is as follows: For each surface phase with \( N_s \) surface sites, it consists of the surface coverages \( \theta_k \) for \( k = 0, N_s - 1 \)

Definition at line 58 of file ImplicitSurfChem.h.

Public Member Functions

 ImplicitSurfChem (vector< InterfaceKinetics * > k, double rtol=1.e-7, double atol=1.e-14, double maxStepSize=0, size_t maxSteps=20000, size_t maxErrTestFails=7)
 Constructor for multiple surfaces. More...
 
void initialize (double t0=0.0)
 Must be called before calling method 'advance'. More...
 
void setMaxStepSize (double maxstep=0.0)
 Set the maximum integration step-size. More...
 
void setTolerances (double rtol=1.e-7, double atol=1.e-14)
 Set the relative and absolute integration tolerances. More...
 
void setMaxSteps (size_t maxsteps=20000)
 Set the maximum number of CVODES integration steps. More...
 
void setMaxErrTestFails (size_t maxErrTestFails=7)
 Set the maximum number of CVODES error test failures. More...
 
void integrate (double t0, double t1)
 Integrate from t0 to t1. The integrator is reinitialized first. More...
 
void integrate0 (double t0, double t1)
 Integrate from t0 to t1 without reinitializing the integrator. More...
 
void solvePseudoSteadyStateProblem (int ifuncOverride=-1, double timeScaleOverride=1.0)
 Solve for the pseudo steady-state of the surface problem. More...
 
size_t neq () const override
 Return the number of equations. More...
 
void eval (double t, double *y, double *ydot, double *p) override
 Evaluate the value of ydot[k] at the current conditions. More...
 
void getState (double *y) override
 Get the current state of the solution vector. More...
 
void getConcSpecies (double *const vecConcSpecies) const
 Get the specifications for the problem from the values in the ThermoPhase objects for all phases. More...
 
void setConcSpecies (const double *const vecConcSpecies)
 Sets the concentrations within phases that are unknowns in the surface problem. More...
 
void setCommonState_TP (double TKelvin, double PresPa)
 Sets the state variable in all thermodynamic phases (surface and surrounding bulk phases) to the input temperature and pressure. More...
 
vector< InterfaceKinetics * > & getObjects ()
 Returns a reference to the vector of pointers to the InterfaceKinetics objects. More...
 
int checkMatch (vector< ThermoPhase * > m_vec, ThermoPhase *thPtr)
 
void setIOFlag (int ioFlag)
 
- Public Member Functions inherited from FuncEval
virtual void evalDae (double t, double *y, double *ydot, double *p, double *residual)
 Evaluate the right-hand-side DAE function. More...
 
virtual void getConstraints (double *constraints)
 Given a vector of length neq(), mark which variables should be considered algebraic constraints. More...
 
int evalNoThrow (double t, double *y, double *ydot)
 Evaluate the right-hand side using return code to indicate status. More...
 
int evalDaeNoThrow (double t, double *y, double *ydot, double *residual)
 Evaluate the right-hand side using return code to indicate status. More...
 
virtual void preconditionerSetup (double t, double *y, double gamma)
 Evaluate the setup processes for the Jacobian preconditioner. More...
 
virtual void preconditionerSolve (double *rhs, double *output)
 Evaluate the linear system Ax=b where A is the preconditioner. More...
 
virtual void updatePreconditioner (double gamma)
 Update the preconditioner based on already computed jacobian values. More...
 
int preconditioner_setup_nothrow (double t, double *y, double gamma)
 Preconditioner setup that doesn't throw an error but returns a CVODES flag. More...
 
int preconditioner_solve_nothrow (double *rhs, double *output)
 Preconditioner solve that doesn't throw an error but returns a CVODES flag. More...
 
virtual void getStateDae (double *y, double *ydot)
 Fill in the vectors y and ydot with the current state of the system. More...
 
virtual size_t nparams () const
 Number of sensitivity parameters. More...
 
void suppressErrors (bool suppress)
 Enable or disable suppression of errors when calling eval() More...
 
bool suppressErrors () const
 Get current state of error suppression. More...
 
string getErrors () const
 Return a string containing the text of any suppressed errors. More...
 
void clearErrors ()
 Clear any previously-stored suppressed errors. More...
 

Protected Member Functions

void updateState (double *y)
 Set the mixture to a state consistent with solution vector y. More...
 

Protected Attributes

vector< SurfPhase * > m_surf
 vector of pointers to surface phases. More...
 
vector< ThermoPhase * > m_bulkPhases
 Vector of pointers to bulk phases. More...
 
vector< InterfaceKinetics * > m_vecKinPtrs
 vector of pointers to InterfaceKinetics objects More...
 
vector< size_t > m_nsp
 Vector of number of species in each Surface Phase. More...
 
vector< size_t > m_specStartIndex
 
size_t m_nv = 0
 Total number of surface species in all surface phases. More...
 
size_t m_numTotalBulkSpecies = 0
 
size_t m_numTotalSpecies = 0
 
vector< vector< int > > pLocVec
 
unique_ptr< Integratorm_integ
 Pointer to the CVODE integrator. More...
 
double m_atol
 
double m_rtol
 
double m_maxstep
 max step size More...
 
size_t m_nmax
 maximum number of steps allowed More...
 
size_t m_maxErrTestFails
 maximum number of error test failures allowed More...
 
vector< double > m_work
 
vector< double > m_concSpecies
 Temporary vector - length num species in the Kinetics object. More...
 
vector< double > m_concSpeciesSave
 
int m_mediumSpeciesStart = -1
 Index into the species vector of the kinetics manager, pointing to the first species from the surrounding medium. More...
 
int m_bulkSpeciesStart = -1
 Index into the species vector of the kinetics manager, pointing to the first species from the condensed phase of the particles. More...
 
int m_surfSpeciesStart = -1
 Index into the species vector of the kinetics manager, pointing to the first species from the surface of the particles. More...
 
unique_ptr< solveSPm_surfSolver
 Pointer to the helper method, Placid, which solves the surface problem. More...
 
bool m_commonTempPressForPhases = true
 If true, a common temperature and pressure for all surface and bulk phases associated with the surface problem is imposed. More...
 
- Protected Attributes inherited from FuncEval
bool m_suppress_errors = false
 
vector< string > m_errors
 Errors occurring during function evaluations. More...
 

Private Attributes

int m_ioFlag = 0
 Controls the amount of printing from this routine and underlying routines. More...
 

Additional Inherited Members

- Public Attributes inherited from FuncEval
vector< double > m_sens_params
 Values for the problem parameters for which sensitivities are computed This is the array which is perturbed and passed back as the fourth argument to eval(). More...
 
vector< double > m_paramScales
 Scaling factors for each sensitivity parameter. More...
 

Constructor & Destructor Documentation

◆ ImplicitSurfChem()

ImplicitSurfChem ( vector< InterfaceKinetics * >  k,
double  rtol = 1.e-7,
double  atol = 1.e-14,
double  maxStepSize = 0,
size_t  maxSteps = 20000,
size_t  maxErrTestFails = 7 
)

Constructor for multiple surfaces.

Parameters
kVector of pointers to InterfaceKinetics objects Each object consists of a surface or an edge containing internal degrees of freedom representing the concentration of surface adsorbates.
rtolThe relative tolerance for the integrator
atolThe absolute tolerance for the integrator
maxStepSizeThe maximum step-size the integrator is allowed to take. If zero, this option is disabled.
maxStepsThe maximum number of time-steps the integrator can take. If not supplied, uses the default value in the CVodesIntegrator (20000).
maxErrTestFailsthe maximum permissible number of error test failures If not supplied, uses the default value in CVODES (7).

Definition at line 18 of file ImplicitSurfChem.cpp.

Member Function Documentation

◆ initialize()

void initialize ( double  t0 = 0.0)

Must be called before calling method 'advance'.

Definition at line 124 of file ImplicitSurfChem.cpp.

◆ setMaxStepSize()

void setMaxStepSize ( double  maxstep = 0.0)

Set the maximum integration step-size.

Note, setting this value to zero disables this option

Definition at line 97 of file ImplicitSurfChem.cpp.

◆ setTolerances()

void setTolerances ( double  rtol = 1.e-7,
double  atol = 1.e-14 
)

Set the relative and absolute integration tolerances.

Definition at line 105 of file ImplicitSurfChem.cpp.

◆ setMaxSteps()

void setMaxSteps ( size_t  maxsteps = 20000)

Set the maximum number of CVODES integration steps.

Definition at line 112 of file ImplicitSurfChem.cpp.

◆ setMaxErrTestFails()

void setMaxErrTestFails ( size_t  maxErrTestFails = 7)

Set the maximum number of CVODES error test failures.

Definition at line 118 of file ImplicitSurfChem.cpp.

◆ integrate()

void integrate ( double  t0,
double  t1 
)

Integrate from t0 to t1. The integrator is reinitialized first.

This routine does a time accurate solve from t = t0 to t = t1. of the surface problem.

Parameters
t0Initial Time -> this is an input
t1Final Time -> This is an input

Definition at line 133 of file ImplicitSurfChem.cpp.

◆ integrate0()

void integrate0 ( double  t0,
double  t1 
)

Integrate from t0 to t1 without reinitializing the integrator.

Use when the coverages have not changed from their values on return from the last call to integrate or integrate0.

Parameters
t0Initial Time -> this is an input
t1Final Time -> This is an input

Definition at line 144 of file ImplicitSurfChem.cpp.

◆ solvePseudoSteadyStateProblem()

void solvePseudoSteadyStateProblem ( int  ifuncOverride = -1,
double  timeScaleOverride = 1.0 
)

Solve for the pseudo steady-state of the surface problem.

Solve for the steady state of the surface problem. This is the same thing as the advanceCoverages() function, but at infinite times.

Note, a direct solve is carried out under the hood here, to reduce the computational time.

Parameters
ifuncOverrideOne of the values defined in Surface Problem Solver Methods. The default is -1, which means that the program will decide.
timeScaleOverrideWhen a pseudo transient is selected this value can be used to override the default time scale for integration which is one. When SFLUX_TRANSIENT is used, this is equal to the time over which the equations are integrated. When SFLUX_INITIALIZE is used, this is equal to the time used in the initial transient algorithm, before the equation system is solved directly.

Definition at line 176 of file ImplicitSurfChem.cpp.

◆ neq()

size_t neq ( ) const
inlineoverridevirtual

Return the number of equations.

Implements FuncEval.

Definition at line 151 of file ImplicitSurfChem.h.

◆ eval()

void eval ( double  t,
double *  y,
double *  ydot,
double *  p 
)
overridevirtual

Evaluate the value of ydot[k] at the current conditions.

Parameters
tTime (seconds)
yVector containing the current solution vector
ydotOutput vector containing the value of the derivative of the surface coverages.
pUnused parameter pass-through parameter vector

Reimplemented from FuncEval.

Definition at line 159 of file ImplicitSurfChem.cpp.

◆ getState()

void getState ( double *  y)
overridevirtual

Get the current state of the solution vector.

Parameters
yValue of the solution vector to be used. On output, this contains the initial value of the solution.

Reimplemented from FuncEval.

Definition at line 88 of file ImplicitSurfChem.cpp.

◆ getConcSpecies()

void getConcSpecies ( double *const  vecConcSpecies) const

Get the specifications for the problem from the values in the ThermoPhase objects for all phases.

  1. concentrations of all species in all phases, m_concSpecies
  2. Temperature and pressure
Parameters
vecConcSpeciesVector of concentrations. The phase concentration vectors are contiguous within the object, in the same order as the unknown vector.

Definition at line 252 of file ImplicitSurfChem.cpp.

◆ setConcSpecies()

void setConcSpecies ( const double *const  vecConcSpecies)

Sets the concentrations within phases that are unknowns in the surface problem.

Fills the local concentration vector for all of the species in all of the phases that are unknowns in the surface problem.

Parameters
vecConcSpeciesVector of concentrations. The phase concentration vectors are contiguous within the object, in the same order as the unknown vector.

Definition at line 268 of file ImplicitSurfChem.cpp.

◆ setCommonState_TP()

void setCommonState_TP ( double  TKelvin,
double  PresPa 
)

Sets the state variable in all thermodynamic phases (surface and surrounding bulk phases) to the input temperature and pressure.

Parameters
TKelvininput temperature (kelvin)
PresPainput pressure in pascal.

Definition at line 284 of file ImplicitSurfChem.cpp.

◆ getObjects()

vector<InterfaceKinetics*>& getObjects ( )
inline

Returns a reference to the vector of pointers to the InterfaceKinetics objects.

This should probably go away in the future, as it opens up the class.

Definition at line 211 of file ImplicitSurfChem.h.

◆ updateState()

void updateState ( double *  y)
protected

Set the mixture to a state consistent with solution vector y.

This function will set the surface site factions in the underlying SurfPhase objects to the current value of the solution vector.

Parameters
yCurrent value of the solution vector. The length is equal to the sum of the number of surface sites in all the surface phases.

Definition at line 150 of file ImplicitSurfChem.cpp.

Member Data Documentation

◆ m_surf

vector<SurfPhase*> m_surf
protected

vector of pointers to surface phases.

Definition at line 234 of file ImplicitSurfChem.h.

◆ m_bulkPhases

vector<ThermoPhase*> m_bulkPhases
protected

Vector of pointers to bulk phases.

Definition at line 237 of file ImplicitSurfChem.h.

◆ m_vecKinPtrs

vector<InterfaceKinetics*> m_vecKinPtrs
protected

vector of pointers to InterfaceKinetics objects

Definition at line 240 of file ImplicitSurfChem.h.

◆ m_nsp

vector<size_t> m_nsp
protected

Vector of number of species in each Surface Phase.

Definition at line 243 of file ImplicitSurfChem.h.

◆ m_nv

size_t m_nv = 0
protected

Total number of surface species in all surface phases.

This is the total number of unknowns in m_mode 0 problem

Definition at line 251 of file ImplicitSurfChem.h.

◆ m_integ

unique_ptr<Integrator> m_integ
protected

Pointer to the CVODE integrator.

Definition at line 258 of file ImplicitSurfChem.h.

◆ m_maxstep

double m_maxstep
protected

max step size

Definition at line 260 of file ImplicitSurfChem.h.

◆ m_nmax

size_t m_nmax
protected

maximum number of steps allowed

Definition at line 261 of file ImplicitSurfChem.h.

◆ m_maxErrTestFails

size_t m_maxErrTestFails
protected

maximum number of error test failures allowed

Definition at line 262 of file ImplicitSurfChem.h.

◆ m_concSpecies

vector<double> m_concSpecies
protected

Temporary vector - length num species in the Kinetics object.

This is the sum of the number of species in each phase included in the kinetics object.

Definition at line 270 of file ImplicitSurfChem.h.

◆ m_mediumSpeciesStart

int m_mediumSpeciesStart = -1
protected

Index into the species vector of the kinetics manager, pointing to the first species from the surrounding medium.

Definition at line 277 of file ImplicitSurfChem.h.

◆ m_bulkSpeciesStart

int m_bulkSpeciesStart = -1
protected

Index into the species vector of the kinetics manager, pointing to the first species from the condensed phase of the particles.

Definition at line 282 of file ImplicitSurfChem.h.

◆ m_surfSpeciesStart

int m_surfSpeciesStart = -1
protected

Index into the species vector of the kinetics manager, pointing to the first species from the surface of the particles.

Definition at line 287 of file ImplicitSurfChem.h.

◆ m_surfSolver

unique_ptr<solveSP> m_surfSolver
protected

Pointer to the helper method, Placid, which solves the surface problem.

Definition at line 291 of file ImplicitSurfChem.h.

◆ m_commonTempPressForPhases

bool m_commonTempPressForPhases = true
protected

If true, a common temperature and pressure for all surface and bulk phases associated with the surface problem is imposed.

Definition at line 295 of file ImplicitSurfChem.h.

◆ m_ioFlag

int m_ioFlag = 0
private

Controls the amount of printing from this routine and underlying routines.

Definition at line 300 of file ImplicitSurfChem.h.


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