16#include "cantera/numerics/eigen_dense.h"
55 double& mcpdTdt = RHS[0];
56 auto dndt = RHS.subspan(m_sidx);
73 for (
size_t n = 0; n <
m_nsp; n++) {
82 for (
auto outlet : m_outlet) {
83 for (
size_t n = 0; n <
m_nsp; n++) {
90 for (
auto inlet : m_inlet) {
93 for (
size_t n = 0; n <
m_nsp; n++) {
97 mcpdTdt -=
m_hk[n] * imw[n] * mdot_spec;
119 for (
int k = 0; k < dnk_dnj.outerSize(); k++) {
120 for (Eigen::SparseMatrix<double>::InnerIterator it(dnk_dnj, k); it; ++it) {
121 trips.emplace_back(it.row() +
offset, it.col() +
offset, it.value());
129 * std::sqrt(std::numeric_limits<double>::epsilon());
131 vector<double> yCurrent(
m_nv);
134 vector<double> lhsPerturbed(
m_nv, 1.0), lhsCurrent(
m_nv, 1.0);
135 vector<double> rhsPerturbed(
m_nv, 0.0), rhsCurrent(
m_nv, 0.0);
136 vector<double> yPerturbed = yCurrent;
138 yPerturbed[0] += deltaTemp;
142 eval(time, lhsPerturbed, rhsPerturbed);
145 eval(time, lhsCurrent, rhsCurrent);
147 for (
size_t j = 0; j <
m_nv; j++) {
148 double ydotPerturbed = rhsPerturbed[j] / lhsPerturbed[j];
149 double ydotCurrent = rhsCurrent[j] / lhsCurrent[j];
150 trips.emplace_back(
static_cast<int>(j +
m_offset),
static_cast<int>(
m_offset),
151 (ydotPerturbed - ydotCurrent) / deltaTemp);
155 Eigen::VectorXd netProductionRates = Eigen::VectorXd::Zero(
m_nsp);
156 Eigen::VectorXd enthalpy = Eigen::VectorXd::Zero(
m_nsp);
157 Eigen::VectorXd specificHeat = Eigen::VectorXd::Zero(
m_nsp);
163 for (
size_t i = 0; i <
m_nsp; i++) {
164 netProductionRates[i] *=
m_vol;
166 double qdot = enthalpy.dot(netProductionRates);
168 Eigen::VectorXd hk_dnkdnj_sums = dnk_dnj.transpose() * enthalpy;
170 for (
size_t j = 0; j <
m_nsp; j++) {
172 (specificHeat[j] * qdot -
m_TotalCp * hk_dnkdnj_sums[j]) * denom);
178 bool includePressureSpecies =
183 for (
auto outlet : m_outlet) {
184 for (
size_t n = 0; n <
m_nsp; n++) {
187 includeComposition, includePressureSpecies);
191 for (
auto inlet : m_inlet) {
192 for (
size_t n = 0; n <
m_nsp; n++) {
195 includeComposition, includePressureSpecies);
203 for (
size_t n = 0; n <
m_nsp; n++) {
206 includeComposition, includePressureSpecies);
213 for (
size_t i = 0; i < m_wall.size(); i++) {
214 int f = 2 *
m_lr[i] - 1;
224 double& f_species, span<double> f_energy)
226 f_species = 1.0 /
m_vol;
227 for (
size_t k = 0; k <
m_nsp; k++) {
233 SparseTriplets& trips,
size_t row,
double coeff)
const
239 trips.emplace_back(row,
m_offset, coeff);
244 SparseTriplets& trips,
size_t row,
size_t k,
double coeff)
const
249 for (
size_t j = 0; j <
m_nsp; j++) {
250 double dYdn = -Yk * mw[j] /
m_mass;
254 trips.emplace_back(row,
m_offset + m_sidx + j, coeff * dYdn);
259 size_t row,
double coeff,
bool includeComposition)
const
263 if (includeComposition) {
267 for (
size_t k = 0; k <
m_nsp; k++) {
268 trips.emplace_back(row,
m_offset + m_sidx + k,
276 if (nm ==
"temperature") {
282 throw CanteraError(
"IdealGasConstPressureReactor::componentIndex",
283 "Component '{}' not found", nm);
289 return "temperature";
290 }
else if (k >= m_sidx && k <
neq()) {
293 throw IndexError(
"IdealGasConstPressureMoleReactor::componentName",
294 "components", k,
m_nv);
300 return 1.5 * m_thermo->
maxTemp();
309 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,...
Header file for base class WallBase.
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.
void addInletEnthalpyJacobian(SparseTriplets &trips, size_t row, double coeff, bool includeComposition=true)
Add Jacobian terms for the inlet enthalpy dependence on upstream state.
virtual void addMassFlowRateJacobian(SparseTriplets &trips, size_t row, double coeff, bool includePressureSpecies=true)
Add Jacobian terms proportional to derivatives of the mass flow rate.
double enthalpy_mass()
specific enthalpy
void addOutletSpeciesMassFlowRateJacobian(SparseTriplets &trips, size_t row, size_t k, double coeff, bool includeComposition=true, bool includePressureSpecies=true)
Add Jacobian terms proportional to derivatives of outletSpeciesMassFlowRate(k).
double massFlowRate()
Mass flow rate (kg/s).
double upperBound(size_t k) const override
Get the upper bound on the k-th component of the local state vector.
void addTemperatureJacobian(SparseTriplets &trips, size_t row, double coeff) const override
Add terms proportional to derivatives with respect to this reactor's temperature.
void eval(double t, span< double > LHS, span< double > RHS) override
Evaluate the reactor governing equations.
double m_TotalCp
Total heat capacity ( ) [J/K].
void addSpeciesMassFractionJacobian(SparseTriplets &trips, size_t row, size_t k, double coeff) const override
Add terms proportional to derivatives with respect to a species mass fraction.
size_t componentIndex(const string &nm) const override
Return the index in the solution vector for this reactor of the component named nm.
void addEnthalpyJacobian(SparseTriplets &trips, size_t row, double coeff, bool includeComposition=true) const override
Add terms proportional to derivatives of this reactor's specific enthalpy.
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 getJacobianElements(SparseTriplets &trips) override
Calculate an approximate Jacobian to accelerate preconditioned solvers.
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].
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.
double massFraction(size_t k) const
Return the mass fraction of a single species.
virtual void setMolesNoTruncate(span< const double > N)
Set the state of the object with moles in [kmol].
size_t speciesIndex(const string &name, bool raise=true) const
Returns the index of a species named 'name' within the Phase object.
span< const double > molecularWeights() const
Return a const reference to the internal vector of molecular weights.
double temperature() const
Temperature (K).
string speciesName(size_t k) const
Name of the species with index k.
virtual double density() const
Density (kg/m^3).
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.
bool isJacobianLocalRow(size_t row) const
Check whether a global Jacobian row belongs to this reactor.
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.
vector< int > m_lr
Vector of length nWalls(), indicating whether this reactor is on the left (0) or right (1) of each wa...
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.
bool m_jac_skip_connector_pressure_composition_dependence
Omit species terms in pressure derivatives, which can add fill-in.
bool m_jac_skip_connector_composition_dependence
Omit flow composition derivatives, which can add dense connector blocks.
Kinetics * m_kin
Pointer to the homogeneous Kinetics object that handles the reactions.
vector< double > m_wdot
Species net molar production rates.
bool m_jac_skip_flow_devices
Omit flow-device terms from the sparse Jacobian.
double m_Qdot
net heat transfer into the reactor, through walls [W]
bool m_jac_skip_walls
Omit wall terms from the sparse Jacobian.
void updateSurfaceProductionRates()
Update m_sdot to reflect current production rates of bulk phase species due to reactions on adjacent ...
vector< double > m_sdot
Total production rate of bulk phase species on surfaces [kmol/s].
void initialize(double t0=0.0) override
Initialize the reactor.
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.
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].
double enthalpy_mass() const
Specific enthalpy. Units: J/kg.
Eigen::SparseMatrix< double > netProductionRates_ddCi()
Calculate derivatives for species net production rates with respect to species concentration at const...
Namespace for the Cantera kernel.
span< double > asSpan(Eigen::DenseBase< Derived > &v)
Convenience wrapper for accessing Eigen vector/array/map data as a span.
offset
Offsets of solution components in the 1D solution array.
const double BigNumber
largest number to compare to inf.
Various templated functions that carry out common vector and polynomial operations (see Templated Arr...