Cantera  3.0.0
Loading...
Searching...
No Matches

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]

Detailed Description

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

Definition at line 44 of file StFlow.h.

Public Member Functions

 StFlow (ThermoPhase *ph=0, size_t nsp=1, size_t points=1)
 Create a new flow domain.
 
 StFlow (shared_ptr< ThermoPhase > th, size_t nsp=1, size_t points=1)
 Delegating constructor.
 
 StFlow (shared_ptr< Solution > sol, const string &id="", size_t points=1)
 Create a new flow domain.
 
string type () const override
 String indicating the domain implemented.
 
string componentName (size_t n) const override
 Name of the nth component. May be overloaded.
 
size_t componentIndex (const string &name) const override
 index of component with name name.
 
virtual bool componentActive (size_t n) const
 Returns true if the specified component is an active part of the solver state.
 
void show (const double *x) override
 Print the solution.
 
shared_ptr< SolutionArrayasArray (const double *soln) const override
 Save the state of this domain as a SolutionArray.
 
void fromArray (SolutionArray &arr, double *soln) override
 Restore the solution for this domain from a SolutionArray.
 
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.
 
void setAxisymmetricFlow ()
 Set flow configuration for axisymmetric counterflow flames, using specified inlet mass fluxes.
 
void setUnstrainedFlow ()
 Set flow configuration for burner-stabilized flames, using specified inlet mass fluxes.
 
virtual string flowType () const
 Return the type of flow domain being represented, either "Free Flame" or "Axisymmetric Stagnation".
 
void solveEnergyEqn (size_t j=npos)
 
virtual size_t getSolvingStage () const
 Get the solving stage (used by IonFlow specialization)
 
virtual void setSolvingStage (const size_t stage)
 Solving stage mode for handling ionized species (used by IonFlow specialization)
 
virtual void solveElectricField (size_t j=npos)
 Set to solve electric field in a point (used by IonFlow specialization)
 
virtual void fixElectricField (size_t j=npos)
 Set to fix voltage in a point (used by IonFlow specialization)
 
virtual bool doElectricField (size_t j) const
 Retrieve flag indicating whether electric field is solved or not (used by IonFlow specialization)
 
void enableRadiation (bool doRadiation)
 Turn radiation on / off.
 
bool radiationEnabled () const
 Returns true if the radiation term in the energy equation is enabled.
 
double radiativeHeatLoss (size_t j) const
 Return radiative heat loss at grid point j.
 
void setBoundaryEmissivities (double e_left, double e_right)
 Set the emissivities for the boundary values.
 
double leftEmissivity () const
 Return emissivity at left boundary.
 
double rightEmissivity () const
 Return emissivity at right boundary.
 
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.
 
void setGas (const double *x, size_t j)
 Set the gas object state to be consistent with the solution at point j.
 
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.
 
double density (size_t j) const
 
virtual bool fixed_mdot ()
 
bool isFree () const
 Retrieve flag indicating whether flow is freely propagating.
 
bool isStrained () const
 Retrieve flag indicating whether flow uses radial momentum.
 
void setViscosityFlag (bool dovisc)
 
void eval (size_t j, double *x, double *r, integer *mask, double rdt) override
 Evaluate the residual function for axisymmetric stagnation flow.
 
virtual void evalRightBoundary (double *x, double *res, int *diag, double rdt)
 Evaluate all residual components at the right boundary.
 
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.
 
size_t leftExcessSpecies () const
 Index of the species on the left boundary with the largest mass fraction.
 
size_t rightExcessSpecies () const
 Index of the species on the right boundary with the largest mass fraction.
 
Problem Specification
void setupGrid (size_t n, const double *z) override
 called to set up initial grid, and after grid refinement
 
void resetBadValues (double *xg) override
 When called, this function should reset "bad" values in the state vector such as negative species concentrations.
 
ThermoPhasephase ()
 
Kineticskinetics ()
 
void setThermo (ThermoPhase &th)
 Set the thermo manager.
 
void setKinetics (shared_ptr< Kinetics > kin) override
 Set the kinetics manager.
 
void setKinetics (Kinetics &kin)
 Set the kinetics manager.
 
void setTransport (shared_ptr< Transport > trans) override
 Set transport model to existing instance.
 
void setTransport (Transport &trans)
 Set transport model to existing instance.
 
void setTransportModel (const string &trans)
 Set the transport model.
 
string transportModel () const
 Retrieve transport model.
 
void enableSoret (bool withSoret)
 Enable thermal diffusion, also known as Soret diffusion.
 
bool withSoret () const
 
void setPressure (double p)
 Set the pressure.
 
double pressure () const
 The current pressure [Pa].
 
void _getInitialSoln (double *x) override
 Write the initial solution estimate into array x.
 
void _finalize (const double *x) override
 In some cases, a domain may need to set parameters that depend on the initial solution estimate.
 
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.
 
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.
 
double T_fixed (size_t j) const
 The fixed temperature value at point j.
 
- Public Member Functions inherited from Domain1D
 Domain1D (size_t nv=1, size_t points=1, double time=0.0)
 Constructor.
 
 Domain1D (const Domain1D &)=delete
 
Domain1Doperator= (const Domain1D &)=delete
 
int domainType ()
 Domain type flag.
 
virtual string type () const
 String indicating the domain implemented.
 
size_t domainIndex ()
 The left-to-right location of this domain.
 
virtual bool isConnector ()
 True if the domain is a connector domain.
 
void setSolution (shared_ptr< Solution > sol)
 Set the solution manager.
 
virtual void setKinetics (shared_ptr< Kinetics > kin)
 Set the kinetics manager.
 
virtual void setTransport (shared_ptr< Transport > trans)
 Set transport model to existing instance.
 
const OneDimcontainer () const
 The container holding this domain.
 
void setContainer (OneDim *c, size_t index)
 Specify the container object for this domain, and the position of this domain in the list.
 
void setBandwidth (int bw=-1)
 Set the Jacobian bandwidth. See the discussion of method bandwidth().
 
size_t bandwidth ()
 Set the Jacobian bandwidth for this domain.
 
virtual void init ()
 Initialize.
 
virtual void setInitialState (double *xlocal=0)
 
virtual void setState (size_t point, const double *state, double *x)
 
virtual void resetBadValues (double *xg)
 When called, this function should reset "bad" values in the state vector such as negative species concentrations.
 
virtual void resize (size_t nv, size_t np)
 Resize the domain to have nv components and np grid points.
 
Refinerrefiner ()
 Return a reference to the grid refiner.
 
size_t nComponents () const
 Number of components at each grid point.
 
void checkComponentIndex (size_t n) const
 Check that the specified component index is in range.
 
void checkComponentArraySize (size_t nn) const
 Check that an array size is at least nComponents().
 
size_t nPoints () const
 Number of grid points in this domain.
 
void checkPointIndex (size_t n) const
 Check that the specified point index is in range.
 
void checkPointArraySize (size_t nn) const
 Check that an array size is at least nPoints().
 
virtual string componentName (size_t n) const
 Name of the nth component. May be overloaded.
 
void setComponentName (size_t n, const string &name)
 
virtual size_t componentIndex (const string &name) const
 index of component with name 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.
 
void setSteadyTolerances (double rtol, double atol, size_t n=npos)
 Set tolerances for steady-state mode.
 
double rtol (size_t n)
 Relative tolerance of the nth component.
 
double atol (size_t n)
 Absolute tolerance of the nth component.
 
double steady_rtol (size_t n)
 Steady relative tolerance of the nth component.
 
double steady_atol (size_t n)
 Steady absolute tolerance of the nth component.
 
double transient_rtol (size_t n)
 Transient relative tolerance of the nth component.
 
double transient_atol (size_t n)
 Transient absolute tolerance of the nth component.
 
double upperBound (size_t n) const
 Upper bound on the nth component.
 
double lowerBound (size_t n) const
 Lower bound on the nth component.
 
void initTimeInteg (double dt, const double *x0)
 Prepare to do time stepping with time step dt.
 
void setSteadyMode ()
 Prepare to solve the steady-state problem.
 
bool steady ()
 True if in steady-state mode.
 
bool transient ()
 True if not in steady-state mode.
 
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.
 
virtual void eval (size_t j, double *x, double *r, integer *mask, double rdt=0.0)
 Evaluate the residual function at point j.
 
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)
 
AnyMap serialize (const double *soln) const
 Save the state of this domain as an AnyMap.
 
virtual shared_ptr< SolutionArrayasArray (const double *soln) const
 Save the state of this domain as a SolutionArray.
 
shared_ptr< SolutionArraytoArray (bool normalize=false) const
 Save the state of this domain to a SolutionArray.
 
void restore (const AnyMap &state, double *soln, int loglevel)
 Restore the solution for this domain from an AnyMap.
 
virtual void fromArray (SolutionArray &arr, double *soln)
 Restore the solution for this domain from a SolutionArray.
 
void fromArray (const shared_ptr< SolutionArray > &arr)
 Restore the solution for this domain from a SolutionArray.
 
shared_ptr< Solutionsolution () const
 Return thermo/kinetics/transport manager used in the domain.
 
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.
 
virtual size_t loc (size_t j=0) const
 Location of the start of the local solution vector in the global solution vector,.
 
size_t firstPoint () const
 The index of the first (that is, left-most) grid point belonging to this domain.
 
size_t lastPoint () const
 The index of the last (that is, right-most) grid point belonging to this domain.
 
void linkLeft (Domain1D *left)
 Set the left neighbor to domain 'left.
 
void linkRight (Domain1D *right)
 Set the right neighbor to domain 'right.'.
 
void append (Domain1D *right)
 Append domain 'right' to this one, and update all links.
 
Domain1Dleft () const
 Return a pointer to the left neighbor.
 
Domain1Dright () const
 Return a pointer to the right neighbor.
 
double prevSoln (size_t n, size_t j) const
 Value of component n at point j in the previous solution.
 
void setID (const string &s)
 Specify an identifying tag for this domain.
 
string id () const
 
virtual void showSolution_s (std::ostream &s, const double *x)
 
virtual void showSolution (const double *x)
 Print the solution.
 
virtual void show (std::ostream &s, const double *x)
 Print the solution.
 
virtual void show (const double *x)
 Print the solution.
 
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 void setupGrid (size_t n, const double *z)
 called to set up initial grid, and after grid refinement
 
virtual void _getInitialSoln (double *x)
 Writes some or all initial solution values into the global solution array, beginning at the location pointed to by x.
 
virtual double initialValue (size_t n, size_t j)
 Initial value of solution component n at grid point j.
 
virtual void _finalize (const double *x)
 In some cases, a domain may need to set parameters that depend on the initial solution estimate.
 
void forceFullUpdate (bool update)
 In some cases, for computational efficiency some properties (such as transport coefficients) may not be updated during Jacobian evaluations.
 
void setData (shared_ptr< vector< double > > &data)
 Set shared data pointer.
 

Public Attributes

double m_zfixed = Undef
 Location of the point where temperature is fixed.
 
double m_tfixed = -1.0
 Temperature at the point used to fix the flame location.
 

Protected Member Functions

AnyMap getMeta () const override
 Retrieve meta data.
 
void setMeta (const AnyMap &state) override
 Retrieve meta data.
 
double wdot (size_t k, size_t j) const
 
void getWdot (double *x, size_t j)
 Write the net production rates at point j into array m_wdot
 
virtual void updateProperties (size_t jg, double *x, size_t jmin, size_t jmax)
 Update the properties (thermo, transport, and diffusion flux).
 
virtual void evalResidual (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax)
 Evaluate the residual function.
 
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.
 
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.
 
virtual void grad_hk (const double *x, size_t j)
 Get the gradient of species specific molar enthalpies.
 
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.
 
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
 
virtual AnyMap getMeta () const
 Retrieve meta data.
 
virtual void setMeta (const AnyMap &meta)
 Retrieve meta data.
 

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
 
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.
 
Array2D m_dhk_dz
 Array of size m_nsp by m_points-1 for saving enthalpy fluxes.
 
Array2D m_wdot
 
size_t m_nsp
 Number of species in the mechanism.
 
ThermoPhasem_thermo = nullptr
 
Kineticsm_kin = nullptr
 
Transportm_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.
 
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
 
vector< double > m_qdotRadiation
 radiative heat loss vector
 
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.
 
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
 
double m_rdt = 0.0
 
size_t m_nv = 0
 
size_t m_points
 Number of grid points.
 
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
 
OneDimm_container = nullptr
 
size_t m_index
 
int m_type = 0
 
size_t m_iloc = 0
 Starting location within the solution vector for unknowns that correspond to this domain.
 
size_t m_jstart = 0
 
Domain1Dm_left = nullptr
 
Domain1Dm_right = nullptr
 
string m_id
 Identity tag for the domain.
 
unique_ptr< Refinerm_refiner
 
vector< string > m_name
 
int m_bw = -1
 
bool m_force_full_update = false
 
shared_ptr< Solutionm_solution
 Composite thermo/kinetics/transport handler.
 

Private Attributes

vector< double > m_ybar
 

Constructor & Destructor Documentation

◆ StFlow() [1/3]

StFlow ( ThermoPhase 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 19 of file StFlow.cpp.

◆ StFlow() [2/3]

StFlow ( shared_ptr< ThermoPhase th,
size_t  nsp = 1,
size_t  points = 1 
)

Delegating constructor.

Definition at line 93 of file StFlow.cpp.

◆ StFlow() [3/3]

StFlow ( shared_ptr< Solution sol,
const string &  id = "",
size_t  points = 1 
)

Create a new flow domain.

Parameters
solSolution object used to evaluate all thermodynamic, kinetic, and transport properties
idname of flow domain
pointsinitial number of grid points

Definition at line 100 of file StFlow.cpp.

◆ ~StFlow()

~StFlow ( )

Definition at line 121 of file StFlow.cpp.

Member Function Documentation

◆ type()

string type ( ) const
overridevirtual

String indicating the domain implemented.

Since
New in Cantera 3.0.
Todo:
Transition back to domainType after Cantera 3.0

Reimplemented from Domain1D.

Definition at line 128 of file StFlow.cpp.

◆ setupGrid()

void setupGrid ( size_t  n,
const double *  z 
)
overridevirtual

called to set up initial grid, and after grid refinement

Reimplemented from Domain1D.

Definition at line 212 of file StFlow.cpp.

◆ resetBadValues()

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

◆ phase()

ThermoPhase & phase ( )
inline

Definition at line 79 of file StFlow.h.

◆ kinetics()

Kinetics & kinetics ( )
inline

Definition at line 83 of file StFlow.h.

◆ setThermo()

void setThermo ( ThermoPhase th)

Set the thermo manager.

Deprecated:
To be removed after Cantera 3.0 (unused)

Definition at line 138 of file StFlow.cpp.

◆ setKinetics() [1/2]

void setKinetics ( shared_ptr< Kinetics kin)
overridevirtual

Set the kinetics manager.

Since
New in Cantera 3.0.

Reimplemented from Domain1D.

Definition at line 143 of file StFlow.cpp.

◆ setKinetics() [2/2]

void setKinetics ( Kinetics kin)

Set the kinetics manager.

Deprecated:
To be removed after Cantera 3.0; replaced by Domain1D::setKinetics()

Definition at line 154 of file StFlow.cpp.

◆ setTransport() [1/2]

void setTransport ( shared_ptr< Transport trans)
overridevirtual

Set transport model to existing instance.

Since
New in Cantera 3.0.

Reimplemented from Domain1D.

Definition at line 161 of file StFlow.cpp.

◆ setTransport() [2/2]

void setTransport ( Transport trans)

Set transport model to existing instance.

Deprecated:
To be removed after Cantera 3.0; replaced by Domain1D::setKinetics()

Definition at line 252 of file StFlow.cpp.

◆ setTransportModel()

void setTransportModel ( const string &  trans)

Set the transport model.

Since
New in Cantera 3.0.

Definition at line 237 of file StFlow.cpp.

◆ transportModel()

string transportModel ( ) const

Retrieve transport model.

Since
New in Cantera 3.0.

Definition at line 248 of file StFlow.cpp.

◆ enableSoret()

void enableSoret ( bool  withSoret)
inline

Enable thermal diffusion, also known as Soret diffusion.

Requires that multicomponent transport properties be enabled to carry out calculations.

Definition at line 119 of file StFlow.h.

◆ withSoret()

bool withSoret ( ) const
inline

Definition at line 122 of file StFlow.h.

◆ setPressure()

void setPressure ( double  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 128 of file StFlow.h.

◆ pressure()

double pressure ( ) const
inline

The current pressure [Pa].

Definition at line 133 of file StFlow.h.

◆ _getInitialSoln()

void _getInitialSoln ( double *  x)
overridevirtual

Write the initial solution estimate into array x.

Reimplemented from Domain1D.

Definition at line 271 of file StFlow.cpp.

◆ _finalize()

void _finalize ( const double *  x)
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 306 of file StFlow.cpp.

◆ setFixedTempProfile()

void setFixedTempProfile ( vector< double > &  zfixed,
vector< double > &  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 145 of file StFlow.h.

◆ setTemperature()

void setTemperature ( size_t  j,
double  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 154 of file StFlow.h.

◆ T_fixed()

double T_fixed ( size_t  j) const
inline

The fixed temperature value at point j.

Definition at line 160 of file StFlow.h.

◆ componentName()

string componentName ( size_t  n) const
overridevirtual

Name of the nth component. May be overloaded.

Reimplemented from Domain1D.

Definition at line 677 of file StFlow.cpp.

◆ componentIndex()

size_t componentIndex ( const string &  name) const
overridevirtual

index of component with name name.

Reimplemented from Domain1D.

Definition at line 699 of file StFlow.cpp.

◆ componentActive()

bool componentActive ( size_t  n) const
virtual

Returns true if the specified component is an active part of the solver state.

Reimplemented in IonFlow.

Definition at line 722 of file StFlow.cpp.

◆ show()

void show ( const double *  x)
overridevirtual

Print the solution.

Reimplemented from Domain1D.

Definition at line 618 of file StFlow.cpp.

◆ asArray()

shared_ptr< SolutionArray > asArray ( const double *  soln) const
overridevirtual

Save the state of this domain as a SolutionArray.

Parameters
solnlocal solution vector for this domain
Todo:
Despite the method's name, data are copied; the intent is to access data directly in future revisions, where a non-const version will be implemented.
Since
New in Cantera 3.0.

Reimplemented from Domain1D.

Definition at line 785 of file StFlow.cpp.

◆ fromArray()

void fromArray ( SolutionArray arr,
double *  soln 
)
overridevirtual

Restore the solution for this domain from a SolutionArray.

Parameters
[in]arrSolutionArray defining the state of this domain
[out]solnValue of the solution vector, local to this domain
Since
New in Cantera 3.0.

Reimplemented from Domain1D.

Definition at line 819 of file StFlow.cpp.

◆ setFreeFlow()

void setFreeFlow ( )
inline

Set flow configuration for freely-propagating flames, using an internal point with a fixed temperature as the condition to determine the inlet mass flux.

Definition at line 181 of file StFlow.h.

◆ setAxisymmetricFlow()

void setAxisymmetricFlow ( )
inline

Set flow configuration for axisymmetric counterflow flames, using specified inlet mass fluxes.

Definition at line 190 of file StFlow.h.

◆ setUnstrainedFlow()

void setUnstrainedFlow ( )
inline

Set flow configuration for burner-stabilized flames, using specified inlet mass fluxes.

Definition at line 199 of file StFlow.h.

◆ flowType()

string flowType ( ) const
virtual

Return the type of flow domain being represented, either "Free Flame" or "Axisymmetric Stagnation".

See also
setFreeFlow setAxisymmetricFlow
Deprecated:
To be removed after Cantera 3.0; replaced by type().

Definition at line 849 of file StFlow.cpp.

◆ solveEnergyEqn()

void solveEnergyEqn ( size_t  j = npos)

Definition at line 917 of file StFlow.cpp.

◆ getSolvingStage()

size_t getSolvingStage ( ) const
virtual

Get the solving stage (used by IonFlow specialization)

Since
New in Cantera 3.0

Reimplemented in IonFlow.

Definition at line 941 of file StFlow.cpp.

◆ setSolvingStage()

void setSolvingStage ( const size_t  stage)
virtual

Solving stage mode for handling ionized species (used by IonFlow specialization)

  • stage=1: the fluxes of charged species are set to zero
  • stage=2: the electric field equation is solved, and the drift flux for ionized species is evaluated

Reimplemented in IonFlow.

Definition at line 947 of file StFlow.cpp.

◆ solveElectricField()

void solveElectricField ( size_t  j = npos)
virtual

Set to solve electric field in a point (used by IonFlow specialization)

Reimplemented in IonFlow.

Definition at line 953 of file StFlow.cpp.

◆ fixElectricField()

void fixElectricField ( size_t  j = npos)
virtual

Set to fix voltage in a point (used by IonFlow specialization)

Reimplemented in IonFlow.

Definition at line 959 of file StFlow.cpp.

◆ doElectricField()

bool doElectricField ( size_t  j) const
virtual

Retrieve flag indicating whether electric field is solved or not (used by IonFlow specialization)

Reimplemented in IonFlow.

Definition at line 965 of file StFlow.cpp.

◆ enableRadiation()

void enableRadiation ( bool  doRadiation)
inline

Turn radiation on / off.

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

◆ radiationEnabled()

bool radiationEnabled ( ) const
inline

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

Definition at line 252 of file StFlow.h.

◆ radiativeHeatLoss()

double radiativeHeatLoss ( size_t  j) const
inline

Return radiative heat loss at grid point j.

Definition at line 257 of file StFlow.h.

◆ setBoundaryEmissivities()

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

◆ leftEmissivity()

double leftEmissivity ( ) const
inline

Return emissivity at left boundary.

Definition at line 270 of file StFlow.h.

◆ rightEmissivity()

double rightEmissivity ( ) const
inline

Return emissivity at right boundary.

Definition at line 273 of file StFlow.h.

◆ fixTemperature()

void fixTemperature ( size_t  j = npos)

Definition at line 985 of file StFlow.cpp.

◆ doEnergy()

bool doEnergy ( size_t  j)
inline

Definition at line 277 of file StFlow.h.

◆ resize()

void resize ( size_t  components,
size_t  points 
)
overridevirtual

Change the grid size. Called after grid refinement.

Reimplemented from Domain1D.

Definition at line 186 of file StFlow.cpp.

◆ setGas()

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

◆ setGasAtMidpoint()

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

◆ density()

double density ( size_t  j) const
inline

Definition at line 291 of file StFlow.h.

◆ fixed_mdot()

bool fixed_mdot ( )
virtual
Deprecated:
To be removed after Cantera 3.0. Superseded by isFree()

Definition at line 300 of file StFlow.cpp.

◆ isFree()

bool isFree ( ) const
inline

Retrieve flag indicating whether flow is freely propagating.

The flow is unstrained and the axial mass flow rate is not specified. For free flame propagation, the axial velocity is determined by the solver.

Since
New in Cantera 3.0

Definition at line 304 of file StFlow.h.

◆ isStrained()

bool isStrained ( ) const
inline

Retrieve flag indicating whether flow uses radial momentum.

If true, radial momentum equation for \( V \) as well as \( d\Lambda/dz = 0 \) are solved; if false, \( \Lambda(z) = 0 \) and \( V(z) = 0 \) by definition.

Since
New in Cantera 3.0

Definition at line 315 of file StFlow.h.

◆ setViscosityFlag()

void setViscosityFlag ( bool  dovisc)
inline

Definition at line 319 of file StFlow.h.

◆ eval()

void eval ( size_t  j,
double *  x,
double *  r,
integer *  mask,
double  rdt 
)
overridevirtual

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

◆ evalRightBoundary()

void evalRightBoundary ( double *  x,
double *  res,
int *  diag,
double  rdt 
)
virtual

Evaluate all residual components at the right boundary.

Definition at line 1009 of file StFlow.cpp.

◆ evalContinuity()

void evalContinuity ( size_t  j,
double *  x,
double *  r,
int *  diag,
double  rdt 
)
virtual

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

Definition at line 1039 of file StFlow.cpp.

◆ leftExcessSpecies()

size_t leftExcessSpecies ( ) const
inline

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

Definition at line 340 of file StFlow.h.

◆ rightExcessSpecies()

size_t rightExcessSpecies ( ) const
inline

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

Definition at line 345 of file StFlow.h.

◆ getMeta()

AnyMap getMeta ( ) const
overrideprotectedvirtual

Retrieve meta data.

Reimplemented from Domain1D.

Definition at line 736 of file StFlow.cpp.

◆ setMeta()

void setMeta ( const AnyMap meta)
overrideprotectedvirtual

Retrieve meta data.

Reimplemented from Domain1D.

Definition at line 861 of file StFlow.cpp.

◆ wdot()

double wdot ( size_t  k,
size_t  j 
) const
inlineprotected

Definition at line 353 of file StFlow.h.

◆ getWdot()

void getWdot ( double *  x,
size_t  j 
)
inlineprotected

Write the net production rates at point j into array m_wdot

Definition at line 358 of file StFlow.h.

◆ updateProperties()

void updateProperties ( size_t  jg,
double *  x,
size_t  jmin,
size_t  jmax 
)
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 380 of file StFlow.cpp.

◆ evalResidual()

void evalResidual ( double *  x,
double *  rsd,
int *  diag,
double  rdt,
size_t  jmin,
size_t  jmax 
)
protectedvirtual

Evaluate the residual function.

This function is called in eval after updateProperties is called.

Reimplemented in IonFlow.

Definition at line 404 of file StFlow.cpp.

◆ updateThermo()

void updateThermo ( const double *  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 377 of file StFlow.h.

◆ T() [1/2]

double T ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 390 of file StFlow.h.

◆ T() [2/2]

double & T ( double *  x,
size_t  j 
)
inlineprotected

Definition at line 393 of file StFlow.h.

◆ T_prev()

double T_prev ( size_t  j) const
inlineprotected

Definition at line 396 of file StFlow.h.

◆ rho_u()

double rho_u ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 400 of file StFlow.h.

◆ u()

double u ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 404 of file StFlow.h.

◆ V()

double V ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 408 of file StFlow.h.

◆ V_prev()

double V_prev ( size_t  j) const
inlineprotected

Definition at line 411 of file StFlow.h.

◆ lambda()

double lambda ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 415 of file StFlow.h.

◆ Y() [1/2]

double Y ( const double *  x,
size_t  k,
size_t  j 
) const
inlineprotected

Definition at line 419 of file StFlow.h.

◆ Y() [2/2]

double & Y ( double *  x,
size_t  k,
size_t  j 
)
inlineprotected

Definition at line 423 of file StFlow.h.

◆ Y_prev()

double Y_prev ( size_t  k,
size_t  j 
) const
inlineprotected

Definition at line 427 of file StFlow.h.

◆ X()

double X ( const double *  x,
size_t  k,
size_t  j 
) const
inlineprotected

Definition at line 431 of file StFlow.h.

◆ flux()

double flux ( size_t  k,
size_t  j 
) const
inlineprotected

Definition at line 435 of file StFlow.h.

◆ dVdz()

double dVdz ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 444 of file StFlow.h.

◆ dYdz()

double dYdz ( const double *  x,
size_t  k,
size_t  j 
) const
inlineprotected

Definition at line 449 of file StFlow.h.

◆ dTdz()

double dTdz ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 454 of file StFlow.h.

◆ shear()

double shear ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 460 of file StFlow.h.

◆ divHeatFlux()

double divHeatFlux ( const double *  x,
size_t  j 
) const
inlineprotected

Definition at line 466 of file StFlow.h.

◆ mindex()

size_t mindex ( size_t  k,
size_t  j,
size_t  m 
)
inlineprotected

Definition at line 472 of file StFlow.h.

◆ updateDiffFluxes()

void updateDiffFluxes ( const double *  x,
size_t  j0,
size_t  j1 
)
protectedvirtual

Update the diffusive mass fluxes.

Reimplemented in IonFlow.

Definition at line 635 of file StFlow.cpp.

◆ grad_hk()

void grad_hk ( const double *  x,
size_t  j 
)
protectedvirtual

Get the gradient of species specific molar enthalpies.

Definition at line 1077 of file StFlow.cpp.

◆ updateTransport()

void updateTransport ( double *  x,
size_t  j0,
size_t  j1 
)
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 588 of file StFlow.cpp.

Member Data Documentation

◆ m_press

double m_press = -1.0
protected

Definition at line 486 of file StFlow.h.

◆ m_dz

vector<double> m_dz
protected

Definition at line 489 of file StFlow.h.

◆ m_rho

vector<double> m_rho
protected

Definition at line 492 of file StFlow.h.

◆ m_wtm

vector<double> m_wtm
protected

Definition at line 493 of file StFlow.h.

◆ m_wt

vector<double> m_wt
protected

Definition at line 496 of file StFlow.h.

◆ m_cp

vector<double> m_cp
protected

Definition at line 497 of file StFlow.h.

◆ m_visc

vector<double> m_visc
protected

Definition at line 500 of file StFlow.h.

◆ m_tcon

vector<double> m_tcon
protected

Definition at line 501 of file StFlow.h.

◆ m_diff

vector<double> m_diff
protected

Definition at line 502 of file StFlow.h.

◆ m_multidiff

vector<double> m_multidiff
protected

Definition at line 503 of file StFlow.h.

◆ m_dthermal

Array2D m_dthermal
protected

Definition at line 504 of file StFlow.h.

◆ m_flux

Array2D m_flux
protected

Definition at line 505 of file StFlow.h.

◆ m_hk

Array2D m_hk
protected

Array of size m_nsp by m_points for saving molar enthalpies.

Definition at line 508 of file StFlow.h.

◆ m_dhk_dz

Array2D m_dhk_dz
protected

Array of size m_nsp by m_points-1 for saving enthalpy fluxes.

Definition at line 511 of file StFlow.h.

◆ m_wdot

Array2D m_wdot
protected

Definition at line 514 of file StFlow.h.

◆ m_nsp

size_t m_nsp
protected

Number of species in the mechanism.

Definition at line 516 of file StFlow.h.

◆ m_thermo

ThermoPhase* m_thermo = nullptr
protected

Definition at line 518 of file StFlow.h.

◆ m_kin

Kinetics* m_kin = nullptr
protected

Definition at line 519 of file StFlow.h.

◆ m_trans

Transport* m_trans = nullptr
protected

Definition at line 520 of file StFlow.h.

◆ m_epsilon_left

double m_epsilon_left = 0.0
protected

Definition at line 523 of file StFlow.h.

◆ m_epsilon_right

double m_epsilon_right = 0.0
protected

Definition at line 524 of file StFlow.h.

◆ m_kRadiating

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

◆ m_do_energy

vector<bool> m_do_energy
protected

Definition at line 531 of file StFlow.h.

◆ m_do_soret

bool m_do_soret = false
protected

Definition at line 532 of file StFlow.h.

◆ m_do_species

vector<bool> m_do_species
protected

Definition at line 533 of file StFlow.h.

◆ m_do_multicomponent

bool m_do_multicomponent = false
protected

Definition at line 534 of file StFlow.h.

◆ m_do_radiation

bool m_do_radiation = false
protected

flag for the radiative heat loss

Definition at line 537 of file StFlow.h.

◆ m_qdotRadiation

vector<double> m_qdotRadiation
protected

radiative heat loss vector

Definition at line 540 of file StFlow.h.

◆ m_fixedtemp

vector<double> m_fixedtemp
protected

Definition at line 543 of file StFlow.h.

◆ m_zfix

vector<double> m_zfix
protected

Definition at line 544 of file StFlow.h.

◆ m_tfix

vector<double> m_tfix
protected

Definition at line 545 of file StFlow.h.

◆ m_kExcessLeft

size_t m_kExcessLeft = 0
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 550 of file StFlow.h.

◆ m_kExcessRight

size_t m_kExcessRight = 0
protected

Definition at line 551 of file StFlow.h.

◆ m_dovisc

bool m_dovisc
protected

Definition at line 553 of file StFlow.h.

◆ m_isFree

bool m_isFree
protected

Definition at line 554 of file StFlow.h.

◆ m_usesLambda

bool m_usesLambda
protected

Definition at line 555 of file StFlow.h.

◆ m_zfixed

double m_zfixed = Undef

Location of the point where temperature is fixed.

Definition at line 563 of file StFlow.h.

◆ m_tfixed

double m_tfixed = -1.0

Temperature at the point used to fix the flame location.

Definition at line 566 of file StFlow.h.

◆ m_ybar

vector<double> m_ybar
private

Definition at line 569 of file StFlow.h.


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