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

Abstract base class for components. More...

#include <Component.hh>

Inheritance diagram for Garfield::Component:
Garfield::ComponentTcadBase< 2 > Garfield::ComponentTcadBase< 3 > Garfield::ComponentAnalyticField Garfield::ComponentChargedRing Garfield::ComponentConstant Garfield::ComponentFieldMap Garfield::ComponentGrid Garfield::ComponentNeBem2d Garfield::ComponentNeBem3d Garfield::ComponentNeBem3dMap Garfield::ComponentParallelPlate Garfield::ComponentTcadBase< N > Garfield::ComponentUser Garfield::ComponentVoxel

Public Member Functions

 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.
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
 Calculate the drift field at given point.
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).
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 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)
 Calculate the weighting field at a given point and for a given electrode.
virtual double WeightingPotential (const double x, const double y, const double z, const std::string &label)
 Calculate the weighting potential at a given point.
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 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 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.

Protected Member Functions

virtual void Reset ()=0
 Reset the component.
virtual void UpdatePeriodicity ()=0
 Verify periodicities.

Protected Attributes

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.

Private Member Functions

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, const unsigned int nV, const bool wfield, const std::string &label)

Detailed Description

Abstract base class for components.

Definition at line 16 of file Component.hh.

Constructor & Destructor Documentation

◆ Component() [1/2]

Garfield::Component::Component ( )
delete

Default constructor.

◆ Component() [2/2]

Garfield::Component::Component ( const std::string & name)

Constructor.

◆ ~Component()

virtual Garfield::Component::~Component ( )
virtualdefault

Destructor.

Member Function Documentation

◆ CellSizeX()

double Garfield::Component::CellSizeX ( )

◆ CellSizeY()

double Garfield::Component::CellSizeY ( )

◆ CellSizeZ()

double Garfield::Component::CellSizeZ ( )

◆ Clear()

virtual void Garfield::Component::Clear ( )
virtual

Reset.

Reimplemented in Garfield::ComponentGrid, and Garfield::ComponentVoxel.

◆ CreateGPUTransferObject()

virtual double Garfield::Component::CreateGPUTransferObject ( ComponentGPU *& comp_gpu)
virtual

Create and initialise GPU Transfer class.

Reimplemented in Garfield::ComponentAnsys123, and Garfield::ComponentFieldMap.

◆ CrossedPlane()

virtual bool Garfield::Component::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 )
virtual

Determine whether the line between two points crosses a plane.

Reimplemented in Garfield::ComponentAnalyticField.

◆ CrossedWire()

virtual bool Garfield::Component::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 )
virtual

Determine whether the line between two points crosses a wire.

Parameters
x0,y0,z0first point [cm].
x1,y1,z1second point [cm]
xc,yc,zcpoint [cm] where the line crosses the wire or the coordinates of the wire centre.
centreflag whether to return the coordinates of the line-wire crossing point or of the wire centre.
rcradius [cm] of the wire.

Reimplemented in Garfield::ComponentAnalyticField, and Garfield::ComponentNeBem2d.

◆ DelayedSignalTimes()

virtual const std::vector< double > & Garfield::Component::DelayedSignalTimes ( const std::string & )
inlinevirtual

Return the time steps at which the delayed weighting potential/field are stored/evaluated.

Reimplemented in Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 85 of file Component.hh.

86 {
87 return m_wdtimes;
88 }
std::vector< double > m_wdtimes
Time steps at which the delayed weighting potentials/fields are stored.
Definition Component.hh:462

◆ DelayedWeightingField()

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

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 in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.

◆ DelayedWeightingPotential()

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

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 in Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.

◆ DelayedWeightingPotentials()

virtual void Garfield::Component::DelayedWeightingPotentials ( const double x,
const double y,
const double z,
const std::string & label,
std::vector< double > & dwp )
virtual

Calculate the delayed weighting potentials at a given point and for a given electrode, for a set of pre-defined times.

Reimplemented in Garfield::ComponentFieldMap, Garfield::ComponentTcad2d, and Garfield::ComponentTcad3d.

◆ DisableDebugging()

void Garfield::Component::DisableDebugging ( )
inline

Switch off debugging messages.

Definition at line 359 of file Component.hh.

359{ m_debug = false; }
bool m_debug
Switch on/off debugging messages.
Definition Component.hh:444

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

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

◆ ElectricField() [3/3]

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

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)

Implemented in Garfield::ComponentAnalyticField, Garfield::ComponentChargedRing, Garfield::ComponentConstant, Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem2d, Garfield::ComponentNeBem3d, Garfield::ComponentNeBem3dMap, Garfield::ComponentParallelPlate, Garfield::ComponentTcad2d, Garfield::ComponentTcad3d, Garfield::ComponentUser, and Garfield::ComponentVoxel.

◆ ElectricPotential()

virtual double Garfield::Component::ElectricPotential ( const double x,
const double y,
const double z )
virtual

Calculate the (drift) electrostatic potential [V] at (x, y, z).

◆ ElectronAttachment()

virtual bool Garfield::Component::ElectronAttachment ( const double ,
const double ,
const double ,
double & eta )
inlinevirtual

Get the electron attachment coefficient.

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 374 of file Component.hh.

375 {
376 eta = 0;
377 return false;
378 }

◆ ElectronMobility()

virtual bool Garfield::Component::ElectronMobility ( const double ,
const double ,
const double ,
double & mu )
inlinevirtual

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 387 of file Component.hh.

388 {
389 mu = 0;
390 return false;
391 }

◆ ElectronTownsend()

virtual bool Garfield::Component::ElectronTownsend ( const double ,
const double ,
const double ,
double & alpha )
inlinevirtual

Get the electron Townsend coefficient.

Reimplemented in Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 400 of file Component.hh.

401 {
402 alpha = 0;
403 return false;
404 }

◆ ElectronVelocity()

virtual bool Garfield::Component::ElectronVelocity ( const double ,
const double ,
const double ,
double & vx,
double & vy,
double & vz )
inlinevirtual

Get the electron drift velocity.

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 412 of file Component.hh.

414 {
415 vx = vy = vz = 0;
416 return false;
417 }

◆ EnableAxialPeriodicityX()

void Garfield::Component::EnableAxialPeriodicityX ( const bool on = true)
inline

Enable axial periodicity in the $x$ direction.

Definition at line 291 of file Component.hh.

291 {
292 m_axiallyPeriodic[0] = on;
294 }
virtual void UpdatePeriodicity()=0
Verify periodicities.
std::array< bool, 3 > m_axiallyPeriodic
Axial periodicity in x, y, z.
Definition Component.hh:451

◆ EnableAxialPeriodicityY()

void Garfield::Component::EnableAxialPeriodicityY ( const bool on = true)
inline

Enable axial periodicity in the $y$ direction.

Definition at line 296 of file Component.hh.

296 {
297 m_axiallyPeriodic[1] = on;
299 }

◆ EnableAxialPeriodicityZ()

void Garfield::Component::EnableAxialPeriodicityZ ( const bool on = true)
inline

Enable axial periodicity in the $z$ direction.

Definition at line 301 of file Component.hh.

301 {
302 m_axiallyPeriodic[2] = on;
304 }

◆ EnableDebugging()

void Garfield::Component::EnableDebugging ( const bool on = true)
inline

Switch on debugging messages.

Definition at line 357 of file Component.hh.

357{ m_debug = on; }

◆ EnableMirrorPeriodicityX()

void Garfield::Component::EnableMirrorPeriodicityX ( const bool on = true)
inline

Enable mirror periodicity in the $x$ direction.

Definition at line 269 of file Component.hh.

269 {
270 m_mirrorPeriodic[0] = on;
272 }
std::array< bool, 3 > m_mirrorPeriodic
Mirror periodicity in x, y, z.
Definition Component.hh:449

◆ EnableMirrorPeriodicityY()

void Garfield::Component::EnableMirrorPeriodicityY ( const bool on = true)
inline

Enable mirror periodicity in the $y$ direction.

Definition at line 274 of file Component.hh.

274 {
275 m_mirrorPeriodic[1] = on;
277 }

◆ EnableMirrorPeriodicityZ()

void Garfield::Component::EnableMirrorPeriodicityZ ( const bool on = true)
inline

Enable mirror periodicity in the $y$ direction.

Definition at line 279 of file Component.hh.

279 {
280 m_mirrorPeriodic[2] = on;
282 }

◆ EnablePeriodicityX()

void Garfield::Component::EnablePeriodicityX ( const bool on = true)
inline

Enable simple periodicity in the $x$ direction.

Definition at line 247 of file Component.hh.

247 {
248 m_periodic[0] = on;
250 }
std::array< bool, 3 > m_periodic
Simple periodicity in x, y, z.
Definition Component.hh:447

◆ EnablePeriodicityY()

void Garfield::Component::EnablePeriodicityY ( const bool on = true)
inline

Enable simple periodicity in the $y$ direction.

Definition at line 252 of file Component.hh.

252 {
253 m_periodic[1] = on;
255 }

◆ EnablePeriodicityZ()

void Garfield::Component::EnablePeriodicityZ ( const bool on = true)
inline

Enable simple periodicity in the $z$ direction.

Definition at line 257 of file Component.hh.

257 {
258 m_periodic[2] = on;
260 }

◆ EnableRotationSymmetryX()

void Garfield::Component::EnableRotationSymmetryX ( const bool on = true)
inline

Enable rotation symmetry around the $x$ axis.

Definition at line 313 of file Component.hh.

313 {
314 m_rotationSymmetric[0] = on;
316 }
std::array< bool, 3 > m_rotationSymmetric
Rotation symmetry around x-axis, y-axis, z-axis.
Definition Component.hh:453

◆ EnableRotationSymmetryY()

void Garfield::Component::EnableRotationSymmetryY ( const bool on = true)
inline

Enable rotation symmetry around the $y$ axis.

Definition at line 318 of file Component.hh.

318 {
319 m_rotationSymmetric[1] = on;
321 }

◆ EnableRotationSymmetryZ()

void Garfield::Component::EnableRotationSymmetryZ ( const bool on = true)
inline

Enable rotation symmetry around the $z$ axis.

Definition at line 323 of file Component.hh.

323 {
324 m_rotationSymmetric[2] = on;
326 }

◆ EnableTriangleSymmetricXY()

void Garfield::Component::EnableTriangleSymmetricXY ( const bool on = true,
const bool oct = 2 )
inline

Enable triangular periodicity in the $xy$ plane.

Definition at line 335 of file Component.hh.

335 {
336 m_triangleSymmetric[0] = on;
338 m_mirrorPeriodic[0] = on;
339 m_mirrorPeriodic[1] = on;
340 }
int m_triangleSymmetricOct
Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).
Definition Component.hh:457
std::array< bool, 3 > m_triangleSymmetric
Triangle symmetric in the xy, xz, and yz plane.
Definition Component.hh:455

◆ EnableTriangleSymmetricXZ()

void Garfield::Component::EnableTriangleSymmetricXZ ( const bool on = true,
const bool oct = 2 )
inline

Enable triangular periodicity in the $xz$ plane.

Definition at line 342 of file Component.hh.

342 {
343 m_triangleSymmetric[1] = on;
345 m_mirrorPeriodic[0] = on;
346 m_mirrorPeriodic[2] = on;
347 }

◆ EnableTriangleSymmetricYZ()

void Garfield::Component::EnableTriangleSymmetricYZ ( const bool on = true,
const bool oct = 2 )
inline

Enable triangular periodicity in the $yz$ plane.

Definition at line 349 of file Component.hh.

349 {
350 m_triangleSymmetric[2] = on;
352 m_mirrorPeriodic[1] = on;
353 m_mirrorPeriodic[2] = on;
354 }

◆ GetBoundingBox()

virtual bool Garfield::Component::GetBoundingBox ( double & xmin,
double & ymin,
double & zmin,
double & xmax,
double & ymax,
double & zmax )
virtual

◆ GetElementaryCell()

virtual bool Garfield::Component::GetElementaryCell ( double & xmin,
double & ymin,
double & zmin,
double & xmax,
double & ymax,
double & zmax )
virtual

◆ GetElementNodes()

virtual bool Garfield::Component::GetElementNodes ( const size_t ,
std::vector< size_t > &  ) const
inlinevirtual

Get the indices of the nodes constituting a given element.

Reimplemented in Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 150 of file Component.hh.

151 {
152 return false;
153 }

◆ GetElementRegion()

virtual bool Garfield::Component::GetElementRegion ( const size_t ,
size_t & ,
bool &  ) const
inlinevirtual

Get the region/material of a mesh element and a flag whether it is associated to an active medium.

Reimplemented in Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 156 of file Component.hh.

157 {
158 return false;
159 }

◆ GetMedium()

◆ GetNode()

virtual bool Garfield::Component::GetNode ( const size_t i,
double & x,
double & y,
double & z ) const
virtual

Get the coordinates of a mesh node.

Reimplemented in Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentTcad2d, and Garfield::ComponentTcad3d.

◆ GetNumberOfElements()

virtual size_t Garfield::Component::GetNumberOfElements ( ) const
inlinevirtual

◆ GetNumberOfNodes()

virtual size_t Garfield::Component::GetNumberOfNodes ( ) const
inlinevirtual

Return the number of mesh nodes.

Reimplemented in Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 161 of file Component.hh.

161{ return 0; }

◆ GetVoltageRange()

◆ HasAttachmentMap()

virtual bool Garfield::Component::HasAttachmentMap ( ) const
inlinevirtual

Does the component have maps of the attachment coefficient?

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 367 of file Component.hh.

367{ return false; }

◆ HasMagneticField()

virtual bool Garfield::Component::HasMagneticField ( ) const
virtual

Does the component have a non-zero magnetic field?

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentUser, and Garfield::ComponentVoxel.

◆ HasMobilityMap()

virtual bool Garfield::Component::HasMobilityMap ( ) const
inlinevirtual

Does the component have maps of the low-field mobility?

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 369 of file Component.hh.

369{ return false; }

◆ HasTownsendMap()

virtual bool Garfield::Component::HasTownsendMap ( ) const
inlinevirtual

Does the component have maps of the Townsend coefficient?

Reimplemented in Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 365 of file Component.hh.

365{ return false; }

◆ HasVelocityMap()

virtual bool Garfield::Component::HasVelocityMap ( ) const
inlinevirtual

Does the component have velocity maps?

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 371 of file Component.hh.

371{ return false; }

◆ HoleAttachment()

virtual bool Garfield::Component::HoleAttachment ( const double ,
const double ,
const double ,
double & eta )
inlinevirtual

Get the hole attachment coefficient.

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 380 of file Component.hh.

381 {
382 eta = 0;
383 return false;
384 }

◆ HoleMobility()

virtual bool Garfield::Component::HoleMobility ( const double ,
const double ,
const double ,
double & mu )
inlinevirtual

Get the hole Mobility coefficient.

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 393 of file Component.hh.

394 {
395 mu = 0;
396 return false;
397 }

◆ HoleTownsend()

virtual bool Garfield::Component::HoleTownsend ( const double ,
const double ,
const double ,
double & alpha )
inlinevirtual

Get the hole Townsend coefficient.

Reimplemented in Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 406 of file Component.hh.

407 {
408 alpha = 0;
409 return false;
410 }

◆ HoleVelocity()

virtual bool Garfield::Component::HoleVelocity ( const double ,
const double ,
const double ,
double & vx,
double & vy,
double & vz )
inlinevirtual

Get the hole drift velocity.

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.

Definition at line 419 of file Component.hh.

421 {
422 vx = vy = vz = 0;
423 return false;
424 }

◆ IntegrateFluxCircle()

double Garfield::Component::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.

Parameters
xc,yccentre of the circle [cm]
rradius [cm]
nInumber of intervals for the integration
Returns
charge enclosed in the circle [fC / cm]

◆ IntegrateFluxLine()

double Garfield::Component::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).

Parameters
x0,y0,z0coordinates of the starting point
x1,y1,z1coordinates of the end point
xp,yp,zpnormal vector
nInumber of intervals for the integration
isigninclude both negative and positive contributions (0) or only contributions with a given polarity (+1,-1)

◆ IntegrateFluxParallelogram() [1/2]

double Garfield::Component::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,
const unsigned int nV,
const bool wfield,
const std::string & label )
private

◆ IntegrateFluxParallelogram() [2/2]

double Garfield::Component::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.

Parameters
x0,y0,z0coordinates of one of the corners [cm]
dx1,dy1,dz1vector to one of the adjacent corners [cm]
dx2,dy2,dz2vector to the other adjacent corner [cm]
nU,nVnumber of integration points in the two directions
Returns
flux [V cm]

◆ IntegrateFluxSphere()

double Garfield::Component::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.

Parameters
xc,yc,zccentre of the sphere [cm]
rradius of the sphere [cm]
nInumber of integration intervals in phi and theta
Returns
charge enclosed in the sphere [fC]

◆ IntegrateWeightingFluxParallelogram()

double Garfield::Component::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.

◆ InTrapRadius()

virtual bool Garfield::Component::InTrapRadius ( const double q0,
const double x0,
const double y0,
const double z0,
double & xw,
double & yw,
double & rw )
virtual

Determine whether a particle is inside the trap radius of a wire.

Parameters
q0charge of the particle [in elementary charges].
x0,y0,z0position [cm] of the particle.
xw,ywcoordinates of the wire (if applicable).
rwradius of the wire (if applicable).

Reimplemented in Garfield::ComponentAnalyticField, and Garfield::ComponentNeBem2d.

◆ Is3d()

virtual bool Garfield::Component::Is3d ( )
inlinevirtual

Does the component have a 3D field (map)?

Reimplemented in Garfield::ComponentFieldMap.

Definition at line 131 of file Component.hh.

131{ return true; }

◆ IsAxiallyPeriodic()

void Garfield::Component::IsAxiallyPeriodic ( bool & perx,
bool & pery,
bool & perz )
inline

Return axial periodicity flags.

Definition at line 306 of file Component.hh.

306 {
307 perx = m_axiallyPeriodic[0];
308 pery = m_axiallyPeriodic[1];
309 perz = m_axiallyPeriodic[2];
310 }

◆ IsMirrorPeriodic()

void Garfield::Component::IsMirrorPeriodic ( bool & perx,
bool & pery,
bool & perz )
inline

Return mirror periodicity flags.

Definition at line 284 of file Component.hh.

284 {
285 perx = m_mirrorPeriodic[0];
286 pery = m_mirrorPeriodic[1];
287 perz = m_mirrorPeriodic[2];
288 }

◆ IsPeriodic()

void Garfield::Component::IsPeriodic ( bool & perx,
bool & pery,
bool & perz )
inline

Return periodicity flags.

Definition at line 262 of file Component.hh.

262 {
263 perx = m_periodic[0];
264 pery = m_periodic[1];
265 perz = m_periodic[2];
266 }

◆ IsReady()

virtual bool Garfield::Component::IsReady ( )
inlinevirtual

Ready for use?

Definition at line 129 of file Component.hh.

129{ return m_ready; }
bool m_ready
Ready for use?
Definition Component.hh:441

◆ IsRotationSymmetric()

void Garfield::Component::IsRotationSymmetric ( bool & rotx,
bool & roty,
bool & rotz )
inline

Return rotation symmetry flags.

Definition at line 328 of file Component.hh.

328 {
329 rotx = m_rotationSymmetric[0];
330 roty = m_rotationSymmetric[1];
331 rotz = m_rotationSymmetric[2];
332 }

◆ MagneticField()

virtual void Garfield::Component::MagneticField ( const double x,
const double y,
const double z,
double & bx,
double & by,
double & bz,
int & status )
virtual

Calculate the magnetic field at a given point.

Parameters
x,y,zcoordinates [cm].
bx,by,bzcomponents of the magnetic field [Tesla].
statusstatus flag.

Reimplemented in Garfield::ComponentGrid, Garfield::ComponentNeBem3dMap, Garfield::ComponentUser, and Garfield::ComponentVoxel.

◆ Reset()

◆ SetGeometry()

virtual void Garfield::Component::SetGeometry ( Geometry * geo)
virtual

Define the geometry.

◆ SetMagneticField()

void Garfield::Component::SetMagneticField ( const double bx,
const double by,
const double bz )

Set a constant magnetic field.

◆ StepSizeHint()

virtual double Garfield::Component::StepSizeHint ( )
inlinevirtual

Reimplemented in Garfield::ComponentAnalyticField.

Definition at line 426 of file Component.hh.

426{ return -1.; }

◆ UpdatePeriodicity()

◆ WeightingField()

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

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 in Garfield::ComponentAnalyticField, Garfield::ComponentConstant, Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem3d, Garfield::ComponentNeBem3dMap, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.

◆ WeightingPotential()

virtual double Garfield::Component::WeightingPotential ( const double x,
const double y,
const double z,
const std::string & label )
virtual

Member Data Documentation

◆ m_axiallyPeriodic

std::array<bool, 3> Garfield::Component::m_axiallyPeriodic = {{false, false, false}}
protected

Axial periodicity in x, y, z.

Definition at line 451 of file Component.hh.

451{{false, false, false}};

◆ m_b0

std::array<double, 3> Garfield::Component::m_b0 = {{0., 0., 0.}}
protected

Constant magnetic field.

Definition at line 439 of file Component.hh.

439{{0., 0., 0.}};

◆ m_className

std::string Garfield::Component::m_className {"Component"}
protected

Class name.

Definition at line 433 of file Component.hh.

433{"Component"};

◆ m_debug

bool Garfield::Component::m_debug {false}
protected

Switch on/off debugging messages.

Definition at line 444 of file Component.hh.

444{false};

◆ m_geometry

Geometry* Garfield::Component::m_geometry {nullptr}
protected

Pointer to the geometry.

Definition at line 436 of file Component.hh.

436{nullptr};

◆ m_mirrorPeriodic

std::array<bool, 3> Garfield::Component::m_mirrorPeriodic = {{false, false, false}}
protected

Mirror periodicity in x, y, z.

Definition at line 449 of file Component.hh.

449{{false, false, false}};

◆ m_outsideCone

bool Garfield::Component::m_outsideCone = false
protected

Definition at line 460 of file Component.hh.

◆ m_periodic

std::array<bool, 3> Garfield::Component::m_periodic = {{false, false, false}}
protected

Simple periodicity in x, y, z.

Definition at line 447 of file Component.hh.

447{{false, false, false}};

◆ m_ready

bool Garfield::Component::m_ready {false}
protected

Ready for use?

Definition at line 441 of file Component.hh.

441{false};

◆ m_rotationSymmetric

std::array<bool, 3> Garfield::Component::m_rotationSymmetric = {{false, false, false}}
protected

Rotation symmetry around x-axis, y-axis, z-axis.

Definition at line 453 of file Component.hh.

453{{false, false, false}};

◆ m_triangleOctRules

const std::array<int, 4> Garfield::Component::m_triangleOctRules = {1, 4, 5, 8}
protected

Octants where |x| >= |y|.

Definition at line 459 of file Component.hh.

459{1, 4, 5, 8};

◆ m_triangleSymmetric

std::array<bool, 3> Garfield::Component::m_triangleSymmetric = {{false, false, false}}
protected

Triangle symmetric in the xy, xz, and yz plane.

Definition at line 455 of file Component.hh.

455{{false, false, false}};

◆ m_triangleSymmetricOct

int Garfield::Component::m_triangleSymmetricOct = 0
protected

Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).

Definition at line 457 of file Component.hh.

◆ m_wdtimes

std::vector<double> Garfield::Component::m_wdtimes
protected

Time steps at which the delayed weighting potentials/fields are stored.

Definition at line 462 of file Component.hh.


The documentation for this class was generated from the following file:
  • /builds/garfield/docs/source/Include/Garfield/Component.hh