Cantera  3.1.0a1

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

#include <IonFlow.h>

Inheritance diagram for IonFlow:
[legend]

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. See Pedersen and Brown [30] for details.

Definition at line 28 of file IonFlow.h.

Public Member Functions

 IonFlow (ThermoPhase *ph=0, size_t nsp=1, size_t points=1)
 
 IonFlow (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...
 
size_t getSolvingStage () const override
 Get the solving stage (used by IonFlow specialization) More...
 
void setSolvingStage (const size_t stage) override
 Solving stage mode for handling ionized species (used by IonFlow specialization) More...
 
void resize (size_t components, size_t points) override
 Resize the domain to have nv components and np grid points. More...
 
bool componentActive (size_t n) const override
 Returns true if the specified component is an active part of the solver state. 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 solveElectricField (size_t j=npos) override
 Set to solve electric field in a point (used by IonFlow specialization) More...
 
void fixElectricField (size_t j=npos) override
 Set to fix voltage in a point (used by IonFlow specialization) More...
 
bool doElectricField (size_t j) const override
 Retrieve flag indicating whether electric field is solved or not (used by IonFlow specialization) More...
 
void setElectronTransport (vector< double > &tfix, vector< double > &diff_e, vector< double > &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...
 
 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...
 
void show (const double *x) override
 Print the solution. More...
 
shared_ptr< SolutionArrayasArray (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)
 
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...
 
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...
 
ThermoPhasephase ()
 
Kineticskinetics ()
 
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
 
Domain1Doperator= (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 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 ()
 Initialize. More...
 
virtual void setInitialState (double *xlocal=0)
 
virtual void setState (size_t point, const double *state, double *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 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< SolutionArraytoArray (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< Solutionsolution () 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...
 
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 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...
 

Protected Member Functions

void evalElectricField (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) override
 Evaluate the electric field equation residual by Gauss's law. More...
 
void evalSpecies (double *x, double *rsd, int *diag, double rdt, size_t jmin, size_t jmax) override
 Evaluate the species equations' residual. More...
 
void updateTransport (double *x, size_t j0, size_t j1) override
 Update the transport properties at grid points in the range from j0 to j1, based on solution x. More...
 
void updateDiffFluxes (const double *x, size_t j0, size_t j1) override
 Update the diffusive mass fluxes. More...
 
void frozenIonMethod (const double *x, size_t j0, size_t j1)
 Solving phase one: the fluxes of charged species are turned off. More...
 
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
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...
 
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 grad_hk (const double *x, size_t j)
 Get the gradient of species specific molar enthalpies. More...
 
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
 
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

vector< bool > m_do_electric_field
 flag for solving electric field or not More...
 
bool m_import_electron_transport = false
 flag for importing transport of electron More...
 
vector< double > m_speciesCharge
 electrical properties More...
 
vector< size_t > m_kCharge
 index of species with charges More...
 
vector< size_t > m_kNeutral
 index of neutral species More...
 
vector< double > m_mobi_e_fix
 coefficients of polynomial fitting of fixed electron transport profile More...
 
vector< double > m_diff_e_fix
 
vector< double > m_mobility
 mobility More...
 
size_t m_stage = 1
 solving stage More...
 
size_t m_kElectron = npos
 index of electron More...
 
- Protected Attributes inherited from StFlow
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...
 
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. 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
 
OneDimm_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
 
Domain1Dm_left = nullptr
 
Domain1Dm_right = nullptr
 
string m_id
 Identity tag for the domain. More...
 
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. More...
 

Additional Inherited Members

- Public Attributes inherited from StFlow
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...
 

Constructor & Destructor Documentation

◆ IonFlow()

IonFlow ( 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 56 of file IonFlow.cpp.

Member Function Documentation

◆ domainType()

string domainType ( ) const
overridevirtual

Domain type flag.

Since
Starting in Cantera 3.1, the return type is a string.

Reimplemented from Domain1D.

Definition at line 74 of file IonFlow.cpp.

◆ getSolvingStage()

size_t getSolvingStage ( ) const
inlineoverridevirtual

Get the solving stage (used by IonFlow specialization)

Since
New in Cantera 3.0

Reimplemented from StFlow.

Definition at line 42 of file IonFlow.h.

◆ setSolvingStage()

void setSolvingStage ( const size_t  stage)
overridevirtual

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 from StFlow.

Definition at line 189 of file IonFlow.cpp.

◆ resize()

void resize ( size_t  nv,
size_t  np 
)
overridevirtual

Resize the domain to have nv components and np grid points.

This method is virtual so that subclasses can perform other actions required to resize the domain.

Reimplemented from Domain1D.

Definition at line 84 of file IonFlow.cpp.

◆ componentActive()

bool componentActive ( size_t  n) const
overridevirtual

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

Reimplemented from StFlow.

Definition at line 91 of file IonFlow.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 312 of file IonFlow.cpp.

◆ solveElectricField()

void solveElectricField ( size_t  j = npos)
overridevirtual

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

Reimplemented from StFlow.

Definition at line 245 of file IonFlow.cpp.

◆ fixElectricField()

void fixElectricField ( size_t  j = npos)
overridevirtual

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

Reimplemented from StFlow.

Definition at line 270 of file IonFlow.cpp.

◆ doElectricField()

bool doElectricField ( size_t  j) const
inlineoverridevirtual

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

Reimplemented from StFlow.

Definition at line 54 of file IonFlow.h.

◆ setElectronTransport()

void setElectronTransport ( vector< double > &  tfix,
vector< double > &  diff_e,
vector< double > &  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).

See Bisetti and El Morsli [2]. If in the future the class GasTransport is improved, this method may be discarded. This method specifies this profile.

Definition at line 295 of file IonFlow.cpp.

◆ evalElectricField()

void evalElectricField ( double *  x,
double *  rsd,
int *  diag,
double  rdt,
size_t  jmin,
size_t  jmax 
)
overrideprotectedvirtual

Evaluate the electric field equation residual by Gauss's law.

Evaluate the electric field equation residual.

The function calculates the electric field equation as:

\[ \frac{dE}{dz} = \frac{e}{\varepsilon_0} \sum (q_k \cdot n_k) \]

and

\[ E = -\frac{dV}{dz} \]

The electric field equation is based on Gauss's law, accounting for charge density and permittivity of free space ( \( \varepsilon_0 \)). The zero electric field is first evaluated and if the solution state is 2, then the alternative form the electric field equation is evaluated.

For argument explanation, see evalContinuity() base class.

Reimplemented from StFlow.

Definition at line 202 of file IonFlow.cpp.

◆ evalSpecies()

void evalSpecies ( double *  x,
double *  rsd,
int *  diag,
double  rdt,
size_t  jmin,
size_t  jmax 
)
overrideprotectedvirtual

Evaluate the species equations' residual.

This function overloads the original species function.

A Neumann boundary for the charged species at the left boundary is added, and the default boundary condition from the overloaded method is left the same for the right boundary.

For argument explanation, see evalContinuity() base class.

Reimplemented from StFlow.

Definition at line 227 of file IonFlow.cpp.

◆ updateTransport()

void updateTransport ( double *  x,
size_t  j0,
size_t  j1 
)
overrideprotectedvirtual

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

Reimplemented from StFlow.

Definition at line 100 of file IonFlow.cpp.

◆ updateDiffFluxes()

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

Update the diffusive mass fluxes.

Reimplemented from StFlow.

Definition at line 117 of file IonFlow.cpp.

◆ frozenIonMethod()

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

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

Definition at line 127 of file IonFlow.cpp.

◆ electricFieldMethod()

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

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

Definition at line 152 of file IonFlow.cpp.

◆ E()

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

electric field

Definition at line 144 of file IonFlow.h.

◆ ND()

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

number density

Definition at line 153 of file IonFlow.h.

◆ rho_e()

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

total charge density

Definition at line 158 of file IonFlow.h.

Member Data Documentation

◆ m_do_electric_field

vector<bool> m_do_electric_field
protected

flag for solving electric field or not

Definition at line 116 of file IonFlow.h.

◆ m_import_electron_transport

bool m_import_electron_transport = false
protected

flag for importing transport of electron

Definition at line 119 of file IonFlow.h.

◆ m_speciesCharge

vector<double> m_speciesCharge
protected

electrical properties

Definition at line 122 of file IonFlow.h.

◆ m_kCharge

vector<size_t> m_kCharge
protected

index of species with charges

Definition at line 125 of file IonFlow.h.

◆ m_kNeutral

vector<size_t> m_kNeutral
protected

index of neutral species

Definition at line 128 of file IonFlow.h.

◆ m_mobi_e_fix

vector<double> m_mobi_e_fix
protected

coefficients of polynomial fitting of fixed electron transport profile

Definition at line 131 of file IonFlow.h.

◆ m_mobility

vector<double> m_mobility
protected

mobility

Definition at line 135 of file IonFlow.h.

◆ m_stage

size_t m_stage = 1
protected

solving stage

Definition at line 138 of file IonFlow.h.

◆ m_kElectron

size_t m_kElectron = npos
protected

index of electron

Definition at line 141 of file IonFlow.h.


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