chombo-discharge
Classes | Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Physics::ItoPlasma::ItoPlasmaAir3LFA Class Reference

Hard-coded implementation of ItoPlasmaPhysics for simplified discharges in STP air. @detail Transport data for the electrons are taken from BOLSIG+. More...

#include <CD_ItoPlasmaAir3LFA.H>

Inheritance diagram for Physics::ItoPlasma::ItoPlasmaAir3LFA:
Inheritance graph
[legend]
Collaboration diagram for Physics::ItoPlasma::ItoPlasmaAir3LFA:
Collaboration graph
[legend]

Classes

class  Electron
 Electron species. More...
 
class  Negative
 
class  PhotonZ
 
class  Positive
 Positive ion species. More...
 

Public Member Functions

 ItoPlasmaAir3LFA ()
 Full constructor. Initializes species and sets up reactions.
 
virtual ~ItoPlasmaAir3LFA ()
 Destructor. Does nothing of importance.
 
virtual void parseRuntimeOptions () noexcept override
 Parse run-time options.
 
virtual Real computeDt (const RealVect a_E, const RealVect a_pos, const Vector< Real > a_densities) const noexcept override
 Compute a physics-baseded time step. More...
 
virtual Real computeAlpha (const RealVect a_E) const override
 Compute Townsend ionization coefficient. More...
 
virtual void updateReactionRates (const RealVect a_E, const Real a_dx, const Real a_kappa) const noexcept override
 Update reaction rates. More...
 
Vector< Real > computeItoMobilities (const Real a_time, const RealVect a_pos, const RealVect a_E) const noexcept override
 Compute the Ito solver mobilities. More...
 
Vector< Real > computeItoDiffusion (const Real a_time, const RealVect a_pos, const RealVect a_E) const noexcept override
 Compute the Ito solver diffusion coefficients. More...
 
void injectParticlesEB (Vector< List< ItoParticle >> &a_outgoingParticles, Vector< List< Photon >> &a_outgoingPhotons, const Vector< List< ItoParticle >> &a_incomingParticles, const Vector< List< Photon >> &a_incomingPhotons, const Vector< FPR > &a_newNumParticles, const Vector< FPR > &a_oldNumParticles, const RealVect &a_E, const RealVect &a_cellCenter, const RealVect &a_cellCentroid, const RealVect &a_bndryCentroid, const RealVect &a_bndryNormal, const Real a_bndryArea, const Real a_dx, const Real a_dt, const bool a_isDielectric, const int a_matIndex) const noexcept override
 Resolve particle and photon injection at the EB. More...
 
- Public Member Functions inherited from Physics::ItoPlasma::ItoPlasmaPhysics
 ItoPlasmaPhysics ()
 Constructor. Does nothing.
 
virtual ~ItoPlasmaPhysics ()
 Destructor. Does nothing.
 
const Vector< RefCountedPtr< ItoSpecies > > & getItoSpecies () const
 Get all particle species. More...
 
const Vector< RefCountedPtr< RtSpecies > > & getRtSpecies () const
 Get all photon species. More...
 
int getNumPlasmaSpecies () const
 Return number of ion equations.
 
int getNumPhotonSpecies () const
 Return number of RTE equations.
 
virtual Real initialSigma (const Real a_time, const RealVect a_pos) const
 Set initial surface charge. Default is 0, override if you want. More...
 
void advanceKMC (Vector< FPR > &a_numParticles, Vector< FPR > &a_numNewPhotons, const Real a_dt, const RealVect a_E, const Real a_dx, const Real a_kappa) const
 Advance particles. More...
 
void reconcileParticles (Vector< List< ItoParticle > * > &a_particles, const Vector< FPR > &a_newNumParticles, const Vector< FPR > &a_oldNumParticles, const RealVect a_cellPos, const RealVect a_centroidPos, const RealVect a_lo, const RealVect a_hi, const RealVect a_bndryCentroid, const RealVect a_bndryNormal, const Real a_dx, const Real a_kappa) const noexcept
 Reconcile the number of particles. More...
 
void reconcilePhotons (Vector< List< Photon > * > &a_newPhotons, const Vector< FPR > &a_numNewPhotons, const RealVect a_cellPos, const RealVect a_centroidPos, const RealVect a_lo, const RealVect a_hi, const RealVect a_bndryCentroid, const RealVect a_bndryNormal, const Real a_dx, const Real a_kappa) const noexcept
 Generate new photons. More...
 
void reconcilePhotoionization (Vector< List< ItoParticle > * > &a_particles, const Vector< List< Photon > * > &a_absorbedPhotons) const noexcept
 Reconcile photoionization reactions. More...
 

Protected Member Functions

virtual Real excitationRates (const Real a_E) const
 Excitation rate for photoionization term.
 
virtual Real sergeyFactor (const Real a_O2frac) const
 Sergey factor for photoionization term.
 
virtual void addTable (const std::string a_tableName, const std::string a_file) noexcept
 Add a table to m_tables. More...
 
virtual void readTables () noexcept
 Read various tabulated data into the lookup tables in this class.
 
virtual void parseDx () noexcept
 Parse the m_dX variable used for time step restriction.
 
virtual void parseTransport () noexcept
 Parse some auxiliary transport settings.
 
- Protected Member Functions inherited from Physics::ItoPlasma::ItoPlasmaPhysics
void defineKMC () noexcept
 Define the KMC solver and state.
 
void parsePPC () noexcept
 Parse the maximum number of particles generated per cell.
 
void parseDebug () noexcept
 Parse the maximum number of particles generated per cell.
 
void parseAlgorithm () noexcept
 Parse reaction algorithm.
 
void removeParticles (List< ItoParticle > &a_particles, const long long a_numToRemove) const
 Remove particles from the input list. More...
 

Protected Attributes

bool m_extrapBC
 Extrapolate BC or not.
 
Real m_deltaX
 User-defined parameter for time step restriction.
 
int m_numPlasmaSpecies
 Number of plasma species.
 
int m_numRtSpecies
 Number of radiative transfer species.
 
Real m_ionMobility
 Ion mobility.
 
Real m_ionDiffCo
 Ion diffusion coefficient.
 
Real m_ionImpactEfficiency
 Ion impact efficiency.
 
Real m_quantumEfficiency
 Quantum efficiency.
 
Real m_N
 Neutral number density.
 
Real m_p
 Gas pressure (SI)
 
Real m_pq
 Quenching pressure (SI)
 
Real m_T
 Gas temperature.
 
Real m_O2frac
 O2 content.
 
Real m_N2frac
 N2 content.
 
Real m_photoIonizationFactor
 Photo-ionization factor.
 
std::map< std::string, LookupTable1D< 2 > > m_tables
 Lookup tables holding mobilities, rate coefficients etc.
 
- Protected Attributes inherited from Physics::ItoPlasma::ItoPlasmaPhysics
Algorithm m_algorithm
 Algorithm to use for KMC advance.
 
ParticlePlacement m_particlePlacement
 Particle placement algorithm.
 
std::string m_className
 Class name. Used for options parsing.
 
bool m_debug
 Turn on/off debugging.
 
KMCState m_kmcState
 Kinetic Monte Carlo state.
 
KMCSolverType m_kmcSolver
 Kinetic Monte Carlo solver.
 
std::vector< std::shared_ptr< const KMCReaction > > m_kmcReactions
 List of reactions for the KMC solver.
 
std::vector< std::shared_ptr< const ItoPlasmaPhotoReaction > > m_photoReactions
 List of photoionization reactions.
 
Vector< RefCountedPtr< ItoSpecies > > m_plasmaSpecies
 List of solver-tracked particle species.
 
Vector< RefCountedPtr< RtSpecies > > m_rtSpecies
 List of solver-tracked photon species.
 
int m_maxNewParticles
 Maximum new number of particles generated by the chemistry advance.
 
int m_maxNewPhotons
 Maximum new number of photons generated by the chemistry advance.
 
int m_Ncrit
 Solver setting for the Cao et. al algorithm. More...
 
int m_NSSA
 Solver setting for the Cao et. al algorithm. More...
 
Real m_SSAlim
 Solver setting for the Cao et. al. algorithm. More...
 
Real m_eps
 Solver setting for the Cao et. al. algorithm. More...
 

Additional Inherited Members

- Protected Types inherited from Physics::ItoPlasma::ItoPlasmaPhysics
enum class  Algorithm { SSA , Tau , Hybrid }
 Enum for switching between KMC algorithms. More...
 
enum class  ParticlePlacement { Random , Centroid }
 Enum for switching between various particle placement algorithms.
 

Detailed Description

Hard-coded implementation of ItoPlasmaPhysics for simplified discharges in STP air. @detail Transport data for the electrons are taken from BOLSIG+.

Member Function Documentation

◆ addTable()

void ItoPlasmaAir3LFA::addTable ( const std::string  a_tableName,
const std::string  a_file 
)
protectedvirtualnoexcept

Add a table to m_tables.

Parameters
[in]a_tableNameTable identifier for m_table
[in]a_fileFile. Must be in column x/y format

◆ computeAlpha()

Real ItoPlasmaAir3LFA::computeAlpha ( const RealVect  a_E) const
overridevirtual

Compute Townsend ionization coefficient.

Parameters
[in]a_EElectric field.

Implements Physics::ItoPlasma::ItoPlasmaPhysics.

◆ computeDt()

Real ItoPlasmaAir3LFA::computeDt ( const RealVect  a_E,
const RealVect  a_pos,
const Vector< Real >  a_densities 
) const
overridevirtualnoexcept

Compute a physics-baseded time step.

Parameters
[in]a_EElectric field.
[in]a_posPosition
[in]a_numParticlesNumber of particles per cell
Returns
Returns m_deltaX/k where k is Townsend ionization coefficient and m_deltaX is user-specified.

Reimplemented from Physics::ItoPlasma::ItoPlasmaPhysics.

◆ computeItoDiffusion()

Vector< Real > ItoPlasmaAir3LFA::computeItoDiffusion ( const Real  a_time,
const RealVect  a_pos,
const RealVect  a_E 
) const
overridevirtualnoexcept

Compute the Ito solver diffusion coefficients.

Parameters
[in]a_timeTime
[in]a_posPosition
[in]a_EElectric field

Implements Physics::ItoPlasma::ItoPlasmaPhysics.

◆ computeItoMobilities()

Vector< Real > ItoPlasmaAir3LFA::computeItoMobilities ( const Real  a_time,
const RealVect  a_pos,
const RealVect  a_E 
) const
overridevirtualnoexcept

Compute the Ito solver mobilities.

Parameters
[in]a_timeTime
[in]a_posPosition
[in]a_EElectric field

Implements Physics::ItoPlasma::ItoPlasmaPhysics.

◆ injectParticlesEB()

void ItoPlasmaAir3LFA::injectParticlesEB ( Vector< List< ItoParticle >> &  a_outgoingParticles,
Vector< List< Photon >> &  a_outgoingPhotons,
const Vector< List< ItoParticle >> &  a_incomingParticles,
const Vector< List< Photon >> &  a_incomingPhotons,
const Vector< FPR > &  a_newNumParticles,
const Vector< FPR > &  a_oldNumParticles,
const RealVect &  a_E,
const RealVect &  a_cellCenter,
const RealVect &  a_cellCentroid,
const RealVect &  a_bndryCentroid,
const RealVect &  a_bndryNormal,
const Real  a_bndryArea,
const Real  a_dx,
const Real  a_dt,
const bool  a_isDielectric,
const int  a_matIndex 
) const
overridevirtualnoexcept

Resolve particle and photon injection at the EB.

Routine is here to handle charge injection, secondary emission etc.

Parameters
[out]a_outgoingParticlesOutgoing plasma species particles.
[out]a_outgoingPhotonsPhotons injected through the EB
[in]a_incomingParticlesParticles that left the computational domain through the EB
[in]a_incomingPhotonsPhotons that left the computational domain through the EB
[in]a_electricFieldElectric field
[in]a_cellCenterPhysical position of the cell center.
[in]a_cellCentroidCell centroid relative to the cell center (not multiplied by dx)
[in]a_bndryCentroidEB face centroid relative to the cell center (not multiplied by dx)
[in]a_bndryNormalCut-cell normal vector.
[in]a_bndryAreaCut-cell boundary area - not multiplied by dx (2D) or dx^2 (3D)
[in]a_dxGrid resolution on this level.
[in]a_dtTime step
[in]a_isDielectricDielectric or electrode.
[in]a_matIndexMaterial index (taken from computationalGeometry)

Implements Physics::ItoPlasma::ItoPlasmaPhysics.

◆ updateReactionRates()

void ItoPlasmaAir3LFA::updateReactionRates ( const RealVect  a_E,
const Real  a_dx,
const Real  a_kappa 
) const
overridevirtualnoexcept

Update reaction rates.

Parameters
[in]a_EElectric field
[in]a_dxGrid resolution
[in]a_kappaCut-cell volume fraction

Implements Physics::ItoPlasma::ItoPlasmaPhysics.


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