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