Cantera  2.5.1
ImplicitSurfChem.cpp
Go to the documentation of this file.
1 /**
2  * @file ImplicitSurfChem.cpp
3  * Definitions for the implicit integration of surface site density equations
4  * (see \ref kineticsmgr and class
5  * \link Cantera::ImplicitSurfChem ImplicitSurfChem\endlink).
6  */
7 
8 // This file is part of Cantera. See License.txt in the top-level directory or
9 // at https://cantera.org/license.txt for license and copyright information.
10 
14 
15 using namespace std;
16 
17 namespace Cantera
18 {
19 
20 ImplicitSurfChem::ImplicitSurfChem(
21  vector<InterfaceKinetics*> k, double rtol, double atol,
22  double maxStepSize, size_t maxSteps,
23  size_t maxErrTestFails) :
24  m_nv(0),
25  m_numTotalBulkSpecies(0),
26  m_numTotalSpecies(0),
27  m_atol(atol),
28  m_rtol(rtol),
29  m_maxstep(maxStepSize),
30  m_nmax(maxSteps),
31  m_maxErrTestFails(maxErrTestFails),
32  m_mediumSpeciesStart(-1),
33  m_bulkSpeciesStart(-1),
34  m_surfSpeciesStart(-1),
35  m_commonTempPressForPhases(true),
36  m_ioFlag(0)
37 {
38  size_t ntmax = 0;
39  size_t kinSpIndex = 0;
40  // Loop over the number of surface kinetics objects
41  for (size_t n = 0; n < k.size(); n++) {
42  InterfaceKinetics* kinPtr = k[n];
43  m_vecKinPtrs.push_back(kinPtr);
44  size_t ns = k[n]->surfacePhaseIndex();
45  if (ns == npos) {
46  throw CanteraError("ImplicitSurfChem::ImplicitSurfChem",
47  "kinetics manager contains no surface phase");
48  }
49  m_surfindex.push_back(ns);
50  m_surf.push_back((SurfPhase*)&k[n]->thermo(ns));
51  size_t nsp = m_surf.back()->nSpecies();
52  m_nsp.push_back(nsp);
53  m_nv += m_nsp.back();
54  size_t nt = k[n]->nTotalSpecies();
55  ntmax = std::max(nt, ntmax);
56  m_specStartIndex.push_back(kinSpIndex);
57  kinSpIndex += nsp;
58  size_t nPhases = kinPtr->nPhases();
59  vector_int pLocTmp(nPhases);
60  size_t imatch = npos;
61  for (size_t ip = 0; ip < nPhases; ip++) {
62  if (ip != ns) {
63  ThermoPhase* thPtr = & kinPtr->thermo(ip);
64  if ((imatch = checkMatch(m_bulkPhases, thPtr)) == npos) {
65  m_bulkPhases.push_back(thPtr);
66  nsp = thPtr->nSpecies();
67  m_numTotalBulkSpecies += nsp;
68  imatch = m_bulkPhases.size() - 1;
69  }
70  pLocTmp[ip] = int(imatch);
71  } else {
72  pLocTmp[ip] = -int(n);
73  }
74  }
75  pLocVec.push_back(pLocTmp);
76  }
77  m_numTotalSpecies = m_nv + m_numTotalBulkSpecies;
78  m_concSpecies.resize(m_numTotalSpecies, 0.0);
79  m_concSpeciesSave.resize(m_numTotalSpecies, 0.0);
80 
81  m_integ.reset(newIntegrator("CVODE"));
82 
83  // use backward differencing, with a full Jacobian computed
84  // numerically, and use a Newton linear iterator
85  m_integ->setMethod(BDF_Method);
86  m_integ->setProblemType(DENSE + NOJAC);
87  m_work.resize(ntmax);
88 }
89 
90 int ImplicitSurfChem::checkMatch(std::vector<ThermoPhase*> m_vec, ThermoPhase* thPtr)
91 {
92  int retn = -1;
93  for (int i = 0; i < (int) m_vec.size(); i++) {
94  ThermoPhase* th = m_vec[i];
95  if (th == thPtr) {
96  return i;
97  }
98  }
99  return retn;
100 }
101 
102 void ImplicitSurfChem::getState(doublereal* c)
103 {
104  size_t loc = 0;
105  for (size_t n = 0; n < m_surf.size(); n++) {
106  m_surf[n]->getCoverages(c + loc);
107  loc += m_nsp[n];
108  }
109 }
110 
112 {
113  m_maxstep = maxstep;
114  if (m_maxstep > 0) {
115  m_integ->setMaxStepSize(m_maxstep);
116  }
117 }
118 
119 void ImplicitSurfChem::setTolerances(double rtol, double atol)
120 {
121  m_rtol = rtol;
122  m_atol = atol;
123  m_integ->setTolerances(m_rtol, m_atol);
124 }
125 
126 void ImplicitSurfChem::setMaxSteps(size_t maxsteps)
127 {
128  m_nmax = maxsteps;
129  m_integ->setMaxSteps(static_cast<int>(m_nmax));
130 }
131 
132 void ImplicitSurfChem::setMaxErrTestFails(size_t maxErrTestFails)
133 {
134  m_maxErrTestFails = maxErrTestFails;
135  m_integ->setMaxErrTestFails(static_cast<int>(m_maxErrTestFails));
136 }
137 
138 void ImplicitSurfChem::initialize(doublereal t0)
139 {
140  this->setTolerances(m_rtol, m_atol);
141  this->setMaxStepSize(m_maxstep);
142  this->setMaxSteps(m_nmax);
144  m_integ->initialize(t0, *this);
145 }
146 
147 void ImplicitSurfChem::integrate(doublereal t0, doublereal t1)
148 {
149  this->initialize(t0);
150  if (fabs(t1 - t0) < m_maxstep || m_maxstep == 0) {
151  // limit max step size on this run to t1 - t0
152  m_integ->setMaxStepSize(t1 - t0);
153  }
154  m_integ->integrate(t1);
155  updateState(m_integ->solution());
156 }
157 
158 void ImplicitSurfChem::integrate0(doublereal t0, doublereal t1)
159 {
160  m_integ->integrate(t1);
161  updateState(m_integ->solution());
162 }
163 
164 void ImplicitSurfChem::updateState(doublereal* c)
165 {
166  size_t loc = 0;
167  for (size_t n = 0; n < m_surf.size(); n++) {
168  m_surf[n]->setCoverages(c + loc);
169  loc += m_nsp[n];
170  }
171 }
172 
173 void ImplicitSurfChem::eval(doublereal time, doublereal* y,
174  doublereal* ydot, doublereal* p)
175 {
176  updateState(y); // synchronize the surface state(s) with y
177  size_t loc = 0;
178  for (size_t n = 0; n < m_surf.size(); n++) {
179  double rs0 = 1.0/m_surf[n]->siteDensity();
180  m_vecKinPtrs[n]->getNetProductionRates(m_work.data());
181  size_t kstart = m_vecKinPtrs[n]->kineticsSpeciesIndex(0,m_surfindex[n]);
182  double sum = 0.0;
183  for (size_t k = 1; k < m_nsp[n]; k++) {
184  ydot[k + loc] = m_work[kstart + k] * rs0 * m_surf[n]->size(k);
185  sum -= ydot[k];
186  }
187  ydot[loc] = sum;
188  loc += m_nsp[n];
189  }
190 }
191 
193  doublereal timeScaleOverride)
194 {
195  int ifunc;
196  // set bulkFunc. We assume that the bulk concentrations are constant.
197  int bulkFunc = BULK_ETCH;
198  // time scale - time over which to integrate equations
199  doublereal time_scale = timeScaleOverride;
200  if (!m_surfSolver) {
201  m_surfSolver.reset(new solveSP(this, bulkFunc));
202  // set ifunc, which sets the algorithm.
203  ifunc = SFLUX_INITIALIZE;
204  } else {
205  ifunc = SFLUX_RESIDUAL;
206  }
207 
208  // Possibly override the ifunc value
209  if (ifuncOverride >= 0) {
210  ifunc = ifuncOverride;
211  }
212 
213  // Get the specifications for the problem from the values
214  // in the ThermoPhase objects for all phases.
215  //
216  // 1) concentrations of all species in all phases, m_concSpecies[]
217  // 2) Temperature and pressure
220  ThermoPhase& tp = ik->thermo(ik->reactionPhaseIndex());
221  doublereal TKelvin = tp.temperature();
222  doublereal PGas = tp.pressure();
223 
224  // Make sure that there is a common temperature and pressure for all
225  // ThermoPhase objects belonging to the interfacial kinetics object, if it
226  // is required by the problem statement.
228  setCommonState_TP(TKelvin, PGas);
229  }
230 
231  doublereal reltol = 1.0E-6;
232  doublereal atol = 1.0E-20;
233 
234  // Install a filter for negative concentrations. One of the few ways solveSS
235  // can fail is if concentrations on input are below zero.
236  bool rset = false;
237  for (size_t k = 0; k < m_nv; k++) {
238  if (m_concSpecies[k] < 0.0) {
239  rset = true;
240  m_concSpecies[k] = 0.0;
241  }
242  }
243  if (rset) {
245  }
246 
247  m_surfSolver->m_ioflag = m_ioFlag;
248 
249  // Save the current solution
250  m_concSpeciesSave = m_concSpecies;
251 
252  int retn = m_surfSolver->solveSurfProb(ifunc, time_scale, TKelvin, PGas,
253  reltol, atol);
254  if (retn != 1) {
255  // reset the concentrations
256  m_concSpecies = m_concSpeciesSave;
257  setConcSpecies(m_concSpeciesSave.data());
258  ifunc = SFLUX_INITIALIZE;
259  retn = m_surfSolver->solveSurfProb(ifunc, time_scale, TKelvin, PGas,
260  reltol, atol);
261  if (retn != 1) {
262  throw CanteraError("ImplicitSurfChem::solvePseudoSteadyStateProblem",
263  "solveSP return an error condition!");
264  }
265  }
266 }
267 
268 void ImplicitSurfChem::getConcSpecies(doublereal* const vecConcSpecies) const
269 {
270  size_t kstart;
271  for (size_t ip = 0; ip < m_surf.size(); ip++) {
272  ThermoPhase* TP_ptr = m_surf[ip];
273  kstart = m_specStartIndex[ip];
274  TP_ptr->getConcentrations(vecConcSpecies + kstart);
275  }
276  kstart = m_nv;
277  for (size_t ip = 0; ip < m_bulkPhases.size(); ip++) {
278  ThermoPhase* TP_ptr = m_bulkPhases[ip];
279  TP_ptr->getConcentrations(vecConcSpecies + kstart);
280  kstart += TP_ptr->nSpecies();
281  }
282 }
283 
284 void ImplicitSurfChem::setConcSpecies(const doublereal* const vecConcSpecies)
285 {
286  size_t kstart;
287  for (size_t ip = 0; ip < m_surf.size(); ip++) {
288  ThermoPhase* TP_ptr = m_surf[ip];
289  kstart = m_specStartIndex[ip];
290  TP_ptr->setConcentrations(vecConcSpecies + kstart);
291  }
292  kstart = m_nv;
293  for (size_t ip = 0; ip < m_bulkPhases.size(); ip++) {
294  ThermoPhase* TP_ptr = m_bulkPhases[ip];
295  TP_ptr->setConcentrations(vecConcSpecies + kstart);
296  kstart += TP_ptr->nSpecies();
297  }
298 }
299 
300 void ImplicitSurfChem::setCommonState_TP(doublereal TKelvin, doublereal PresPa)
301 {
302  for (size_t ip = 0; ip < m_surf.size(); ip++) {
303  ThermoPhase* TP_ptr = m_surf[ip];
304  TP_ptr->setState_TP(TKelvin, PresPa);
305  }
306  for (size_t ip = 0; ip < m_bulkPhases.size(); ip++) {
307  ThermoPhase* TP_ptr = m_bulkPhases[ip];
308  TP_ptr->setState_TP(TKelvin, PresPa);
309  }
310 }
311 
312 }
Declarations for the implicit integration of surface site density equations (see Kinetics Managers an...
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
Base class for exceptions thrown by Cantera classes.
Definition: ctexceptions.h:61
void integrate0(doublereal t0, doublereal t1)
Integrate from t0 to t1 without reinitializing the integrator.
void integrate(doublereal t0, doublereal t1)
Integrate from t0 to t1. The integrator is reinitialized first.
void updateState(doublereal *y)
Set the mixture to a state consistent with solution vector y.
std::vector< ThermoPhase * > m_bulkPhases
Vector of pointers to bulk phases.
void setConcSpecies(const doublereal *const vecConcSpecies)
Sets the concentrations within phases that are unknowns in the surface problem.
size_t m_nv
Total number of surface species in all surface phases.
doublereal m_maxstep
max step size
int m_ioFlag
Controls the amount of printing from this routine and underlying routines.
virtual void getState(doublereal *y)
Get the current state of the solution vector.
std::unique_ptr< solveSP > m_surfSolver
Pointer to the helper method, Placid, which solves the surface problem.
std::vector< size_t > m_surfindex
index of the surface phase in each InterfaceKinetics object
bool m_commonTempPressForPhases
If true, a common temperature and pressure for all surface and bulk phases associated with the surfac...
virtual void setTolerances(double rtol=1.e-7, double atol=1.e-14)
virtual void setMaxStepSize(double maxstep=0.0)
std::vector< SurfPhase * > m_surf
vector of pointers to surface phases.
virtual void setMaxSteps(size_t maxsteps=20000)
vector_fp m_concSpecies
Temporary vector - length num species in the Kinetics object.
void getConcSpecies(doublereal *const vecConcSpecies) const
virtual void initialize(doublereal t0=0.0)
size_t m_maxErrTestFails
maximum number of error test failures allowed
void solvePseudoSteadyStateProblem(int ifuncOverride=-1, doublereal timeScaleOverride=1.0)
Solve for the pseudo steady-state of the surface problem.
std::vector< size_t > m_nsp
Vector of number of species in each Surface Phase.
virtual void setMaxErrTestFails(size_t maxErrTestFails=7)
void setCommonState_TP(doublereal TKelvin, doublereal PresPa)
Sets the state variable in all thermodynamic phases (surface and surrounding bulk phases) to the inpu...
virtual void eval(doublereal t, doublereal *y, doublereal *ydot, doublereal *p)
Evaluate the value of ydot[k] at the current conditions.
std::unique_ptr< Integrator > m_integ
Pointer to the CVODE integrator.
size_t m_nmax
maximum number of steps allowed
std::vector< InterfaceKinetics * > m_vecKinPtrs
vector of pointers to InterfaceKinetics objects
A kinetics manager for heterogeneous reaction mechanisms.
size_t reactionPhaseIndex() const
Phase where the reactions occur.
Definition: Kinetics.h:214
thermo_t & thermo(size_t n=0)
This method returns a reference to the nth ThermoPhase object defined in this kinetics mechanism.
Definition: Kinetics.h:227
size_t nPhases() const
The number of phases participating in the reaction mechanism.
Definition: Kinetics.h:167
void getConcentrations(double *const c) const
Get the species concentrations (kmol/m^3).
Definition: Phase.cpp:625
size_t nSpecies() const
Returns the number of species in the phase.
Definition: Phase.h:285
virtual void setConcentrations(const double *const conc)
Set the concentrations to the specified values within the phase.
Definition: Phase.cpp:630
doublereal temperature() const
Temperature (K).
Definition: Phase.h:667
virtual double pressure() const
Return the thermodynamic pressure (Pa).
Definition: Phase.h:679
A simple thermodynamic model for a surface phase, assuming an ideal solution model.
Definition: SurfPhase.h:143
Base class for a phase with thermodynamic properties.
Definition: ThermoPhase.h:102
virtual void setState_TP(doublereal t, doublereal p)
Set the temperature (K) and pressure (Pa)
Method to solve a pseudo steady state surface problem.
Definition: solveSP.h:139
const size_t npos
index returned by functions to indicate "no position"
Definition: ct_defs.h:188
std::vector< int > vector_int
Vector of ints.
Definition: ct_defs.h:182
const int BULK_ETCH
Etching of a bulk phase is to be expected.
Definition: solveSP.h:56
const int SFLUX_RESIDUAL
Need to solve the surface problem in order to calculate the surface fluxes of gas-phase species.
Definition: solveSP.h:28
const int SFLUX_INITIALIZE
This assumes that the initial guess supplied to the routine is far from the correct one.
Definition: solveSP.h:21
Namespace for the Cantera kernel.
Definition: AnyMap.cpp:264
@ BDF_Method
Backward Differentiation.
Definition: Integrator.h:33
Header file for implicit surface problem solver (see Chemical Kinetics and class solveSP).