Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation 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:

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)
 
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.
 
- 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 MediumGetMedium (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 $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 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.
 
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.
 

Detailed Description

Component for parallel-plate geometries.

Definition at line 13 of file ComponentParallelPlate.hh.

Constructor & Destructor Documentation

◆ ComponentParallelPlate()

Garfield::ComponentParallelPlate::ComponentParallelPlate ( )

Constructor.

Definition at line 14 of file ComponentParallelPlate.cc.

14: Component("ParallelPlate") {}
Component()=delete
Default constructor.

◆ ~ComponentParallelPlate()

Garfield::ComponentParallelPlate::~ComponentParallelPlate ( )
inline

Destructor.

Definition at line 18 of file ComponentParallelPlate.hh.

18{}

Member Function Documentation

◆ AddPixel()

void Garfield::ComponentParallelPlate::AddPixel ( double  x,
double  y,
double  lx,
double  ly,
const std::string &  label 
)

Add a pixel electrode.

Parameters
x,zposition of the center of the electrode in the xy-plane.
lxwidth in the along $x$.
lywidth in the along $y$.
labelgive name using a string.

Definition at line 586 of file ComponentParallelPlate.cc.

588 {
589 const auto it = std::find(m_readout.cbegin(), m_readout.cend(), label);
590 if (it == m_readout.end() && m_readout.size() > 0) {
591 std::cerr << m_className << "::AddPixel:\n"
592 << "Note that the label " << label << " is already in use.\n";
593 }
594 Electrode pixel;
595 pixel.label = label;
596 pixel.ind = structureelectrode::Pixel;
597 pixel.xpos = x;
598 pixel.ypos = y;
599 pixel.lx = lx_input;
600 pixel.ly = ly_input;
601
602 m_readout.push_back(label);
603 m_readout_p.push_back(std::move(pixel));
604 std::cout << m_className << "::AddPixel: Added pixel electrode.\n";
605}
std::string m_className
Class name.
Definition: Component.hh:329

◆ AddPlane()

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.

626 {
627 const auto it = std::find(m_readout.cbegin(), m_readout.cend(), label);
628 if (it == m_readout.end() && m_readout.size() > 0) {
629 std::cerr << m_className << "::AddPlane:\n"
630 << "Note that the label " << label << " is already in use.\n";
631 }
632 Electrode plate;
633 plate.label = label;
634 plate.ind = structureelectrode::Plane;
635
636 if (!anode) plate.flip = -1;
637
638 m_readout.push_back(label);
639 m_readout_p.push_back(std::move(plate));
640
641 std::cout << m_className << "::AddPlane: Added plane electrode.\n";
642}

◆ AddStrip()

void Garfield::ComponentParallelPlate::AddStrip ( double  x,
double  lx,
const std::string &  label 
)

Add strip electrode.

Definition at line 607 of file ComponentParallelPlate.cc.

608 {
609 const auto it = std::find(m_readout.cbegin(), m_readout.cend(), label);
610 if (it == m_readout.end() && m_readout.size() > 0) {
611 std::cerr << m_className << "::AddStrip:\n"
612 << "Note that the label " << label << " is already in use.\n";
613 }
614 Electrode strip;
615 strip.label = label;
616 strip.ind = structureelectrode::Strip;
617 strip.xpos = x;
618 strip.lx = lx_input;
619
620 m_readout.push_back(label);
621 m_readout_p.push_back(std::move(strip));
622
623 std::cout << m_className << "::AddStrip: Added strip electrode.\n";
624}

◆ DelayedWeightingField()

void Garfield::ComponentParallelPlate::DelayedWeightingField ( const double  x,
const double  y,
const double  z,
const double  t,
double &  wx,
double &  wy,
double &  wz,
const std::string &  label 
)
overridevirtual

Calculate the delayed weighting field at a given point and time and for a given electrode.

Parameters
x,y,zcoordinates [cm].
ttime [ns].
wx,wy,wzcomponents of the weighting field [1/cm].
labelname of the electrode

Reimplemented from Garfield::Component.

Definition at line 542 of file ComponentParallelPlate.cc.

544 {
545 wx = 0.;
546 wy = 0.;
547 wz = 0.;
548
549 if (m_sigma == 0) {
550 if (m_debug) {
551 std::cout << m_className << "::DelayedWeightingField:\n"
552 << " Conductivity is set to zero.\n";
553 }
554 return;
555 }
556
557 for (const auto& electrode : m_readout_p) {
558 if (electrode.label == label) {
559 wx = electrode.flip *
560 IntegrateDelayedField(electrode, fieldcomponent::xcomp, x, y, z, t);
561 wy = electrode.flip *
562 IntegrateDelayedField(electrode, fieldcomponent::ycomp, x, y, z, t);
563 wz = electrode.flip *
564 IntegrateDelayedField(electrode, fieldcomponent::zcomp, x, y, z, t);
565 }
566 }
567}
bool m_debug
Switch on/off debugging messages.
Definition: Component.hh:341

◆ DelayedWeightingPotential()

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

Calculate the delayed weighting potential at a given point and time and for a given electrode.

Parameters
x,y,zcoordinates [cm].
ttime [ns].
labelname of the electrode

Reimplemented from Garfield::Component.

Definition at line 515 of file ComponentParallelPlate.cc.

517 {
518 if (m_sigma == 0) {
519 if (m_debug) {
520 std::cout << m_className << "::DelayedWeightingPotential:\n"
521 << " Conductivity is set to zero.\n";
522 }
523 return 0.;
524 }
525
526 double ret = 0.;
527
528 for (const auto& electrode : m_readout_p) {
529 if (electrode.label == label) {
530 if (!electrode.m_usegrid) {
531 ret +=
532 electrode.flip * IntegrateDelayedPotential(electrode, x, y, z, t);
533 } else {
534 ret += FindDelayedWeightingPotentialInGrid(electrode, x, y, z, t);
535 }
536 }
537 }
538
539 return ret;
540}

◆ ElectricField() [1/2]

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.

Definition at line 429 of file ComponentParallelPlate.cc.

432 {
433 ex = ey = 0.;
434
435 if (z > 0.) {
436 ez = m_ezg;
437 } else {
438 ez = m_ezb;
439 }
440
441 if (m_sigma == 0) {
442 v = -m_eps * m_V * (m_g - z) / (m_b + m_eps * m_g);
443 } else {
444 v = -m_eps * m_V * (m_g - z) / (m_eps * m_g);
445 }
446
447 m = m_geometry ? m_geometry->GetMedium(x, y, z) : m_medium;
448 if (!m) {
449 if (m_debug) {
450 std::cout << m_className << "::ElectricField: No medium at (" << x << ", "
451 << y << ", " << z << ").\n";
452 }
453 status = -6;
454 return;
455 }
456
457 if (z > 0) {
458 status = 0;
459 } else {
460 status = -5;
461 }
462}
Geometry * m_geometry
Pointer to the geometry.
Definition: Component.hh:332
virtual Medium * GetMedium(const double x, const double y, const double z, const bool tesselated=false) const =0
Retrieve the medium at a given point.

◆ ElectricField() [2/2]

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.

Definition at line 399 of file ComponentParallelPlate.cc.

402 {
403 ex = ey = 0.;
404
405 if (z < 0) {
406 ez = m_ezb;
407 } else {
408 ez = m_ezb;
409 }
410
411 m = m_geometry ? m_geometry->GetMedium(x, y, z) : m_medium;
412
413 if (!m) {
414 if (m_debug) {
415 std::cout << m_className << "::ElectricField: No medium at (" << x << ", "
416 << y << ", " << z << ").\n";
417 }
418 status = -6;
419 return;
420 }
421
422 if (z > 0) {
423 status = 0;
424 } else {
425 status = -5;
426 }
427}

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

Definition at line 40 of file ComponentParallelPlate.cc.

42 {
43 // If a geometry is present, try to get the bounding box from there.
44 if (m_geometry) {
45 if (m_geometry->GetBoundingBox(x0, y0, z0, x1, y1, z1)) return true;
46 }
47 x0 = -std::numeric_limits<double>::infinity();
48 y0 = -std::numeric_limits<double>::infinity();
49 x1 = +std::numeric_limits<double>::infinity();
50 y1 = +std::numeric_limits<double>::infinity();
51 // TODO: check!
52 z0 = 0.;
53 z1 = m_g;
54 return true;
55}
virtual bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)=0
Get the bounding box (envelope of the geometry).

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

Definition at line 644 of file ComponentParallelPlate.cc.

645 {
646 if (m_geometry) {
647 return m_geometry->GetMedium(x, y, z);
648 } else if (m_medium) {
649 return m_medium;
650 }
651 return nullptr;
652}

◆ GetVoltageRange()

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

Calculate the voltage range [V].

Implements Garfield::Component.

Definition at line 464 of file ComponentParallelPlate.cc.

464 {
465 if (m_V == 0) return false;
466
467 if (m_V < 0) {
468 vmin = m_V;
469 vmax = 0;
470 } else {
471 vmin = 0;
472 vmax = m_V;
473 }
474 return true;
475}

◆ SetMedium()

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

Definition at line 69 of file ComponentParallelPlate.hh.

69{ m_medium = medium; }

◆ Setup()

void Garfield::ComponentParallelPlate::Setup ( double  g,
double  b,
double  eps,
double  v,
double  sigma = 0. 
)

Define the geometry.

Parameters
gsize of the gap along positive $z$.
bthickness of the resistive layer along negative $z$.
epsrelative permittivity of the resistive layer.
sigmaconductivity of the resistive layer (must be larger then zero, otherwise do not pass it in the function).
Vapplied potential difference between the parallel plates.

Definition at line 16 of file ComponentParallelPlate.cc.

17 {
18 // TODO: can g, b be negative?
19 m_g = g;
20 m_b = b;
21 if (eps < 1.) {
22 std::cerr << m_className << "::Setup: Epsilon must be >= 1.\n";
23 return;
24 }
25 m_eps = eps;
26 m_V = V;
27 // TODO: can sigma be negative?
28 m_sigma = sigma;
29 if (sigma == 0) {
30 m_ezg = -m_eps * m_V / (m_b + m_eps * m_g);
31 m_ezb = -m_V / (m_b + m_eps * m_g);
32 } else {
33 // For large times the resistive layer will act as a perfect conductor.
34 m_ezg = -m_V / m_g;
35 m_ezb = 0.;
36 }
37 std::cout << m_className << "::Setup: Geometry set.\n";
38}

◆ SetWeightingPotentialGrid()

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.

659 {
660 for (auto& electrode : m_readout_p) {
661 if (electrode.label == label) {
662 electrode.gridXSteps = xsteps;
663 electrode.gridYSteps = ysteps;
664 electrode.gridZSteps = zsteps;
665 electrode.gridTSteps = tsteps;
666
667 if (xsteps == 0) electrode.gridXSteps = 1;
668 if (ysteps == 0) electrode.gridYSteps = 1;
669
670 electrode.gridX0 = xmin;
671 electrode.gridY0 = ymin;
672 electrode.gridZ0 = zmin;
673 electrode.gridT0 = tmin;
674
675 electrode.gridXStepSize = (xmax - xmin) / xsteps;
676 electrode.gridYStepSize = (ymax - ymin) / ysteps;
677 electrode.gridZStepSize = (zmax - zmin) / zsteps;
678 electrode.gridTStepSize = (tmax - tmin) / tsteps;
679
680 std::vector<double> nhz(zsteps, 0);
681 std::vector<std::vector<double>> nhy(ysteps, nhz);
682 std::vector<std::vector<std::vector<double>>> nhx(xsteps, nhy);
683 electrode.gridPromptV = nhx;
684
685 std::vector<double> nht(tsteps, 0);
686 std::vector<std::vector<double>> nhzd(zsteps, nht);
687 std::vector<std::vector<std::vector<double>>> nhyd(ysteps, nhzd);
688 std::vector<std::vector<std::vector<std::vector<double>>>> nhxd(xsteps,
689 nhyd);
690 electrode.gridDelayedV = nhxd;
691
692 for (int ix = 0; ix < xsteps; ix++) {
693 for (int iy = 0; iy < xsteps; iy++) {
694 for (int iz = 0; iz < xsteps; iz++) {
695 if (iz * zsteps + zmin >= 0)
696 electrode.gridPromptV[ix][iy][iz] =
697 electrode.flip * IntegratePromptPotential(
698 electrode, ix * xsteps + xmin,
699 iy * ysteps + ymin, iz * zsteps + zmin);
700
701 for (int it = 0; it < tsteps; it++) {
702 if (iz * zsteps + zmin >= 0)
703 electrode.gridDelayedV[ix][iy][iz][it] =
704 electrode.flip * IntegrateDelayedPotential(
705 electrode, ix * xsteps + xmin,
706 iy * ysteps + ymin, iz * zsteps + zmin,
707 it * tsteps + tmin);
708 }
709 }
710 }
711 }
712
713 electrode.m_usegrid = true;
714 }
715 }
716}

Referenced by SetWeightingPotentialGrids().

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

722 {
723 for (const auto& electrode : m_readout_p) {
724 SetWeightingPotentialGrid(electrode.label, xmin, xmax, xsteps, ymin, ymax,
725 ysteps, zmin, zmax, zsteps, tmin, tmax, tsteps);
726 }
727}
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)

◆ WeightingField()

void Garfield::ComponentParallelPlate::WeightingField ( const double  x,
const double  y,
const double  z,
double &  wx,
double &  wy,
double &  wz,
const std::string &  label 
)
overridevirtual

Calculate the weighting field at a given point and for a given electrode.

Parameters
x,y,zcoordinates [cm].
wx,wy,wzcomponents of the weighting field [1/cm].
labelname of the electrode

Reimplemented from Garfield::Component.

Definition at line 477 of file ComponentParallelPlate.cc.

480 {
481 wx = 0;
482 wy = 0;
483 wz = 0;
484
485 for (const auto& electrode : m_readout_p) {
486 if (electrode.label == label) {
487 wx = electrode.flip *
488 IntegrateField(electrode, fieldcomponent::xcomp, x, y, z);
489 wy = electrode.flip *
490 IntegrateField(electrode, fieldcomponent::ycomp, x, y, z);
491 wz = electrode.flip *
492 IntegrateField(electrode, fieldcomponent::zcomp, x, y, z);
493 }
494 }
495}

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

Definition at line 497 of file ComponentParallelPlate.cc.

500 {
501 double ret = 0.;
502
503 for (const auto& electrode : m_readout_p) {
504 if (electrode.label == label) {
505 if (!electrode.m_usegrid) {
506 ret += electrode.flip * IntegratePromptPotential(electrode, x, y, z);
507 } else {
508 ret += FindWeightingPotentialInGrid(electrode, x, y, z);
509 }
510 }
511 }
512 return ret;
513}

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