chombo-discharge
|
Base class for interaction between Kinetic Monte Carlo and Ito-based plasma solvers. More...
#include <CD_ItoKMCPhysics.H>
Public Member Functions | |
ItoKMCPhysics () noexcept | |
Constructor. Does nothing. | |
virtual | ~ItoKMCPhysics () noexcept |
Destructor. Does nothing. | |
void | defineKMC () const noexcept |
Define the KMC solver and state. | |
void | killKMC () const noexcept |
Kill the KMC solver. | |
virtual Real | computeDt (const RealVect a_E, const RealVect a_pos, const Vector< FPR > a_numParticles) const noexcept |
Compute a time step to use. More... | |
virtual Real | computeAlpha (const Real a_E, const RealVect a_x) const =0 |
Compute Townsend ionization coefficient. More... | |
virtual Real | computeEta (const Real a_E, const RealVect a_x) const =0 |
Compute Townsend attachment coefficient. More... | |
const Vector< RefCountedPtr< ItoSpecies > > & | getItoSpecies () const |
Get all particle drift-diffusion species. More... | |
const Vector< RefCountedPtr< CdrSpecies > > & | getCdrSpecies () const |
Get all fluid drift-diffusion species. More... | |
const Vector< RefCountedPtr< RtSpecies > > & | getRtSpecies () const |
Get all photon species. More... | |
virtual Vector< std::string > | getPlotVariableNames () const noexcept |
Get number of plot variables. | |
virtual Vector< Real > | getPlotVariables (const RealVect a_E, const RealVect a_pos, const Vector< Real > &a_phi, const Vector< RealVect > &a_gradPhi, const Real a_dx, const Real a_kappa) const noexcept |
Get plot variables. More... | |
virtual int | getNumberOfPlotVariables () const noexcept |
Get number of plot variables. | |
int | getNumItoSpecies () const |
Return number of Ito solvers. | |
int | getNumCdrSpecies () const |
Return number of CDR solvers. | |
int | getNumPlasmaSpecies () const |
Return total number of plasma species. | |
int | getNumPhotonSpecies () const |
Return number of RTE solvers. | |
virtual bool | needGradients () const noexcept |
Return true/false if physics model needs species gradients. | |
const std::map< int, std::pair< SpeciesType, int > > & | getSpeciesMap () const noexcept |
Get the internal mapping between plasma species and Ito solvers. | |
virtual void | parseRuntimeOptions () noexcept |
Parse run-time options. | |
virtual Real | initialSigma (const Real a_time, const RealVect a_pos) const |
Set initial surface charge. Default is 0, override if you want. More... | |
virtual Vector< Real > | computeMobilities (const Real a_time, const RealVect a_pos, const RealVect a_E) const noexcept=0 |
Compute the Ito solver mobilities. More... | |
virtual Vector< Real > | computeDiffusionCoefficients (const Real a_time, const RealVect a_pos, const RealVect a_E) const noexcept=0 |
Compute the Ito solver diffusion coefficients. More... | |
virtual void | secondaryEmissionEB (Vector< List< ItoParticle >> &a_secondaryParticles, Vector< Real > &a_cdrFluxes, Vector< List< Photon >> &a_secondaryPhotons, const Vector< List< ItoParticle >> &a_primaryParticles, const Vector< Real > &a_cdrFluxesExtrap, const Vector< List< Photon >> &a_primaryPhotons, const RealVect &a_E, const RealVect &a_physicalCellCenter, 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=0 |
Resolve secondary emission at the EB. More... | |
void | advanceKMC (Vector< FPR > &a_numParticles, Vector< FPR > &a_numNewPhotons, const Vector< Real > &a_phi, const Vector< RealVect > &a_gradPhi, const Real a_dt, const RealVect a_E, const RealVect a_pos, 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_itoParticles, Vector< List< PointParticle > * > &a_cdrParticles, const Vector< List< Photon > * > &a_absorbedPhotons) const noexcept |
Reconcile photoionization reactions. More... | |
Protected Types | |
enum class | Algorithm { SSA , TauPlain , TauMidpoint , HybridPlain , HybridMidpoint } |
Enum for switching between KMC algorithms. More... | |
enum class | ParticlePlacement { Random , Centroid } |
Enum for switching between various particle placement algorithms. | |
Protected Member Functions | |
void | define () noexcept |
Define method – defines all the internal machinery. | |
void | defineSpeciesMap () noexcept |
Build internal representation of how we distinguish the Ito and CDR solvers. More... | |
void | definePhotoPathways () noexcept |
Define pathways for photo-reactions. | |
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. | |
virtual void | updateReactionRates (std::vector< std::shared_ptr< const KMCReaction >> &a_kmcReactions, const RealVect a_E, const RealVect a_pos, const Vector< Real > &a_phi, const Vector< RealVect > &a_gradPhi, const Real a_dx, const Real a_kappa) const noexcept=0 |
Update reaction rates. More... | |
void | removeParticles (List< ItoParticle > &a_particles, const long long a_numToRemove) const |
Remove particles from the input list. More... | |
Protected Attributes | |
Algorithm | m_algorithm |
Algorithm to use for KMC advance. | |
ParticlePlacement | m_particlePlacement |
Particle placement algorithm. | |
std::map< int, std::pair< SpeciesType, int > > | m_speciesMap |
Map for associating a plasma species with an Ito solver or CDR solver. | |
std::string | m_className |
Class name. Used for options parsing. | |
bool | m_debug |
Turn on/off debugging. | |
bool | m_isDefined |
Is defined or not. | |
std::vector< KMCReaction > | m_kmcReactions |
List of reactions for the KMC solver. | |
std::vector< ItoKMCPhotoReaction > | m_photoReactions |
List of photoionization reactions. | |
std::map< int, std::pair< std::discrete_distribution< int >, std::map< int, int > > > | m_photoPathways |
Random number generators for photoionization pathways. More... | |
ItoKMCSurfaceReactionSet | m_surfaceReactions |
Surface reactions. | |
Vector< RefCountedPtr< ItoSpecies > > | m_itoSpecies |
List of solver-tracked particle drift-diffusion species. | |
Vector< RefCountedPtr< CdrSpecies > > | m_cdrSpecies |
List of solver-tracked fluid drift-diffusion 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... | |
Static Protected Attributes | |
static thread_local bool | m_hasKMCSolver |
Is the KMC solver defined or not. | |
static thread_local KMCSolverType | m_kmcSolver |
Kinetic Monte Carlo solver used in advanceReactionNetwork. | |
static thread_local KMCState | m_kmcState |
KMC state used in advanceReactionNetwork. | |
static thread_local std::vector< std::shared_ptr< const KMCReaction > > | m_kmcReactionsThreadLocal |
KMC reactions used in advanceReactionNetowkr. More... | |
Base class for interaction between Kinetic Monte Carlo and Ito-based plasma solvers.
When using this class, the user should populate m_kmcReactions with the appropriate reactions AND implement routines for computing the mobility/diffusion coefficeints.
|
strongprotected |
Enum for switching between KMC algorithms.
'SSA' is the Gillespie algorithm, 'Tau' is tau-leaping and 'Hybrid' is the Cao et. al. algorithm.
|
inline |
Advance particles.
[in,out] | a_numParticles | Number of physical particles |
[out] | a_numNewPhotons | Number of new physical photons to generate (of each type) |
[in] | a_phi | Plasma species densities. |
[in] | a_gradPhi | Plasma species density gradients. |
[in] | a_dt | Time step |
[in] | a_E | Electric field |
[in] | a_pos | Physical position |
[in] | a_dx | Grid resolution |
[in] | a_kappa | Cut-cell volume fraction. |
|
pure virtual |
Compute Townsend ionization coefficient.
[in] | a_E | Electric field magnitude |
[in] | a_x | Physical coordinate |
Implemented in Physics::ItoKMC::ItoKMCJSON.
|
pure virtualnoexcept |
Compute the Ito solver diffusion coefficients.
[in] | a_time | Time |
[in] | a_pos | Position |
[in] | a_E | Electric field |
Implemented in Physics::ItoKMC::ItoKMCJSON.
|
inlinevirtualnoexcept |
Compute a time step to use.
[in] | a_E | Electric field. |
[in] | a_pos | Position |
[in] | a_numParticles | Number of reactive particles of each species |
Reimplemented in Physics::ItoKMC::ItoKMCJSON.
|
pure virtual |
Compute Townsend attachment coefficient.
[in] | a_E | Electric field magnitude |
[in] | a_x | Physical coordinate |
Implemented in Physics::ItoKMC::ItoKMCJSON.
|
pure virtualnoexcept |
Compute the Ito solver mobilities.
[in] | a_time | Time |
[in] | a_pos | Position |
[in] | a_E | Electric field |
Implemented in Physics::ItoKMC::ItoKMCJSON.
|
inlineprotectednoexcept |
Build internal representation of how we distinguish the Ito and CDR solvers.
|
inline |
Get all fluid drift-diffusion species.
|
inline |
Get all particle drift-diffusion species.
|
virtualnoexcept |
Get plot variables.
[in] | a_E | Electric field |
[in] | a_pos | Physical position |
[in] | a_phi | Plasma species densities |
[in] | a_gradPhi | Density gradients for plasma species. |
[in] | a_dx | Grid resolution |
[in] | a_kappa | Cut-cell volume fraction |
Reimplemented in Physics::ItoKMC::ItoKMCJSON.
|
inline |
Get all photon species.
|
inlinevirtual |
Set initial surface charge. Default is 0, override if you want.
[in] | a_time | Simulation time |
[in] | a_pos | Physical coordinate |
|
inlinenoexcept |
Reconcile the number of particles.
This will add/remove particles and potentially also adjust the particle weights.
[in,out] | a_particles | Computational particles |
[in] | a_newNumParticles | New number of particles (i.e., after the KMC advance) |
[in] | a_oldNumParticles | Previous number of particles (i.e., before the KMC advance) |
[in] | a_cellPos | Cell center position |
[in] | a_centroidPos | Cell centroid position |
[in] | a_loCorner | Low corner of minimum box enclosing the cut-cell |
[in] | a_hiCorner | High corner of minimum box enclosing the cut-cell |
[in] | a_bndryCentroid | Cut-cell boundary centroid |
[in] | a_bndryNormal | Cut-cell normal (pointing into the domain) |
[in] | a_dx | Grid resolution |
[in] | a_kappa | Cut-cell volume fraction. |
|
inlinenoexcept |
Reconcile photoionization reactions.
[in,out] | a_itoPrticles | Particle products placed in Ito solvers. |
[in,out] | a_cdrParticles | Particle products placed in CDR solvers. |
[in] | a_absorbedPhotons | Photons absorbed on the mesh. |
This runs through the photo-reactions and associates photo-ionization products.
|
inlinenoexcept |
Generate new photons.
This will add photons
[in] | a_newPhotons | New photons |
[in] | a_numNewPhotons | Number of physical photons to be generated. |
[in] | a_cellPos | Cell center position |
[in] | a_centroidPos | Cell centroid position |
[in] | a_loCorner | Low corner of minimum box enclosing the cut-cell |
[in] | a_hiCorner | High corner of minimum box enclosing the cut-cell |
[in] | a_bndryCentroid | Cut-cell boundary centroid |
[in] | a_bndryNormal | Cut-cell normal (pointing into the domain) |
[in] | a_dx | Grid resolution |
[in] | a_kappa | Cut-cell volume fraction. |
|
inlineprotected |
Remove particles from the input list.
This will remove weight from the input particles if we can. Otherwise we remove full particles.
[in,out] | a_particles | List of (super-)particles to remove from. |
[in] | a_numToRemove | Number of physical particles to remove from the input list |
|
pure virtualnoexcept |
Resolve secondary emission at the EB.
Routine is here to handle charge injection, secondary emission etc.
[out] | a_secondaryParticles | Outgoing plasma species particles. |
[out] | a_cdrFluxes | CDR fluxes for CDR species. |
[out] | a_secondaryPhotons | Photons injected through the EB. |
[in] | a_primaryParticles | Particles that left the computational domain through the EB. |
[in] | a_cdrFluxesExtrap | Extrapolated CDR fluxes. |
[in] | a_primaryPhotons | Photons that left the computational domain through the EB. |
[in] | a_newNumParticles | Total number of particles in the cut-cell AFTER the transport step. |
[in] | a_oldNumParticles | Total number of particles in the cut-cell BEFORE the transport step. |
[in] | a_electricField | Electric field. |
[in] | a_physicalCellCenter | Physical position of the cell center. |
[in] | a_cellCentroid | Cell centroid relative to the cell center (not multiplied by dx). |
[in] | a_bndryCentroid | EB face centroid relative to the cell center (not multiplied by dx). |
[in] | a_bndryNormal | Cut-cell normal vector. |
[in] | a_bndryArea | Cut-cell boundary area - not multiplied by dx (2D) or dx^2 (3D). |
[in] | a_dx | Grid resolution on this level. |
[in] | a_dt | Time step. |
[in] | a_isDielectric | Dielectric or electrode. |
[in] | a_matIndex | Material index (taken from computationalGeometry). |
Implemented in Physics::ItoKMC::ItoKMCJSON.
|
protectedpure virtualnoexcept |
Update reaction rates.
[out] | a_kmcReactions | Reaction rates to be set. |
[in] | a_E | Electric field |
[in] | a_pos | Physical position |
[in] | a_phi | Plasma species densities |
[in] | a_gradPhi | Density gradients for plasma species. |
[in] | a_dx | Grid resolution |
[in] | a_kappa | Cut-cell volume fraction |
Implemented in Physics::ItoKMC::ItoKMCJSON.
|
protected |
Solver setting for the Cao et. al. algorithm.
Equal to the maximum permitted change in the relative propensity for non-critical reactions
|
staticprotected |
KMC reactions used in advanceReactionNetowkr.
|
protected |
Solver setting for the Cao et. al algorithm.
Determines critical reactions. A reaction is critical if it is m_Ncrit firings away from depleting a reactant.
|
protected |
Solver setting for the Cao et. al algorithm.
Maximum number of SSA steps to run when switching into SSA-based advancement for non-critical reactions.
|
protected |
Random number generators for photoionization pathways.
The first index is the photon index, i.e. entry in m_photonSpecies. The second index in the map is an RNG generator for selecting photo-reactions, and a map which associates the returned reaction from the RNG generator with an index in m_photoReactions.
|
protected |
Solver setting for the Cao et. al. algorithm.
Equal to the maximum permitted change in the relative propensity for non-critical reactions