chombo-discharge
Public Member Functions | Protected Types | Protected Member Functions | Protected Attributes | List of all members
Physics::ItoPlasma::ItoPlasmaPhysics Class Referenceabstract

Base class for interaction between Kinetic Monte Carlo and Ito-based plasma solvers. More...

#include <CD_ItoPlasmaPhysics.H>

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

Public Member Functions

 ItoPlasmaPhysics ()
 Constructor. Does nothing.
 
virtual ~ItoPlasmaPhysics ()
 Destructor. Does nothing.
 
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 RealVect a_E) const =0
 Compute Townsend ionization coefficient. More...
 
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 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 > computeItoMobilities (const Real a_time, const RealVect a_pos, const RealVect a_E) const noexcept=0
 Compute the Ito solver mobilities. More...
 
virtual Vector< Real > computeItoDiffusion (const Real a_time, const RealVect a_pos, const RealVect a_E) const noexcept=0
 Compute the Ito solver diffusion coefficients. More...
 
virtual 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=0
 Resolve particle and photon injection at the EB. 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 Types

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.
 

Protected Member Functions

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.
 
virtual void updateReactionRates (const RealVect a_E, 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::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...
 

Detailed Description

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.

Member Enumeration Documentation

◆ Algorithm

Enum for switching between KMC algorithms.

'SSA' is the Gillespie algorithm, 'Tau' is tau-leaping and 'Hybrid' is the Cao et. al. algorithm.

Member Function Documentation

◆ advanceKMC()

void ItoPlasmaPhysics::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
inline

Advance particles.

Parameters
[in,out]a_numParticlesNumber of physical particles
[out]a_numNewPhotonsNumber of new physical photons to generate (of each type)
[in]a_dtTime step
[in]a_EElectric field
[in]a_dxGrid resolution
[in]a_kappaCut-cell volume fraction.

◆ computeAlpha()

virtual Real Physics::ItoPlasma::ItoPlasmaPhysics::computeAlpha ( const RealVect  a_E) const
pure virtual

Compute Townsend ionization coefficient.

Parameters
[in]a_EElectric field.
Returns
Should return the Townsend ionization coefficient.

Implemented in Physics::ItoPlasma::ItoPlasmaAir3LFA.

◆ computeDt()

Real ItoPlasmaPhysics::computeDt ( const RealVect  a_E,
const RealVect  a_pos,
const Vector< FPR a_numParticles 
) const
inlinevirtualnoexcept

Compute a time step to use.

Parameters
[in]a_EElectric field.
[in]a_posPosition
[in]a_numParticlesNumber of reactive particles of each species
Returns
Default returns infinity but user can implement a physics based time step.

Reimplemented in Physics::ItoPlasma::ItoPlasmaAir3LFA.

◆ computeItoDiffusion()

virtual Vector<Real> Physics::ItoPlasma::ItoPlasmaPhysics::computeItoDiffusion ( const Real  a_time,
const RealVect  a_pos,
const RealVect  a_E 
) const
pure virtualnoexcept

Compute the Ito solver diffusion coefficients.

Parameters
[in]a_timeTime
[in]a_posPosition
[in]a_EElectric field
Returns
Must return a vector of non-negative diffusion coefficients for the plasma species

Implemented in Physics::ItoPlasma::ItoPlasmaAir3LFA.

◆ computeItoMobilities()

virtual Vector<Real> Physics::ItoPlasma::ItoPlasmaPhysics::computeItoMobilities ( const Real  a_time,
const RealVect  a_pos,
const RealVect  a_E 
) const
pure virtualnoexcept

Compute the Ito solver mobilities.

Parameters
[in]a_timeTime
[in]a_posPosition
[in]a_EElectric field
Returns
Must return a vector of non-negative mobility coefficients for the plasma species

Implemented in Physics::ItoPlasma::ItoPlasmaAir3LFA.

◆ getItoSpecies()

const Vector< RefCountedPtr< ItoSpecies > > & ItoPlasmaPhysics::getItoSpecies ( ) const
inline

Get all particle species.

Returns
m_plasmaSpecies

◆ getRtSpecies()

const Vector< RefCountedPtr< RtSpecies > > & ItoPlasmaPhysics::getRtSpecies ( ) const
inline

Get all photon species.

Returns
m_rtSpecies

◆ initialSigma()

Real ItoPlasmaPhysics::initialSigma ( const Real  a_time,
const RealVect  a_pos 
) const
inlinevirtual

Set initial surface charge. Default is 0, override if you want.

Parameters
[in]a_timeSimulation time
[in]a_posPhysical coordinate

◆ injectParticlesEB()

virtual void Physics::ItoPlasma::ItoPlasmaPhysics::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
pure virtualnoexcept

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_newNumParticlesTotal number of particles in the cut-cell AFTER the transport step
[in]a_oldNumParticlesTotal number of particles in the cut-cell BEFORE the transport step
[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)

Implemented in Physics::ItoPlasma::ItoPlasmaAir3LFA.

◆ reconcileParticles()

void ItoPlasmaPhysics::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
inlinenoexcept

Reconcile the number of particles.

This will add/remove particles and potentially also adjust the particle weights.

Parameters
[in,out]a_particlesComputational particles
[in]a_newNumParticlesNew number of particles (i.e., after the KMC advance)
[in]a_oldNumParticlesPrevious number of particles (i.e., before the KMC advance)
[in]a_cellPosCell center position
[in]a_centroidPosCell centroid position
[in]a_loCornerLow corner of minimum box enclosing the cut-cell
[in]a_hiCornerHigh corner of minimum box enclosing the cut-cell
[in]a_bndryCentroidCut-cell boundary centroid
[in]a_bndryNormalCut-cell normal (pointing into the domain)
[in]a_dxGrid resolution
[in]a_kappaCut-cell volume fraction.
Note
Public because this is called by ItoPlasmaStepper

◆ reconcilePhotoionization()

void ItoPlasmaPhysics::reconcilePhotoionization ( Vector< List< ItoParticle > * > &  a_particles,
const Vector< List< Photon > * > &  a_absorbedPhotons 
) const
inlinenoexcept

Reconcile photoionization reactions.

Parameters
[in,out]a_particlesParticle products.
[in]a_absorbedPhotonsPhotons absorbed on the mesh.

This runs through the photo-reactions and associates photo-ionization products.

◆ reconcilePhotons()

void ItoPlasmaPhysics::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
inlinenoexcept

Generate new photons.

This will add photons

Parameters
[in]a_newPhotonsNew photons
[in]a_numNewPhotonsNumber of physical photons to be generated.
[in]a_cellPosCell center position
[in]a_centroidPosCell centroid position
[in]a_loCornerLow corner of minimum box enclosing the cut-cell
[in]a_hiCornerHigh corner of minimum box enclosing the cut-cell
[in]a_bndryCentroidCut-cell boundary centroid
[in]a_bndryNormalCut-cell normal (pointing into the domain)
[in]a_dxGrid resolution
[in]a_kappaCut-cell volume fraction.
Note
Public because this is called by ItoPlasmaStepper

◆ removeParticles()

void ItoPlasmaPhysics::removeParticles ( List< ItoParticle > &  a_particles,
const long long  a_numToRemove 
) const
inlineprotected

Remove particles from the input list.

This will remove weight from the input particles if we can. Otherwise we remove full particles.

Parameters
[in,out]a_particlesList of (super-)particles to remove from.
[in]a_numToRemoveNumber of physical particles to remove from the input list

◆ updateReactionRates()

virtual void Physics::ItoPlasma::ItoPlasmaPhysics::updateReactionRates ( const RealVect  a_E,
const Real  a_dx,
const Real  a_kappa 
) const
protectedpure virtualnoexcept

Update reaction rates.

Parameters
[in]a_EElectric field
[in]a_dxGrid resolution
[in]a_kappaCut-cell volume fraction
Note
Must be implemented by the user.

Implemented in Physics::ItoPlasma::ItoPlasmaAir3LFA.

Member Data Documentation

◆ m_eps

Real Physics::ItoPlasma::ItoPlasmaPhysics::m_eps
protected

Solver setting for the Cao et. al. algorithm.

Equal to the maximum permitted change in the relative propensity for non-critical reactions

◆ m_Ncrit

int Physics::ItoPlasma::ItoPlasmaPhysics::m_Ncrit
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.

◆ m_NSSA

int Physics::ItoPlasma::ItoPlasmaPhysics::m_NSSA
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.

◆ m_SSAlim

Real Physics::ItoPlasma::ItoPlasmaPhysics::m_SSAlim
protected

Solver setting for the Cao et. al. algorithm.

Equal to the maximum permitted change in the relative propensity for non-critical reactions


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