Cantera 2.6.0
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
IonFlow Class Reference

This class models the ion transportation in a flame. More...

#include <IonFlow.h>

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

Public Member Functions

 IonFlow (IdealGasPhase *ph=0, size_t nsp=1, size_t points=1)
 
virtual void setSolvingStage (const size_t phase)
 set the solving stage More...
 
virtual void resize (size_t components, size_t points)
 Change the grid size. Called after grid refinement. More...
 
virtual bool componentActive (size_t n) const
 Returns true if the specified component is an active part of the solver state. More...
 
virtual void _finalize (const double *x)
 In some cases, a domain may need to set parameters that depend on the initial solution estimate. More...
 
void solveElectricField (size_t j=npos)
 set to solve electric field on a point More...
 
void fixElectricField (size_t j=npos)
 set to fix voltage on a point More...
 
bool doElectricField (size_t j)
 
void setElectronTransport (vector_fp &tfix, vector_fp &diff_e, vector_fp &mobi_e)
 Sometimes it is desired to carry out the simulation using a specified electron transport profile, rather than assuming it as a constant (0.4). More...
 
- Public Member Functions inherited from StFlow
 StFlow (ThermoPhase *ph=0, size_t nsp=1, size_t points=1)
 Create a new flow domain. More...
 
 StFlow (shared_ptr< ThermoPhase > th, size_t nsp=1, size_t points=1)
 Delegating constructor. More...
 
virtual std::string componentName (size_t n) const
 Name of the nth component. May be overloaded. More...
 
virtual size_t componentIndex (const std::string &name) const
 index of component with name name. More...
 
virtual 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 AnyMap serialize (const double *soln) const
 Save the state of this domain as an AnyMap. More...
 
virtual void restore (const AnyMap &state, double *soln, int loglevel)
 Restore the solution for this domain from an AnyMap. More...
 
void setFreeFlow ()
 Set flow configuration for freely-propagating flames, using an internal point with a fixed temperature as the condition to determine the inlet mass flux. More...
 
void setAxisymmetricFlow ()
 Set flow configuration for axisymmetric counterflow or burner-stabilized flames, using specified inlet mass fluxes. More...
 
virtual std::string flowType () const
 Return the type of flow domain being represented, either "Free Flame" or "Axisymmetric Stagnation". More...
 
void solveEnergyEqn (size_t j=npos)
 
void enableRadiation (bool doRadiation)
 Turn radiation on / off. More...
 
bool radiationEnabled () const
 Returns true if the radiation term in the energy equation is enabled. More...
 
double radiativeHeatLoss (size_t j) const
 Return radiative heat loss at grid point j. More...
 
void setBoundaryEmissivities (double e_left, double e_right)
 Set the emissivities for the boundary values. More...
 
double leftEmissivity () const
 Return emissivitiy at left boundary. More...
 
double rightEmissivity () const
 Return emissivitiy at right boundary. More...
 
void fixTemperature (size_t j=npos)
 
bool doEnergy (size_t j)
 
void setGas (const doublereal *x, size_t j)
 Set the gas object state to be consistent with the solution at point j. More...
 
void setGasAtMidpoint (const doublereal *x, size_t j)
 Set the gas state to be consistent with the solution at the midpoint between j and j + 1. More...
 
doublereal density (size_t j) const
 
virtual bool fixed_mdot ()
 
void setViscosityFlag (bool dovisc)
 
virtual void eval (size_t j, doublereal *x, doublereal *r, integer *mask, doublereal rdt)
 
virtual void evalRightBoundary (double *x, double *res, int *diag, double rdt)
 Evaluate all residual components at the right boundary. More...
 
virtual void evalContinuity (size_t j, double *x, double *r, int *diag, double rdt)
 Evaluate the residual corresponding to the continuity equation at all interior grid points. More...
 
size_t leftExcessSpecies () const
 Index of the species on the left boundary with the largest mass fraction. More...
 
size_t rightExcessSpecies () const
 Index of the species on the right boundary with the largest mass fraction. More...
 
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)
 
ThermoPhasephase ()
 
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 enableSoret (bool withSoret)
 Enable thermal diffusion, also known as Soret diffusion. More...
 
bool withSoret () const
 
void setPressure (doublereal p)
 Set the pressure. More...
 
doublereal pressure () const
 The current pressure [Pa]. More...
 
virtual void _getInitialSoln (double *x)
 Write the initial solution estimate into array x. More...
 
void setFixedTempProfile (vector_fp &zfixed, vector_fp &tfixed)
 Sometimes it is desired to carry out the simulation using a specified temperature profile, rather than computing it by solving the energy equation. More...
 
void setTemperature (size_t j, doublereal t)
 
doublereal T_fixed (size_t j) const
 The fixed temperature value at point j. More...
 
- Public Member Functions inherited from Domain1D
 Domain1D (size_t nv=1, size_t points=1, double time=0.0)
 Constructor. More...
 
 Domain1D (const Domain1D &)=delete
 
Domain1Doperator= (const Domain1D &)=delete
 
int domainType ()
 Domain type flag. More...
 
size_t domainIndex ()
 The left-to-right location of this domain. More...
 
bool isConnector ()
 True if the domain is a connector domain. More...
 
const 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 setBounds (size_t n, doublereal lower, doublereal upper)
 
void setTransientTolerances (doublereal rtol, doublereal atol, size_t n=npos)
 Set tolerances for time-stepping mode. More...
 
void setSteadyTolerances (doublereal rtol, doublereal atol, size_t n=npos)
 Set tolerances for steady-state mode. More...
 
doublereal rtol (size_t n)
 Relative tolerance of the nth component. More...
 
doublereal atol (size_t n)
 Absolute tolerance of the nth component. More...
 
double steady_rtol (size_t n)
 Steady relative tolerance of the nth component. More...
 
double steady_atol (size_t n)
 Steady absolute tolerance of the nth component. More...
 
double transient_rtol (size_t n)
 Transient relative tolerance of the nth component. More...
 
double transient_atol (size_t n)
 Transient absolute tolerance of the nth component. More...
 
doublereal upperBound (size_t n) const
 Upper bound on the nth component. More...
 
doublereal lowerBound (size_t n) const
 Lower bound on the nth component. More...
 
void initTimeInteg (doublereal dt, const doublereal *x0)
 Prepare to do time stepping with time step dt. More...
 
void setSteadyMode ()
 Prepare to solve the steady-state problem. More...
 
bool steady ()
 True if in steady-state mode. More...
 
bool transient ()
 True if not in steady-state mode. More...
 
void needJacUpdate ()
 
size_t index (size_t n, size_t j) const
 
doublereal value (const doublereal *x, size_t n, size_t j) const
 
virtual void setJac (MultiJac *jac)
 
size_t size () const
 
void locate ()
 Find the index of the first grid point in this domain, and the start of its variables in the global solution vector. More...
 
virtual size_t loc (size_t j=0) const
 Location of the start of the local solution vector in the global solution vector,. More...
 
size_t firstPoint () const
 The index of the first (that is, left-most) grid point belonging to this domain. More...
 
size_t lastPoint () const
 The index of the last (that is, right-most) grid point belonging to this domain. More...
 
void linkLeft (Domain1D *left)
 Set the left neighbor to domain 'left. More...
 
void linkRight (Domain1D *right)
 Set the right neighbor to domain 'right.'. More...
 
void append (Domain1D *right)
 Append domain 'right' to this one, and update all links. More...
 
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
 
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 (such as transport coefficients) may not be updated during Jacobian evaluations. More...
 

Protected Member Functions

virtual void evalResidual (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax)
 
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...
 
virtual void updateDiffFluxes (const double *x, size_t j0, size_t j1)
 Update the diffusive mass fluxes. More...
 
virtual void frozenIonMethod (const double *x, size_t j0, size_t j1)
 Solving phase one: the fluxes of charged species are turned off. More...
 
virtual void electricFieldMethod (const double *x, size_t j0, size_t j1)
 Solving phase two: the electric field equation is added coupled by the electrical drift. More...
 
double E (const double *x, size_t j) const
 electric field More...
 
double dEdz (const double *x, size_t j) const
 
double ND (const double *x, size_t k, size_t j) const
 number density More...
 
double rho_e (double *x, size_t j) const
 total charge density More...
 
- Protected Member Functions inherited from StFlow
doublereal wdot (size_t k, size_t j) const
 
void getWdot (doublereal *x, size_t j)
 Write the net production rates at point j into array m_wdot More...
 
virtual void updateProperties (size_t jg, double *x, size_t jmin, size_t jmax)
 Update the properties (thermo, transport, and diffusion flux). More...
 
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)
 
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
 
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

std::vector< bool > m_do_electric_field
 flag for solving electric field or not More...
 
bool m_import_electron_transport
 flag for importing transport of electron More...
 
vector_fp m_speciesCharge
 electrical properties More...
 
std::vector< size_t > m_kCharge
 index of species with charges More...
 
std::vector< size_t > m_kNeutral
 index of neutral species More...
 
vector_fp m_mobi_e_fix
 coefficients of polynomial fitting of fixed electron transport profile More...
 
vector_fp m_diff_e_fix
 
vector_fp m_mobility
 mobility More...
 
size_t m_stage
 solving stage More...
 
size_t m_kElectron
 index of electron More...
 
- Protected Attributes inherited from StFlow
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
 
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::unique_ptr< Refinerm_refiner
 
std::vector< std::string > m_name
 
int m_bw
 
bool m_force_full_update
 

Additional Inherited Members

- Public Attributes inherited from StFlow
double m_zfixed
 Location of the point where temperature is fixed. More...
 
double m_tfixed
 Temperature at the point used to fix the flame location. More...
 

Detailed Description

This class models the ion transportation in a flame.

There are three stages of the simulation.

The first stage turns off the diffusion of ions due to the fast diffusion rate of electron without internal electric forces (ambi- polar diffusion effect).

The second stage evaluates drift flux from electric field calculated from Poisson's equation, which is solved together with other equations. Poisson's equation is coupled because the total charge densities depends on the species' concentration. Reference: Pederson, Timothy, and R. C. Brown. "Simulation of electric field effects in premixed methane flames." Combustion and Flames 94.4(1993): 433-448.

Definition at line 31 of file IonFlow.h.

Constructor & Destructor Documentation

◆ IonFlow()

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

Definition at line 20 of file IonFlow.cpp.

Member Function Documentation

◆ setSolvingStage()

void setSolvingStage ( const size_t  phase)
virtual

set the solving stage

Definition at line 169 of file IonFlow.cpp.

◆ resize()

void resize ( size_t  components,
size_t  points 
)
virtual

Change the grid size. Called after grid refinement.

Reimplemented from StFlow.

Definition at line 61 of file IonFlow.cpp.

◆ componentActive()

bool componentActive ( size_t  n) const
virtual

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

Reimplemented from StFlow.

Definition at line 68 of file IonFlow.cpp.

◆ _finalize()

void _finalize ( const double *  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 StFlow.

Definition at line 283 of file IonFlow.cpp.

◆ solveElectricField()

void solveElectricField ( size_t  j = npos)

set to solve electric field on a point

Definition at line 216 of file IonFlow.cpp.

References Cantera::npos.

◆ fixElectricField()

void fixElectricField ( size_t  j = npos)

set to fix voltage on a point

Definition at line 241 of file IonFlow.cpp.

References Cantera::npos.

◆ doElectricField()

bool doElectricField ( size_t  j)
inline

Definition at line 46 of file IonFlow.h.

◆ setElectronTransport()

void setElectronTransport ( vector_fp tfix,
vector_fp diff_e,
vector_fp mobi_e 
)

Sometimes it is desired to carry out the simulation using a specified electron transport profile, rather than assuming it as a constant (0.4).

Reference: Bisetti, Fabrizio, and Mbark El Morsli. "Calculation and analysis of the mobility and diffusion coefficient of thermal electrons in methane/air premixed flames." Combustion and flame 159.12 (2012): 3518-3521. If in the future the class GasTranport is improved, this method may be discard. This method specifies this profile.

Definition at line 266 of file IonFlow.cpp.

References Cantera::polyfit().

◆ evalResidual()

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

This function overloads the original function. The residual function of electric field is added.

Reimplemented from StFlow.

Definition at line 181 of file IonFlow.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 from StFlow.

Definition at line 77 of file IonFlow.cpp.

References Cantera::poly5().

◆ updateDiffFluxes()

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

Update the diffusive mass fluxes.

Reimplemented from StFlow.

Definition at line 92 of file IonFlow.cpp.

◆ frozenIonMethod()

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

Solving phase one: the fluxes of charged species are turned off.

Definition at line 102 of file IonFlow.cpp.

◆ electricFieldMethod()

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

Solving phase two: the electric field equation is added coupled by the electrical drift.

Definition at line 129 of file IonFlow.cpp.

◆ E()

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

electric field

Definition at line 108 of file IonFlow.h.

◆ dEdz()

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

Definition at line 112 of file IonFlow.h.

◆ ND()

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

number density

Definition at line 117 of file IonFlow.h.

References Cantera::Avogadro.

Referenced by IonFlow::rho_e().

◆ rho_e()

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

total charge density

Definition at line 122 of file IonFlow.h.

References Cantera::ElectronCharge, IonFlow::m_kCharge, IonFlow::m_speciesCharge, and IonFlow::ND().

Member Data Documentation

◆ m_do_electric_field

std::vector<bool> m_do_electric_field
protected

flag for solving electric field or not

Definition at line 80 of file IonFlow.h.

◆ m_import_electron_transport

bool m_import_electron_transport
protected

flag for importing transport of electron

Definition at line 83 of file IonFlow.h.

◆ m_speciesCharge

vector_fp m_speciesCharge
protected

electrical properties

Definition at line 86 of file IonFlow.h.

Referenced by IonFlow::rho_e().

◆ m_kCharge

std::vector<size_t> m_kCharge
protected

index of species with charges

Definition at line 89 of file IonFlow.h.

Referenced by IonFlow::rho_e().

◆ m_kNeutral

std::vector<size_t> m_kNeutral
protected

index of neutral species

Definition at line 92 of file IonFlow.h.

◆ m_mobi_e_fix

vector_fp m_mobi_e_fix
protected

coefficients of polynomial fitting of fixed electron transport profile

Definition at line 95 of file IonFlow.h.

◆ m_diff_e_fix

vector_fp m_diff_e_fix
protected

Definition at line 96 of file IonFlow.h.

◆ m_mobility

vector_fp m_mobility
protected

mobility

Definition at line 99 of file IonFlow.h.

◆ m_stage

size_t m_stage
protected

solving stage

Definition at line 102 of file IonFlow.h.

◆ m_kElectron

size_t m_kElectron
protected

index of electron

Definition at line 105 of file IonFlow.h.


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