Cantera  2.0
Public Member Functions | Public Attributes | Private Member Functions | Private Attributes | List of all members
solveSP Class Reference

Method to solve a pseudo steady state surface problem. More...

#include <solveSP.h>

Collaboration diagram for solveSP:
[legend]

Public Member Functions

 solveSP (ImplicitSurfChem *surfChemPtr, int bulkFunc=BULK_ETCH)
 Constructor for the object.
 
 ~solveSP ()
 Destructor. Deletes the integrator.
 
int solveSurfProb (int ifunc, doublereal time_scale, doublereal TKelvin, doublereal PGas, doublereal reltol, doublereal abstol)
 Main routine that actually calculates the pseudo steady state of the surface problem.
 

Public Attributes

int m_ioflag
 

Private Member Functions

 solveSP (const solveSP &right)
 Unimplemented private copy constructor.
 
solveSPoperator= (const solveSP &right)
 Unimplemented private assignment operator.
 
void print_header (int ioflag, int ifunc, doublereal time_scale, int damping, doublereal reltol, doublereal abstol, doublereal TKelvin, doublereal PGas, doublereal netProdRate[], doublereal XMolKinSpecies[])
 Printing routine that gets called at the start of every invocation.
 
void printIteration (int ioflag, doublereal damp, int label_d, int label_t, doublereal inv_t, doublereal t_real, size_t iter, doublereal update_norm, doublereal resid_norm, doublereal netProdRate[], doublereal CSolnSP[], doublereal resid[], doublereal XMolSolnSP[], doublereal wtSpecies[], size_t dim, bool do_time)
 Printing routine that gets called after every iteration.
 
void printFinal (int ioflag, doublereal damp, int label_d, int label_t, doublereal inv_t, doublereal t_real, size_t iter, doublereal update_norm, doublereal resid_norm, doublereal netProdRateKinSpecies[], const doublereal CSolnSP[], const doublereal resid[], doublereal XMolSolnSP[], const doublereal wtSpecies[], const doublereal wtRes[], size_t dim, bool do_time, doublereal TKelvin, doublereal PGas)
 Print a summary of the solution.
 
doublereal calc_t (doublereal netProdRateSolnSP[], doublereal XMolSolnSP[], int *label, int *label_old, doublereal *label_factor, int ioflag)
 Calculate a conservative delta T to use in a pseudo-steady state algorithm.
 
void calcWeights (doublereal wtSpecies[], doublereal wtResid[], const Array2D &Jac, const doublereal CSolnSP[], const doublereal abstol, const doublereal reltol)
 Calculate the solution and residual weights.
 
void updateState (const doublereal *cSurfSpec)
 Update the surface states of the surface phases.
 
void updateMFSolnSP (doublereal *XMolSolnSP)
 Update mole fraction vector consisting of unknowns in surface problem.
 
void updateMFKinSpecies (doublereal *XMolKinSp, int isp)
 Update the mole fraction vector for a specific kinetic species vector corresponding to one InterfaceKinetics object.
 
void evalSurfLarge (const doublereal *CSolnSP)
 Update the vector that keeps track of the largest species in each surface phase.
 
void fun_eval (doublereal *resid, const doublereal *CSolnSP, const doublereal *CSolnOldSP, const bool do_time, const doublereal deltaT)
 Main Function evaluation.
 
void resjac_eval (std::vector< doublereal * > &JacCol, doublereal *resid, doublereal *CSolnSP, const doublereal *CSolnSPOld, const bool do_time, const doublereal deltaT)
 Main routine that calculates the current residual and Jacobian.
 

Private Attributes

ImplicitSurfChemm_SurfChemPtr
 Pointer to the manager of the implicit surface chemistry problem.
 
std::vector< InterfaceKinetics * > & m_objects
 Vector of interface kinetics objects.
 
size_t m_neq
 Total number of equations to solve in the implicit problem.
 
int m_bulkFunc
 This variable determines how the bulk phases are to be handled.
 
size_t m_numSurfPhases
 Number of surface phases in the surface problem.
 
size_t m_numTotSurfSpecies
 Total number of surface species in all surface phases.
 
std::vector< size_t > m_indexKinObjSurfPhase
 Mapping between the surface phases and the InterfaceKinetics objects.
 
std::vector< size_t > m_nSpeciesSurfPhase
 Vector of length number of surface phases containing the number of surface species in each phase.
 
std::vector< SurfPhase * > m_ptrsSurfPhase
 Vector of surface phase pointers.
 
std::vector< size_t > m_eqnIndexStartSolnPhase
 Index of the start of the unknowns for each solution phase.
 
std::vector< size_t > m_kinObjPhaseIDSurfPhase
 Phase ID in the InterfaceKinetics object of the surface phase.
 
size_t m_numBulkPhasesSS
 Total number of volumetric condensed phases included in the steady state problem handled by this routine.
 
std::vector< size_t > m_numBulkSpecies
 Vector of number of species in the m_numBulkPhases phases.
 
size_t m_numTotBulkSpeciesSS
 Total number of species in all bulk phases.
 
std::vector< ThermoPhase * > m_bulkPhasePtrs
 Vector of bulk phase pointers, length is equal to m_numBulkPhases.
 
std::vector< size_t > m_kinSpecIndex
 Index between the equation index and the position in the kinetic species array for the appropriate kinetics operator.
 
std::vector< size_t > m_kinObjIndex
 Index between the equation index and the index of the InterfaceKinetics object.
 
std::vector< size_t > m_spSurfLarge
 Vector containing the indices of the largest species in each surface phase.
 
doublereal m_atol
 m_atol is the absolute tolerance in real units.
 
doublereal m_rtol
 m_rtol is the relative error tolerance.
 
doublereal m_maxstep
 maximum value of the time step
 
size_t m_maxTotSpecies
 Maximum number of species in any single kinetics operator -> also maxed wrt the total # of solution species.
 
vector_fp m_netProductionRatesSave
 Temporary vector with length equal to max m_maxTotSpecies.
 
vector_fp m_numEqn1
 Temporary vector with length equal to max m_maxTotSpecies.
 
vector_fp m_numEqn2
 Temporary vector with length equal to max m_maxTotSpecies.
 
vector_fp m_CSolnSave
 Temporary vector with length equal to max m_maxTotSpecies.
 
vector_fp m_CSolnSP
 Solution vector.
 
vector_fp m_CSolnSPInit
 Saved initial solution vector.
 
vector_fp m_CSolnSPOld
 Saved solution vector at the old time step.
 
vector_fp m_wtResid
 Weights for the residual norm calculation.
 
vector_fp m_wtSpecies
 Weights for the species concentrations norm calculation.
 
vector_fp m_resid
 Residual for the surface problem.
 
vector_fp m_XMolKinSpecies
 Vector of mole fractions.
 
vector_int m_ipiv
 pivots
 
std::vector< doublereal * > m_JacCol
 Vector of pointers to the top of the columns of the jacobians.
 
Array2D m_Jac
 Jacobian.
 

Detailed Description

Method to solve a pseudo steady state surface problem.

The following class handles solving the surface problem. The calculation uses Newton's method to obtain the surface fractions of the surface and bulk species by requiring that the surface species production rate = 0 and that the either the bulk fractions are proportional to their production rates or they are constants.

Currently, the bulk mole fractions are treated as constants. Implementation of their being added to the unknown solution vector is delayed.

Lets introduce the unknown vector for the "surface problem". The surface problem is defined as the evaluation of the surface site fractions for multiple surface phases. The unknown vector will consist of the vector of surface concentrations for each species in each surface vector. Species are grouped first by their surface phases

C_i_j = Concentration of the ith species in the jth surface phase Nj = number of surface species in the jth surface phase

The unknown solution vector is defined as follows:

         kindexSP

C_0_0 0 C_1_0 1 C_2_0 2 . . . ... C_N0-1_0 N0-1 C_0_1 N0 C_1_1 N0+1 C_2_1 N0+2 . . . ... C_N1-1_1 NO+N1-1

Note there are a couple of different types of species indices floating around in the formulation of this object.

kindexSP This is the species index in the contiguous vector of unknowns for the surface problem.

Note, in the future, BULK_DEPOSITION systems will be added, and the solveSP unknown vector will get more complicated. It will include the mole fraction and growth rates of specified bulk phases

Indices which relate to individual kinetics objects use the suffix KSI (kinetics species index).

Solution Method

This routine is typically used within a residual calculation in a large code. It's typically invoked millions of times for large calculations, and it must work every time. Therefore, requirements demand that it be robust but also efficient.

The solution methodology is largely determined by the ifunc<> parameter, that is input to the solution object. This parameter may have the following 4 values:

1: SFLUX_INITIALIZE = This assumes that the initial guess supplied to the routine is far from the correct one. Substantial work plus transient time-stepping is to be expected to find a solution.

2: SFLUX_RESIDUAL = Need to solve the surface problem in order to calculate the surface fluxes of gas-phase species. (Can expect a moderate change in the solution vector -> try to solve the system by direct methods with no damping first -> then, try time-stepping if the first method fails) A "time_scale" supplied here is used in the algorithm to determine when to shut off time-stepping.

3: SFLUX_JACOBIAN = Calculation of the surface problem is due to the need for a numerical jacobian for the gas-problem. The solution is expected to be very close to the initial guess, and extra accuracy is needed because solution variables have been delta'd from nominal values to create jacobian entries.

4: SFLUX_TRANSIENT = The transient calculation is performed here for an amount of time specified by "time_scale". It is not guaranteed to be time-accurate - just stable and fairly fast. The solution after del_t time is returned, whether it's converged to a steady state or not. This is a poor man's time stepping algorithm.

Pseudo time stepping algorithm: The time step is determined from sdot[], so so that the time step doesn't ever change the value of a variable by more than 100%.

This algorithm does use a damped Newton's method to relax the equations. Damping is based on a "delta damping" technique. The solution unknowns are not allowed to vary too much between iterations.

EXTRA_ACCURACY:A constant that is the ratio of the required update norm in this Newton iteration compared to that in the nonlinear solver. A value of 0.1 is used so surface species are safely overconverged.

Functions called:

ct_dgetrf – First half of LAPACK direct solve of a full Matrix

ct_dgetrs – Second half of LAPACK direct solve of a full matrix. Returns solution vector in the right-hand-side vector, resid.


Definition at line 197 of file solveSP.h.

Constructor & Destructor Documentation

solveSP ( ImplicitSurfChem surfChemPtr,
int  bulkFunc = BULK_ETCH 
)
~solveSP ( )

Destructor. Deletes the integrator.

Definition at line 147 of file solveSP.cpp.

solveSP ( const solveSP right)
private

Unimplemented private copy constructor.

Member Function Documentation

solveSP& operator= ( const solveSP right)
private

Unimplemented private assignment operator.

int solveSurfProb ( int  ifunc,
doublereal  time_scale,
doublereal  TKelvin,
doublereal  PGas,
doublereal  reltol,
doublereal  abstol 
)

Main routine that actually calculates the pseudo steady state of the surface problem.

The actual converged solution is returned as part of the internal state of the InterfaceKinetics objects.

Parameters
ifuncDetermines the type of solution algorithm to be used. Possible values are SFLUX_INITIALIZE , SFLUX_RESIDUAL SFLUX_JACOBIAN SFLUX_TRANSIENT .
time_scaleTime over which to integrate the surface equations, where applicable
TKelvinTemperature (kelvin)
PGasPressure (pascals)
reltolRelative tolerance to use
abstolabsolute tolerance.
Returns
Returns 1 if the surface problem is successfully solved. Returns -1 if the surface problem wasn't solved successfully. Note the actual converged solution is returned as part of the internal state of the InterfaceKinetics objects.

Definition at line 158 of file solveSP.cpp.

References solveSP::calc_t(), solveSP::calcWeights(), DATA_PTR, solveSP::evalSurfLarge(), solveSP::fun_eval(), Phase::getConcentrations(), Cantera::int2str(), solveSP::m_CSolnSP, solveSP::m_CSolnSPInit, solveSP::m_CSolnSPOld, solveSP::m_ipiv, solveSP::m_Jac, solveSP::m_JacCol, solveSP::m_kinSpecIndex, solveSP::m_neq, solveSP::m_netProductionRatesSave, solveSP::m_nSpeciesSurfPhase, solveSP::m_numEqn1, solveSP::m_numSurfPhases, solveSP::m_ptrsSurfPhase, solveSP::m_resid, solveSP::m_wtResid, solveSP::m_wtSpecies, solveSP::m_XMolKinSpecies, ckr::max(), solveSP::print_header(), solveSP::printFinal(), solveSP::printIteration(), solveSP::resjac_eval(), SFLUX_INITIALIZE, and solveSP::updateState().

Referenced by ImplicitSurfChem::solvePseudoSteadyStateProblem().

void print_header ( int  ioflag,
int  ifunc,
doublereal  time_scale,
int  damping,
doublereal  reltol,
doublereal  abstol,
doublereal  TKelvin,
doublereal  PGas,
doublereal  netProdRate[],
doublereal  XMolKinSpecies[] 
)
private

Printing routine that gets called at the start of every invocation.

Definition at line 901 of file solveSP.cpp.

References solveSP::m_bulkFunc, and SFLUX_INITIALIZE.

Referenced by solveSP::solveSurfProb().

void printIteration ( int  ioflag,
doublereal  damp,
int  label_d,
int  label_t,
doublereal  inv_t,
doublereal  t_real,
size_t  iter,
doublereal  update_norm,
doublereal  resid_norm,
doublereal  netProdRate[],
doublereal  CSolnSP[],
doublereal  resid[],
doublereal  XMolSolnSP[],
doublereal  wtSpecies[],
size_t  dim,
bool  do_time 
)
private

Printing routine that gets called after every iteration.

Definition at line 953 of file solveSP.cpp.

References Cantera::int2str(), Kinetics::kineticsSpeciesName(), solveSP::m_kinObjIndex, solveSP::m_kinSpecIndex, and solveSP::m_objects.

Referenced by solveSP::solveSurfProb().

void printFinal ( int  ioflag,
doublereal  damp,
int  label_d,
int  label_t,
doublereal  inv_t,
doublereal  t_real,
size_t  iter,
doublereal  update_norm,
doublereal  resid_norm,
doublereal  netProdRateKinSpecies[],
const doublereal  CSolnSP[],
const doublereal  resid[],
doublereal  XMolSolnSP[],
const doublereal  wtSpecies[],
const doublereal  wtRes[],
size_t  dim,
bool  do_time,
doublereal  TKelvin,
doublereal  PGas 
)
private

Print a summary of the solution.

Definition at line 1002 of file solveSP.cpp.

References Cantera::int2str(), Kinetics::kineticsSpeciesName(), solveSP::m_kinObjIndex, solveSP::m_kinSpecIndex, and solveSP::m_objects.

Referenced by solveSP::solveSurfProb().

doublereal calc_t ( doublereal  netProdRateSolnSP[],
doublereal  XMolSolnSP[],
int *  label,
int *  label_old,
doublereal *  label_factor,
int  ioflag 
)
private

Calculate a conservative delta T to use in a pseudo-steady state algorithm.

This routine calculates a pretty conservative 1/del_t based on MAX_i(sdot_i/(X_i*SDen0)). This probably guarantees diagonal dominance.

Small surface fractions are allowed to intervene in the del_t determination, no matter how small. This may be changed. Now minimum changed to 1.0e-12,

Maximum time step set to time_scale.

Parameters
netProdRateSolnSPOutput variable. Net production rate of all of the species in the solution vector.
XMolSolnSPoutput variable. Mole fraction of all of the species in the solution vector
labelOutput variable. Pointer to the value of the species index (kindexSP) that is controlling the time step
label_oldOutput variable. Pointer to the value of the species index (kindexSP) that controlled the time step at the previous iteration
label_factorOutput variable. Pointer to the current factor that is used to indicate the same species is controlling the time step.
ioflagLevel of the output requested.
Returns
Returns the 1. / delta T to be used on the next step

Definition at line 839 of file solveSP.cpp.

References DATA_PTR, InterfaceKinetics::getNetProductionRates(), Kinetics::kineticsSpeciesIndex(), solveSP::m_nSpeciesSurfPhase, solveSP::m_numEqn1, solveSP::m_numSurfPhases, solveSP::m_objects, Phase::molarDensity(), Kinetics::surfacePhaseIndex(), Kinetics::thermo(), and solveSP::updateMFSolnSP().

Referenced by solveSP::solveSurfProb().

void calcWeights ( doublereal  wtSpecies[],
doublereal  wtResid[],
const Array2D Jac,
const doublereal  CSolnSP[],
const doublereal  abstol,
const doublereal  reltol 
)
private

Calculate the solution and residual weights.

Parameters
wtSpeciesWeights to use for the soln unknowns. These are in concentration units
wtResidWeights to sue for the residual unknowns.
JacJacobian. Row sum scaling is used for the Jacobian
CSolnSPSolution vector for the surface problem
abstolAbsolute error tolerance
reltolRelative error tolerance

Definition at line 786 of file solveSP.cpp.

References solveSP::m_bulkFunc, solveSP::m_bulkPhasePtrs, solveSP::m_neq, solveSP::m_nSpeciesSurfPhase, solveSP::m_numBulkPhasesSS, solveSP::m_numBulkSpecies, solveSP::m_numSurfPhases, and solveSP::m_ptrsSurfPhase.

Referenced by solveSP::solveSurfProb().

void updateState ( const doublereal *  cSurfSpec)
private

Update the surface states of the surface phases.

Definition at line 450 of file solveSP.cpp.

References solveSP::m_nSpeciesSurfPhase, solveSP::m_numSurfPhases, and solveSP::m_ptrsSurfPhase.

Referenced by solveSP::fun_eval(), and solveSP::solveSurfProb().

void updateMFSolnSP ( doublereal *  XMolSolnSP)
private

Update mole fraction vector consisting of unknowns in surface problem.

Parameters
XMolSolnSPVector of mole fractions for the unknowns in the surface problem.

Definition at line 462 of file solveSP.cpp.

References solveSP::m_eqnIndexStartSolnPhase, solveSP::m_numSurfPhases, and solveSP::m_ptrsSurfPhase.

Referenced by solveSP::calc_t().

void updateMFKinSpecies ( doublereal *  XMolKinSp,
int  isp 
)
private

Update the mole fraction vector for a specific kinetic species vector corresponding to one InterfaceKinetics object.

Parameters
XMolKinSpMole fraction vector corresponding to a particular kinetic species for a single InterfaceKinetics Object This is a vector over all the species in all of the phases in the InterfaceKinetics object
ispID of the InterfaceKinetics Object.

Definition at line 474 of file solveSP.cpp.

References Phase::getMoleFractions(), Kinetics::kineticsSpeciesIndex(), solveSP::m_objects, Kinetics::nPhases(), and Kinetics::thermo().

void evalSurfLarge ( const doublereal *  CSolnSP)
private

Update the vector that keeps track of the largest species in each surface phase.

Parameters
CsolnSPVector of the current values of the surface concentrations in all of the surface species.

Definition at line 489 of file solveSP.cpp.

References solveSP::m_nSpeciesSurfPhase, solveSP::m_numSurfPhases, and solveSP::m_spSurfLarge.

Referenced by solveSP::solveSurfProb().

void fun_eval ( doublereal *  resid,
const doublereal *  CSolnSP,
const doublereal *  CSolnOldSP,
const bool  do_time,
const doublereal  deltaT 
)
private

Main Function evaluation.

Parameters
residoutput Vector of residuals, length = m_neq
CSolnSPVector of species concentrations, unknowns in the problem, length = m_neq
CSolnSPOldOld Vector of species concentrations, unknowns in the problem, length = m_neq
do_timeCalculate a time dependent residual
deltaTDelta time for time dependent problem.

Definition at line 516 of file solveSP.cpp.

References DATA_PTR, InterfaceKinetics::getNetProductionRates(), Kinetics::kineticsSpeciesIndex(), solveSP::m_bulkFunc, solveSP::m_bulkPhasePtrs, solveSP::m_indexKinObjSurfPhase, solveSP::m_netProductionRatesSave, solveSP::m_nSpeciesSurfPhase, solveSP::m_numBulkPhasesSS, solveSP::m_numEqn1, solveSP::m_numSurfPhases, solveSP::m_numTotSurfSpecies, solveSP::m_objects, solveSP::m_ptrsSurfPhase, solveSP::m_spSurfLarge, Kinetics::surfacePhaseIndex(), and solveSP::updateState().

Referenced by solveSP::resjac_eval(), and solveSP::solveSurfProb().

void resjac_eval ( std::vector< doublereal * > &  JacCol,
doublereal *  resid,
doublereal *  CSolnSP,
const doublereal *  CSolnSPOld,
const bool  do_time,
const doublereal  deltaT 
)
private

Main routine that calculates the current residual and Jacobian.

Parameters
JacColVector of pointers to the tops of columns of the Jacobian to be evaluated.
residoutput Vector of residuals, length = m_neq
CSolnSPVector of species concentrations, unknowns in the problem, length = m_neq. These are tweaked in order to derive the columns of the jacobian.
CSolnSPOldOld Vector of species concentrations, unknowns in the problem, length = m_neq
do_timeCalculate a time dependent residual
deltaTDelta time for time dependent problem.

Definition at line 637 of file solveSP.cpp.

References DATA_PTR, solveSP::fun_eval(), solveSP::m_bulkFunc, solveSP::m_bulkPhasePtrs, solveSP::m_neq, solveSP::m_nSpeciesSurfPhase, solveSP::m_numBulkPhasesSS, solveSP::m_numBulkSpecies, solveSP::m_numEqn2, solveSP::m_numSurfPhases, solveSP::m_ptrsSurfPhase, and ckr::max().

Referenced by solveSP::solveSurfProb().

Member Data Documentation

ImplicitSurfChem* m_SurfChemPtr
private

Pointer to the manager of the implicit surface chemistry problem.

This object actually calls the current object. Thus, we are providing a loop-back functionality here.

Definition at line 405 of file solveSP.h.

std::vector<InterfaceKinetics*>& m_objects
private

Vector of interface kinetics objects.

Each of these is associated with one and only one surface phase.

Definition at line 411 of file solveSP.h.

Referenced by solveSP::calc_t(), solveSP::fun_eval(), solveSP::printFinal(), solveSP::printIteration(), solveSP::solveSP(), and solveSP::updateMFKinSpecies().

size_t m_neq
private

Total number of equations to solve in the implicit problem.

Note, this can be zero, and frequently is

Definition at line 417 of file solveSP.h.

Referenced by solveSP::calcWeights(), solveSP::resjac_eval(), solveSP::solveSP(), and solveSP::solveSurfProb().

int m_bulkFunc
private

This variable determines how the bulk phases are to be handled.

= BULK_ETCH (default) The concentrations of the bulk phases are considered constant, just as the gas phase is. They are not part of the solution vector. = BULK_DEPOSITION = We solve here for the composition of the bulk phases by calculating a growth rate. The equations for the species in the bulk phases are unknowns in this calculation.

Definition at line 430 of file solveSP.h.

Referenced by solveSP::calcWeights(), solveSP::fun_eval(), solveSP::print_header(), and solveSP::resjac_eval().

size_t m_numSurfPhases
private

Number of surface phases in the surface problem.

This number is equal to the number of InterfaceKinetics objects in the problem. (until further noted)

Definition at line 437 of file solveSP.h.

Referenced by solveSP::calc_t(), solveSP::calcWeights(), solveSP::evalSurfLarge(), solveSP::fun_eval(), solveSP::resjac_eval(), solveSP::solveSP(), solveSP::solveSurfProb(), solveSP::updateMFSolnSP(), and solveSP::updateState().

size_t m_numTotSurfSpecies
private

Total number of surface species in all surface phases.

This is also the number of equations to solve for m_mode=0 system It's equal to the sum of the number of species in each of the m_numSurfPhases.

Definition at line 445 of file solveSP.h.

Referenced by solveSP::fun_eval(), and solveSP::solveSP().

std::vector<size_t> m_indexKinObjSurfPhase
private

Mapping between the surface phases and the InterfaceKinetics objects.

Currently this is defined to be a 1-1 mapping (and probably assumed in some places) m_surfKinObjID[i] = i

Definition at line 453 of file solveSP.h.

Referenced by solveSP::fun_eval(), and solveSP::solveSP().

std::vector<size_t> m_nSpeciesSurfPhase
private

Vector of length number of surface phases containing the number of surface species in each phase.

Length is equal to the number of surface phases, m_numSurfPhases

Definition at line 460 of file solveSP.h.

Referenced by solveSP::calc_t(), solveSP::calcWeights(), solveSP::evalSurfLarge(), solveSP::fun_eval(), solveSP::resjac_eval(), solveSP::solveSP(), solveSP::solveSurfProb(), and solveSP::updateState().

std::vector<SurfPhase*> m_ptrsSurfPhase
private

Vector of surface phase pointers.

This is created during the constructor Length is equal to the number of surface phases, m_numSurfPhases

Definition at line 467 of file solveSP.h.

Referenced by solveSP::calcWeights(), solveSP::fun_eval(), solveSP::resjac_eval(), solveSP::solveSP(), solveSP::solveSurfProb(), solveSP::updateMFSolnSP(), and solveSP::updateState().

std::vector<size_t> m_eqnIndexStartSolnPhase
private

Index of the start of the unknowns for each solution phase.

  i_eqn = m_eqnIndexStartPhase[isp]

isp is the phase id in the list of phases solved by the surface problem.

i_eqn is the equation number of the first unknown in the solution vector corresponding to isp'th phase.

Definition at line 479 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::updateMFSolnSP().

std::vector<size_t> m_kinObjPhaseIDSurfPhase
private

Phase ID in the InterfaceKinetics object of the surface phase.

For each surface phase, this lists the PhaseId of the surface phase in the corresponding InterfaceKinetics object

Length is equal to m_numSurfPhases

Definition at line 488 of file solveSP.h.

Referenced by solveSP::solveSP().

size_t m_numBulkPhasesSS
private

Total number of volumetric condensed phases included in the steady state problem handled by this routine.

This is equal to or less than the total number of volumetric phases in all of the InterfaceKinetics objects. We usually do not include bulk phases. Bulk phases are only included in the calculation when their domain isn't included in the underlying continuum model conservation equation system.

This is equal to 0, for the time being

Definition at line 501 of file solveSP.h.

Referenced by solveSP::calcWeights(), solveSP::fun_eval(), solveSP::resjac_eval(), and solveSP::solveSP().

std::vector<size_t> m_numBulkSpecies
private

Vector of number of species in the m_numBulkPhases phases.

Length is number of bulk phases

Definition at line 507 of file solveSP.h.

Referenced by solveSP::calcWeights(), and solveSP::resjac_eval().

size_t m_numTotBulkSpeciesSS
private

Total number of species in all bulk phases.

This is also the number of bulk equations to solve when bulk equation solving is turned on.

Definition at line 518 of file solveSP.h.

Referenced by solveSP::solveSP().

std::vector<ThermoPhase*> m_bulkPhasePtrs
private

Vector of bulk phase pointers, length is equal to m_numBulkPhases.

Definition at line 524 of file solveSP.h.

Referenced by solveSP::calcWeights(), solveSP::fun_eval(), and solveSP::resjac_eval().

std::vector<size_t> m_kinSpecIndex
private

Index between the equation index and the position in the kinetic species array for the appropriate kinetics operator.

Length = m_neq.

ksp = m_kinSpecIndex[ieq] ksp is the kinetic species index for the ieq'th equation.

Definition at line 535 of file solveSP.h.

Referenced by solveSP::printFinal(), solveSP::printIteration(), solveSP::solveSP(), and solveSP::solveSurfProb().

std::vector<size_t> m_kinObjIndex
private

Index between the equation index and the index of the InterfaceKinetics object.

Length m_neq

Definition at line 542 of file solveSP.h.

Referenced by solveSP::printFinal(), solveSP::printIteration(), and solveSP::solveSP().

std::vector<size_t> m_spSurfLarge
private

Vector containing the indices of the largest species in each surface phase.

      k = m_spSurfLarge[i]

where k is the local species index, i.e., it varies from 0 num species in phase-1 i is the surface phase index in the problem

length is equal to m_numSurfPhases

Definition at line 555 of file solveSP.h.

Referenced by solveSP::evalSurfLarge(), solveSP::fun_eval(), and solveSP::solveSP().

doublereal m_atol
private

m_atol is the absolute tolerance in real units.

units are (kmol/m2)

Definition at line 561 of file solveSP.h.

doublereal m_rtol
private

m_rtol is the relative error tolerance.

Definition at line 564 of file solveSP.h.

doublereal m_maxstep
private

maximum value of the time step

units = seconds

Definition at line 570 of file solveSP.h.

size_t m_maxTotSpecies
private

Maximum number of species in any single kinetics operator -> also maxed wrt the total # of solution species.

Definition at line 574 of file solveSP.h.

Referenced by solveSP::solveSP().

vector_fp m_netProductionRatesSave
private

Temporary vector with length equal to max m_maxTotSpecies.

Definition at line 577 of file solveSP.h.

Referenced by solveSP::fun_eval(), solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_numEqn1
private

Temporary vector with length equal to max m_maxTotSpecies.

Definition at line 580 of file solveSP.h.

Referenced by solveSP::calc_t(), solveSP::fun_eval(), solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_numEqn2
private

Temporary vector with length equal to max m_maxTotSpecies.

Definition at line 583 of file solveSP.h.

Referenced by solveSP::resjac_eval(), and solveSP::solveSP().

vector_fp m_CSolnSave
private

Temporary vector with length equal to max m_maxTotSpecies.

Definition at line 586 of file solveSP.h.

Referenced by solveSP::solveSP().

vector_fp m_CSolnSP
private

Solution vector.

length MAX(1, m_neq)

Definition at line 592 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_CSolnSPInit
private

Saved initial solution vector.

length MAX(1, m_neq)

Definition at line 598 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_CSolnSPOld
private

Saved solution vector at the old time step.

length MAX(1, m_neq)

Definition at line 604 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_wtResid
private

Weights for the residual norm calculation.

length MAX(1, m_neq)

Definition at line 610 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_wtSpecies
private

Weights for the species concentrations norm calculation.

length MAX(1, m_neq)

Definition at line 616 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_resid
private

Residual for the surface problem.

The residual vector of length "dim" that, that has the value of "sdot" for surface species. The residuals for the bulk species are a function of the sdots for all species in the bulk phase. The last residual of each phase enforces {Sum(fractions) = 1}. After linear solve (dgetrf_ & dgetrs_), resid holds the update vector.

length MAX(1, m_neq)

Definition at line 629 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_fp m_XMolKinSpecies
private

Vector of mole fractions.

*length m_maxTotSpecies

Definition at line 635 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

vector_int m_ipiv
private

pivots

length MAX(1, m_neq)

Definition at line 641 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

std::vector<doublereal*> m_JacCol
private

Vector of pointers to the top of the columns of the jacobians.

The "dim" by "dim" computed Jacobian matrix for the local Newton's method.

Definition at line 649 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().

Array2D m_Jac
private

Jacobian.

m_neq by m_neq computed Jacobian matrix for the local Newton's method.

Definition at line 656 of file solveSP.h.

Referenced by solveSP::solveSP(), and solveSP::solveSurfProb().


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