Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ComponentParallelPlate Class Reference

Component for parallel-plate geometries. More...

#include <ComponentParallelPlate.hh>

Inheritance diagram for Garfield::ComponentParallelPlate:
Garfield::Component

Classes

struct  Electrode
 Structure that captures the information of the electrodes under study. More...

Public Member Functions

 ComponentParallelPlate ()
 Constructor.
 ~ComponentParallelPlate ()
 Destructor.
void Setup (const unsigned int N, std::vector< double > eps, std::vector< double > d, const double V, std::vector< int > sigmaIndex={})
 Define the geometry.
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status) override
 Calculate the drift field at given point.
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).
double WeightingPotential (const double x, const double y, const double z, const std::string &label) override
 Calculate the weighting potential at a given point.
bool GetVoltageRange (double &vmin, double &vmax) override
 Calculate the voltage range [V].
void AddPixel (double x, double z, double lx, double lz, const std::string &label, bool fromAnode=true)
 Add a pixel electrode.
void AddStrip (double z, double lz, const std::string &label, bool fromAnode=true)
 Add strip electrode.
void AddPlane (const std::string &label, bool fromAnode=true)
 Add plane electrode, if you want to read the signal from the cathode set the second argument to false.
void SetMedium (Medium *medium)
 Setting the medium.
void SetWeightingPotentialGrid (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 std::string &label)
 Calculate time-dependent weighting potential on a grid.
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)
 This will calculate all electrodes time-dependent weighting potential on the specified grid.
void LoadWeightingPotentialGrid (const std::string &label)
 This will load a previously calculated grid of time-dependent weighting potential values.
MediumGetMedium (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.
bool getLayer (const double y, int &m, double &epsM)
void getPermittivityFromLayer (int m, double &eps)
void getZBoundFromLayer (int m, double &zbottom, double &ztop)
int NumberOfLayers ()
void IndexOfGasGaps (std::vector< int > &indexGasGap)
void SetIntegrationPrecision (const double eps)
void SetIntegrationUpperbound (const double p)
void DisablePotentialCalculationOutsideGasGap ()
std::array< double, 3 > ElectricField (const double x, const double y, const double z)
 Calculate the drift field [V/cm] at (x, y, z).
Public Member Functions inherited from Garfield::Component
 Component ()=delete
 Default constructor.
 Component (const std::string &name)
 Constructor.
virtual ~Component ()=default
 Destructor.
virtual void SetGeometry (Geometry *geo)
 Define the geometry.
virtual void Clear ()
 Reset.
std::array< double, 3 > ElectricField (const double x, const double y, const double z)
 Calculate the drift field [V/cm] at (x, y, z).
virtual double ElectricPotential (const double x, const double y, const double z)
 Calculate the (drift) electrostatic potential [V] at (x, y, z).
virtual void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
 Calculate the weighting field at a given point and for a given electrode.
virtual const std::vector< double > & DelayedSignalTimes (const std::string &)
 Return the time steps at which the delayed weighting potential/field are stored/evaluated.
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)
 Calculate the delayed weighting field at a given point and time and for a given electrode.
virtual double DelayedWeightingPotential (const double x, const double y, const double z, const double t, const std::string &label)
 Calculate the delayed weighting potential at a given point and time and for a given electrode.
virtual void DelayedWeightingPotentials (const double x, const double y, const double z, const std::string &label, std::vector< double > &dwp)
 Calculate the delayed weighting potentials at a given point and for a given electrode, for a set of pre-defined times.
virtual void MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
 Calculate the magnetic field at a given point.
void SetMagneticField (const double bx, const double by, const double bz)
 Set a constant magnetic field.
virtual bool IsReady ()
 Ready for use?
virtual bool Is3d ()
 Does the component have a 3D field (map)?
virtual bool GetElementaryCell (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
 Get the coordinates of the elementary cell.
double CellSizeX ()
double CellSizeY ()
double CellSizeZ ()
virtual size_t GetNumberOfElements () const
 Return the number of mesh elements.
virtual bool GetElementNodes (const size_t, std::vector< size_t > &) const
 Get the indices of the nodes constituting a given element.
virtual bool GetElementRegion (const size_t, size_t &, bool &) const
 Get the region/material of a mesh element and a flag whether it is associated to an active medium.
virtual size_t GetNumberOfNodes () const
 Return the number of mesh nodes.
virtual bool GetNode (const size_t i, double &x, double &y, double &z) const
 Get the coordinates of a mesh node.
double IntegrateFluxCircle (const double xc, const double yc, const double r, const unsigned int nI=50)
 Integrate the normal component of the electric field over a circle.
double IntegrateFluxSphere (const double xc, const double yc, const double zc, const double r, const unsigned int nI=20)
 Integrate the normal component of the electric field over a sphere.
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)
 Integrate the normal component of the electric field over a parallelogram.
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)
 Integrate the normal component of the weighting field over a parallelogram.
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)
 Integrate the electric field flux through a line from (x0,y0,z0) to (x1,y1,z1) along a direction (xp,yp,zp).
virtual bool CrossedWire (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)
 Determine whether the line between two points crosses a wire.
virtual bool InTrapRadius (const double q0, const double x0, const double y0, const double z0, double &xw, double &yw, double &rw)
 Determine whether a particle is inside the trap radius of a wire.
virtual bool CrossedPlane (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc)
 Determine whether the line between two points crosses a plane.
void EnablePeriodicityX (const bool on=true)
 Enable simple periodicity in the $x$ direction.
void EnablePeriodicityY (const bool on=true)
 Enable simple periodicity in the $y$ direction.
void EnablePeriodicityZ (const bool on=true)
 Enable simple periodicity in the $z$ direction.
void IsPeriodic (bool &perx, bool &pery, bool &perz)
 Return periodicity flags.
void EnableMirrorPeriodicityX (const bool on=true)
 Enable mirror periodicity in the $x$ direction.
void EnableMirrorPeriodicityY (const bool on=true)
 Enable mirror periodicity in the $y$ direction.
void EnableMirrorPeriodicityZ (const bool on=true)
 Enable mirror periodicity in the $y$ direction.
void IsMirrorPeriodic (bool &perx, bool &pery, bool &perz)
 Return mirror periodicity flags.
void EnableAxialPeriodicityX (const bool on=true)
 Enable axial periodicity in the $x$ direction.
void EnableAxialPeriodicityY (const bool on=true)
 Enable axial periodicity in the $y$ direction.
void EnableAxialPeriodicityZ (const bool on=true)
 Enable axial periodicity in the $z$ direction.
void IsAxiallyPeriodic (bool &perx, bool &pery, bool &perz)
 Return axial periodicity flags.
void EnableRotationSymmetryX (const bool on=true)
 Enable rotation symmetry around the $x$ axis.
void EnableRotationSymmetryY (const bool on=true)
 Enable rotation symmetry around the $y$ axis.
void EnableRotationSymmetryZ (const bool on=true)
 Enable rotation symmetry around the $z$ axis.
void IsRotationSymmetric (bool &rotx, bool &roty, bool &rotz)
 Return rotation symmetry flags.
void EnableTriangleSymmetricXY (const bool on=true, const bool oct=2)
 Enable triangular periodicity in the $xy$ plane.
void EnableTriangleSymmetricXZ (const bool on=true, const bool oct=2)
 Enable triangular periodicity in the $xz$ plane.
void EnableTriangleSymmetricYZ (const bool on=true, const bool oct=2)
 Enable triangular periodicity in the $yz$ plane.
void EnableDebugging (const bool on=true)
 Switch on debugging messages.
void DisableDebugging ()
 Switch off debugging messages.
virtual bool HasMagneticField () const
 Does the component have a non-zero magnetic field?
virtual bool HasTownsendMap () const
 Does the component have maps of the Townsend coefficient?
virtual bool HasAttachmentMap () const
 Does the component have maps of the attachment coefficient?
virtual bool HasMobilityMap () const
 Does the component have maps of the low-field mobility?
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 ElectronMobility (const double, const double, const double, double &mu)
virtual bool HoleMobility (const double, const double, const double, double &mu)
 Get the hole Mobility coefficient.
virtual bool ElectronTownsend (const double, const double, const double, double &alpha)
 Get the electron Townsend coefficient.
virtual bool HoleTownsend (const double, const double, const double, double &alpha)
 Get the hole Townsend 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 double StepSizeHint ()
virtual double CreateGPUTransferObject (ComponentGPU *&comp_gpu)
 Create and initialise GPU Transfer class.

Private Types

enum  structureelectrode { NotSet = -1 , Plane , Strip , Pixel }
 Possible readout groups. More...

Private Member Functions

double IntegratePromptPotential (const Electrode &el, const double x, const double y, const double z)
void CalculateDynamicalWeightingPotential (const Electrode &el)
double FindWeightingPotentialInGrid (Electrode &el, const double x, const double y, const double z)
bool Nsigma (int N, std::vector< std::vector< int > > &sigmaMatrix)
bool Ntheta (int N, std::vector< std::vector< int > > &thetaMatrix, std::vector< std::vector< int > > &sigmaMatrix)
void constructGeometryMatrices (const int N)
void constructGeometryFunction (const int N, const std::vector< double > &d)
void setHIntegrand ()
void setwpPixelIntegrand ()
void setwpStripIntegrand ()
double constWEFieldLayer (const int indexLayer)
double wpPlane (const double z)
double constEFieldLayer (const int indexLayer)
bool decToBinary (int n, std::vector< int > &binaryNum)
void LayerUpdate (const double z, const int im, const double epsM)
void UpdatePeriodicity () override
 Verify periodicities.
void Reset () override
 Reset the component.

Private Attributes

double m_precision = 1.e-12
double m_V = 0.
 Voltage difference between the parallel plates.
bool m_getPotentialInPlate = true
int m_N = 0
 Number of layers.
double m_upperBoundIntegration = 30
std::vector< double > m_eps
 relative permittivity of each layer
std::vector< double > m_epsHolder
std::vector< double > m_d
 thickness of each layer
std::vector< double > m_z
std::vector< bool > m_conductive
 Flag whether a layer is conductive.
TF2 m_hIntegrand
TF1 m_wpStripIntegral
 Weighting potential integrand for strips.
TF2 m_wpPixelIntegral
 Weighting potential integrand for pixels.
std::vector< std::vector< std::vector< int > > > m_sigmaMatrix
std::vector< std::vector< std::vector< int > > > m_thetaMatrix
std::vector< std::vector< double > > m_cMatrix
 c-matrixl.
std::vector< std::vector< double > > m_vMatrix
 v-matrixl.
std::vector< std::vector< double > > m_gMatrix
 g-matrixl.
std::vector< std::vector< double > > m_wMatrix
 w-matrixl.
int m_currentLayer = 0
 Index of the current layer.
double m_currentPosition = -1
Mediumm_medium = nullptr
std::vector< std::string > m_readout
std::vector< Electrodem_readout_p

Static Private Attributes

static constexpr double m_Vw = 1.

Additional Inherited Members

Protected Attributes inherited from Garfield::Component
std::string m_className {"Component"}
 Class name.
Geometrym_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.
std::array< bool, 3 > m_triangleSymmetric = {{false, false, false}}
 Triangle symmetric in the xy, xz, and yz plane.
int m_triangleSymmetricOct = 0
 Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).
const std::array< int, 4 > m_triangleOctRules = {1, 4, 5, 8}
 Octants where |x| >= |y|.
bool m_outsideCone = false
std::vector< double > m_wdtimes
 Time steps at which the delayed weighting potentials/fields are stored.

Detailed Description

Component for parallel-plate geometries.

Definition at line 19 of file ComponentParallelPlate.hh.

Member Enumeration Documentation

◆ structureelectrode

Constructor & Destructor Documentation

◆ ComponentParallelPlate()

Garfield::ComponentParallelPlate::ComponentParallelPlate ( )

Constructor.

◆ ~ComponentParallelPlate()

Garfield::ComponentParallelPlate::~ComponentParallelPlate ( )
inline

Destructor.

Definition at line 24 of file ComponentParallelPlate.hh.

24{}

Member Function Documentation

◆ AddPixel()

void Garfield::ComponentParallelPlate::AddPixel ( double x,
double z,
double lx,
double lz,
const std::string & label,
bool fromAnode = true )

Add a pixel electrode.

Parameters
x,zposition of the center of the electrode in the xz-plane.
lxwidth in the along $x$.
lzwidth in the along $z$.
labelgive name using a string.
fromAnodeis $true$ is the electrode is the andode and $false$ if it is the cathode.

◆ AddPlane()

void Garfield::ComponentParallelPlate::AddPlane ( const std::string & label,
bool fromAnode = true )

Add plane electrode, if you want to read the signal from the cathode set the second argument to false.

◆ AddStrip()

void Garfield::ComponentParallelPlate::AddStrip ( double z,
double lz,
const std::string & label,
bool fromAnode = true )

Add strip electrode.

◆ CalculateDynamicalWeightingPotential()

void Garfield::ComponentParallelPlate::CalculateDynamicalWeightingPotential ( const Electrode & el)
private

◆ constEFieldLayer()

double Garfield::ComponentParallelPlate::constEFieldLayer ( const int indexLayer)
inlineprivate

Definition at line 278 of file ComponentParallelPlate.hh.

278 {
279 if (m_conductive[indexLayer]) return 0.;
280 double invEz = 0;
281 for (int i = 1; i <= m_N - 1; i++) {
282 // TODO!
283 if (m_conductive[indexLayer]) continue;
284 invEz -= m_d[i - 1] / m_epsHolder[i - 1];
285 }
286 return m_V / (m_epsHolder[indexLayer - 1] * invEz);
287 }
std::vector< bool > m_conductive
Flag whether a layer is conductive.
double m_V
Voltage difference between the parallel plates.
std::vector< double > m_d
thickness of each layer

◆ constructGeometryFunction()

void Garfield::ComponentParallelPlate::constructGeometryFunction ( const int N,
const std::vector< double > & d )
private

◆ constructGeometryMatrices()

void Garfield::ComponentParallelPlate::constructGeometryMatrices ( const int N)
private

◆ constWEFieldLayer()

double Garfield::ComponentParallelPlate::constWEFieldLayer ( const int indexLayer)
inlineprivate

Definition at line 256 of file ComponentParallelPlate.hh.

256 {
257 double invEz = 0;
258 for (int i = 1; i <= m_N - 1; i++) {
259 invEz += m_d[i - 1] / m_epsHolder[i - 1];
260 }
261 return 1 / (m_epsHolder[indexLayer - 1] * invEz);
262 }

◆ decToBinary()

bool Garfield::ComponentParallelPlate::decToBinary ( int n,
std::vector< int > & binaryNum )
private

◆ DisablePotentialCalculationOutsideGasGap()

void Garfield::ComponentParallelPlate::DisablePotentialCalculationOutsideGasGap ( )
inline

Definition at line 155 of file ComponentParallelPlate.hh.

◆ ElectricField() [1/3]

std::array< double, 3 > Garfield::Component::ElectricField ( const double x,
const double y,
const double z )

Calculate the drift field [V/cm] at (x, y, z).

◆ ElectricField() [2/3]

void Garfield::ComponentParallelPlate::ElectricField ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
double & v,
Medium *& m,
int & status )
overridevirtual

Calculate the drift field [V/cm] and potential [V] at (x, y, z).

Implements Garfield::Component.

◆ ElectricField() [3/3]

void Garfield::ComponentParallelPlate::ElectricField ( const double x,
const double y,
const double z,
double & ex,
double & ey,
double & ez,
Medium *& m,
int & status )
overridevirtual

Calculate the drift field at given point.

Parameters
x,y,zcoordinates [cm].
ex,ey,ezcomponents of the electric field [V/cm].
mpointer to the medium at this location.
statusstatus 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.

◆ FindWeightingPotentialInGrid()

double Garfield::ComponentParallelPlate::FindWeightingPotentialInGrid ( Electrode & el,
const double x,
const double y,
const double z )
private

◆ GetBoundingBox()

bool Garfield::ComponentParallelPlate::GetBoundingBox ( double & xmin,
double & ymin,
double & zmin,
double & xmax,
double & ymax,
double & zmax )
overridevirtual

Get the bounding box coordinates.

Reimplemented from Garfield::Component.

◆ getLayer()

bool Garfield::ComponentParallelPlate::getLayer ( const double y,
int & m,
double & epsM )
inline

Definition at line 119 of file ComponentParallelPlate.hh.

119 {
120 m = -1;
121 if (y < m_z[0]) return false;
122 for (int i = 1; i < m_N; i++) {
123 if (y <= m_z[i]) {
124 m = i;
125 break;
126 }
127 }
128 if (m == -1) return false;
129 epsM = m_epsHolder[m - 1];
130 return true;
131 }

◆ GetMedium()

Medium * Garfield::ComponentParallelPlate::GetMedium ( const double x,
const double y,
const double z )
overridevirtual

Get the medium at a given location (x, y, z).

Reimplemented from Garfield::Component.

◆ getPermittivityFromLayer()

void Garfield::ComponentParallelPlate::getPermittivityFromLayer ( int m,
double & eps )
inline

Definition at line 133 of file ComponentParallelPlate.hh.

133 {
134 eps = m_epsHolder.at(m - 1);
135 }

◆ GetVoltageRange()

bool Garfield::ComponentParallelPlate::GetVoltageRange ( double & vmin,
double & vmax )
overridevirtual

Calculate the voltage range [V].

Implements Garfield::Component.

◆ getZBoundFromLayer()

void Garfield::ComponentParallelPlate::getZBoundFromLayer ( int m,
double & zbottom,
double & ztop )
inline

Definition at line 137 of file ComponentParallelPlate.hh.

137 {
138 ztop = m_z.at(m);
139 zbottom = m_z.at(m - 1);
140 }

◆ IndexOfGasGaps()

void Garfield::ComponentParallelPlate::IndexOfGasGaps ( std::vector< int > & indexGasGap)
inline

Definition at line 144 of file ComponentParallelPlate.hh.

144 {
145 indexGasGap = {};
146 for (int i = 1; i < m_N; i++) {
147 if (!m_conductive[i]) indexGasGap.push_back(i);
148 }
149 }

◆ IntegratePromptPotential()

double Garfield::ComponentParallelPlate::IntegratePromptPotential ( const Electrode & el,
const double x,
const double y,
const double z )
private

◆ LayerUpdate()

void Garfield::ComponentParallelPlate::LayerUpdate ( const double z,
const int im,
const double epsM )
inlineprivate

Definition at line 293 of file ComponentParallelPlate.hh.

293 {
294 if (z == m_currentPosition) return;
295
297
298 if (im != m_currentLayer) {
299 m_currentLayer = im;
300 for (int i = 0; i < im - 1; i++) m_eps[i] = m_epsHolder[i];
301 m_eps[im - 1] = epsM;
302 m_eps[im] = epsM;
303 for (int i = im + 1; i < m_N; i++) m_eps[i] = m_epsHolder[i - 1];
304 }
305
306 double diff1 = m_z[im] - z;
307 double diff2 = z - m_z[im - 1];
308
309 std::vector<double> d(m_N, 0.);
310 for (int i = 0; i < im - 1; i++) d[i] = m_d[i];
311 d[im - 1] = diff2;
312 d[im] = diff1;
313 for (int i = im + 1; i < m_N; i++) d[i] = m_d[i - 1];
314 // TODO::Construct c and g matrices only for im != m_currentLayer.
316 };
int m_currentLayer
Index of the current layer.
std::vector< double > m_eps
relative permittivity of each layer
void constructGeometryFunction(const int N, const std::vector< double > &d)

◆ LoadWeightingPotentialGrid()

void Garfield::ComponentParallelPlate::LoadWeightingPotentialGrid ( const std::string & label)
inline

This will load a previously calculated grid of time-dependent weighting potential values.

Definition at line 97 of file ComponentParallelPlate.hh.

97 {
98 for (auto &electrode : m_readout_p) {
99 if (electrode.label != label) continue;
100 if (electrode.grid.LoadWeightingField(label + "map", "xyz", true)) {
101 std::cout << m_className << "::LoadWeightingPotentialGrid: "
102 << "Weighting potential set for " << label << ".\n";
103 electrode.m_usegrid = true;
104 return;
105 }
106 }
107 std::cerr << m_className
108 << "::LoadWeightingPotentialGrid: Could not find file for "
109 << label << ".\n";
110 }
std::string m_className
Class name.
Definition Component.hh:433

◆ Nsigma()

bool Garfield::ComponentParallelPlate::Nsigma ( int N,
std::vector< std::vector< int > > & sigmaMatrix )
private

◆ Ntheta()

bool Garfield::ComponentParallelPlate::Ntheta ( int N,
std::vector< std::vector< int > > & thetaMatrix,
std::vector< std::vector< int > > & sigmaMatrix )
private

◆ NumberOfLayers()

int Garfield::ComponentParallelPlate::NumberOfLayers ( )
inline

Definition at line 142 of file ComponentParallelPlate.hh.

142{ return m_N - 1; }

◆ Reset()

void Garfield::ComponentParallelPlate::Reset ( )
overrideprivatevirtual

Reset the component.

Implements Garfield::Component.

◆ setHIntegrand()

void Garfield::ComponentParallelPlate::setHIntegrand ( )
private

◆ SetIntegrationPrecision()

void Garfield::ComponentParallelPlate::SetIntegrationPrecision ( const double eps)
inline

Definition at line 151 of file ComponentParallelPlate.hh.

◆ SetIntegrationUpperbound()

void Garfield::ComponentParallelPlate::SetIntegrationUpperbound ( const double p)
inline

◆ SetMedium()

void Garfield::ComponentParallelPlate::SetMedium ( Medium * medium)
inline

Setting the medium.

Definition at line 72 of file ComponentParallelPlate.hh.

◆ Setup()

void Garfield::ComponentParallelPlate::Setup ( const unsigned int N,
std::vector< double > eps,
std::vector< double > d,
const double V,
std::vector< int > sigmaIndex = {} )

Define the geometry.

Parameters
Namount of layers in the geometry, this includes the gas gaps $y$.
dthickness of the layers starting from the bottom to the top layer along $y$.
epsrelative permittivities of the layers starting from the bottom to the top layer along $y$ . Here, the gas gaps having a value of 1.
sigmaIndexIndices of the resistive layers (optional).
Vapplied potential difference between the parallel plates.

◆ SetWeightingPotentialGrid()

void Garfield::ComponentParallelPlate::SetWeightingPotentialGrid ( 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 std::string & label )

Calculate time-dependent weighting potential on a grid.

Parameters
xmin,ymin,zminminimum value of the interval in the $x$-, $y$- and $z$-direction.
xmax,ymax,zmaxmaximum value of the interval in the $x$-, $y$- and $z$-direction.
xsteps,ysteps,zstepsmumber of grid nodes in the $x$-, $y$- and $z$-direction.
labelgive name using a string.

◆ 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 )

This will calculate all electrodes time-dependent weighting potential on the specified grid.

◆ setwpPixelIntegrand()

void Garfield::ComponentParallelPlate::setwpPixelIntegrand ( )
private

◆ setwpStripIntegrand()

void Garfield::ComponentParallelPlate::setwpStripIntegrand ( )
private

◆ UpdatePeriodicity()

void Garfield::ComponentParallelPlate::UpdatePeriodicity ( )
overrideprivatevirtual

Verify periodicities.

Implements Garfield::Component.

◆ WeightingPotential()

double Garfield::ComponentParallelPlate::WeightingPotential ( const double x,
const double y,
const double z,
const std::string & label )
overridevirtual

Calculate the weighting potential at a given point.

Parameters
x,y,zcoordinates [cm].
labelname of the electrode.
Returns
weighting potential [dimensionless].

Reimplemented from Garfield::Component.

◆ wpPlane()

double Garfield::ComponentParallelPlate::wpPlane ( const double z)
inlineprivate

Definition at line 265 of file ComponentParallelPlate.hh.

265 {
266 int im = -1;
267 double epsM = -1;
268 if (!getLayer(z, im, epsM)) return 0.;
269 double v = 1 - (z - m_z[im - 1]) * constWEFieldLayer(im);
270 for (int i = 1; i <= im - 1; i++) {
271 v -= m_d[i - 1] * constWEFieldLayer(i);
272 }
273
274 return v;
275 }
double constWEFieldLayer(const int indexLayer)
bool getLayer(const double y, int &m, double &epsM)

Member Data Documentation

◆ m_cMatrix

std::vector<std::vector<double> > Garfield::ComponentParallelPlate::m_cMatrix
private

c-matrixl.

Definition at line 191 of file ComponentParallelPlate.hh.

◆ m_conductive

std::vector<bool> Garfield::ComponentParallelPlate::m_conductive
private

Flag whether a layer is conductive.

Definition at line 177 of file ComponentParallelPlate.hh.

◆ m_currentLayer

int Garfield::ComponentParallelPlate::m_currentLayer = 0
private

Index of the current layer.

Definition at line 196 of file ComponentParallelPlate.hh.

◆ m_currentPosition

double Garfield::ComponentParallelPlate::m_currentPosition = -1
private

Definition at line 197 of file ComponentParallelPlate.hh.

◆ m_d

std::vector<double> Garfield::ComponentParallelPlate::m_d
private

thickness of each layer

Definition at line 173 of file ComponentParallelPlate.hh.

◆ m_eps

std::vector<double> Garfield::ComponentParallelPlate::m_eps
private

relative permittivity of each layer

Definition at line 171 of file ComponentParallelPlate.hh.

◆ m_epsHolder

std::vector<double> Garfield::ComponentParallelPlate::m_epsHolder
private

Definition at line 172 of file ComponentParallelPlate.hh.

◆ m_getPotentialInPlate

bool Garfield::ComponentParallelPlate::m_getPotentialInPlate = true
private

Definition at line 165 of file ComponentParallelPlate.hh.

◆ m_gMatrix

std::vector<std::vector<double> > Garfield::ComponentParallelPlate::m_gMatrix
private

g-matrixl.

Definition at line 193 of file ComponentParallelPlate.hh.

◆ m_hIntegrand

TF2 Garfield::ComponentParallelPlate::m_hIntegrand
private

Definition at line 179 of file ComponentParallelPlate.hh.

◆ m_medium

Medium* Garfield::ComponentParallelPlate::m_medium = nullptr
private

Definition at line 199 of file ComponentParallelPlate.hh.

◆ m_N

int Garfield::ComponentParallelPlate::m_N = 0
private

Number of layers.

Definition at line 167 of file ComponentParallelPlate.hh.

◆ m_precision

double Garfield::ComponentParallelPlate::m_precision = 1.e-12
private

Definition at line 160 of file ComponentParallelPlate.hh.

◆ m_readout

std::vector<std::string> Garfield::ComponentParallelPlate::m_readout
private

Definition at line 217 of file ComponentParallelPlate.hh.

◆ m_readout_p

std::vector<Electrode> Garfield::ComponentParallelPlate::m_readout_p
private

Definition at line 218 of file ComponentParallelPlate.hh.

◆ m_sigmaMatrix

std::vector<std::vector<std::vector<int> > > Garfield::ComponentParallelPlate::m_sigmaMatrix
private

Definition at line 184 of file ComponentParallelPlate.hh.

◆ m_thetaMatrix

std::vector<std::vector<std::vector<int> > > Garfield::ComponentParallelPlate::m_thetaMatrix
private

Definition at line 187 of file ComponentParallelPlate.hh.

◆ m_upperBoundIntegration

double Garfield::ComponentParallelPlate::m_upperBoundIntegration = 30
private

Definition at line 169 of file ComponentParallelPlate.hh.

◆ m_V

double Garfield::ComponentParallelPlate::m_V = 0.
private

Voltage difference between the parallel plates.

Definition at line 163 of file ComponentParallelPlate.hh.

◆ m_vMatrix

std::vector<std::vector<double> > Garfield::ComponentParallelPlate::m_vMatrix
private

v-matrixl.

Definition at line 192 of file ComponentParallelPlate.hh.

◆ m_Vw

double Garfield::ComponentParallelPlate::m_Vw = 1.
staticconstexprprivate

Definition at line 161 of file ComponentParallelPlate.hh.

◆ m_wMatrix

std::vector<std::vector<double> > Garfield::ComponentParallelPlate::m_wMatrix
private

w-matrixl.

Definition at line 194 of file ComponentParallelPlate.hh.

◆ m_wpPixelIntegral

TF2 Garfield::ComponentParallelPlate::m_wpPixelIntegral
private

Weighting potential integrand for pixels.

Definition at line 182 of file ComponentParallelPlate.hh.

◆ m_wpStripIntegral

TF1 Garfield::ComponentParallelPlate::m_wpStripIntegral
private

Weighting potential integrand for strips.

Definition at line 181 of file ComponentParallelPlate.hh.

◆ m_z

std::vector<double> Garfield::ComponentParallelPlate::m_z
private

Definition at line 174 of file ComponentParallelPlate.hh.


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