23 throw CanteraError(
"IdealGasConstPressureMoleReactor::getState",
24 "Error: reactor is empty.");
26 m_thermo->restoreState(m_state);
30 y[0] = m_thermo->temperature();
39 if (m_thermo->type() !=
"ideal-gas") {
40 throw CanteraError(
"IdealGasConstPressureMoleReactor::initialize",
41 "Incompatible phase type '{}' provided", m_thermo->type());
53 m_thermo->setMolesNoTruncate(y + m_sidx);
62 double& mcpdTdt = RHS[0];
63 double* dndt = RHS + m_sidx;
67 m_thermo->restoreState(m_state);
69 m_thermo->getPartialMolarEnthalpies(&
m_hk[0]);
70 const vector<double>& imw = m_thermo->inverseMolecularWeights();
82 for (
size_t n = 0; n <
m_nsp; n++) {
91 for (
auto outlet : m_outlet) {
92 for (
size_t n = 0; n <
m_nsp; n++) {
94 dndt[n] -=
outlet->outletSpeciesMassFlowRate(n) * imw[n];
99 for (
auto inlet : m_inlet) {
100 double mdot =
inlet->massFlowRate();
101 mcpdTdt +=
inlet->enthalpy_mass() * mdot;
102 for (
size_t n = 0; n <
m_nsp; n++) {
103 double mdot_spec =
inlet->outletSpeciesMassFlowRate(n);
105 dndt[n] +=
inlet->outletSpeciesMassFlowRate(n) * imw[n];
106 mcpdTdt -=
m_hk[n] * imw[n] * mdot_spec;
111 LHS[0] =
m_mass * m_thermo->cp_mass();
120 throw CanteraError(
"IdealGasConstPressureMoleReactor::jacobian",
121 "Reactor must be initialized first.");
127 Eigen::SparseMatrix<double> dnk_dnj =
m_kin->netProductionRates_ddCi();
129 size_t ssize = m_nv - m_sidx;
132 if (!m_surfaces.empty()) {
133 vector<Eigen::Triplet<double>> species_trips(dnk_dnj.nonZeros());
134 for (
int k = 0; k < dnk_dnj.outerSize(); k++) {
135 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
136 species_trips.emplace_back(
static_cast<int>(it.row()),
137 static_cast<int>(it.col()), it.value());
141 dnk_dnj.resize(ssize, ssize);
142 dnk_dnj.setFromTriplets(species_trips.begin(), species_trips.end());
145 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(ssize);
147 m_kin->getNetProductionRates(netProductionRates.data());
149 for (
auto &S: m_surfaces) {
150 auto curr_kin = S->kinetics();
151 vector<double> prod_rates(curr_kin->nTotalSpecies());
152 curr_kin->getNetProductionRates(prod_rates.data());
153 for (
size_t i = 0; i < curr_kin->nTotalSpecies(); i++) {
155 size_t row =
speciesIndex(curr_kin->kineticsSpeciesName(i));
156 netProductionRates[row] += prod_rates[i];
162 double molarVol = m_thermo->molarVolume();
166 for (
int k = 0; k < dnk_dnj.outerSize(); k++) {
167 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
169 if (
static_cast<size_t>(it.row()) <
m_nsp) {
170 it.valueRef() = it.value() + netProductionRates[it.row()] * molarVol;
172 m_jac_trips.emplace_back(
static_cast<int>(it.row() + m_sidx),
173 static_cast<int>(it.col() + m_sidx), it.value());
179 double deltaTemp = m_thermo->temperature()
180 * std::sqrt(std::numeric_limits<double>::epsilon());
182 vector<double> yCurrent(m_nv);
185 vector<double> lhsPerturbed(m_nv, 1.0), lhsCurrent(m_nv, 1.0);
186 vector<double> rhsPerturbed(m_nv, 0.0), rhsCurrent(m_nv, 0.0);
187 vector<double> yPerturbed = yCurrent;
189 yPerturbed[0] += deltaTemp;
192 double time = (
m_net !=
nullptr) ?
m_net->time() : 0.0;
193 eval(time, lhsPerturbed.data(), rhsPerturbed.data());
196 eval(time, lhsCurrent.data(), rhsCurrent.data());
198 for (
size_t j = 0; j < m_nv; j++) {
199 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
200 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
202 (ydotPerturbed - ydotCurrent) / deltaTemp);
206 Eigen::VectorXd enthalpy = Eigen::VectorXd::Zero(ssize);
207 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(ssize);
209 m_thermo->getPartialMolarCp(specificHeat.data());
210 m_thermo->getPartialMolarEnthalpies(enthalpy.data());
212 for (
size_t i = 0; i <
m_nsp; i++) {
213 netProductionRates[i] *=
m_vol;
215 double qdot = enthalpy.dot(netProductionRates);
218 double* moles = yCurrent.data() + m_sidx;
219 for (
size_t i = 0; i < ssize; i++) {
220 NCp += moles[i] * specificHeat[i];
222 double denom = 1 / (NCp * NCp);
223 Eigen::VectorXd hk_dnkdnj_sums = dnk_dnj.transpose() * enthalpy;
225 for (
size_t j = 0; j < ssize; j++) {
226 m_jac_trips.emplace_back(0,
static_cast<int>(j + m_sidx),
227 (specificHeat[j] * qdot - NCp * hk_dnkdnj_sums[j]) * denom);
231 Eigen::SparseMatrix<double> jac(m_nv, m_nv);
238 if (nm ==
"temperature") {
244 throw CanteraError(
"IdealGasConstPressureReactor::componentIndex",
245 "Component '{}' not found", nm);
251 return "temperature";
252 }
else if (k >= m_sidx && k <
neq()) {
254 if (k < m_thermo->nSpecies()) {
255 return m_thermo->speciesName(k);
257 k -= m_thermo->nSpecies();
259 for (
auto& S : m_surfaces) {
261 if (k < th->nSpecies()) {
268 throw IndexError(
"IdealGasConstPressureMoleReactor::componentName",
269 "components", k, m_nv);
275 return 1.5 * m_thermo->maxTemp();
284 return 0.5 * m_thermo->minTemp();
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.
void initialize(double t0=0.0) override
Initialize the reactor.
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.
Eigen::SparseMatrix< double > jacobian() override
Calculate an approximate Jacobian to accelerate preconditioned solvers.
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.
vector< double > m_hk
Species molar enthalpies.
An array index is out of range.
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.
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...
void updateSurfaceState(double *y) override
Update the state of SurfPhase objects attached to this reactor.
size_t nSpecies() const
Returns the number of species in the phase.
string speciesName(size_t k) const
Name of the species with index k.
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.
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.
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
vector< double > m_wdot
Species net molar production rates.
virtual void evalWalls(double t)
Evaluate terms related to Walls.
size_t neq()
Number of equations (state variables) for this reactor.
double m_Qdot
net heat transfer into the reactor, through walls [W]
vector< Eigen::Triplet< double > > m_jac_trips
Vector of triplets representing the jacobian.
vector< double > m_sdot
Production rates of gas phase species on surfaces [kmol/s].
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...
virtual void updateConnected(bool updatePressure)
Update the state information needed by connected reactors, flow devices, and reactor walls.
Base class for a phase with thermodynamic properties.
Namespace for the Cantera kernel.
const double BigNumber
largest number to compare to inf.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...