Cantera  4.0.0a1
Loading...
Searching...
No Matches
IdealGasConstPressureMoleReactor.cpp
Go to the documentation of this file.
1//! @file IdealGasConstPressureMoleReactor.cpp A constant pressure
2//! zero-dimensional 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
11#include "cantera/zeroD/Wall.h"
16#include "cantera/numerics/eigen_dense.h"
17#include <limits>
18
19namespace Cantera
20{
21
23{
24 // get mass for calculations
25 m_mass = m_thermo->density() * m_vol;
26 // set the first component to the temperature
27 y[0] = m_thermo->temperature();
28 // get moles of species in remaining state
29 getMoles(y.subspan(m_sidx));
30}
31
33{
35 m_hk.resize(m_nsp, 0.0);
36}
37
39{
40 // the components of y are: [0] the temperature, [1...K+1) are the
41 // moles of each species, and [K+1...] are the moles of surface
42 // species on each wall.
43 setMassFromMoles(y.subspan(m_sidx));
44 m_thermo->setMolesNoTruncate(y.subspan(m_sidx, m_nsp));
45 m_thermo->setState_TP(y[0], m_pressure);
46 m_vol = m_mass / m_thermo->density();
48 m_TotalCp = m_mass * m_thermo->cp_mass();
49 updateConnected(false);
50}
51
52void IdealGasConstPressureMoleReactor::eval(double time, span<double> LHS,
53 span<double> RHS)
54{
55 double& mcpdTdt = RHS[0]; // m * c_p * dT/dt
56 auto dndt = RHS.subspan(m_sidx); // kmol per s
57
58 evalWalls(time);
60 auto imw = m_thermo->inverseMolecularWeights();
61
62 if (m_chem) {
63 m_kin->getNetProductionRates(m_wdot); // "omega dot"
64 }
65
66 // external heat transfer
67 mcpdTdt += m_Qdot;
68
69 if (m_energy) {
70 mcpdTdt += m_thermo->intrinsicHeating() * m_vol;
71 }
72
73 for (size_t n = 0; n < m_nsp; n++) {
74 // heat release from gas phase and surface reactions
75 mcpdTdt -= m_wdot[n] * m_hk[n] * m_vol;
76 mcpdTdt -= m_sdot[n] * m_hk[n];
77 // production in gas phase and from surfaces
78 dndt[n] = m_wdot[n] * m_vol + m_sdot[n];
79 }
80
81 // add terms for outlets
82 for (auto outlet : m_outlet) {
83 for (size_t n = 0; n < m_nsp; n++) {
84 // flow of species out of system
85 dndt[n] -= outlet->outletSpeciesMassFlowRate(n) * imw[n];
86 }
87 }
88
89 // add terms for inlets
90 for (auto inlet : m_inlet) {
91 double mdot = inlet->massFlowRate();
92 mcpdTdt += inlet->enthalpy_mass() * mdot;
93 for (size_t n = 0; n < m_nsp; n++) {
94 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
95 // flow of species into system and dilution by other species
96 dndt[n] += inlet->outletSpeciesMassFlowRate(n) * imw[n];
97 mcpdTdt -= m_hk[n] * imw[n] * mdot_spec;
98 }
99 }
100
101 if (m_energy) {
102 LHS[0] = m_TotalCp;
103 } else {
104 RHS[0] = 0.0;
105 }
106}
107
109{
110 // dnk_dnj represents d(dot(n_k)) / d (n_j) but is first assigned as
111 // d (dot(omega)) / d c_j, it is later transformed appropriately.
112 Eigen::SparseMatrix<double> dnk_dnj = m_kin->netProductionRates_ddCi();
113
114 // add species to species derivatives elements to the jacobian
115 // calculate ROP derivatives, excluding the terms -n_i / (V * N) dc_i/dn_j
116 // as it substantially reduces matrix sparsity
117 size_t offset = static_cast<int>(m_offset + m_sidx);
118 // double molarVol = m_thermo->molarVolume();
119 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
120 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
121 trips.emplace_back(it.row() + offset, it.col() + offset, it.value());
122 }
123 }
124
125 // Temperature Derivatives
126 if (m_energy) {
127 // getting perturbed state for finite difference
128 double deltaTemp = m_thermo->temperature()
129 * std::sqrt(std::numeric_limits<double>::epsilon());
130 // get current state
131 vector<double> yCurrent(m_nv);
132 getState(yCurrent);
133 // finite difference temperature derivatives
134 vector<double> lhsPerturbed(m_nv, 1.0), lhsCurrent(m_nv, 1.0);
135 vector<double> rhsPerturbed(m_nv, 0.0), rhsCurrent(m_nv, 0.0);
136 vector<double> yPerturbed = yCurrent;
137 // perturb temperature
138 yPerturbed[0] += deltaTemp;
139 // getting perturbed state
140 updateState(yPerturbed);
141 double time = (m_net != nullptr) ? m_net->time() : 0.0;
142 eval(time, lhsPerturbed, rhsPerturbed);
143 // reset and get original state
144 updateState(yCurrent);
145 eval(time, lhsCurrent, rhsCurrent);
146 // d ydot_j/dT
147 for (size_t j = 0; j < m_nv; j++) {
148 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
149 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
150 trips.emplace_back(static_cast<int>(j + m_offset), static_cast<int>(m_offset),
151 (ydotPerturbed - ydotCurrent) / deltaTemp);
152 }
153 // d T_dot/dnj
154 // allocate vectors for whole system
155 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(m_nsp);
156 Eigen::VectorXd enthalpy = Eigen::VectorXd::Zero(m_nsp);
157 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(m_nsp);
158 // gas phase
159 m_thermo->getPartialMolarCp(asSpan(specificHeat));
160 m_thermo->getPartialMolarEnthalpies(asSpan(enthalpy));
161 m_kin->getNetProductionRates(asSpan(netProductionRates));
162 // scale production rates by the volume for gas species
163 for (size_t i = 0; i < m_nsp; i++) {
164 netProductionRates[i] *= m_vol;
165 }
166 double qdot = enthalpy.dot(netProductionRates);
167 double denom = 1 / (m_TotalCp * m_TotalCp);
168 Eigen::VectorXd hk_dnkdnj_sums = dnk_dnj.transpose() * enthalpy;
169 // Add derivatives to jac by spanning columns
170 for (size_t j = 0; j < m_nsp; j++) {
171 trips.emplace_back(m_offset, static_cast<int>(j + m_offset + m_sidx),
172 (specificHeat[j] * qdot - m_TotalCp * hk_dnkdnj_sums[j]) * denom);
173 }
174 }
175
176 // Add cross-reactor terms due to flow devices and walls.
177 bool includeComposition = !m_jac_skip_connector_composition_dependence;
178 bool includePressureSpecies =
180
182 auto imw = m_thermo->inverseMolecularWeights();
183 for (auto outlet : m_outlet) {
184 for (size_t n = 0; n < m_nsp; n++) {
186 m_offset + m_sidx + n, n, -imw[n],
187 includeComposition, includePressureSpecies);
188 }
189 }
190
191 for (auto inlet : m_inlet) {
192 for (size_t n = 0; n < m_nsp; n++) {
194 m_offset + m_sidx + n, n, imw[n],
195 includeComposition, includePressureSpecies);
196 }
197 if (m_energy) {
199 inlet->enthalpy_mass() / m_TotalCp, includePressureSpecies);
200 // d(h_in)/d(upstream_state): temperature and (optionally) composition.
202 1.0 / m_TotalCp, includeComposition);
203 for (size_t n = 0; n < m_nsp; n++) {
205 -m_hk[n] * imw[n] / m_TotalCp,
206 includeComposition, includePressureSpecies);
207 }
208 }
209 }
210 }
211
212 if (!m_jac_skip_walls && m_energy) {
213 for (size_t i = 0; i < m_wall.size(); i++) {
214 int f = 2 * m_lr[i] - 1;
215 // Connector preconditioner terms include numerator derivatives for
216 // wall heat transfer. Derivatives of m_TotalCp are omitted here to avoid
217 // dense temperature/composition fill-in.
218 m_wall[i]->addHeatRateJacobian(trips, m_offset, f / m_TotalCp);
219 }
220 }
221}
222
224 double& f_species, span<double> f_energy)
225{
226 f_species = 1.0 / m_vol;
227 for (size_t k = 0; k < m_nsp; k++) {
228 f_energy[k] = - m_hk[k] / m_TotalCp;
229 }
230}
231
233 SparseTriplets& trips, size_t row, double coeff) const
234{
235 if (!isJacobianLocalRow(row)) {
236 // Local temperature-column entries are already captured by the finite
237 // difference temperature column in getJacobianElements. Cross-reactor
238 // temperature entries are not, so they are added here.
239 trips.emplace_back(row, m_offset, coeff);
240 }
241}
242
244 SparseTriplets& trips, size_t row, size_t k, double coeff) const
245{
246 auto mw = m_thermo->molecularWeights();
247 double Yk = m_thermo->massFraction(k);
248 // Convert flow-carried composition derivatives from dY_k/dy to dY_k/dn_j.
249 for (size_t j = 0; j < m_nsp; j++) {
250 double dYdn = -Yk * mw[j] / m_mass;
251 if (j == k) {
252 dYdn += mw[k] / m_mass;
253 }
254 trips.emplace_back(row, m_offset + m_sidx + j, coeff * dYdn);
255 }
256}
257
259 size_t row, double coeff, bool includeComposition) const
260{
261 // d(h_mass)/dT = cp_mass
262 addTemperatureJacobian(trips, row, coeff * m_thermo->cp_mass());
263 if (includeComposition) {
264 // d(h_mass)/d(n_k) = (h_k_molar - h_mass * W_k) / m_mass
265 double h_mass = m_thermo->enthalpy_mass();
266 auto mw = m_thermo->molecularWeights();
267 for (size_t k = 0; k < m_nsp; k++) {
268 trips.emplace_back(row, m_offset + m_sidx + k,
269 coeff * (m_hk[k] - h_mass * mw[k]) / m_mass);
270 }
271 }
272}
273
275{
276 if (nm == "temperature") {
277 return 0;
278 }
279 try {
280 return m_thermo->speciesIndex(nm) + m_sidx;
281 } catch (const CanteraError&) {
282 throw CanteraError("IdealGasConstPressureReactor::componentIndex",
283 "Component '{}' not found", nm);
284 }
285}
286
288 if (k == 0) {
289 return "temperature";
290 } else if (k >= m_sidx && k < neq()) {
291 return m_thermo->speciesName(k - m_sidx);
292 }
293 throw IndexError("IdealGasConstPressureMoleReactor::componentName",
294 "components", k, m_nv);
295}
296
298 if (k == 0) {
299 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
300 return 1.5 * m_thermo->maxTemp();
301 } else {
302 return BigNumber; // moles of a bulk or surface species
303 }
304}
305
307 if (k == 0) {
308 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
309 return 0.5 * m_thermo->minTemp();
310 } else {
312 }
313}
314
315}
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,...
Header file for base class WallBase.
Base class for exceptions thrown by Cantera classes.
double lowerBound(size_t k) const override
Get the lower bound on the k-th component of the local state vector.
double outletSpeciesMassFlowRate(size_t k)
Mass flow rate (kg/s) of outlet species k.
void addInletEnthalpyJacobian(SparseTriplets &trips, size_t row, double coeff, bool includeComposition=true)
Add Jacobian terms for the inlet enthalpy dependence on upstream state.
virtual void addMassFlowRateJacobian(SparseTriplets &trips, size_t row, double coeff, bool includePressureSpecies=true)
Add Jacobian terms proportional to derivatives of the mass flow rate.
double enthalpy_mass()
specific enthalpy
void addOutletSpeciesMassFlowRateJacobian(SparseTriplets &trips, size_t row, size_t k, double coeff, bool includeComposition=true, bool includePressureSpecies=true)
Add Jacobian terms proportional to derivatives of outletSpeciesMassFlowRate(k).
double massFlowRate()
Mass flow rate (kg/s).
Definition FlowDevice.h:37
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
void addTemperatureJacobian(SparseTriplets &trips, size_t row, double coeff) const override
Add terms proportional to derivatives with respect to this reactor's temperature.
void eval(double t, span< double > LHS, span< double > RHS) override
Evaluate the reactor governing equations.
void addSpeciesMassFractionJacobian(SparseTriplets &trips, size_t row, size_t k, double coeff) const override
Add terms proportional to derivatives with respect to a species mass fraction.
size_t componentIndex(const string &nm) const override
Return the index in the solution vector for this reactor of the component named nm.
void addEnthalpyJacobian(SparseTriplets &trips, size_t row, double coeff, bool includeComposition=true) const override
Add terms proportional to derivatives of this reactor's specific enthalpy.
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 getJacobianElements(SparseTriplets &trips) override
Calculate an approximate Jacobian to accelerate preconditioned solvers.
void initialize(double t0=0.0) override
Initialize the reactor.
void updateState(span< const double > y) override
Set the state of the reactor to correspond to the state vector y.
vector< double > m_hk
Species molar enthalpies.
void getState(span< double > y) override
Get the current state of the reactor.
void getJacobianScalingFactors(double &f_species, span< double > f_energy) override
Get scaling factors for the Jacobian matrix terms proportional to .
An array index is out of range.
virtual void getNetProductionRates(span< double > wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:447
void getMoles(span< double > y)
Get moles of the system from mass fractions stored by thermo object.
void setMassFromMoles(span< const double > y)
Set internal mass variable based on moles given.
span< const double > inverseMolecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition Phase.cpp:416
double massFraction(size_t k) const
Return the mass fraction of a single species.
Definition Phase.cpp:473
virtual void setMolesNoTruncate(span< const double > N)
Set the state of the object with moles in [kmol].
Definition Phase.cpp:550
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
Definition Phase.cpp:127
span< const double > molecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition Phase.cpp:411
double temperature() const
Temperature (K).
Definition Phase.h:586
string speciesName(size_t k) const
Name of the species with index k.
Definition Phase.cpp:143
virtual double density() const
Density (kg/m^3).
Definition Phase.h:611
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.
bool isJacobianLocalRow(size_t row) const
Check whether a global Jacobian row belongs to this reactor.
size_t neq()
Number of equations (state variables) for this reactor.
size_t m_nv
Number of state variables for this reactor.
FlowDevice & inlet(size_t n=0)
Return a reference to the n-th inlet FlowDevice connected to this reactor.
vector< int > m_lr
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wa...
double m_vol
Current volume of the reactor [m^3].
size_t m_offset
Offset into global ReactorNet state vector.
double m_mass
Current mass of the reactor [kg].
size_t m_nsp
Number of homogeneous species in the mixture.
virtual void updateConnected(bool updatePressure)
Update state information needed by connected reactors, flow devices, and walls.
double time()
Current value of the simulation time [s], for reactor networks that are solved in the time domain.
void evalWalls(double t) override
Evaluate terms related to Walls.
Definition Reactor.cpp:209
bool m_jac_skip_connector_pressure_composition_dependence
Omit species terms in pressure derivatives, which can add fill-in.
Definition Reactor.h:189
bool m_jac_skip_connector_composition_dependence
Omit flow composition derivatives, which can add dense connector blocks.
Definition Reactor.h:186
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition Reactor.h:166
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:170
bool m_jac_skip_flow_devices
Omit flow-device terms from the sparse Jacobian.
Definition Reactor.h:180
double m_Qdot
net heat transfer into the reactor, through walls [W]
Definition Reactor.h:169
bool m_jac_skip_walls
Omit wall terms from the sparse Jacobian.
Definition Reactor.h:183
void updateSurfaceProductionRates()
Update m_sdot to reflect current production rates of bulk phase species due to reactions on adjacent ...
Definition Reactor.cpp:306
vector< double > m_sdot
Total production rate of bulk phase species on surfaces [kmol/s].
Definition Reactor.h:174
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:78
virtual void getPartialMolarEnthalpies(span< double > hbar) const
Returns an array of partial molar enthalpies for the species in the mixture.
virtual void setState_TP(double t, double p)
Set the temperature (K) and pressure (Pa)
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(span< double > cpbar) const
Return an array of partial molar heat capacities 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.
virtual double intrinsicHeating()
Intrinsic volumetric heating rate [W/m³].
double cp_mass() const
Specific heat at constant pressure and composition [J/kg/K].
double enthalpy_mass() const
Specific enthalpy. Units: J/kg.
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:609
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
Definition eigen_dense.h:46
offset
Offsets of solution components in the 1D solution array.
Definition Flow1D.h:25
const double BigNumber
largest number to compare to inf.
Definition ct_defs.h:163
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...