Cantera  3.2.0a5
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 if (nm == "temperature") {
47 return 0;
48 }
49 if (nm == "volume") {
50 return 1;
51 }
52 try {
53 return speciesIndex(nm) + m_sidx;
54 } catch (const CanteraError&) {
55 throw CanteraError("IdealGasMoleReactor::componentIndex",
56 "Component '{}' not found", nm);
57 }
58}
59
61{
62 if (k == 0) {
63 return "temperature";
64 } else {
66 }
67}
68
70{
71 if (m_thermo->type() != "ideal-gas") {
72 throw CanteraError("IdealGasMoleReactor::initialize",
73 "Incompatible phase type '{}' provided", m_thermo->type());
74 }
76 m_uk.resize(m_nsp, 0.0);
77}
78
79double IdealGasMoleReactor::upperBound(size_t k) const {
80 if (k == 0) {
81 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
82 return 1.5 * m_thermo->maxTemp();
83 } else {
85 }
86}
87
88double IdealGasMoleReactor::lowerBound(size_t k) const {
89 if (k == 0) {
90 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
91 return 0.5 * m_thermo->minTemp();
92 } else {
94 }
95}
96
98{
99 if (nSurfs() != 0) {
100 throw CanteraError("IdealGasMoleReactor::steadyConstraints",
101 "Steady state solver cannot currently be used with IdealGasMoleReactor"
102 " when reactor surfaces are present.\n"
103 "See https://github.com/Cantera/enhancements/issues/234");
104 }
105 if (energyEnabled()) {
106 return {1}; // volume
107 } else {
108 return {0, 1}; // temperature and volume
109 }
110}
111
113{
114 // the components of y are: [0] the temperature, [1] the volume, [2...K+1) are the
115 // moles of each species, and [K+1...] are the moles of surface
116 // species on each wall.
118 m_vol = y[1];
119 // set state
120 m_thermo->setMolesNoTruncate(y + m_sidx);
121 m_thermo->setState_TD(y[0], m_mass / m_vol);
122 updateConnected(true);
124}
125
126void IdealGasMoleReactor::eval(double time, double* LHS, double* RHS)
127{
128 double& mcvdTdt = RHS[0]; // m * c_v * dT/dt
129 double* dndt = RHS + m_sidx; // kmol per s
130
131 evalWalls(time);
132
133 m_thermo->restoreState(m_state);
134
135 m_thermo->getPartialMolarIntEnergies(&m_uk[0]);
136 const vector<double>& imw = m_thermo->inverseMolecularWeights();
137
138 if (m_chem) {
139 m_kin->getNetProductionRates(&m_wdot[0]); // "omega dot"
140 }
141
142 // evaluate surfaces
143 evalSurfaces(LHS + m_nsp + m_sidx, RHS + m_nsp + m_sidx, m_sdot.data());
144
145 // external heat transfer and compression work
146 mcvdTdt += - m_pressure * m_vdot + m_Qdot;
147
148 for (size_t n = 0; n < m_nsp; n++) {
149 // heat release from gas phase and surface reactions
150 mcvdTdt -= m_wdot[n] * m_uk[n] * m_vol;
151 mcvdTdt -= m_sdot[n] * m_uk[n];
152 // production in gas phase and from surfaces
153 dndt[n] = (m_wdot[n] * m_vol + m_sdot[n]);
154 }
155
156 // add terms for outlets
157 for (auto outlet : m_outlet) {
158 for (size_t n = 0; n < m_nsp; n++) {
159 // flow of species into system and dilution by other species
160 dndt[n] -= outlet->outletSpeciesMassFlowRate(n) * imw[n];
161 }
162 double mdot = outlet->massFlowRate();
163 mcvdTdt -= mdot * m_pressure * m_vol / m_mass; // flow work
164 }
165
166 // add terms for inlets
167 for (auto inlet : m_inlet) {
168 double mdot = inlet->massFlowRate();
169 mcvdTdt += inlet->enthalpy_mass() * mdot;
170 for (size_t n = 0; n < m_nsp; n++) {
171 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
172 // flow of species into system and dilution by other species
173 dndt[n] += inlet->outletSpeciesMassFlowRate(n) * imw[n];
174 mcvdTdt -= m_uk[n] * imw[n] * mdot_spec;
175 }
176 }
177
178 RHS[1] = m_vdot;
179 if (m_energy) {
180 LHS[0] = m_mass * m_thermo->cv_mass();
181 } else {
182 RHS[0] = 0;
183 }
184}
185
186Eigen::SparseMatrix<double> IdealGasMoleReactor::jacobian()
187{
188 if (m_nv == 0) {
189 throw CanteraError("IdealGasMoleReactor::jacobian",
190 "Reactor must be initialized first.");
191 }
192 // clear former jacobian elements
193 m_jac_trips.clear();
194 // dnk_dnj represents d(dot(n_k)) / d (n_j) but is first assigned as
195 // d (dot(omega)) / d c_j, it is later transformed appropriately.
196 Eigen::SparseMatrix<double> dnk_dnj = m_kin->netProductionRates_ddCi();
197 // species size that accounts for surface species
198 size_t ssize = m_nv - m_sidx;
199 // map derivatives from the surface chemistry jacobian
200 // to the reactor jacobian
201 if (!m_surfaces.empty()) {
202 vector<Eigen::Triplet<double>> species_trips;
203 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
204 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
205 species_trips.emplace_back(static_cast<int>(it.row()),
206 static_cast<int>(it.col()), it.value());
207 }
208 }
209 addSurfaceJacobian(species_trips);
210 dnk_dnj.resize(ssize, ssize);
211 dnk_dnj.setFromTriplets(species_trips.begin(), species_trips.end());
212 }
213 // add species to species derivatives elements to the jacobian
214 // calculate ROP derivatives, excluding the terms -n_i / (V * N) dc_i/dn_j
215 // as it substantially reduces matrix sparsity
216 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
217 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
218 m_jac_trips.emplace_back(static_cast<int>(it.row() + m_sidx),
219 static_cast<int>(it.col() + m_sidx), it.value());
220 }
221 }
222 // Temperature Derivatives
223 if (m_energy) {
224 // getting perturbed state for finite difference
225 double deltaTemp = m_thermo->temperature()
226 * std::sqrt(std::numeric_limits<double>::epsilon());
227 // finite difference temperature derivatives
228 vector<double> lhsPerturbed(m_nv, 1.0), lhsCurrent(m_nv, 1.0);
229 vector<double> rhsPerturbed(m_nv, 0.0), rhsCurrent(m_nv, 0.0);
230 vector<double> yCurrent(m_nv);
231 getState(yCurrent.data());
232 vector<double> yPerturbed = yCurrent;
233 // perturb temperature
234 yPerturbed[0] += deltaTemp;
235 // getting perturbed state
236 updateState(yPerturbed.data());
237 double time = (m_net != nullptr) ? m_net->time() : 0.0;
238 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
239 // reset and get original state
240 updateState(yCurrent.data());
241 eval(time, lhsCurrent.data(), rhsCurrent.data());
242 // d ydot_j/dT
243 for (size_t j = 0; j < m_nv; j++) {
244 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
245 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
246 m_jac_trips.emplace_back(static_cast<int>(j), 0,
247 (ydotPerturbed - ydotCurrent) / deltaTemp);
248 }
249 // d T_dot/dnj
250 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(ssize);
251 Eigen::VectorXd internal_energy = Eigen::VectorXd::Zero(ssize);
252 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(ssize);
253 // getting species data
254 m_thermo->getPartialMolarIntEnergies(internal_energy.data());
255 m_kin->getNetProductionRates(netProductionRates.data());
256 m_thermo->getPartialMolarCp(specificHeat.data());
257 // convert Cp to Cv for ideal gas as Cp - Cv = R
258 for (size_t i = 0; i < m_nsp; i++) {
259 specificHeat[i] -= GasConstant;
260 netProductionRates[i] *= m_vol;
261 }
262 // scale net production rates by volume to get molar rate
263 double qdot = internal_energy.dot(netProductionRates);
264 // find the sum of n_i and cp_i
265 double NCv = 0.0;
266 double* moles = yCurrent.data() + m_sidx;
267 for (size_t i = 0; i < ssize; i++) {
268 NCv += moles[i] * specificHeat[i];
269 }
270 // make denominator beforehand
271 double denom = 1 / (NCv * NCv);
272 Eigen::VectorXd uk_dnkdnj_sums = dnk_dnj.transpose() * internal_energy;
273 // add derivatives to jacobian
274 for (size_t j = 0; j < ssize; j++) {
275 m_jac_trips.emplace_back(0, static_cast<int>(j + m_sidx),
276 (specificHeat[j] * qdot - NCv * uk_dnkdnj_sums[j]) * denom);
277 }
278 }
279 // convert triplets to sparse matrix
280 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
281 jac.setFromTriplets(m_jac_trips.begin(), m_jac_trips.end());
282 return jac;
283}
284
285}
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:489
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:323
void setState_TD(double t, double rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition Phase.cpp:435
double temperature() const
Temperature (K).
Definition Phase.h:629
const vector< double > & inverseMolecularWeights() const
Return a const reference to the internal vector of molecular weights.
Definition Phase.cpp:465
virtual double density() const
Density (kg/m^3).
Definition Phase.h:654
virtual void setMolesNoTruncate(const double *const N)
Set the state of the object with moles in [kmol].
Definition Phase.cpp:593
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:286
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:297
virtual void evalWalls(double t)
Evaluate terms related to Walls.
Definition Reactor.cpp:285
bool energyEnabled() const override
Returns true if solution of the energy equation is enabled.
Definition Reactor.h:88
double m_Qdot
net heat transfer into the reactor, through walls [W]
Definition Reactor.h:290
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
Definition Reactor.h:307
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
Definition Reactor.h:295
double m_vdot
net rate of volume change from moving walls [m^3/s]
Definition Reactor.h:288
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:452
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
Definition Reactor.cpp:190
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:641
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:120
Namespace for the Cantera kernel.
Definition AnyMap.cpp:595
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...