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