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
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 // set the first component to the temperature
26 y[0] = m_thermo->temperature();
27 // get moles of species in remaining state
28 getMoles(y.subspan(m_sidx));
29}
30
32{
33 if (m_thermo->type() != "ideal-gas" && m_thermo->type() != "plasma") {
34 throw CanteraError("IdealGasConstPressureMoleReactor::initialize",
35 "Incompatible phase type '{}' provided", m_thermo->type());
36 }
38 m_hk.resize(m_nsp, 0.0);
39}
40
42{
43 // the components of y are: [0] the temperature, [1...K+1) are the
44 // moles of each species, and [K+1...] are the moles of surface
45 // species on each wall.
46 setMassFromMoles(y.subspan(m_sidx));
47 m_thermo->setMolesNoTruncate(y.subspan(m_sidx, m_nsp));
48 m_thermo->setState_TP(y[0], m_pressure);
49 m_vol = m_mass / m_thermo->density();
51 m_TotalCp = m_mass * m_thermo->cp_mass();
52 updateConnected(false);
53}
54
55void IdealGasConstPressureMoleReactor::eval(double time, span<double> LHS,
56 span<double> RHS)
57{
58 double& mcpdTdt = RHS[0]; // m * c_p * dT/dt
59 auto dndt = RHS.subspan(m_sidx); // kmol per s
60
61 evalWalls(time);
63 auto imw = m_thermo->inverseMolecularWeights();
64
65 if (m_chem) {
66 m_kin->getNetProductionRates(m_wdot); // "omega dot"
67 }
68
69 // external heat transfer
70 mcpdTdt += m_Qdot;
71
72 if (m_energy) {
73 mcpdTdt += m_thermo->intrinsicHeating() * m_vol;
74 }
75
76 for (size_t n = 0; n < m_nsp; n++) {
77 // heat release from gas phase and surface reactions
78 mcpdTdt -= m_wdot[n] * m_hk[n] * m_vol;
79 mcpdTdt -= m_sdot[n] * m_hk[n];
80 // production in gas phase and from surfaces
81 dndt[n] = m_wdot[n] * m_vol + m_sdot[n];
82 }
83
84 // add terms for outlets
85 for (auto outlet : m_outlet) {
86 for (size_t n = 0; n < m_nsp; n++) {
87 // flow of species out of system
88 dndt[n] -= outlet->outletSpeciesMassFlowRate(n) * imw[n];
89 }
90 }
91
92 // add terms for inlets
93 for (auto inlet : m_inlet) {
94 double mdot = inlet->massFlowRate();
95 mcpdTdt += inlet->enthalpy_mass() * mdot;
96 for (size_t n = 0; n < m_nsp; n++) {
97 double mdot_spec = inlet->outletSpeciesMassFlowRate(n);
98 // flow of species into system and dilution by other species
99 dndt[n] += inlet->outletSpeciesMassFlowRate(n) * imw[n];
100 mcpdTdt -= m_hk[n] * imw[n] * mdot_spec;
101 }
102 }
103
104 if (m_energy) {
105 LHS[0] = m_TotalCp;
106 } else {
107 RHS[0] = 0.0;
108 }
109}
110
112 vector<Eigen::Triplet<double>>& trips)
113{
114 // dnk_dnj represents d(dot(n_k)) / d (n_j) but is first assigned as
115 // d (dot(omega)) / d c_j, it is later transformed appropriately.
116 Eigen::SparseMatrix<double> dnk_dnj = m_kin->netProductionRates_ddCi();
117
118 // add species to species derivatives elements to the jacobian
119 // calculate ROP derivatives, excluding the terms -n_i / (V * N) dc_i/dn_j
120 // as it substantially reduces matrix sparsity
121 size_t offset = static_cast<int>(m_offset + m_sidx);
122 // double molarVol = m_thermo->molarVolume();
123 for (int k = 0; k < dnk_dnj.outerSize(); k++) {
124 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
125 trips.emplace_back(it.row() + offset, it.col() + offset, it.value());
126 }
127 }
128
129 // Temperature Derivatives
130 if (m_energy) {
131 // getting perturbed state for finite difference
132 double deltaTemp = m_thermo->temperature()
133 * std::sqrt(std::numeric_limits<double>::epsilon());
134 // get current state
135 vector<double> yCurrent(m_nv);
136 getState(yCurrent);
137 // finite difference temperature derivatives
138 vector<double> lhsPerturbed(m_nv, 1.0), lhsCurrent(m_nv, 1.0);
139 vector<double> rhsPerturbed(m_nv, 0.0), rhsCurrent(m_nv, 0.0);
140 vector<double> yPerturbed = yCurrent;
141 // perturb temperature
142 yPerturbed[0] += deltaTemp;
143 // getting perturbed state
144 updateState(yPerturbed);
145 double time = (m_net != nullptr) ? m_net->time() : 0.0;
146 eval(time, lhsPerturbed, rhsPerturbed);
147 // reset and get original state
148 updateState(yCurrent);
149 eval(time, lhsCurrent, rhsCurrent);
150 // d ydot_j/dT
151 for (size_t j = 0; j < m_nv; j++) {
152 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
153 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
154 trips.emplace_back(static_cast<int>(j + m_offset), static_cast<int>(m_offset),
155 (ydotPerturbed - ydotCurrent) / deltaTemp);
156 }
157 // d T_dot/dnj
158 // allocate vectors for whole system
159 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(m_nsp);
160 Eigen::VectorXd enthalpy = Eigen::VectorXd::Zero(m_nsp);
161 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(m_nsp);
162 // gas phase
163 m_thermo->getPartialMolarCp(asSpan(specificHeat));
164 m_thermo->getPartialMolarEnthalpies(asSpan(enthalpy));
165 m_kin->getNetProductionRates(asSpan(netProductionRates));
166 // scale production rates by the volume for gas species
167 for (size_t i = 0; i < m_nsp; i++) {
168 netProductionRates[i] *= m_vol;
169 }
170 double qdot = enthalpy.dot(netProductionRates);
171 double denom = 1 / (m_TotalCp * m_TotalCp);
172 Eigen::VectorXd hk_dnkdnj_sums = dnk_dnj.transpose() * enthalpy;
173 // Add derivatives to jac by spanning columns
174 for (size_t j = 0; j < m_nsp; j++) {
175 trips.emplace_back(m_offset, static_cast<int>(j + m_offset + m_sidx),
176 (specificHeat[j] * qdot - m_TotalCp * hk_dnkdnj_sums[j]) * denom);
177 }
178 }
179}
180
182 double& f_species, span<double> f_energy)
183{
184 f_species = 1.0 / m_vol;
185 for (size_t k = 0; k < m_nsp; k++) {
186 f_energy[k] = - m_hk[k] / m_TotalCp;
187 }
188}
189
191{
192 if (nm == "temperature") {
193 return 0;
194 }
195 try {
196 return m_thermo->speciesIndex(nm) + m_sidx;
197 } catch (const CanteraError&) {
198 throw CanteraError("IdealGasConstPressureReactor::componentIndex",
199 "Component '{}' not found", nm);
200 }
201}
202
204 if (k == 0) {
205 return "temperature";
206 } else if (k >= m_sidx && k < neq()) {
207 return m_thermo->speciesName(k - m_sidx);
208 }
209 throw IndexError("IdealGasConstPressureMoleReactor::componentName",
210 "components", k, m_nv);
211}
212
214 if (k == 0) {
215 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
216 return 1.5 * m_thermo->maxTemp();
217 } else {
218 return BigNumber; // moles of a bulk or surface species
219 }
220}
221
223 if (k == 0) {
224 //@todo: Revise pending resolution of https://github.com/Cantera/enhancements/issues/229
225 return 0.5 * m_thermo->minTemp();
226 } else {
228 }
229}
230
231}
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 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.
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.
size_t componentIndex(const string &nm) const override
Return the index in the solution vector for this reactor of the component named nm.
void getJacobianElements(vector< Eigen::Triplet< double > > &trips) override
Calculate an approximate Jacobian to accelerate preconditioned solvers.
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 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: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
double temperature() const
Temperature (K).
Definition Phase.h:585
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: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 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.
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
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
void initialize(double t0=0.0) override
Initialize the reactor.
Definition Reactor.cpp:63
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.
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.
virtual double intrinsicHeating()
Intrinsic volumetric heating rate [W/m³].
double cp_mass() const
Specific heat at constant pressure and composition [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: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...