Cantera  2.3.0
Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
StFlow Class Referenceabstract

This class represents 1D flow domains that satisfy the one-dimensional similarity solution for chemically-reacting, axisymmetric flows. More...

#include <StFlow.h>

Inheritance diagram for StFlow:
[legend]
Collaboration diagram for StFlow:
[legend]

Public Member Functions

 StFlow (IdealGasPhase *ph=0, size_t nsp=1, size_t points=1)
 Create a new flow domain. 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
 
virtual void showSolution (const doublereal *x)
 Print the solution. More...
 
virtual XML_Nodesave (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 std::string flowType ()
 
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...
 
void setBoundaryEmissivities (doublereal e_left, doublereal e_right)
 Set the emissivities for the boundary values. More...
 
void fixTemperature (size_t j=npos)
 
bool doEnergy (size_t j)
 
void resize (size_t components, size_t points)
 Change the grid size. Called after grid refinement. More...
 
virtual void setFixedPoint (int j0, doublereal t0)
 
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 (doublereal *x, doublereal *res, integer *diag, doublereal rdt)=0
 Evaluate all residual components at the right boundary. More...
 
virtual void evalContinuity (size_t j, doublereal *x, doublereal *r, integer *diag, doublereal rdt)=0
 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)
 
thermo_tphase ()
 
Kineticskinetics ()
 
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 setTransport (Transport &trans, bool withSoret)
 set the transport manager More...
 
void enableSoret (bool withSoret)
 
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...
 
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 OneDimcontainer () 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)
 
Refinerrefiner ()
 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 setComponentType (size_t n, int ctype)
 
size_t componentIndex (const std::string &name) const
 index of component with name name. More...
 
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 ()
 
void evalss (doublereal *x, doublereal *r, integer *mask)
 
virtual doublereal residual (doublereal *x, size_t n, size_t j)
 
int timeDerivativeFlag (size_t n)
 
void setAlgebraic (size_t n)
 
virtual void update (doublereal *x)
 
doublereal time () const
 
void incrementTime (doublereal dt)
 
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 (i.e., left-most) grid point belonging to this domain. More...
 
size_t lastPoint () const
 The index of the last (i.e., 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...
 
Domain1Dleft () const
 Return a pointer to the left neighbor. More...
 
Domain1Dright () 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
 
void setDesc (const std::string &s)
 Specify descriptive text for this domain. More...
 
const std::string & desc ()
 
virtual void getTransientMask (integer *mask)
 
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_fpgrid ()
 
const vector_fpgrid () 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 (e.g. 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...
 
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)
 
void updateDiffFluxes (const doublereal *x, size_t j0, size_t j1)
 Update the diffusive mass fluxes. More...
 
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
 
IdealGasPhasem_thermo
 
Kineticsm_kin
 
Transportm_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
 
doublereal m_time
 
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
 
OneDimm_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
 
Domain1Dm_left
 
Domain1Dm_right
 
std::string m_id
 Identity tag for the domain. More...
 
std::string m_desc
 
std::unique_ptr< Refinerm_refiner
 
vector_int m_td
 
std::vector< std::string > m_name
 
int m_bw
 
bool m_force_full_update
 

Private Attributes

vector_fp m_ybar
 

Detailed Description

This class represents 1D flow domains that satisfy the one-dimensional similarity solution for chemically-reacting, axisymmetric flows.

Definition at line 35 of file StFlow.h.

Constructor & Destructor Documentation

◆ StFlow()

StFlow ( IdealGasPhase ph = 0,
size_t  nsp = 1,
size_t  points = 1 
)

Create a new flow domain.

Parameters
phObject representing the gas phase. This object will be used to evaluate all thermodynamic, kinetic, and transport properties.
nspNumber of species.
pointsInitial number of grid points

Definition at line 16 of file StFlow.cpp.

Member Function Documentation

◆ setupGrid()

void setupGrid ( size_t  n,
const doublereal *  z 
)
virtual

called to set up initial grid, and after grid refinement

Reimplemented from Domain1D.

Definition at line 120 of file StFlow.cpp.

References StFlow::resize().

◆ resetBadValues()

void resetBadValues ( double *  xg)
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 135 of file StFlow.cpp.

References Phase::getMassFractions(), Domain1D::loc(), and Phase::setMassFractions().

◆ setThermo()

void setThermo ( IdealGasPhase th)
inline

Set the thermo manager.

Note that the flow equations assume the ideal gas equation.

Definition at line 68 of file StFlow.h.

◆ setKinetics()

void setKinetics ( Kinetics kin)
inline

Set the kinetics manager. The kinetics manager must.

Definition at line 73 of file StFlow.h.

◆ setTransport() [1/2]

void setTransport ( Transport trans)

set the transport manager

Definition at line 166 of file StFlow.cpp.

References Array2D::resize(), and Transport::transportType().

◆ setTransport() [2/2]

void setTransport ( Transport trans,
bool  withSoret 
)

set the transport manager

Deprecated:
The withSoret argument is deprecated and unused. Use the form of setTransport with signature setTransport(Transport& trans). To be removed after Cantera 2.3.

Definition at line 144 of file StFlow.cpp.

References Array2D::resize(), Transport::transportType(), and Cantera::warn_deprecated().

◆ setPressure()

void setPressure ( doublereal  p)
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 95 of file StFlow.h.

Referenced by StFlow::restore().

◆ pressure()

doublereal pressure ( ) const
inline

The current pressure [Pa].

Definition at line 100 of file StFlow.h.

◆ _getInitialSoln()

void _getInitialSoln ( double *  x)
virtual

Write the initial solution estimate into array x.

Reimplemented from Domain1D.

Definition at line 189 of file StFlow.cpp.

References Phase::getMassFractions(), and Phase::temperature().

◆ _finalize()

void _finalize ( const doublereal *  x)
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 FreeFlame.

Definition at line 217 of file StFlow.cpp.

References Cantera::linearInterp().

Referenced by FreeFlame::_finalize().

◆ setFixedTempProfile()

void setFixedTempProfile ( vector_fp zfixed,
vector_fp tfixed 
)
inline

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 112 of file StFlow.h.

◆ setTemperature()

void setTemperature ( size_t  j,
doublereal  t 
)
inline

Set the temperature fixed point at grid point j, and disable the energy equation so that the solution will be held to this value.

Definition at line 121 of file StFlow.h.

◆ T_fixed()

doublereal T_fixed ( size_t  j) const
inline

The fixed temperature value at point j.

Definition at line 127 of file StFlow.h.

Referenced by StFlow::eval().

◆ componentName()

string componentName ( size_t  n) const
virtual

Name of the nth component. May be overloaded.

Reimplemented from Domain1D.

Definition at line 556 of file StFlow.cpp.

References Phase::speciesName().

◆ showSolution()

void showSolution ( const doublereal *  x)
virtual

Print the solution.

Reimplemented from Domain1D.

Definition at line 497 of file StFlow.cpp.

References StFlow::m_do_radiation, StFlow::m_qdotRadiation, Domain1D::showSolution(), and Cantera::writelog().

◆ save()

XML_Node & save ( XML_Node o,
const doublereal *const  sol 
)
virtual

Save the current solution for this domain into an XML_Node.

Parameters
oXML_Node to save the solution to.
solCurrent value of the solution vector. The object will pick out which part of the solution vector pertains to this object.

Reimplemented from Domain1D.

Reimplemented in FreeFlame.

Definition at line 758 of file StFlow.cpp.

References XML_Node::addAttribute(), Cantera::addFloat(), Cantera::addFloatArray(), Cantera::addNamedFloatArray(), Cantera::addString(), Array2D::getRow(), Domain1D::loc(), StFlow::m_do_radiation, StFlow::m_qdotRadiation, Array2D::nColumns(), Domain1D::nPoints(), Domain1D::refiner(), Domain1D::save(), and Phase::speciesName().

Referenced by FreeFlame::save().

◆ restore()

void restore ( const XML_Node dom,
doublereal *  soln,
int  loglevel 
)
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.

Parameters
domXML_Node for this domain
solnCurrent value of the solution vector, local to this object.
loglevel0 to suppress all output; 1 to show warnings; 2 for verbose output

Reimplemented from Domain1D.

Reimplemented in FreeFlame.

Definition at line 596 of file StFlow.cpp.

References XML_Node::getChildren(), Cantera::getFloat(), Phase::nSpecies(), Domain1D::restore(), StFlow::setPressure(), XML_Node::value(), and Cantera::writelog().

Referenced by FreeFlame::restore().

◆ enableRadiation()

void enableRadiation ( bool  doRadiation)
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 166 of file StFlow.h.

References StFlow::m_do_radiation.

◆ radiationEnabled()

bool radiationEnabled ( ) const
inline

Returns true if the radiation term in the energy equation is enabled.

Definition at line 171 of file StFlow.h.

References StFlow::m_do_radiation.

◆ setBoundaryEmissivities()

void setBoundaryEmissivities ( doublereal  e_left,
doublereal  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 840 of file StFlow.cpp.

◆ resize()

void resize ( size_t  components,
size_t  points 
)
virtual

Change the grid size. Called after grid refinement.

Reimplemented from Domain1D.

Definition at line 96 of file StFlow.cpp.

References StFlow::m_qdotRadiation, Array2D::resize(), and Domain1D::resize().

Referenced by StFlow::setupGrid().

◆ setGas()

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 197 of file StFlow.cpp.

References Phase::setMassFractions_NoNorm(), IdealGasPhase::setPressure(), and Phase::setTemperature().

Referenced by StFlow::eval(), ReactingSurf1D::eval(), StFlow::getWdot(), and StFlow::updateThermo().

◆ setGasAtMidpoint()

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 205 of file StFlow.cpp.

References Phase::setMassFractions_NoNorm(), IdealGasPhase::setPressure(), and Phase::setTemperature().

Referenced by StFlow::updateTransport().

◆ eval()

void eval ( size_t  j,
doublereal *  x,
doublereal *  r,
integer *  mask,
doublereal  rdt 
)
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 235 of file StFlow.cpp.

References IdealGasPhase::cp_R_ref(), IdealGasPhase::enthalpy_RT_ref(), StFlow::evalContinuity(), StFlow::evalRightBoundary(), Domain1D::firstPoint(), Cantera::GasConstant, StFlow::getWdot(), Domain1D::lastPoint(), StFlow::leftExcessSpecies(), Domain1D::loc(), StFlow::m_do_radiation, StFlow::m_kExcessLeft, StFlow::m_kRadiating, StFlow::m_qdotRadiation, Cantera::npos, Cantera::OneAtm, StFlow::setGas(), Cantera::StefanBoltz, StFlow::T_fixed(), StFlow::updateDiffFluxes(), StFlow::updateThermo(), and StFlow::updateTransport().

◆ evalRightBoundary()

virtual void evalRightBoundary ( doublereal *  x,
doublereal *  res,
integer *  diag,
doublereal  rdt 
)
pure virtual

Evaluate all residual components at the right boundary.

Implemented in FreeFlame, and AxiStagnFlow.

Referenced by StFlow::eval().

◆ evalContinuity()

virtual void evalContinuity ( size_t  j,
doublereal *  x,
doublereal *  r,
integer *  diag,
doublereal  rdt 
)
pure virtual

Evaluate the residual corresponding to the continuity equation at all interior grid points.

Implemented in FreeFlame, and AxiStagnFlow.

Referenced by StFlow::eval().

◆ leftExcessSpecies()

size_t leftExcessSpecies ( ) const
inline

Index of the species on the left boundary with the largest mass fraction.

Definition at line 232 of file StFlow.h.

References StFlow::m_kExcessLeft.

Referenced by StFlow::eval().

◆ rightExcessSpecies()

size_t rightExcessSpecies ( ) const
inline

Index of the species on the right boundary with the largest mass fraction.

Definition at line 237 of file StFlow.h.

Referenced by Outlet1D::eval(), OutletRes1D::eval(), ReactingSurf1D::eval(), AxiStagnFlow::evalRightBoundary(), and FreeFlame::evalRightBoundary().

◆ getWdot()

void getWdot ( doublereal *  x,
size_t  j 
)
inlineprotected

Write the net production rates at point j into array m_wdot

Definition at line 247 of file StFlow.h.

References Kinetics::getNetProductionRates(), and StFlow::setGas().

Referenced by StFlow::eval().

◆ updateThermo()

void updateThermo ( const doublereal *  x,
size_t  j0,
size_t  j1 
)
inlineprotected

Update the thermodynamic properties from point j0 to point j1 (inclusive), based on solution x.

Definition at line 256 of file StFlow.h.

References ThermoPhase::cp_mass(), Phase::density(), Phase::meanMolecularWeight(), and StFlow::setGas().

Referenced by StFlow::eval().

◆ updateDiffFluxes()

void updateDiffFluxes ( const doublereal *  x,
size_t  j0,
size_t  j1 
)
protected

Update the diffusive mass fluxes.

Definition at line 514 of file StFlow.cpp.

Referenced by StFlow::eval().

◆ updateTransport()

void updateTransport ( doublereal *  x,
size_t  j0,
size_t  j1 
)
protected

Update the transport properties at grid points in the range from j0 to j1, based on solution x.

Definition at line 467 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::eval().

Member Data Documentation

◆ m_kRadiating

std::vector<size_t> m_kRadiating
protected

Indices within the ThermoPhase of the radiating species.

First index is for CO2, second is for H2O.

Definition at line 396 of file StFlow.h.

Referenced by StFlow::eval().

◆ m_do_radiation

bool m_do_radiation
protected

flag for the radiative heat loss

Definition at line 405 of file StFlow.h.

Referenced by StFlow::enableRadiation(), StFlow::eval(), StFlow::radiationEnabled(), StFlow::save(), and StFlow::showSolution().

◆ m_qdotRadiation

vector_fp m_qdotRadiation
protected

radiative heat loss vector

Definition at line 408 of file StFlow.h.

Referenced by StFlow::eval(), StFlow::resize(), StFlow::save(), and StFlow::showSolution().

◆ m_kExcessLeft

size_t m_kExcessLeft
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 418 of file StFlow.h.

Referenced by StFlow::eval(), and StFlow::leftExcessSpecies().


The documentation for this class was generated from the following files: