Cantera  3.2.0a4
Loading...
Searching...
No Matches
IdealGasMoleReactor.cpp
Go to the documentation of this file.
1//! @file IdealGasMoleReactor.cpp A constant volume zero-dimensional
2//! reactor with moles as the state
3
4// This file is part of Cantera. See License.txt in the top-level directory or
5// at https://cantera.org/license.txt for license and copyright information.
6
15#include <limits>
16
17namespace Cantera
18{
19
21{
22 if (m_thermo == 0) {
23 throw CanteraError("IdealGasMoleReactor::getState",
24 "Error: reactor is empty.");
25 }
26 m_thermo->restoreState(m_state);
27
28 // get mass for calculations
29 m_mass = m_thermo->density() * m_vol;
30
31 // set the first component to the temperature
32 y[0] = m_thermo->temperature();
33
34 // set the second component to the volume
35 y[1] = m_vol;
36
37 // get moles of species in remaining state
38 getMoles(y + m_sidx);
39 // set the remaining components to the surface species moles on
40 // the walls
42}
43
44size_t IdealGasMoleReactor::componentIndex(const string& nm) const
45{
46 size_t k = speciesIndex(nm);
47 if (k != npos) {
48 return k + m_sidx;
49 } else if (nm == "temperature") {
50 return 0;
51 } else if (nm == "volume") {
52 return 1;
53 } else {
54 return npos;
55 }
56}
57
59{
60 if (k == 0) {
61 return "temperature";
62 } else {
64 }
65}
66
68{
69 if (m_thermo->type() != "ideal-gas") {
70 throw CanteraError("IdealGasMoleReactor::initialize",
71 "Incompatible phase type '{}' provided", m_thermo->type());
72 }
74 m_uk.resize(m_nsp, 0.0);
75}
76
77double IdealGasMoleReactor::upperBound(size_t k) const {
78 if (k == 0) {
79 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
80 return 1.5 * m_thermo->maxTemp();
81 } else {
83 }
84}
85
86double IdealGasMoleReactor::lowerBound(size_t k) const {
87 if (k == 0) {
88 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
89 return 0.5 * m_thermo->minTemp();
90 } else {
92 }
93}
94
96{
97 if (nSurfs() != 0) {
98 throw CanteraError("IdealGasMoleReactor::steadyConstraints",
99 "Steady state solver cannot currently be used with IdealGasMoleReactor"
100 " when reactor surfaces are present.\n"
101 "See https://github.com/Cantera/enhancements/issues/234");
102 }
103 if (energyEnabled()) {
104 return {1}; // volume
105 } else {
106 return {0, 1}; // temperature and volume
107 }
108}
109
111{
112 // the components of y are: [0] the temperature, [1] the volume, [2...K+1) are the
113 // moles of each species, and [K+1...] are the moles of surface
114 // species on each wall.
116 m_vol = y[1];
117 // set state
118 m_thermo->setMolesNoTruncate(y + m_sidx);
119 m_thermo->setState_TD(y[0], m_mass / m_vol);
120 updateConnected(true);
122}
123
124void IdealGasMoleReactor::eval(double time, double* LHS, double* RHS)
125{
126 double& mcvdTdt = RHS[0]; // m * c_v * dT/dt
127 double* dndt = RHS + m_sidx; // kmol per s
128
129 evalWalls(time);
130
131 m_thermo->restoreState(m_state);
132
133 m_thermo->getPartialMolarIntEnergies(&m_uk[0]);
134 const vector<double>& imw = m_thermo->inverseMolecularWeights();
135
136 if (m_chem) {
137 m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
138 }
139
140 // evaluate surfaces
141 evalSurfaces(LHS + m_nsp + m_sidx, RHS + m_nsp + m_sidx, m_sdot.data());
142
143 // external heat transfer and compression work
144 mcvdTdt += - m_pressure * m_vdot + m_Qdot;
145
146 for (size_t n = 0; n < m_nsp; n++) {
147 // heat release from gas phase and surface reactions
148 mcvdTdt -= m_wdot[n] * m_uk[n] * m_vol;
149 mcvdTdt -= m_sdot[n] * m_uk[n];
150 // production in gas phase and from surfaces
151 dndt[n] = (m_wdot[n] * m_vol + m_sdot[n]);
152 }
153
154 // add terms for outlets
155 for (auto outlet : m_outlet) {
156 for (size_t n = 0; n < m_nsp; n++) {
157 // flow of species into system and dilution by other species
158 dndt[n] -= outlet->outletSpeciesMassFlowRate(n) * imw[n];
159 }
160 double mdot = outlet->massFlowRate();
161 mcvdTdt -= mdot * m_pressure * m_vol / m_mass; // flow work
162 }
163
164 // add terms for inlets
165 for (auto inlet : m_inlet) {
166 double mdot = inlet->massFlowRate();
167 mcvdTdt += inlet->enthalpy_mass() * mdot;
168 for (size_t n = 0; n < m_nsp; n++) {
169 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
170 // flow of species into system and dilution by other species
171 dndt[n] += inlet->outletSpeciesMassFlowRate(n) * imw[n];
172 mcvdTdt -= m_uk[n] * imw[n] * mdot_spec;
173 }
174 }
175
176 RHS[1] = m_vdot;
177 if (m_energy) {
178 LHS[0] = m_mass * m_thermo->cv_mass();
179 } else {
180 RHS[0] = 0;
181 }
182}
183
184Eigen::SparseMatrix<double> IdealGasMoleReactor::jacobian()
185{
186 if (m_nv == 0) {
187 throw CanteraError("IdealGasMoleReactor::jacobian",
188 "Reactor must be initialized first.");
189 }
190 // clear former jacobian elements
191 m_jac_trips.clear();
192 // dnk_dnj represents d(dot(n_k)) / d (n_j) but is first assigned as
193 // d (dot(omega)) / d c_j, it is later transformed appropriately.
194 Eigen::SparseMatrix<double> dnk_dnj = m_kin->netProductionRates_ddCi();
195 // species size that accounts for surface species
196 size_t ssize = m_nv - m_sidx;
197 // map derivatives from the surface chemistry jacobian
198 // to the reactor jacobian
199 if (!m_surfaces.empty()) {
200 vector<Eigen::Triplet<double>> species_trips;
201 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
202 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
203 species_trips.emplace_back(static_cast<int>(it.row()),
204 static_cast<int>(it.col()), it.value());
205 }
206 }
207 addSurfaceJacobian(species_trips);
208 dnk_dnj.resize(ssize, ssize);
209 dnk_dnj.setFromTriplets(species_trips.begin(), species_trips.end());
210 }
211 // add species to species derivatives elements to the jacobian
212 // calculate ROP derivatives, excluding the terms -n_i / (V * N) dc_i/dn_j
213 // as it substantially reduces matrix sparsity
214 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
215 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
216 m_jac_trips.emplace_back(static_cast<int>(it.row() + m_sidx),
217 static_cast<int>(it.col() + m_sidx), it.value());
218 }
219 }
220 // Temperature Derivatives
221 if (m_energy) {
222 // getting perturbed state for finite difference
223 double deltaTemp = m_thermo->temperature()
224 * std::sqrt(std::numeric_limits<double>::epsilon());
225 // finite difference temperature derivatives
226 vector<double> lhsPerturbed(m_nv, 1.0), lhsCurrent(m_nv, 1.0);
227 vector<double> rhsPerturbed(m_nv, 0.0), rhsCurrent(m_nv, 0.0);
228 vector<double> yCurrent(m_nv);
229 getState(yCurrent.data());
230 vector<double> yPerturbed = yCurrent;
231 // perturb temperature
232 yPerturbed[0] += deltaTemp;
233 // getting perturbed state
234 updateState(yPerturbed.data());
235 double time = (m_net != nullptr) ? m_net->time() : 0.0;
236 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
237 // reset and get original state
238 updateState(yCurrent.data());
239 eval(time, lhsCurrent.data(), rhsCurrent.data());
240 // d ydot_j/dT
241 for (size_t j = 0; j < m_nv; j++) {
242 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
243 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
244 m_jac_trips.emplace_back(static_cast<int>(j), 0,
245 (ydotPerturbed - ydotCurrent) / deltaTemp);
246 }
247 // d T_dot/dnj
248 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(ssize);
249 Eigen::VectorXd internal_energy = Eigen::VectorXd::Zero(ssize);
250 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(ssize);
251 // getting species data
252 m_thermo->getPartialMolarIntEnergies(internal_energy.data());
253 m_kin->getNetProductionRates(netProductionRates.data());
254 m_thermo->getPartialMolarCp(specificHeat.data());
255 // convert Cp to Cv for ideal gas as Cp - Cv = R
256 for (size_t i = 0; i < m_nsp; i++) {
257 specificHeat[i] -= GasConstant;
258 netProductionRates[i] *= m_vol;
259 }
260 // scale net production rates by volume to get molar rate
261 double qdot = internal_energy.dot(netProductionRates);
262 // find the sum of n_i and cp_i
263 double NCv = 0.0;
264 double* moles = yCurrent.data() + m_sidx;
265 for (size_t i = 0; i < ssize; i++) {
266 NCv += moles[i] * specificHeat[i];
267 }
268 // make denominator beforehand
269 double denom = 1 / (NCv * NCv);
270 Eigen::VectorXd uk_dnkdnj_sums = dnk_dnj.transpose() * internal_energy;
271 // add derivatives to jacobian
272 for (size_t j = 0; j < ssize; j++) {
273 m_jac_trips.emplace_back(0, static_cast<int>(j + m_sidx),
274 (specificHeat[j] * qdot - NCv * uk_dnkdnj_sums[j]) * denom);
275 }
276 }
277 // convert triplets to sparse matrix
278 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
279 jac.setFromTriplets(m_jac_trips.begin(), m_jac_trips.end());
280 return jac;
281}
282
283}
Base class for kinetics managers and also contains the kineticsmgr module documentation (see Kinetics...
Header file for class ReactorSurface.
Header for a simple thermodynamics model of a surface phase derived from ThermoPhase,...
Header file for class ThermoPhase, the base class for phases with thermodynamic properties,...
Base class for exceptions thrown by Cantera classes.
double outletSpeciesMassFlowRate(size_t k)
Mass flow rate (kg/s) of outlet species k.
double enthalpy_mass()
specific enthalpy
double massFlowRate()
Mass flow rate (kg/s).
Definition FlowDevice.h:36
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
vector< size_t > steadyConstraints() const override
Get the indices of equations that are algebraic constraints when solving the steady-state problem.
void eval(double t, double *LHS, double *RHS) override
Evaluate the reactor governing equations.
size_t componentIndex(const string &nm) const override
Return the index in the solution vector for this reactor of the component named nm.
Eigen::SparseMatrix< double > jacobian() override
Calculate an approximate Jacobian to accelerate preconditioned solvers.
vector< double > m_uk
Species molar internal energies.
void getState(double *y) override
Get the the current state of the reactor.
double lowerBound(size_t k) const override
Get the lower bound on the k-th component of the local state vector.
string componentName(size_t k) override
Return the name of the solution component with index i.
void updateState(double *y) override
Set the state of the reactor to correspond to the state vector y.
void initialize(double t0=0.0) override
Initialize the reactor.
virtual void getNetProductionRates(double *wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:428
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
void evalSurfaces(double *LHS, double *RHS, double *sdot) override
Evaluate terms related to surface reactions.
void getSurfaceInitialConditions(double *y) override
Get initial conditions for SurfPhase objects attached to this reactor.
void getMoles(double *y)
Get moles of the system from mass fractions stored by thermo object.
const size_t m_sidx
const value for the species start index
Definition MoleReactor.h:64
void setMassFromMoles(double *y)
Set internal mass variable based on moles given.
virtual void addSurfaceJacobian(vector< Eigen::Triplet< double > > &triplets)
For each surface in the reactor, update vector of triplets with all relevant surface jacobian derivat...
double lowerBound(size_t k) const override
Get the lower bound on the k-th component of the local state vector.
string componentName(size_t k) override
Return the name of the solution component with index i.
void initialize(double t0=0.0) override
Initialize the reactor.
void updateSurfaceState(double *y) override
Update the state of SurfPhase objects attached to this reactor.
void restoreState(const vector< double > &state)
Restore a state saved on a previous call to saveState.
Definition Phase.cpp:286
void setState_TD(double t, double rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition Phase.cpp:398
double temperature() const
Temperature (K).
Definition Phase.h:587
const vector< double > & inverseMolecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition Phase.cpp:428
virtual double density() const
Density (kg/m^3).
Definition Phase.h:612
virtual void setMolesNoTruncate(const double *const N)
Set the state of the object with moles in [kmol].
Definition Phase.cpp:556
FlowDevice & outlet(size_t n=0)
Return a reference to the n-th outlet FlowDevice connected to this reactor.
double m_pressure
Current pressure in the reactor [Pa].
ReactorNet * m_net
The ReactorNet that this reactor is part of.
virtual size_t nSurfs() const
Return the number of surfaces in a reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
double m_vol
Current volume of the reactor [m^3].
double m_mass
Current mass of the reactor [kg].
size_t m_nsp
Number of homogeneous species in the mixture.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition Reactor.h:293
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:304
virtual void evalWalls(double t)
Evaluate terms related to Walls.
Definition Reactor.cpp:297
bool energyEnabled() const
Returns true if solution of the energy equation is enabled.
Definition Reactor.h:95
double m_Qdot
net heat transfer into the reactor, through walls [W]
Definition Reactor.h:297
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
Definition Reactor.h:314
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
Definition Reactor.h:302
double m_vdot
net rate of volume change from moving walls [m^3/s]
Definition Reactor.h:295
virtual size_t speciesIndex(const string &nm) const
Return the index in the solution vector for this reactor of the species named nm, in either the homog...
Definition Reactor.cpp:464
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
Definition Reactor.cpp:202
virtual double minTemp(size_t k=npos) const
Minimum temperature for which the thermodynamic data for the species or phase are valid.
virtual void getPartialMolarCp(double *cpbar) const
Return an array of partial molar heat capacities for the species in the mixture.
string type() const override
String indicating the thermodynamic model implemented.
virtual void getPartialMolarIntEnergies(double *ubar) const
Return an array of partial molar internal energies for the species in the mixture.
virtual double maxTemp(size_t k=npos) const
Maximum temperature for which the thermodynamic data for the species are valid.
double cv_mass() const
Specific heat at constant volume. Units: J/kg/K.
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:580
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
const size_t npos
index returned by functions to indicate "no position"
Definition ct_defs.h:180
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...