1#ifndef GARFIELD_AVALANCHEGRIDSPACECHARGE_HH
2#define GARFIELD_AVALANCHEGRIDSPACECHARGE_HH
67 const int nof_approx = 1.) {
90 void Set2dGrid(
double zmin,
double zmax,
int zsteps,
double rmax,
int rsteps);
97 void AddElectron(
double x,
double y,
double z,
double t = 0,
int n = 1);
122 [[nodiscard]]
const std::vector<std::pair<double, long>>
184 bool SnapTo2dGrid(
double x,
double y,
double z,
long n = 1,
int gasLayer = 0);
194 double nNegIon,
int iz,
int ir,
int gasGap);
198 int ir,
double stepZ,
double stepR,
int gasGap);
203 const std::string &fieldOption,
int gasGap);
213 double &eFieldZ,
double &eFieldR);
219 double &eFieldZ,
double &eFieldR);
225 double &eFieldZ,
double &eFieldR);
229 double &drift,
double &dSigmaL,
double &dSigmaT,
230 double &wv,
double &wr,
double &alphaPT,
231 double &etaPT,
int gasGap);
235 double &yg,
double &zg,
int gasGap);
307 std::vector<std::vector<int>>
310 std::vector<std::vector<GridNode>>
m_grid;
316 1500, 800, 400, 200, 100,
std::vector< std::vector< Point > > m_vElectrons
Electrons to transfer onto grid.
std::vector< double > m_zGrid
Grid points of z-coordinate.
void EnableStickyAnode(const bool option=true)
Enable sticky anode (default on)
bool m_isgridset
Distance between the grid points.
bool AddFieldFromChargeAt(int iz, int ir, int fz, int fr, double N, double &eFieldZ, double &eFieldR)
void Set2dGrid(double zmin, double zmax, int zsteps, double rmax, int rsteps)
void EnableMC(const bool option=true)
Enable Monte Carlo (gain/diffusion) up to a certain total number of electrons e.g.
void ExportGrid(const std::string &filename)
Export the current grid to a txt file (electron, ions numbers and field magnitude) Take care: only wh...
void AddExtraElectron(double y, int n=1)
After calling AddElectron, add more electrons on the same transversal line (y-freedom).
std::vector< std::pair< double, long > > m_vNElectronEvolution
std::vector< double > m_vYPointInGasGap
Example point (y-coord) in each gas gap.
std::vector< std::vector< int > > m_zGasGapBoundaries
[k] -> {izLeft, ..., izRight}
std::vector< double > m_ezBkg
Uniform background field in z direction, can be negative.
std::vector< int > m_vSaturatedGaps
Which gas gaps are saturated if saturation is on.
void SetFieldCalculation(const std::string &option="coulomb", const int nof_approx=1.)
Sets the different options to calculate the space charge field coulomb: free field approximation (rel...
std::vector< double > m_vEElliptic
std::vector< double > m_rGrid
Distance between the grid points.
bool m_bRatesAvailable
Flag if temporal rates are available to the simulation.
void EnableDiffusion(const bool option=true)
Enable diffusion (default off)
double m_time0
Initial time.
std::vector< std::vector< GridNode > > m_grid
grid with nodes on it
void GetEllipticIntegrals(double x, double &K, double &E)
AvalancheGridSpaceCharge()
Default constructor.
bool m_bWrAvailable
Flag if bulk drift velocity is available to the simulation.
void SetStopAtK(bool option=true)
Stop the avalanche if K % field is reached.
int m_rSteps
Number of grid points.
bool m_bPreparedImportAvalanche
std::vector< std::vector< double > > m_vCoNGasLayer
Coordinates of center of electron number.
std::vector< double > m_vKElliptic
void SetSensor(Sensor *sensor)
Set the sensor (and determine if it includes a parallel-plate component).
long m_nTotElectron
Total amount of electrons at time step.
void EnableSpaceChargeEffect(const bool option=true)
Enable space charge calculations (default on)
long m_nTotPosIons
total amount of charge created
std::vector< int > m_vIndexGasGaps
Which layer indices are gas layers.
void DiffuseTimeStep(double dx, long nElectron, double nPosIon, double nNegIon, int iz, int ir, int gasGap)
std::string m_sFieldOption
int GetGasGapNumber(int layerIndex)
bool m_run
Tracking if the charges are still in the drift gap.
long ReachedKPercent() const
Returns if 100 * K % background charge has been reached.
bool SnapTo2dGrid(double x, double y, double z, long n=1, int gasLayer=0)
void GetSwarmParameters(double MagEField, double &alpha, double &eta, double &drift, double &dSigmaL, double &dSigmaT, double &wv, double &wr, double &alphaPT, double &etaPT, int gasGap)
void EnableAdaptiveTimeStepping(const bool option=true)
Enable adaptive time stepping (default on)
int m_zSteps
Number of grid points.
void GetFreeChargedRing(double zi, double ri, double zf, double rf, double &eFieldZ, double &eFieldR)
~AvalancheGridSpaceCharge()=default
Destructor.
void ImportEllipticIntegralValues(const std::string &filename)
void AddElectrons(AvalancheMicroscopic *avmc)
Import electron (no ions) data from AvalancheMicroscopic class to axi-symmetric grid.
void DistributeCharges(long nElectron, double nPosIon, double nNegIon, int iz, int ir, double stepZ, double stepR, int gasGap)
std::vector< long > m_vGroupSizes
same values as Lippmann & Riegler
bool AddFieldFromChargeAt(int iz, int ir, double zf, double rf, double N, double &eFieldZ, double &eFieldR)
double GetMeanDistance()
Return current mean distance of the electrons on the grid.
void SetNCrit(const long NCrit=1e8)
Set electron saturation value applied for each gas gap if space charge effect is turned off (default ...
long GetAvalancheSize() const
Returns the total positive charge in the gap's.
std::vector< double > m_vXElliptic
void GetLocalField(int iz, int ir, double &eFieldZ, double &eFieldR, const std::string &fieldOption, int gasGap)
void GetFreeChargedRing(int iz, int ir, int fz, int fr, double &eFieldZ, double &eFieldR)
void EnableTOF(const bool option=true)
Enable to use TOF swarm parameters (default on)
AvalancheGridSpaceCharge(Sensor *sensor)
Constructor.
void Reset()
Reset the charges.
void SetK(float option=0.95)
Set the streamer-inception criterion constant K in the interval (0, ∞) s.t.
bool m_bUseTOF
Flag if TOF parameters should be used, else Magboltz drift and SST spatial coefficients.
void AddElectron(double x, double y, double z, double t=0, int n=1)
Set n electrons onto the grid.
ComponentParallelPlate * m_pp
void PrepareElectronsFromMicroscopicAvalanche()
void StartGridAvalanche(double dtime=-1)
Starts the simulation with the imported electrons for a time step dt (dt = -1 until there are no elec...
void GetGlobalCoordinates(double r, double z, double phi, double &xg, double &yg, double &zg, int gasGap)
const std::vector< std::pair< double, long > > GetElectronEvolution() const
Returns the total electron number evolution.
void EnableDebugging(const bool option=true)
Enable/disable debugging (log messages) (default off)
Calculate electron drift lines and avalanches using microscopic tracking.
Component for parallel-plate geometries.
long nElectronHolder
at t+dt
double townsend
townsend at this node 1/cm
double attachmentPT
Attachment rate from TOF experiment 1/ns -> 1/cm.
int gasGapIndex
Gas gap index: -1 if not gas gap; starts with 0, 1, ...
double velocity
Magnitude of velocity of the node (not negative) cm/ns.
double eFieldR
Space-charge electric field in R direction (can be negative)
double nPosIonHolder
at t+dt
double nNegIonHolder
at t+dt
double townsendPT
Ionization rate from TOF experiment 1/ns -> 1/cm.
double nPosIon
pos ion on node (smeared values allowed)
double dSigmaT
Diffusion transverse to E (radial, phi dir is net 0).
double eFieldZ
Space-charge electric field in Z direction (can be negative)
double Wr
bulk drift cm/ns
double Wv
flux drift cm/ns
int layerIndex
LayerIndex in ParallelPlate convention != gas gap index.
double nNegIon
neg ion on node (smeared values allowed)
long nElectron
electrons on node
double dSigmaL
Diffusion along E.
double attachment
attachment at this node 1/cm