Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
|
Component for parallel-plate geometries. More...
#include <ComponentParallelPlate.hh>
Public Member Functions | |
ComponentParallelPlate () | |
Constructor. | |
~ComponentParallelPlate () | |
Destructor. | |
void | Setup (double g, double b, double eps, double v, double sigma=0.) |
void | ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override |
void | ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status) override |
Calculate the drift field [V/cm] and potential [V] at (x, y, z). | |
void | WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) override |
double | WeightingPotential (const double x, const double y, const double z, const std::string &label) override |
double | DelayedWeightingPotential (const double x, const double y, const double z, const double t, const std::string &label) override |
void | DelayedWeightingField (const double x, const double y, const double z, const double t, double &wx, double &wy, double &wz, const std::string &label) override |
bool | GetVoltageRange (double &vmin, double &vmax) override |
Calculate the voltage range [V]. | |
void | AddPixel (double x, double y, double lx, double ly, const std::string &label) |
void | AddStrip (double x, double lx, const std::string &label) |
Add strip electrode. | |
void | AddPlane (const std::string &label, bool anode=true) |
void | SetMedium (Medium *medium) |
void | SetWeightingPotentialGrid (const std::string &label, const double xmin, const double xmax, const double xsteps, const double ymin, const double ymax, const double ysteps, const double zmin, const double zmax, const double zsteps, const double tmin, const double tmax, const double tsteps) |
void | SetWeightingPotentialGrids (const double xmin, const double xmax, const double xsteps, const double ymin, const double ymax, const double ysteps, const double zmin, const double zmax, const double zsteps, const double tmin, const double tmax, const double tsteps) |
Medium * | GetMedium (const double x, const double y, const double z) override |
Get the medium at a given location (x, y, z). | |
bool | GetBoundingBox (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) override |
Get the bounding box coordinates. | |
Public Member Functions inherited from Garfield::Component | |
Component ()=delete | |
Default constructor. | |
Component (const std::string &name) | |
Constructor. | |
virtual | ~Component () |
Destructor. | |
virtual void | SetGeometry (Geometry *geo) |
Define the geometry. | |
virtual void | Clear () |
Reset. | |
virtual Medium * | GetMedium (const double x, const double y, const double z) |
Get the medium at a given location (x, y, z). | |
virtual void | ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)=0 |
virtual void | ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status)=0 |
Calculate the drift field [V/cm] and potential [V] at (x, y, z). | |
virtual bool | GetVoltageRange (double &vmin, double &vmax)=0 |
Calculate the voltage range [V]. | |
virtual void | WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label) |
virtual double | WeightingPotential (const double x, const double y, const double z, const std::string &label) |
virtual void | DelayedWeightingField (const double x, const double y, const double z, const double t, double &wx, double &wy, double &wz, const std::string &label) |
virtual double | DelayedWeightingPotential (const double x, const double y, const double z, const double t, const std::string &label) |
virtual void | MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status) |
void | SetMagneticField (const double bx, const double by, const double bz) |
Set a constant magnetic field. | |
virtual bool | IsReady () |
Ready for use? | |
virtual bool | GetBoundingBox (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) |
Get the bounding box coordinates. | |
virtual bool | GetElementaryCell (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax) |
Get the coordinates of the elementary cell. | |
double | IntegrateFluxCircle (const double xc, const double yc, const double r, const unsigned int nI=50) |
double | IntegrateFluxSphere (const double xc, const double yc, const double zc, const double r, const unsigned int nI=20) |
double | IntegrateFluxParallelogram (const double x0, const double y0, const double z0, const double dx1, const double dy1, const double dz1, const double dx2, const double dy2, const double dz2, const unsigned int nU=20, const unsigned int nV=20) |
double | IntegrateWeightingFluxParallelogram (const std::string &label, const double x0, const double y0, const double z0, const double dx1, const double dy1, const double dz1, const double dx2, const double dy2, const double dz2, const unsigned int nU=20, const unsigned int nV=20) |
double | IntegrateFluxLine (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const double xp, const double yp, const double zp, const unsigned int nI, const int isign=0) |
virtual bool | IsWireCrossed (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc, const bool centre, double &rc) |
virtual bool | IsInTrapRadius (const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw) |
void | EnablePeriodicityX (const bool on=true) |
Enable simple periodicity in the direction. | |
void | EnablePeriodicityY (const bool on=true) |
Enable simple periodicity in the direction. | |
void | EnablePeriodicityZ (const bool on=true) |
Enable simple periodicity in the direction. | |
void | IsPeriodic (bool &perx, bool &pery, bool &perz) |
Return periodicity flags. | |
void | EnableMirrorPeriodicityX (const bool on=true) |
Enable mirror periodicity in the direction. | |
void | EnableMirrorPeriodicityY (const bool on=true) |
Enable mirror periodicity in the direction. | |
void | EnableMirrorPeriodicityZ (const bool on=true) |
Enable mirror periodicity in the direction. | |
void | IsMirrorPeriodic (bool &perx, bool &pery, bool &perz) |
Return mirror periodicity flags. | |
void | EnableAxialPeriodicityX (const bool on=true) |
Enable axial periodicity in the direction. | |
void | EnableAxialPeriodicityY (const bool on=true) |
Enable axial periodicity in the direction. | |
void | EnableAxialPeriodicityZ (const bool on=true) |
Enable axial periodicity in the direction. | |
void | IsAxiallyPeriodic (bool &perx, bool &pery, bool &perz) |
Return axial periodicity flags. | |
void | EnableRotationSymmetryX (const bool on=true) |
Enable rotation symmetry around the axis. | |
void | EnableRotationSymmetryY (const bool on=true) |
Enable rotation symmetry around the axis. | |
void | EnableRotationSymmetryZ (const bool on=true) |
Enable rotation symmetry around the axis. | |
void | IsRotationSymmetric (bool &rotx, bool &roty, bool &rotz) |
Return rotation symmetry flags. | |
void | EnableDebugging () |
Switch on debugging messages. | |
void | DisableDebugging () |
Switch off debugging messages. | |
virtual bool | HasAttachmentMap () const |
Does the component have attachment maps? | |
virtual bool | HasVelocityMap () const |
Does the component have velocity maps? | |
virtual bool | ElectronAttachment (const double, const double, const double, double &eta) |
Get the electron attachment coefficient. | |
virtual bool | HoleAttachment (const double, const double, const double, double &eta) |
Get the hole attachment coefficient. | |
virtual bool | ElectronVelocity (const double, const double, const double, double &vx, double &vy, double &vz) |
Get the electron drift velocity. | |
virtual bool | HoleVelocity (const double, const double, const double, double &vx, double &vy, double &vz) |
Get the hole drift velocity. | |
virtual bool | GetElectronLifetime (const double, const double, const double, double &etau) |
virtual bool | GetHoleLifetime (const double, const double, const double, double &htau) |
Additional Inherited Members | |
virtual void | Reset ()=0 |
Reset the component. | |
virtual void | UpdatePeriodicity ()=0 |
Verify periodicities. | |
Protected Attributes inherited from Garfield::Component | |
std::string | m_className = "Component" |
Class name. | |
Geometry * | m_geometry = nullptr |
Pointer to the geometry. | |
std::array< double, 3 > | m_b0 = {{0., 0., 0.}} |
Constant magnetic field. | |
bool | m_ready = false |
Ready for use? | |
bool | m_debug = false |
Switch on/off debugging messages. | |
std::array< bool, 3 > | m_periodic = {{false, false, false}} |
Simple periodicity in x, y, z. | |
std::array< bool, 3 > | m_mirrorPeriodic = {{false, false, false}} |
Mirror periodicity in x, y, z. | |
std::array< bool, 3 > | m_axiallyPeriodic = {{false, false, false}} |
Axial periodicity in x, y, z. | |
std::array< bool, 3 > | m_rotationSymmetric = {{false, false, false}} |
Rotation symmetry around x-axis, y-axis, z-axis. | |
Component for parallel-plate geometries.
Definition at line 13 of file ComponentParallelPlate.hh.
Garfield::ComponentParallelPlate::ComponentParallelPlate | ( | ) |
|
inline |
void Garfield::ComponentParallelPlate::AddPixel | ( | double | x, |
double | y, | ||
double | lx, | ||
double | ly, | ||
const std::string & | label | ||
) |
Add a pixel electrode.
x,z | position of the center of the electrode in the xy-plane. |
lx | width in the along . |
ly | width in the along . |
label | give name using a string. |
Definition at line 586 of file ComponentParallelPlate.cc.
void Garfield::ComponentParallelPlate::AddPlane | ( | const std::string & | label, |
bool | anode = true |
||
) |
Add plane electrode, if you want to read the signal from the cathode set the second argument to false.
Definition at line 626 of file ComponentParallelPlate.cc.
void Garfield::ComponentParallelPlate::AddStrip | ( | double | x, |
double | lx, | ||
const std::string & | label | ||
) |
Add strip electrode.
Definition at line 607 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the delayed weighting field at a given point and time and for a given electrode.
x,y,z | coordinates [cm]. |
t | time [ns]. |
wx,wy,wz | components of the weighting field [1/cm]. |
label | name of the electrode |
Reimplemented from Garfield::Component.
Definition at line 542 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the delayed weighting potential at a given point and time and for a given electrode.
x,y,z | coordinates [cm]. |
t | time [ns]. |
label | name of the electrode |
Reimplemented from Garfield::Component.
Definition at line 515 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the drift field [V/cm] and potential [V] at (x, y, z).
Implements Garfield::Component.
Definition at line 429 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the drift field at given point.
x,y,z | coordinates [cm]. |
ex,ey,ez | components of the electric field [V/cm]. |
m | pointer to the medium at this location. |
status | status flag |
Status flags:
0: Inside an active medium > 0: Inside a wire of type X -4 ... -1: On the side of a plane where no wires are -5: Inside the mesh but not in an active medium -6: Outside the mesh -10: Unknown potential type (should not occur) other: Other cases (should not occur)
Implements Garfield::Component.
Definition at line 399 of file ComponentParallelPlate.cc.
|
overridevirtual |
Get the bounding box coordinates.
Reimplemented from Garfield::Component.
Definition at line 40 of file ComponentParallelPlate.cc.
|
overridevirtual |
Get the medium at a given location (x, y, z).
Reimplemented from Garfield::Component.
Definition at line 644 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the voltage range [V].
Implements Garfield::Component.
Definition at line 464 of file ComponentParallelPlate.cc.
|
inline |
Definition at line 69 of file ComponentParallelPlate.hh.
void Garfield::ComponentParallelPlate::Setup | ( | double | g, |
double | b, | ||
double | eps, | ||
double | v, | ||
double | sigma = 0. |
||
) |
Define the geometry.
g | size of the gap along positive . |
b | thickness of the resistive layer along negative . |
eps | relative permittivity of the resistive layer. |
sigma | conductivity of the resistive layer (must be larger then zero, otherwise do not pass it in the function). |
V | applied potential difference between the parallel plates. |
Definition at line 16 of file ComponentParallelPlate.cc.
void Garfield::ComponentParallelPlate::SetWeightingPotentialGrid | ( | const std::string & | label, |
const double | xmin, | ||
const double | xmax, | ||
const double | xsteps, | ||
const double | ymin, | ||
const double | ymax, | ||
const double | ysteps, | ||
const double | zmin, | ||
const double | zmax, | ||
const double | zsteps, | ||
const double | tmin, | ||
const double | tmax, | ||
const double | tsteps | ||
) |
Definition at line 654 of file ComponentParallelPlate.cc.
Referenced by SetWeightingPotentialGrids().
void Garfield::ComponentParallelPlate::SetWeightingPotentialGrids | ( | const double | xmin, |
const double | xmax, | ||
const double | xsteps, | ||
const double | ymin, | ||
const double | ymax, | ||
const double | ysteps, | ||
const double | zmin, | ||
const double | zmax, | ||
const double | zsteps, | ||
const double | tmin, | ||
const double | tmax, | ||
const double | tsteps | ||
) |
Definition at line 718 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the weighting field at a given point and for a given electrode.
x,y,z | coordinates [cm]. |
wx,wy,wz | components of the weighting field [1/cm]. |
label | name of the electrode |
Reimplemented from Garfield::Component.
Definition at line 477 of file ComponentParallelPlate.cc.
|
overridevirtual |
Calculate the weighting potential at a given point.
x,y,z | coordinates [cm]. |
label | name of the electrode. |
Reimplemented from Garfield::Component.
Definition at line 497 of file ComponentParallelPlate.cc.