Cantera
3.1.0a1
|
This class represents 1D flow domains that satisfy the one-dimensional similarity solution for chemically-reacting, axisymmetric flows. More...
#include <StFlow.h>
This class represents 1D flow domains that satisfy the one-dimensional similarity solution for chemically-reacting, axisymmetric flows.
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... | |
StFlow (shared_ptr< Solution > sol, const string &id="", size_t points=1) | |
Create a new flow domain. More... | |
string | domainType () const override |
Domain type flag. More... | |
string | componentName (size_t n) const override |
Name of the nth component. May be overloaded. More... | |
size_t | componentIndex (const string &name) const override |
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... | |
void | show (const double *x) override |
Print the solution. More... | |
shared_ptr< SolutionArray > | asArray (const double *soln) const override |
Save the state of this domain as a SolutionArray. More... | |
void | fromArray (SolutionArray &arr, double *soln) override |
Restore the solution for this domain from a SolutionArray. 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 flames, using specified inlet mass fluxes. More... | |
void | setUnstrainedFlow () |
Set flow configuration for burner-stabilized flames, using specified inlet mass fluxes. More... | |
void | solveEnergyEqn (size_t j=npos) |
virtual size_t | getSolvingStage () const |
Get the solving stage (used by IonFlow specialization) More... | |
virtual void | setSolvingStage (const size_t stage) |
Solving stage mode for handling ionized species (used by IonFlow specialization) More... | |
virtual void | solveElectricField (size_t j=npos) |
Set to solve electric field in a point (used by IonFlow specialization) More... | |
virtual void | fixElectricField (size_t j=npos) |
Set to fix voltage in a point (used by IonFlow specialization) More... | |
virtual bool | doElectricField (size_t j) const |
Retrieve flag indicating whether electric field is solved or not (used by IonFlow specialization) More... | |
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 emissivity at left boundary. More... | |
double | rightEmissivity () const |
Return emissivity at right boundary. More... | |
void | fixTemperature (size_t j=npos) |
bool | doEnergy (size_t j) |
void | resize (size_t components, size_t points) override |
Change the grid size. Called after grid refinement. More... | |
void | setGas (const double *x, size_t j) |
Set the gas object state to be consistent with the solution at point j. More... | |
void | setGasAtMidpoint (const double *x, size_t j) |
Set the gas state to be consistent with the solution at the midpoint between j and j + 1. More... | |
double | density (size_t j) const |
bool | isFree () const |
Retrieve flag indicating whether flow is freely propagating. More... | |
bool | isStrained () const |
Retrieve flag indicating whether flow uses radial momentum. More... | |
void | setViscosityFlag (bool dovisc) |
void | eval (size_t jGlobal, double *xGlobal, double *rsdGlobal, integer *diagGlobal, double rdt) override |
Evaluate the residual functions for axisymmetric stagnation flow. 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 | |
void | setupGrid (size_t n, const double *z) override |
called to set up initial grid, and after grid refinement More... | |
void | resetBadValues (double *xg) override |
When called, this function should reset "bad" values in the state vector such as negative species concentrations. More... | |
ThermoPhase & | phase () |
Kinetics & | kinetics () |
void | setKinetics (shared_ptr< Kinetics > kin) override |
Set the kinetics manager. More... | |
void | setTransport (shared_ptr< Transport > trans) override |
Set transport model to existing instance. More... | |
void | setTransportModel (const string &trans) |
Set the transport model. More... | |
string | transportModel () const |
Retrieve transport model. More... | |
void | enableSoret (bool withSoret) |
Enable thermal diffusion, also known as Soret diffusion. More... | |
bool | withSoret () const |
void | setPressure (double p) |
Set the pressure. More... | |
double | pressure () const |
The current pressure [Pa]. More... | |
void | _getInitialSoln (double *x) override |
Write the initial solution estimate into array x. More... | |
void | _finalize (const double *x) override |
In some cases, a domain may need to set parameters that depend on the initial solution estimate. More... | |
void | setFixedTempProfile (vector< double > &zfixed, vector< double > &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, double t) |
Set the temperature fixed point at grid point j, and disable the energy equation so that the solution will be held to this value. More... | |
double | 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 |
string | type () const |
String indicating the domain implemented. More... | |
size_t | domainIndex () |
The left-to-right location of this domain. More... | |
virtual bool | isConnector () |
True if the domain is a connector domain. More... | |
void | setSolution (shared_ptr< Solution > sol) |
Set the solution manager. 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 () |
Initialize. More... | |
virtual void | setInitialState (double *xlocal=0) |
virtual void | setState (size_t point, const double *state, double *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 string &name) |
void | setBounds (size_t n, double lower, double upper) |
void | setTransientTolerances (double rtol, double atol, size_t n=npos) |
Set tolerances for time-stepping mode. More... | |
void | setSteadyTolerances (double rtol, double atol, size_t n=npos) |
Set tolerances for steady-state mode. More... | |
double | rtol (size_t n) |
Relative tolerance of the nth component. More... | |
double | 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... | |
double | upperBound (size_t n) const |
Upper bound on the nth component. More... | |
double | lowerBound (size_t n) const |
Lower bound on the nth component. More... | |
void | initTimeInteg (double dt, const double *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 () |
Set this if something has changed in the governing equations (for example, the value of a constant has been changed, so that the last-computed Jacobian is no longer valid. More... | |
size_t | index (size_t n, size_t j) const |
double | value (const double *x, size_t n, size_t j) const |
virtual void | setJac (MultiJac *jac) |
shared_ptr< SolutionArray > | toArray (bool normalize=false) const |
Save the state of this domain to a SolutionArray. More... | |
void | fromArray (const shared_ptr< SolutionArray > &arr) |
Restore the solution for this domain from a SolutionArray. More... | |
shared_ptr< Solution > | solution () const |
Return thermo/kinetics/transport manager used in the domain. More... | |
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 string &s) |
Specify an identifying tag for this domain. More... | |
string | id () const |
virtual void | show (std::ostream &s, const double *x) |
Print the solution. More... | |
double | z (size_t jlocal) const |
double | zmin () const |
double | zmax () const |
void | setProfile (const string &name, double *values, double *soln) |
vector< double > & | grid () |
const vector< double > & | grid () const |
double | grid (size_t point) const |
virtual double | 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... | |
void | setData (shared_ptr< vector< double >> &data) |
Set shared data pointer. More... | |
Public Attributes | |
double | m_zfixed = Undef |
Location of the point where temperature is fixed. More... | |
double | m_tfixed = -1.0 |
Temperature at the point used to fix the flame location. More... | |
Protected Member Functions | |
AnyMap | getMeta () const override |
Retrieve meta data. More... | |
void | setMeta (const AnyMap &state) override |
Retrieve meta data. More... | |
double | wdot (size_t k, size_t j) const |
virtual void | updateProperties (size_t jg, double *x, size_t jmin, size_t jmax) |
Update the properties (thermo, transport, and diffusion flux). More... | |
void | computeRadiation (double *x, size_t jmin, size_t jmax) |
Computes the radiative heat loss vector over points jmin to jmax and stores the data in the qdotRadiation variable. More... | |
virtual void | evalContinuity (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the continuity equation residual. More... | |
virtual void | evalMomentum (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the momentum equation residual. More... | |
virtual void | evalLambda (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the lambda equation residual. More... | |
virtual void | evalEnergy (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the energy equation residual. More... | |
virtual void | evalSpecies (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the species equations' residuals. More... | |
virtual void | evalElectricField (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) |
Evaluate the electric field equation residual to be zero everywhere. More... | |
void | updateThermo (const double *x, size_t j0, size_t j1) |
Update the thermodynamic properties from point j0 to point j1 (inclusive), based on solution x. More... | |
double | shear (const double *x, size_t j) const |
double | divHeatFlux (const double *x, size_t j) const |
size_t | mindex (size_t k, size_t j, size_t m) |
virtual void | updateDiffFluxes (const double *x, size_t j0, size_t j1) |
Update the diffusive mass fluxes. More... | |
virtual void | grad_hk (const double *x, size_t j) |
Get the gradient of species specific molar enthalpies. More... | |
virtual void | updateTransport (double *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 | |
double | T (const double *x, size_t j) const |
double & | T (double *x, size_t j) |
double | T_prev (size_t j) const |
double | rho_u (const double *x, size_t j) const |
double | u (const double *x, size_t j) const |
double | V (const double *x, size_t j) const |
double | V_prev (size_t j) const |
double | lambda (const double *x, size_t j) const |
double | Y (const double *x, size_t k, size_t j) const |
double & | Y (double *x, size_t k, size_t j) |
double | Y_prev (size_t k, size_t j) const |
double | X (const double *x, size_t k, size_t j) const |
double | flux (size_t k, size_t j) const |
convective spatial derivatives. | |
These use upwind differencing, assuming u(z) is negative | |
double | dVdz (const double *x, size_t j) const |
double | dYdz (const double *x, size_t k, size_t j) const |
double | dTdz (const double *x, size_t j) const |
Protected Attributes | |
double | m_press = -1.0 |
vector< double > | m_dz |
vector< double > | m_rho |
vector< double > | m_wtm |
vector< double > | m_wt |
vector< double > | m_cp |
vector< double > | m_visc |
vector< double > | m_tcon |
vector< double > | m_diff |
Array of size m_nsp by m_points for saving density times diffusion coefficient times species molar mass divided by mean molecular weight. More... | |
vector< double > | m_multidiff |
Array2D | m_dthermal |
Array2D | m_flux |
Array2D | m_hk |
Array of size m_nsp by m_points for saving molar enthalpies. More... | |
Array2D | m_dhk_dz |
Array of size m_nsp by m_points-1 for saving enthalpy fluxes. More... | |
Array2D | m_wdot |
size_t | m_nsp |
Number of species in the mechanism. More... | |
ThermoPhase * | m_thermo = nullptr |
Kinetics * | m_kin = nullptr |
Transport * | m_trans = nullptr |
double | m_epsilon_left = 0.0 |
double | m_epsilon_right = 0.0 |
vector< size_t > | m_kRadiating |
Indices within the ThermoPhase of the radiating species. More... | |
vector< bool > | m_do_energy |
bool | m_do_soret = false |
vector< bool > | m_do_species |
bool | m_do_multicomponent = false |
bool | m_do_radiation = false |
flag for the radiative heat loss More... | |
vector< double > | m_qdotRadiation |
radiative heat loss vector More... | |
vector< double > | m_fixedtemp |
vector< double > | m_zfix |
vector< double > | m_tfix |
size_t | m_kExcessLeft = 0 |
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 = 0 |
bool | m_dovisc |
bool | m_isFree |
bool | m_usesLambda |
Protected Attributes inherited from Domain1D | |
shared_ptr< vector< double > > | m_state |
data pointer shared from OneDim More... | |
double | m_rdt = 0.0 |
size_t | m_nv = 0 |
size_t | m_points |
Number of grid points. More... | |
vector< double > | m_slast |
vector< double > | m_max |
vector< double > | m_min |
vector< double > | m_rtol_ss |
vector< double > | m_rtol_ts |
vector< double > | m_atol_ss |
vector< double > | m_atol_ts |
vector< double > | m_z |
OneDim * | m_container = nullptr |
size_t | m_index |
size_t | m_iloc = 0 |
Starting location within the solution vector for unknowns that correspond to this domain. More... | |
size_t | m_jstart = 0 |
Domain1D * | m_left = nullptr |
Domain1D * | m_right = nullptr |
string | m_id |
Identity tag for the domain. More... | |
unique_ptr< Refiner > | m_refiner |
vector< string > | m_name |
int | m_bw = -1 |
bool | m_force_full_update = false |
shared_ptr< Solution > | m_solution |
Composite thermo/kinetics/transport handler. More... | |
Private Attributes | |
vector< double > | m_ybar |
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 19 of file StFlow.cpp.
StFlow | ( | shared_ptr< ThermoPhase > | th, |
size_t | nsp = 1 , |
||
size_t | points = 1 |
||
) |
Delegating constructor.
Definition at line 92 of file StFlow.cpp.
Create a new flow domain.
sol | Solution object used to evaluate all thermodynamic, kinetic, and transport properties |
id | name of flow domain |
points | initial number of grid points |
Definition at line 99 of file StFlow.cpp.
|
overridevirtual |
Domain type flag.
string
. Reimplemented from Domain1D.
Definition at line 124 of file StFlow.cpp.
|
overridevirtual |
called to set up initial grid, and after grid refinement
Reimplemented from Domain1D.
Definition at line 186 of file StFlow.cpp.
|
overridevirtual |
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 201 of file StFlow.cpp.
|
overridevirtual |
Set the kinetics manager.
Reimplemented from Domain1D.
Definition at line 134 of file StFlow.cpp.
|
overridevirtual |
Set transport model to existing instance.
Reimplemented from Domain1D.
Definition at line 140 of file StFlow.cpp.
void setTransportModel | ( | const string & | trans | ) |
string transportModel | ( | ) | const |
|
inline |
|
inline |
|
inline |
|
overridevirtual |
Write the initial solution estimate into array x.
Reimplemented from Domain1D.
Definition at line 220 of file StFlow.cpp.
|
overridevirtual |
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.
Definition at line 249 of file StFlow.cpp.
|
inline |
|
inline |
|
inline |
|
overridevirtual |
Name of the nth component. May be overloaded.
Reimplemented from Domain1D.
Definition at line 700 of file StFlow.cpp.
|
overridevirtual |
index of component with name name.
Reimplemented from Domain1D.
Definition at line 722 of file StFlow.cpp.
|
virtual |
Returns true if the specified component is an active part of the solver state.
Reimplemented in IonFlow.
Definition at line 745 of file StFlow.cpp.
|
overridevirtual |
|
overridevirtual |
Save the state of this domain as a SolutionArray.
soln | local solution vector for this domain |
Reimplemented from Domain1D.
Definition at line 808 of file StFlow.cpp.
|
overridevirtual |
Restore the solution for this domain from a SolutionArray.
[in] | arr | SolutionArray defining the state of this domain |
[out] | soln | Value of the solution vector, local to this domain |
Reimplemented from Domain1D.
Definition at line 842 of file StFlow.cpp.
|
inline |
|
inline |
|
inline |
|
virtual |
Get the solving stage (used by IonFlow specialization)
Reimplemented in IonFlow.
Definition at line 952 of file StFlow.cpp.
|
virtual |
Solving stage mode for handling ionized species (used by IonFlow specialization)
stage=1
: the fluxes of charged species are set to zerostage=2
: the electric field equation is solved, and the drift flux for ionized species is evaluated Reimplemented in IonFlow.
Definition at line 958 of file StFlow.cpp.
|
virtual |
Set to solve electric field in a point (used by IonFlow specialization)
Reimplemented in IonFlow.
Definition at line 964 of file StFlow.cpp.
|
virtual |
Set to fix voltage in a point (used by IonFlow specialization)
Reimplemented in IonFlow.
Definition at line 970 of file StFlow.cpp.
|
virtual |
Retrieve flag indicating whether electric field is solved or not (used by IonFlow specialization)
Reimplemented in IonFlow.
Definition at line 976 of file StFlow.cpp.
|
inline |
|
inline |
|
inline |
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 982 of file StFlow.cpp.
|
inline |
|
inline |
|
overridevirtual |
Change the grid size. Called after grid refinement.
Reimplemented from Domain1D.
Definition at line 160 of file StFlow.cpp.
void setGas | ( | const double * | x, |
size_t | j | ||
) |
Set the gas object state to be consistent with the solution at point j.
Definition at line 229 of file StFlow.cpp.
void setGasAtMidpoint | ( | const double * | 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 237 of file StFlow.cpp.
|
inline |
|
inline |
|
overridevirtual |
Evaluate the residual functions for axisymmetric stagnation flow.
If jGlobal == 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.
These residuals at all the boundary grid points are evaluated using a default boundary condition that may be modified by a boundary object that is attached to the domain. The boundary object connected will modify these equations by subtracting the boundary object's values for V, T, mdot, etc. As a result, these residual equations will force the solution variables to the values of the connected boundary object.
jGlobal | Global grid point at which to update the residual | |
[in] | xGlobal | Global state vector |
[out] | rsdGlobal | Global residual vector |
[out] | diagGlobal | Global boolean mask indicating whether each solution component has a time derivative (1) or not (0). |
[in] | rdt | Reciprocal of the timestep (rdt=0 implies steady-state.) |
Reimplemented from Domain1D.
Definition at line 296 of file StFlow.cpp.
|
inline |
|
inline |
|
overrideprotectedvirtual |
|
overrideprotectedvirtual |
|
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 334 of file StFlow.cpp.
|
protected |
Computes the radiative heat loss vector over points jmin to jmax and stores the data in the qdotRadiation variable.
The simple radiation model used was established by Liu and Rogg [21]. This model considers the radiation of CO2 and H2O.
This model uses the optically thin limit and the gray-gas approximation to simply calculate a volume specified heat flux out of the Planck absorption coefficients, the boundary emissivities and the temperature. Polynomial lines calculate the species Planck coefficients for H2O and CO2. The data for the lines are taken from the RADCAL program [8]. The coefficients for the polynomials are taken from TNF Workshop material.
Definition at line 358 of file StFlow.cpp.
|
protectedvirtual |
Evaluate the continuity equation residual.
This function calculates the residual of the continuity equation
\[ \frac{d(\rho u)}{dz} + 2\rho V = 0 \]
Axisymmetric flame: The continuity equation propagates information from right-to-left. The \( \rho u \) at point 0 is dependent on \( \rho u \) at point 1, but not on \( \dot{m} \) from the inlet.
Freely-propagating flame: The continuity equation propagates information away from a fixed temperature point that is set in the domain.
Unstrained flame: A specified mass flux; the main example being burner-stabilized flames.
The default boundary condition for the continuity equation is ( \( u = 0 \)) at the left and right boundary.
[in] | x | Local domain state vector, includes variables like temperature, density, etc. |
[out] | rsd | Local domain residual vector that stores the continuity equation residuals. |
[out] | diag | Local domain diagonal matrix that controls whether an entry has a time-derivative (used by the solver). |
[in] | rdt | Reciprocal of the timestep. |
[in] | jmin | The index for the starting point in the local domain grid. |
[in] | jmax | The index for the ending point in the local domain grid. |
Definition at line 405 of file StFlow.cpp.
|
protectedvirtual |
Evaluate the momentum equation residual.
The function calculates the radial momentum equation defined as
\[ \rho u \frac{dV}{dz} + \rho V^2 = \frac{d}{dz}\left( \mu \frac{dV}{dz} \right) - \Lambda \]
The radial momentum equation is used for axisymmetric flows, and incorporates terms for time and spatial variations of radial velocity ( \( V \)). The default boundary condition is zero radial velocity ( \( V \)) at the left and right boundary.
For argument explanation, see evalContinuity().
Definition at line 462 of file StFlow.cpp.
|
protectedvirtual |
Evaluate the lambda equation residual.
The function calculates the lambda equation as
\[ \frac{d\Lambda}{dz} = 0 \]
The lambda equation serves as an eigenvalue that allows the momentum equation and continuity equations to be simultaneously satisfied in axisymmetric flows. The lambda equation propagates information from left-to-right. The default boundary condition is \( \Lambda = 0 \) at the left and zero flux at the right boundary.
For argument explanation, see evalContinuity().
Definition at line 494 of file StFlow.cpp.
|
protectedvirtual |
Evaluate the energy equation residual.
The function calculates the energy equation:
\[ \rho c_p u \frac{dT}{dz} = \frac{d}{dz}\left( \lambda \frac{dT}{dz} \right) - \sum_k h_kW_k\dot{\omega}_k - \sum_k j_k \frac{dh_k}{dz} \]
The energy equation includes contributions from chemical reactions and diffusion. Default is zero temperature ( \( T \)) at the left and right boundaries. These boundary values are updated by the specific boundary object connected to the domain.
For argument explanation, see evalContinuity().
Definition at line 522 of file StFlow.cpp.
|
protectedvirtual |
Evaluate the species equations' residuals.
The function calculates the species equations as
\[ \rho u \frac{dY_k}{dz} + \frac{dj_k}{dz} = W_k\dot{\omega}_k \]
The species equations include terms for temporal and spatial variations of species mass fractions ( \( Y_k \)). The default boundary condition is zero flux for species at the left and right boundary.
For argument explanation, see evalContinuity().
Reimplemented in IonFlow.
Definition at line 560 of file StFlow.cpp.
|
protectedvirtual |
Evaluate the electric field equation residual to be zero everywhere.
The electric field equation is implemented in the IonFlow class. The default boundary condition is zero electric field ( \( E \)) at the boundary, and \( E \) is zero within the domain.
For argument explanation, see evalContinuity().
Reimplemented in IonFlow.
Definition at line 599 of file StFlow.cpp.
|
inlineprotected |
Update the thermodynamic properties from point j0 to point j1 (inclusive), based on solution x.
The gas state is set to be consistent with the solution at the points from j0 to j1.
Properties that are computed and cached are:
|
protectedvirtual |
Update the diffusive mass fluxes.
Reimplemented in IonFlow.
Definition at line 660 of file StFlow.cpp.
|
protectedvirtual |
Get the gradient of species specific molar enthalpies.
Definition at line 1020 of file StFlow.cpp.
|
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 608 of file StFlow.cpp.
|
protected |
|
protected |
|
protected |
|
protected |
Indices within the ThermoPhase of the radiating species.
First index is for CO2, second is for H2O.
|
protected |
|
protected |
|
protected |
double m_zfixed = Undef |
double m_tfixed = -1.0 |