Cantera 2.6.0
|
This class represents 1D flow domains that satisfy the one-dimensional similarity solution for chemically-reacting, axisymmetric flows. More...
#include <StFlow.h>
Public Member Functions | |
StFlow (ThermoPhase *ph=0, size_t nsp=1, size_t points=1) | |
Create a new flow domain. More... | |
StFlow (shared_ptr< ThermoPhase > th, size_t nsp=1, size_t points=1) | |
Delegating constructor. More... | |
virtual std::string | componentName (size_t n) const |
Name of the nth component. May be overloaded. More... | |
virtual size_t | componentIndex (const std::string &name) const |
index of component with name name. More... | |
virtual bool | componentActive (size_t n) const |
Returns true if the specified component is an active part of the solver state. More... | |
virtual void | showSolution (const doublereal *x) |
Print the solution. More... | |
virtual XML_Node & | save (XML_Node &o, const doublereal *const sol) |
Save the current solution for this domain into an XML_Node. More... | |
virtual void | restore (const XML_Node &dom, doublereal *soln, int loglevel) |
Restore the solution for this domain from an XML_Node. More... | |
virtual AnyMap | serialize (const double *soln) const |
Save the state of this domain as an AnyMap. More... | |
virtual void | restore (const AnyMap &state, double *soln, int loglevel) |
Restore the solution for this domain from an AnyMap. More... | |
void | setFreeFlow () |
Set flow configuration for freely-propagating flames, using an internal point with a fixed temperature as the condition to determine the inlet mass flux. More... | |
void | setAxisymmetricFlow () |
Set flow configuration for axisymmetric counterflow or burner-stabilized flames, using specified inlet mass fluxes. More... | |
virtual std::string | flowType () const |
Return the type of flow domain being represented, either "Free Flame" or "Axisymmetric Stagnation". More... | |
void | solveEnergyEqn (size_t j=npos) |
void | enableRadiation (bool doRadiation) |
Turn radiation on / off. More... | |
bool | radiationEnabled () const |
Returns true if the radiation term in the energy equation is enabled. More... | |
double | radiativeHeatLoss (size_t j) const |
Return radiative heat loss at grid point j. More... | |
void | setBoundaryEmissivities (double e_left, double e_right) |
Set the emissivities for the boundary values. More... | |
double | leftEmissivity () const |
Return emissivitiy at left boundary. More... | |
double | rightEmissivity () const |
Return emissivitiy at right boundary. More... | |
void | fixTemperature (size_t j=npos) |
bool | doEnergy (size_t j) |
virtual void | resize (size_t components, size_t points) |
Change the grid size. Called after grid refinement. More... | |
void | setGas (const doublereal *x, size_t j) |
Set the gas object state to be consistent with the solution at point j. More... | |
void | setGasAtMidpoint (const doublereal *x, size_t j) |
Set the gas state to be consistent with the solution at the midpoint between j and j + 1. More... | |
doublereal | density (size_t j) const |
virtual bool | fixed_mdot () |
void | setViscosityFlag (bool dovisc) |
virtual void | eval (size_t j, doublereal *x, doublereal *r, integer *mask, doublereal rdt) |
virtual void | evalRightBoundary (double *x, double *res, int *diag, double rdt) |
Evaluate all residual components at the right boundary. More... | |
virtual void | evalContinuity (size_t j, double *x, double *r, int *diag, double rdt) |
Evaluate the residual corresponding to the continuity equation at all interior grid points. More... | |
size_t | leftExcessSpecies () const |
Index of the species on the left boundary with the largest mass fraction. More... | |
size_t | rightExcessSpecies () const |
Index of the species on the right boundary with the largest mass fraction. More... | |
Problem Specification | |
virtual void | setupGrid (size_t n, const doublereal *z) |
called to set up initial grid, and after grid refinement More... | |
virtual void | resetBadValues (double *xg) |
ThermoPhase & | phase () |
Kinetics & | kinetics () |
void | setThermo (IdealGasPhase &th) |
Set the thermo manager. More... | |
void | setKinetics (Kinetics &kin) |
Set the kinetics manager. The kinetics manager must. More... | |
void | setTransport (Transport &trans) |
set the transport manager More... | |
void | enableSoret (bool withSoret) |
Enable thermal diffusion, also known as Soret diffusion. More... | |
bool | withSoret () const |
void | setPressure (doublereal p) |
Set the pressure. More... | |
doublereal | pressure () const |
The current pressure [Pa]. More... | |
virtual void | _getInitialSoln (double *x) |
Write the initial solution estimate into array x. More... | |
virtual void | _finalize (const doublereal *x) |
In some cases, a domain may need to set parameters that depend on the initial solution estimate. More... | |
void | setFixedTempProfile (vector_fp &zfixed, vector_fp &tfixed) |
Sometimes it is desired to carry out the simulation using a specified temperature profile, rather than computing it by solving the energy equation. More... | |
void | setTemperature (size_t j, doublereal t) |
doublereal | T_fixed (size_t j) const |
The fixed temperature value at point j. More... | |
Public Member Functions inherited from Domain1D | |
Domain1D (size_t nv=1, size_t points=1, double time=0.0) | |
Constructor. More... | |
Domain1D (const Domain1D &)=delete | |
Domain1D & | operator= (const Domain1D &)=delete |
int | domainType () |
Domain type flag. More... | |
size_t | domainIndex () |
The left-to-right location of this domain. More... | |
bool | isConnector () |
True if the domain is a connector domain. More... | |
const OneDim & | container () const |
The container holding this domain. More... | |
void | setContainer (OneDim *c, size_t index) |
Specify the container object for this domain, and the position of this domain in the list. More... | |
void | setBandwidth (int bw=-1) |
Set the Jacobian bandwidth. See the discussion of method bandwidth(). More... | |
size_t | bandwidth () |
Set the Jacobian bandwidth for this domain. More... | |
virtual void | init () |
virtual void | setInitialState (doublereal *xlocal=0) |
virtual void | setState (size_t point, const doublereal *state, doublereal *x) |
Refiner & | refiner () |
Return a reference to the grid refiner. More... | |
size_t | nComponents () const |
Number of components at each grid point. More... | |
void | checkComponentIndex (size_t n) const |
Check that the specified component index is in range. More... | |
void | checkComponentArraySize (size_t nn) const |
Check that an array size is at least nComponents(). More... | |
size_t | nPoints () const |
Number of grid points in this domain. More... | |
void | checkPointIndex (size_t n) const |
Check that the specified point index is in range. More... | |
void | checkPointArraySize (size_t nn) const |
Check that an array size is at least nPoints(). More... | |
void | setComponentName (size_t n, const std::string &name) |
void | setBounds (size_t n, doublereal lower, doublereal upper) |
void | setTransientTolerances (doublereal rtol, doublereal atol, size_t n=npos) |
Set tolerances for time-stepping mode. More... | |
void | setSteadyTolerances (doublereal rtol, doublereal atol, size_t n=npos) |
Set tolerances for steady-state mode. More... | |
doublereal | rtol (size_t n) |
Relative tolerance of the nth component. More... | |
doublereal | atol (size_t n) |
Absolute tolerance of the nth component. More... | |
double | steady_rtol (size_t n) |
Steady relative tolerance of the nth component. More... | |
double | steady_atol (size_t n) |
Steady absolute tolerance of the nth component. More... | |
double | transient_rtol (size_t n) |
Transient relative tolerance of the nth component. More... | |
double | transient_atol (size_t n) |
Transient absolute tolerance of the nth component. More... | |
doublereal | upperBound (size_t n) const |
Upper bound on the nth component. More... | |
doublereal | lowerBound (size_t n) const |
Lower bound on the nth component. More... | |
void | initTimeInteg (doublereal dt, const doublereal *x0) |
Prepare to do time stepping with time step dt. More... | |
void | setSteadyMode () |
Prepare to solve the steady-state problem. More... | |
bool | steady () |
True if in steady-state mode. More... | |
bool | transient () |
True if not in steady-state mode. More... | |
void | needJacUpdate () |
size_t | index (size_t n, size_t j) const |
doublereal | value (const doublereal *x, size_t n, size_t j) const |
virtual void | setJac (MultiJac *jac) |
size_t | size () const |
void | locate () |
Find the index of the first grid point in this domain, and the start of its variables in the global solution vector. More... | |
virtual size_t | loc (size_t j=0) const |
Location of the start of the local solution vector in the global solution vector,. More... | |
size_t | firstPoint () const |
The index of the first (that is, left-most) grid point belonging to this domain. More... | |
size_t | lastPoint () const |
The index of the last (that is, right-most) grid point belonging to this domain. More... | |
void | linkLeft (Domain1D *left) |
Set the left neighbor to domain 'left. More... | |
void | linkRight (Domain1D *right) |
Set the right neighbor to domain 'right.'. More... | |
void | append (Domain1D *right) |
Append domain 'right' to this one, and update all links. More... | |
Domain1D * | left () const |
Return a pointer to the left neighbor. More... | |
Domain1D * | right () const |
Return a pointer to the right neighbor. More... | |
double | prevSoln (size_t n, size_t j) const |
Value of component n at point j in the previous solution. More... | |
void | setID (const std::string &s) |
Specify an identifying tag for this domain. More... | |
std::string | id () const |
virtual void | showSolution_s (std::ostream &s, const doublereal *x) |
doublereal | z (size_t jlocal) const |
doublereal | zmin () const |
doublereal | zmax () const |
void | setProfile (const std::string &name, double *values, double *soln) |
vector_fp & | grid () |
const vector_fp & | grid () const |
doublereal | grid (size_t point) const |
virtual doublereal | initialValue (size_t n, size_t j) |
Initial value of solution component n at grid point j. More... | |
void | forceFullUpdate (bool update) |
In some cases, for computational efficiency some properties (such as transport coefficients) may not be updated during Jacobian evaluations. More... | |
Public Attributes | |
double | m_zfixed |
Location of the point where temperature is fixed. More... | |
double | m_tfixed |
Temperature at the point used to fix the flame location. More... | |
Protected Member Functions | |
doublereal | wdot (size_t k, size_t j) const |
void | getWdot (doublereal *x, size_t j) |
Write the net production rates at point j into array m_wdot More... | |
virtual void | updateProperties (size_t jg, double *x, size_t jmin, size_t jmax) |
Update the properties (thermo, transport, and diffusion flux). More... | |
virtual void | evalResidual (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the residual function. More... | |
void | updateThermo (const doublereal *x, size_t j0, size_t j1) |
Update the thermodynamic properties from point j0 to point j1 (inclusive), based on solution x. More... | |
doublereal | shear (const doublereal *x, size_t j) const |
doublereal | divHeatFlux (const doublereal *x, size_t j) const |
size_t | mindex (size_t k, size_t j, size_t m) |
virtual void | updateDiffFluxes (const doublereal *x, size_t j0, size_t j1) |
Update the diffusive mass fluxes. More... | |
virtual void | updateTransport (doublereal *x, size_t j0, size_t j1) |
Update the transport properties at grid points in the range from j0 to j1 , based on solution x . More... | |
Solution components | |
doublereal | T (const doublereal *x, size_t j) const |
doublereal & | T (doublereal *x, size_t j) |
doublereal | T_prev (size_t j) const |
doublereal | rho_u (const doublereal *x, size_t j) const |
doublereal | u (const doublereal *x, size_t j) const |
doublereal | V (const doublereal *x, size_t j) const |
doublereal | V_prev (size_t j) const |
doublereal | lambda (const doublereal *x, size_t j) const |
doublereal | Y (const doublereal *x, size_t k, size_t j) const |
doublereal & | Y (doublereal *x, size_t k, size_t j) |
doublereal | Y_prev (size_t k, size_t j) const |
doublereal | X (const doublereal *x, size_t k, size_t j) const |
doublereal | flux (size_t k, size_t j) const |
convective spatial derivatives. | |
These use upwind differencing, assuming u(z) is negative | |
doublereal | dVdz (const doublereal *x, size_t j) const |
doublereal | dYdz (const doublereal *x, size_t k, size_t j) const |
doublereal | dTdz (const doublereal *x, size_t j) const |
Protected Attributes | |
doublereal | m_press |
vector_fp | m_dz |
vector_fp | m_rho |
vector_fp | m_wtm |
vector_fp | m_wt |
vector_fp | m_cp |
vector_fp | m_visc |
vector_fp | m_tcon |
vector_fp | m_diff |
vector_fp | m_multidiff |
Array2D | m_dthermal |
Array2D | m_flux |
Array2D | m_wdot |
size_t | m_nsp |
IdealGasPhase * | m_thermo |
Kinetics * | m_kin |
Transport * | m_trans |
doublereal | m_epsilon_left |
doublereal | m_epsilon_right |
std::vector< size_t > | m_kRadiating |
Indices within the ThermoPhase of the radiating species. More... | |
std::vector< bool > | m_do_energy |
bool | m_do_soret |
std::vector< bool > | m_do_species |
bool | m_do_multicomponent |
bool | m_do_radiation |
flag for the radiative heat loss More... | |
vector_fp | m_qdotRadiation |
radiative heat loss vector More... | |
vector_fp | m_fixedtemp |
vector_fp | m_zfix |
vector_fp | m_tfix |
size_t | m_kExcessLeft |
Index of species with a large mass fraction at each boundary, for which the mass fraction may be calculated as 1 minus the sum of the other mass fractions. More... | |
size_t | m_kExcessRight |
bool | m_dovisc |
Protected Attributes inherited from Domain1D | |
doublereal | m_rdt |
size_t | m_nv |
size_t | m_points |
vector_fp | m_slast |
vector_fp | m_max |
vector_fp | m_min |
vector_fp | m_rtol_ss |
vector_fp | m_rtol_ts |
vector_fp | m_atol_ss |
vector_fp | m_atol_ts |
vector_fp | m_z |
OneDim * | m_container |
size_t | m_index |
int | m_type |
size_t | m_iloc |
Starting location within the solution vector for unknowns that correspond to this domain. More... | |
size_t | m_jstart |
Domain1D * | m_left |
Domain1D * | m_right |
std::string | m_id |
Identity tag for the domain. More... | |
std::unique_ptr< Refiner > | m_refiner |
std::vector< std::string > | m_name |
int | m_bw |
bool | m_force_full_update |
Private Attributes | |
vector_fp | m_ybar |
This class represents 1D flow domains that satisfy the one-dimensional similarity solution for chemically-reacting, axisymmetric flows.
StFlow | ( | ThermoPhase * | ph = 0 , |
size_t | nsp = 1 , |
||
size_t | points = 1 |
||
) |
Create a new flow domain.
ph | Object representing the gas phase. This object will be used to evaluate all thermodynamic, kinetic, and transport properties. |
nsp | Number of species. |
points | Initial number of grid points |
Definition at line 20 of file StFlow.cpp.
References StFlow::m_kRadiating, StFlow::m_qdotRadiation, ThermoPhase::maxTemp(), Phase::molecularWeights(), Cantera::npos, Phase::nSpecies(), Array2D::resize(), Domain1D::resize(), StFlow::setupGrid(), Phase::speciesIndex(), and ThermoPhase::type().
|
inline |
|
virtual |
called to set up initial grid, and after grid refinement
Reimplemented from Domain1D.
Definition at line 131 of file StFlow.cpp.
References StFlow::resize().
Referenced by StFlow::restore(), and StFlow::StFlow().
|
virtual |
When called, this function should reset "bad" values in the state vector such as negative species concentrations. This function may be called after a failed solution attempt.
Reimplemented from Domain1D.
Definition at line 146 of file StFlow.cpp.
References Phase::getMassFractions(), Domain1D::loc(), and Phase::setMassFractions().
|
inline |
|
inline |
|
inline |
void setTransport | ( | Transport & | trans | ) |
set the transport manager
Definition at line 156 of file StFlow.cpp.
References Array2D::resize(), and Transport::transportType().
|
inline |
|
inline |
Set the pressure.
Since the flow equations are for the limit of small Mach number, the pressure is very nearly constant throughout the flow.
Definition at line 98 of file StFlow.h.
Referenced by StFlow::restore().
|
inline |
|
virtual |
Write the initial solution estimate into array x.
Reimplemented from Domain1D.
Definition at line 168 of file StFlow.cpp.
References Phase::getMassFractions(), and Phase::temperature().
|
virtual |
In some cases, a domain may need to set parameters that depend on the initial solution estimate.
In such cases, the parameters may be set in method _finalize. This method is called just before the Newton solver is called, and the x array is guaranteed to be the local solution vector for this domain that will be used as the initial guess. If no such parameters need to be set, then method _finalize does not need to be overloaded.
Reimplemented from Domain1D.
Reimplemented in IonFlow.
Definition at line 196 of file StFlow.cpp.
References Domain1D::domainType(), Cantera::linearInterp(), StFlow::m_tfixed, StFlow::m_zfixed, and Cantera::Undef.
Sometimes it is desired to carry out the simulation using a specified temperature profile, rather than computing it by solving the energy equation.
This method specifies this profile.
Definition at line 115 of file StFlow.h.
Referenced by StFlow::restore().
|
inline |
|
inline |
The fixed temperature value at point j.
Definition at line 130 of file StFlow.h.
Referenced by StFlow::evalResidual(), and StFlow::evalRightBoundary().
|
virtual |
Name of the nth component. May be overloaded.
Reimplemented from Domain1D.
Definition at line 578 of file StFlow.cpp.
References Phase::speciesName().
Referenced by StFlow::componentIndex(), StFlow::restore(), and StFlow::serialize().
|
virtual |
index of component with name name.
Reimplemented from Domain1D.
Definition at line 600 of file StFlow.cpp.
References StFlow::componentName().
|
virtual |
Returns true if the specified component is an active part of the solver state.
Reimplemented in IonFlow.
Definition at line 623 of file StFlow.cpp.
Referenced by StFlow::restore(), and StFlow::serialize().
|
virtual |
Print the solution.
Reimplemented from Domain1D.
Definition at line 519 of file StFlow.cpp.
References StFlow::m_do_radiation, StFlow::m_qdotRadiation, Domain1D::showSolution(), and Cantera::writelog().
Save the current solution for this domain into an XML_Node.
o | XML_Node to save the solution to. |
sol | Current value of the solution vector. The object will pick out which part of the solution vector pertains to this object. |
Reimplemented from Domain1D.
Definition at line 804 of file StFlow.cpp.
References XML_Node::addAttribute(), Cantera::addFloat(), Cantera::addFloatArray(), Cantera::addNamedFloatArray(), StFlow::flowType(), Array2D::getRow(), Domain1D::loc(), StFlow::m_do_radiation, StFlow::m_qdotRadiation, StFlow::m_tfixed, StFlow::m_zfixed, Array2D::nColumns(), Domain1D::nPoints(), Domain1D::refiner(), Domain1D::save(), Phase::speciesName(), and Cantera::Undef.
|
virtual |
Restore the solution for this domain from an XML_Node.
Base class version of the general Domain1D restore function. Derived classes should call the base class method in addition to restoring their own data.
dom | XML_Node for this domain |
soln | Current value of the solution vector, local to this object. |
loglevel | 0 to suppress all output; 1 to show warnings; 2 for verbose output |
Reimplemented from Domain1D.
Definition at line 637 of file StFlow.cpp.
References XML_Node::child(), Cantera::debuglog(), Domain1D::domainType(), XML_Node::getChildren(), Cantera::getFloat(), Cantera::getFloatArray(), Cantera::getOptionalFloat(), XML_Node::hasChild(), StFlow::m_tfixed, StFlow::m_zfixed, Domain1D::nPoints(), Cantera::npos, Phase::nSpecies(), Domain1D::refiner(), Domain1D::restore(), Refiner::setCriteria(), StFlow::setFixedTempProfile(), Refiner::setGridMin(), StFlow::setPressure(), StFlow::setupGrid(), Phase::speciesIndex(), Phase::speciesName(), XML_Node::value(), Cantera::warn_user(), and Cantera::writelog().
|
virtual |
Save the state of this domain as an AnyMap.
soln | local solution vector for this domain |
Reimplemented from Domain1D.
Definition at line 863 of file StFlow.cpp.
References AnyValue::asString(), StFlow::componentActive(), StFlow::componentName(), AnyValue::empty(), StFlow::flowType(), AnyBase::getMetadata(), ThermoPhase::input(), StFlow::m_do_radiation, StFlow::m_qdotRadiation, StFlow::m_tfixed, StFlow::m_zfixed, Phase::name(), Domain1D::nComponents(), Domain1D::nPoints(), Domain1D::serialize(), Phase::speciesName(), and Cantera::Undef.
|
virtual |
Restore the solution for this domain from an AnyMap.
[in] | state | AnyMap defining the state of this domain |
[out] | soln | Value of the solution vector, local to this domain |
[in] | loglevel | 0 to suppress all output; 1 to show warnings; 2 for verbose output |
Reimplemented from Domain1D.
Definition at line 925 of file StFlow.cpp.
References AnyValue::asBool(), AnyValue::asVector(), StFlow::componentActive(), StFlow::componentName(), AnyMap::getDouble(), AnyMap::hasKey(), AnyValue::isScalar(), StFlow::m_do_radiation, StFlow::m_tfixed, StFlow::m_zfixed, Domain1D::nComponents(), Domain1D::nPoints(), Phase::nSpecies(), Domain1D::restore(), StFlow::setupGrid(), and Cantera::warn_user().
|
inline |
|
inline |
|
inlinevirtual |
Return the type of flow domain being represented, either "Free Flame" or "Axisymmetric Stagnation".
Definition at line 182 of file StFlow.h.
Referenced by StFlow::save(), and StFlow::serialize().
void solveEnergyEqn | ( | size_t | j = npos | ) |
Definition at line 999 of file StFlow.cpp.
|
inline |
Turn radiation on / off.
The simple radiation model used was established by Y. Liu and B. Rogg [Y. Liu and B. Rogg, Modelling of thermally radiating diffusion flames with detailed chemistry and transport, EUROTHERM Seminars, 17:114-127, 1991]. This model considers the radiation of CO2 and H2O.
Definition at line 201 of file StFlow.h.
References StFlow::m_do_radiation.
|
inline |
Returns true
if the radiation term in the energy equation is enabled.
Definition at line 206 of file StFlow.h.
References StFlow::m_do_radiation.
|
inline |
Return radiative heat loss at grid point j.
Definition at line 211 of file StFlow.h.
References StFlow::m_qdotRadiation.
void setBoundaryEmissivities | ( | double | e_left, |
double | e_right | ||
) |
Set the emissivities for the boundary values.
Reads the emissivities for the left and right boundary values in the radiative term and writes them into the variables, which are used for the calculation.
Definition at line 1023 of file StFlow.cpp.
|
inline |
|
inline |
void fixTemperature | ( | size_t | j = npos | ) |
Definition at line 1037 of file StFlow.cpp.
|
virtual |
Change the grid size. Called after grid refinement.
Reimplemented from Domain1D.
Reimplemented in IonFlow.
Definition at line 107 of file StFlow.cpp.
References StFlow::m_qdotRadiation, Array2D::resize(), and Domain1D::resize().
Referenced by StFlow::setupGrid().
void setGas | ( | const doublereal * | x, |
size_t | j | ||
) |
Set the gas object state to be consistent with the solution at point j.
Definition at line 176 of file StFlow.cpp.
References Phase::setMassFractions_NoNorm(), IdealGasPhase::setPressure(), and Phase::setTemperature().
Referenced by ReactingSurf1D::eval(), StFlow::evalResidual(), StFlow::getWdot(), and StFlow::updateThermo().
void setGasAtMidpoint | ( | const doublereal * | x, |
size_t | j | ||
) |
Set the gas state to be consistent with the solution at the midpoint between j and j + 1.
Definition at line 184 of file StFlow.cpp.
References Phase::setMassFractions_NoNorm(), IdealGasPhase::setPressure(), and Phase::setTemperature().
Referenced by StFlow::updateTransport().
|
virtual |
Evaluate the residual function for axisymmetric stagnation flow. If j == npos, the residual function is evaluated at all grid points. Otherwise, the residual function is only evaluated at grid points j-1, j, and j+1. This option is used to efficiently evaluate the Jacobian numerically.
Reimplemented from Domain1D.
Definition at line 243 of file StFlow.cpp.
References StFlow::evalResidual(), Domain1D::firstPoint(), Domain1D::lastPoint(), Domain1D::loc(), Cantera::npos, and StFlow::updateProperties().
|
virtual |
Evaluate all residual components at the right boundary.
Definition at line 1061 of file StFlow.cpp.
References Domain1D::domainType(), StFlow::rightExcessSpecies(), and StFlow::T_fixed().
Referenced by StFlow::evalResidual().
|
virtual |
Evaluate the residual corresponding to the continuity equation at all interior grid points.
Definition at line 1092 of file StFlow.cpp.
References Domain1D::domainType(), StFlow::m_tfixed, and StFlow::m_zfixed.
Referenced by StFlow::evalResidual().
|
inline |
Index of the species on the left boundary with the largest mass fraction.
Definition at line 276 of file StFlow.h.
References StFlow::m_kExcessLeft.
Referenced by StFlow::evalResidual().
|
inline |
Index of the species on the right boundary with the largest mass fraction.
Definition at line 281 of file StFlow.h.
Referenced by Outlet1D::eval(), OutletRes1D::eval(), ReactingSurf1D::eval(), and StFlow::evalRightBoundary().
|
inlineprotected |
|
inlineprotected |
Write the net production rates at point j
into array m_wdot
Definition at line 291 of file StFlow.h.
References Kinetics::getNetProductionRates(), and StFlow::setGas().
Referenced by StFlow::evalResidual().
|
protectedvirtual |
Update the properties (thermo, transport, and diffusion flux).
This function is called in eval after the points which need to be updated are defined.
Definition at line 271 of file StFlow.cpp.
References StFlow::m_kExcessLeft, Cantera::npos, StFlow::updateDiffFluxes(), StFlow::updateThermo(), and StFlow::updateTransport().
Referenced by StFlow::eval().
|
protectedvirtual |
Evaluate the residual function.
This function is called in eval after updateProperties is called.
Reimplemented in IonFlow.
Definition at line 295 of file StFlow.cpp.
References IdealGasPhase::cp_R_ref(), IdealGasPhase::enthalpy_RT_ref(), StFlow::evalContinuity(), StFlow::evalRightBoundary(), Cantera::GasConstant, StFlow::getWdot(), StFlow::leftExcessSpecies(), StFlow::m_do_radiation, StFlow::m_kRadiating, StFlow::m_qdotRadiation, Cantera::npos, Cantera::OneAtm, StFlow::setGas(), Cantera::StefanBoltz, and StFlow::T_fixed().
Referenced by StFlow::eval().
|
inlineprotected |
Update the thermodynamic properties from point j0 to point j1 (inclusive), based on solution x.
Definition at line 310 of file StFlow.h.
References ThermoPhase::cp_mass(), Phase::density(), Phase::meanMolecularWeight(), and StFlow::setGas().
Referenced by StFlow::updateProperties().
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
inlineprotected |
|
protectedvirtual |
Update the diffusive mass fluxes.
Reimplemented in IonFlow.
Definition at line 536 of file StFlow.cpp.
Referenced by StFlow::updateProperties().
|
protectedvirtual |
Update the transport properties at grid points in the range from j0
to j1
, based on solution x
.
Reimplemented in IonFlow.
Definition at line 489 of file StFlow.cpp.
References Phase::density(), Transport::getMixDiffCoeffs(), Transport::getMultiDiffCoeffs(), Transport::getThermalDiffCoeffs(), Phase::meanMolecularWeight(), Array2D::ptrColumn(), StFlow::setGasAtMidpoint(), Transport::thermalConductivity(), and Transport::viscosity().
Referenced by StFlow::updateProperties().
|
protected |
|
protected |
Indices within the ThermoPhase of the radiating species.
First index is for CO2, second is for H2O.
Definition at line 450 of file StFlow.h.
Referenced by StFlow::evalResidual(), and StFlow::StFlow().
|
protected |
flag for the radiative heat loss
Definition at line 459 of file StFlow.h.
Referenced by StFlow::enableRadiation(), StFlow::evalResidual(), StFlow::radiationEnabled(), StFlow::restore(), StFlow::save(), StFlow::serialize(), and StFlow::showSolution().
|
protected |
radiative heat loss vector
Definition at line 462 of file StFlow.h.
Referenced by StFlow::evalResidual(), StFlow::radiativeHeatLoss(), StFlow::resize(), StFlow::save(), StFlow::serialize(), StFlow::showSolution(), and StFlow::StFlow().
|
protected |
Index of species with a large mass fraction at each boundary, for which the mass fraction may be calculated as 1 minus the sum of the other mass fractions.
Definition at line 472 of file StFlow.h.
Referenced by StFlow::leftExcessSpecies(), and StFlow::updateProperties().
double m_zfixed |
Location of the point where temperature is fixed.
Definition at line 483 of file StFlow.h.
Referenced by StFlow::_finalize(), StFlow::evalContinuity(), StFlow::restore(), StFlow::save(), and StFlow::serialize().
double m_tfixed |
Temperature at the point used to fix the flame location.
Definition at line 486 of file StFlow.h.
Referenced by StFlow::_finalize(), StFlow::evalContinuity(), StFlow::restore(), StFlow::save(), and StFlow::serialize().