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