Cantera  4.0.0a1
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 "cantera/numerics/eigen_dense.h"
16#include <limits>
17
18namespace Cantera
19{
20
22{
23 // get mass for calculations
24 m_mass = m_thermo->density() * m_vol;
25
26 // set the first component to the temperature
27 y[0] = m_thermo->temperature();
28
29 // set the second component to the volume
30 y[1] = m_vol;
31
32 // get moles of species in remaining state
33 getMoles(y.subspan(m_sidx));
34}
35
36size_t IdealGasMoleReactor::componentIndex(const string& nm) const
37{
38 if (nm == "temperature") {
39 return 0;
40 }
41 if (nm == "volume") {
42 return 1;
43 }
44 try {
45 return m_thermo->speciesIndex(nm) + m_sidx;
46 } catch (const CanteraError&) {
47 throw CanteraError("IdealGasMoleReactor::componentIndex",
48 "Component '{}' not found", nm);
49 }
50}
51
53{
54 if (k == 0) {
55 return "temperature";
56 } else {
58 }
59}
60
62{
63 if (m_thermo->type() != "ideal-gas" && m_thermo->type() != "plasma") {
64 throw CanteraError("IdealGasMoleReactor::initialize",
65 "Incompatible phase type '{}' provided", m_thermo->type());
66 }
68 m_uk.resize(m_nsp, 0.0);
69}
70
71double IdealGasMoleReactor::upperBound(size_t k) const {
72 if (k == 0) {
73 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
74 return 1.5 * m_thermo->maxTemp();
75 } else {
77 }
78}
79
80double IdealGasMoleReactor::lowerBound(size_t k) const {
81 if (k == 0) {
82 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
83 return 0.5 * m_thermo->minTemp();
84 } else {
86 }
87}
88
90{
93 if (energyEnabled()) {
94 return {1}; // volume
95 } else {
96 return {0, 1}; // temperature and volume
97 }
98}
99
100void IdealGasMoleReactor::updateState(span<const double> y)
101{
102 // the components of y are: [0] the temperature, [1] the volume, [2...K+1) are the
103 // moles of each species, and [K+1...] are the moles of surface
104 // species on each wall.
105 setMassFromMoles(y.subspan(m_sidx));
106 m_vol = y[1];
107 // set state
108 m_thermo->setMolesNoTruncate(y.subspan(m_sidx, m_nsp));
109 m_thermo->setState_TD(y[0], m_mass / m_vol);
111 m_TotalCv = m_mass * m_thermo->cv_mass();
112 updateConnected(true);
113}
114
115void IdealGasMoleReactor::eval(double time, span<double> LHS, span<double> RHS)
116{
117 double& mcvdTdt = RHS[0]; // m * c_v * dT/dt
118 auto dndt = RHS.subspan(m_sidx); // kmol per s
119
120 evalWalls(time);
122 auto imw = m_thermo->inverseMolecularWeights();
123
124 if (m_chem) {
125 m_kin->getNetProductionRates(m_wdot); // "omega dot"
126 }
127
128 // external heat transfer and compression work
129 mcvdTdt += - m_pressure * m_vdot + m_Qdot;
130
131 if (m_energy) {
132 mcvdTdt += m_thermo->intrinsicHeating() * m_vol;
133 }
134
135 for (size_t n = 0; n < m_nsp; n++) {
136 // heat release from gas phase and surface reactions
137 mcvdTdt -= m_wdot[n] * m_uk[n] * m_vol;
138 mcvdTdt -= m_sdot[n] * m_uk[n];
139 // production in gas phase and from surfaces
140 dndt[n] = (m_wdot[n] * m_vol + m_sdot[n]);
141 }
142
143 // add terms for outlets
144 for (auto outlet : m_outlet) {
145 for (size_t n = 0; n < m_nsp; n++) {
146 // flow of species into system and dilution by other species
147 dndt[n] -= outlet->outletSpeciesMassFlowRate(n) * imw[n];
148 }
149 double mdot = outlet->massFlowRate();
150 mcvdTdt -= mdot * m_pressure * m_vol / m_mass; // flow work
151 }
152
153 // add terms for inlets
154 for (auto inlet : m_inlet) {
155 double mdot = inlet->massFlowRate();
156 mcvdTdt += inlet->enthalpy_mass() * mdot;
157 for (size_t n = 0; n < m_nsp; n++) {
158 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
159 // flow of species into system and dilution by other species
160 dndt[n] += inlet->outletSpeciesMassFlowRate(n) * imw[n];
161 mcvdTdt -= m_uk[n] * imw[n] * mdot_spec;
162 }
163 }
164
165 RHS[1] = m_vdot;
166 if (m_energy) {
167 LHS[0] = m_TotalCv;
168 } else {
169 RHS[0] = 0;
170 }
171}
172
173void IdealGasMoleReactor::evalSteady(double t, span<double> LHS, span<double> RHS)
174{
175 eval(t, LHS, RHS);
176 if (!energyEnabled()) {
177 RHS[0] = m_thermo->temperature() - m_initialTemperature;
178 }
179 RHS[1] = m_vol - m_initialVolume;
180}
181
182void IdealGasMoleReactor::getJacobianElements(vector<Eigen::Triplet<double>>& trips)
183{
184 // dnk_dnj represents d(dot(n_k)) / d (n_j) but is first assigned as
185 // d (dot(omega)) / d c_j, it is later transformed appropriately.
186 Eigen::SparseMatrix<double> dnk_dnj = m_kin->netProductionRates_ddCi();
187
188 // add species to species derivatives elements to the jacobian
189 // calculate ROP derivatives, excluding the terms -n_i / (V * N) dc_i/dn_j
190 // as it substantially reduces matrix sparsity
191 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
192 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
193 trips.emplace_back(static_cast<int>(it.row() + m_offset + m_sidx),
194 static_cast<int>(it.col() + m_offset + m_sidx), it.value());
195 }
196 }
197
198 // Temperature Derivatives
199 if (m_energy) {
200 // getting perturbed state for finite difference
201 double deltaTemp = m_thermo->temperature()
202 * std::sqrt(std::numeric_limits<double>::epsilon());
203 // finite difference temperature derivatives
204 vector<double> lhsPerturbed(m_nv, 1.0), lhsCurrent(m_nv, 1.0);
205 vector<double> rhsPerturbed(m_nv, 0.0), rhsCurrent(m_nv, 0.0);
206 vector<double> yCurrent(m_nv);
207 getState(yCurrent);
208 vector<double> yPerturbed = yCurrent;
209 // perturb temperature
210 yPerturbed[0] += deltaTemp;
211 // getting perturbed state
212 updateState(yPerturbed);
213 double time = (m_net != nullptr) ? m_net->time() : 0.0;
214 eval(time, lhsPerturbed, rhsPerturbed);
215 // reset and get original state
216 updateState(yCurrent);
217 eval(time, lhsCurrent, rhsCurrent);
218 // d ydot_j/dT
219 for (size_t j = 0; j < m_nv; j++) {
220 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
221 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
222 trips.emplace_back(static_cast<int>(j + m_offset), m_offset,
223 (ydotPerturbed - ydotCurrent) / deltaTemp);
224 }
225 // d T_dot/dnj
226 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(m_nsp);
227 Eigen::VectorXd internal_energy = Eigen::VectorXd::Zero(m_nsp);
228 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(m_nsp);
229 // getting species data
230 m_thermo->getPartialMolarIntEnergies(asSpan(internal_energy));
231 m_kin->getNetProductionRates(asSpan(netProductionRates));
232 m_thermo->getPartialMolarCp(asSpan(specificHeat));
233 // convert Cp to Cv for ideal gas as Cp - Cv = R
234 for (size_t i = 0; i < m_nsp; i++) {
235 specificHeat[i] -= GasConstant;
236 netProductionRates[i] *= m_vol;
237 }
238 // scale net production rates by volume to get molar rate
239 double qdot = internal_energy.dot(netProductionRates);
240 double denom = 1 / (m_TotalCv * m_TotalCv);
241 Eigen::VectorXd uk_dnkdnj_sums = dnk_dnj.transpose() * internal_energy;
242 // add derivatives to jacobian
243 for (size_t j = 0; j < m_nsp; j++) {
244 trips.emplace_back(m_offset, static_cast<int>(j + m_offset + m_sidx),
245 (specificHeat[j] * qdot - m_TotalCv * uk_dnkdnj_sums[j]) * denom);
246 }
247 }
248}
249
251 double& f_species, span<double> f_energy)
252{
253 f_species = 1.0 / m_vol;
254 for (size_t k = 0; k < m_nsp; k++) {
255 f_energy[k] = - m_uk[k] / m_TotalCv;
256 }
257}
258
259}
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.
void eval(double t, span< double > LHS, span< double > RHS) override
Evaluate the reactor governing equations.
void evalSteady(double t, span< double > LHS, span< double > RHS) override
Evaluate the governing equations with modifications for the steady-state solver.
size_t componentIndex(const string &nm) const override
Return the index in the solution vector for this reactor of the component named nm.
vector< double > m_uk
Species molar internal energies.
void getJacobianElements(vector< Eigen::Triplet< double > > &trips) override
Calculate an approximate Jacobian to accelerate preconditioned solvers.
vector< size_t > initializeSteady() override
Initialize the reactor before solving a steady-state problem.
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.
double m_initialTemperature
Initial temperature [K]; used for steady-state calculations.
double m_initialVolume
Initial volume [m³]; used for steady-state calculations.
void updateState(span< const double > y) override
Set the state of the reactor to correspond to the state vector y.
void getState(span< double > y) override
Get the current state of the reactor.
double m_TotalCv
Total heat capacity ( ) [J/K].
void getJacobianScalingFactors(double &f_species, span< double > f_energy) override
Get scaling factors for the Jacobian matrix terms proportional to .
virtual void getNetProductionRates(span< double > wdot)
Species net production rates [kmol/m^3/s or kmol/m^2/s].
Definition Kinetics.cpp:447
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
void getMoles(span< 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:49
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 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:406
virtual void setMolesNoTruncate(span< const double > N)
Set the state of the object with moles in [kmol].
Definition Phase.cpp:540
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
void setState_TD(double t, double rho)
Set the internally stored temperature (K) and density (kg/m^3)
Definition Phase.cpp:375
double temperature() const
Temperature (K).
Definition Phase.h:585
virtual double density() const
Density (kg/m^3).
Definition Phase.h:610
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.
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.
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:194
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
Definition Reactor.h:147
vector< double > m_wdot
Species net molar production rates.
Definition Reactor.h:151
bool energyEnabled() const override
Returns true if solution of the energy equation is enabled.
Definition Reactor.h:79
double m_Qdot
net heat transfer into the reactor, through walls [W]
Definition Reactor.h:150
void updateSurfaceProductionRates()
Update m_sdot to reflect current production rates of bulk phase species due to reactions on adjacent ...
Definition Reactor.cpp:291
vector< double > m_sdot
Total production rate of bulk phase species on surfaces [kmol/s].
Definition Reactor.h:155
double m_vdot
net rate of volume change from moving walls [m^3/s]
Definition Reactor.h:149
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:63
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.
string type() const override
String indicating the thermodynamic model implemented.
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 and composition [J/kg/K].
virtual void getPartialMolarIntEnergies(span< double > ubar) const
Return an array of partial molar internal energies for the species in the mixture.
virtual double intrinsicHeating()
Intrinsic volumetric heating rate [W/m³].
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Definition Kinetics.cpp:609
const double GasConstant
Universal Gas Constant [J/kmol/K].
Definition ct_defs.h:123
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
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...