9#ifndef CT_PLASMAPHASE_H
10#define CT_PLASMAPHASE_H
13#include "cantera/numerics/eigen_sparse.h"
89 explicit PlasmaPhase(
const std::string& inputFile=
"",
90 const std::string&
id=
"");
92 virtual std::string
type()
const {
117 const double* distrb,
135 return m_isotropicShapeFactor;
182 virtual bool addSpecies(shared_ptr<Species> spec);
237 double m_isotropicShapeFactor;
ThermoPhase object for the ideal gas equation of state - workhorse for Cantera (see Thermodynamic Pro...
A map of string keys to values whose type can vary at runtime.
Class IdealGasPhase represents low-density gases that obey the ideal gas equation of state.
Base class for a phase with plasma properties.
void checkElectronEnergyDistribution() const
Check the electron energy distribution.
double meanElectronEnergy() const
Mean electron energy [eV].
size_t m_nPoints
Number of points of electron energy levels.
virtual bool addSpecies(shared_ptr< Species > spec)
virtual void setParameters(const AnyMap &phaseNode, const AnyMap &rootNode=AnyMap())
Set equation of state parameters from an AnyMap phase description.
virtual void getParameters(AnyMap &phaseNode) const
Store the parameters of a ThermoPhase object such that an identical one could be reconstructed using ...
void normalizeElectronEnergyDistribution()
Electron energy distribution norm.
void setQuadratureMethod(const std::string &method)
Set numerical quadrature method for intergating electron energy distribution function.
std::string m_quadratureMethod
Numerical quadrature method for electron energy distribution.
virtual void setElectronTemperature(double Te)
Set the internally stored electron temperature of the phase (K).
PlasmaPhase(const std::string &inputFile="", const std::string &id="")
Construct and initialize a PlasmaPhase object directly from an input file.
size_t nElectronEnergyLevels() const
Number of electron levels.
Eigen::ArrayXd m_electronEnergyDist
Normalized electron energy distribution vector [-] Length: m_nPoints.
Eigen::ArrayXd m_electronEnergyLevels
electron energy levels [ev]. Length: m_nPoints
void setDiscretizedElectronEnergyDist(const double *levels, const double *distrb, size_t length)
Set discretized electron energy distribution.
void checkElectronEnergyLevels() const
Check the electron energy levels.
void updateElectronTemperatureFromEnergyDist()
Update electron temperature (K) From energy distribution.
void setElectronEnergyDistributionType(const std::string &type)
Set electron energy distribution type.
void updateElectronEnergyDistribution()
Update electron energy distribution.
std::string quadratureMethod() const
Numerical quadrature method. Method: m_quadratureMethod.
virtual double electronTemperature() const
Electron Temperature (K)
virtual void updateThermo() const
virtual void initThermo()
Initialize the ThermoPhase object after all species have been set up.
double m_electronTemp
Electron temperature [K].
void getElectronEnergyDistribution(double *distrb) const
Get electron energy distribution.
void setElectronEnergyLevels(const double *levels, size_t length)
Set electron energy levels.
bool m_do_normalizeElectronEnergyDist
Flag of normalizing electron energy distribution.
virtual std::string type() const
String indicating the thermodynamic model implemented.
bool normalizeElectronEnergyDistEnabled() const
Flag of automatically normalize electron energy distribution.
std::string m_distributionType
Electron energy distribution type.
void setIsotropicElectronEnergyDistribution()
Set isotropic electron energy distribution.
double isotropicShapeFactor() const
The shape factor of isotropic electron energy distribution.
void getElectronEnergyLevels(double *levels) const
Get electron energy levels.
std::string electronEnergyDistributionType() const
Get electron energy distribution type.
void setMeanElectronEnergy(double energy)
Set mean electron energy [eV].
size_t m_electronSpeciesIndex
Index of electron species.
void setIsotropicShapeFactor(double x)
Set the shape factor of isotropic electron energy distribution.
void enableNormalizeElectronEnergyDist(bool enable)
Set flag of automatically normalize electron energy distribution Flag: m_do_normalizeElectronEnergyDi...
Namespace for the Cantera kernel.
const double Boltzmann
Boltzmann constant [J/K].
const double ElectronCharge
Elementary charge [C].