![]() |
Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
|
Abstract base class for components. More...
#include <Component.hh>
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 Medium * | GetMedium (const double x, const double y, const double z) |
Get the medium at a given location (x, y, z). | |
virtual void | ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)=0 |
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 ![]() | |
void | EnablePeriodicityY (const bool on=true) |
Enable simple periodicity in the ![]() | |
void | EnablePeriodicityZ (const bool on=true) |
Enable simple periodicity in the ![]() | |
void | IsPeriodic (bool &perx, bool &pery, bool &perz) |
Return periodicity flags. | |
void | EnableMirrorPeriodicityX (const bool on=true) |
Enable mirror periodicity in the ![]() | |
void | EnableMirrorPeriodicityY (const bool on=true) |
Enable mirror periodicity in the ![]() | |
void | EnableMirrorPeriodicityZ (const bool on=true) |
Enable mirror periodicity in the ![]() | |
void | IsMirrorPeriodic (bool &perx, bool &pery, bool &perz) |
Return mirror periodicity flags. | |
void | EnableAxialPeriodicityX (const bool on=true) |
Enable axial periodicity in the ![]() | |
void | EnableAxialPeriodicityY (const bool on=true) |
Enable axial periodicity in the ![]() | |
void | EnableAxialPeriodicityZ (const bool on=true) |
Enable axial periodicity in the ![]() | |
void | IsAxiallyPeriodic (bool &perx, bool &pery, bool &perz) |
Return axial periodicity flags. | |
void | EnableRotationSymmetryX (const bool on=true) |
Enable rotation symmetry around the ![]() | |
void | EnableRotationSymmetryY (const bool on=true) |
Enable rotation symmetry around the ![]() | |
void | EnableRotationSymmetryZ (const bool on=true) |
Enable rotation symmetry around the ![]() | |
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 ![]() | |
void | EnableTriangleSymmetricXZ (const bool on=true, const bool oct=2) |
Enable triangular periodicity in the ![]() | |
void | EnableTriangleSymmetricYZ (const bool on=true, const bool oct=2) |
Enable triangular periodicity in the ![]() | |
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. | |
Geometry * | m_geometry {nullptr} |
Pointer to the geometry. | |
std::array< double, 3 > | m_b0 = {{0., 0., 0.}} |
Constant magnetic field. | |
bool | m_ready {false} |
Ready for use? | |
bool | m_debug {false} |
Switch on/off debugging messages. | |
std::array< bool, 3 > | m_periodic = {{false, false, false}} |
Simple periodicity in x, y, z. | |
std::array< bool, 3 > | m_mirrorPeriodic = {{false, false, false}} |
Mirror periodicity in x, y, z. | |
std::array< bool, 3 > | m_axiallyPeriodic = {{false, false, false}} |
Axial periodicity in x, y, z. | |
std::array< bool, 3 > | m_rotationSymmetric = {{false, false, false}} |
Rotation symmetry around x-axis, y-axis, z-axis. | |
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) |
Abstract base class for components.
Definition at line 16 of file Component.hh.
|
delete |
Default constructor.
Garfield::Component::Component | ( | const std::string & | name | ) |
Constructor.
|
virtualdefault |
Destructor.
double Garfield::Component::CellSizeX | ( | ) |
double Garfield::Component::CellSizeY | ( | ) |
double Garfield::Component::CellSizeZ | ( | ) |
|
virtual |
Reset.
Reimplemented in Garfield::ComponentGrid, and Garfield::ComponentVoxel.
|
virtual |
Create and initialise GPU Transfer class.
Reimplemented in Garfield::ComponentAnsys123, and Garfield::ComponentFieldMap.
|
virtual |
Determine whether the line between two points crosses a plane.
Reimplemented in Garfield::ComponentAnalyticField.
|
virtual |
Determine whether the line between two points crosses a wire.
x0,y0,z0 | first point [cm]. |
x1,y1,z1 | second point [cm] |
xc,yc,zc | point [cm] where the line crosses the wire or the coordinates of the wire centre. |
centre | flag whether to return the coordinates of the line-wire crossing point or of the wire centre. |
rc | radius [cm] of the wire. |
Reimplemented in Garfield::ComponentAnalyticField, and Garfield::ComponentNeBem2d.
|
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.
|
virtual |
Calculate the delayed weighting field at a given point and time and for a given electrode.
x,y,z | coordinates [cm]. |
t | time [ns]. |
wx,wy,wz | components of the weighting field [1/cm]. |
label | name of the electrode |
Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
virtual |
Calculate the delayed weighting potential at a given point and time and for a given electrode.
x,y,z | coordinates [cm]. |
t | time [ns]. |
label | name of the electrode |
Reimplemented in Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
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.
|
inline |
Switch off debugging messages.
Definition at line 359 of file Component.hh.
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).
|
pure virtual |
Calculate the drift field [V/cm] and potential [V] at (x, y, z).
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.
|
pure virtual |
Calculate the drift field at given point.
x,y,z | coordinates [cm]. |
ex,ey,ez | components of the electric field [V/cm]. |
m | pointer to the medium at this location. |
status | status flag |
Status flags:
0: Inside an active medium > 0: Inside a wire of type X -4 ... -1: On the side of a plane where no wires are -5: Inside the mesh but not in an active medium -6: Outside the mesh -10: Unknown potential type (should not occur) other: Other cases (should not occur)
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.
|
virtual |
Calculate the (drift) electrostatic potential [V] at (x, y, z).
|
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.
|
inlinevirtual |
Reimplemented in Garfield::ComponentGrid, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.
Definition at line 387 of file Component.hh.
|
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.
|
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.
|
inline |
Enable axial periodicity in the
Definition at line 291 of file Component.hh.
|
inline |
Enable axial periodicity in the
Definition at line 296 of file Component.hh.
|
inline |
Enable axial periodicity in the
Definition at line 301 of file Component.hh.
|
inline |
|
inline |
Enable mirror periodicity in the
Definition at line 269 of file Component.hh.
|
inline |
Enable mirror periodicity in the
Definition at line 274 of file Component.hh.
|
inline |
Enable mirror periodicity in the
Definition at line 279 of file Component.hh.
|
inline |
Enable simple periodicity in the
Definition at line 247 of file Component.hh.
|
inline |
Enable simple periodicity in the
Definition at line 252 of file Component.hh.
|
inline |
Enable simple periodicity in the
Definition at line 257 of file Component.hh.
|
inline |
Enable rotation symmetry around the
Definition at line 313 of file Component.hh.
|
inline |
Enable rotation symmetry around the
Definition at line 318 of file Component.hh.
|
inline |
Enable rotation symmetry around the
Definition at line 323 of file Component.hh.
|
inline |
Enable triangular periodicity in the
Definition at line 335 of file Component.hh.
|
inline |
Enable triangular periodicity in the
Definition at line 342 of file Component.hh.
|
inline |
Enable triangular periodicity in the
Definition at line 349 of file Component.hh.
|
virtual |
Get the bounding box coordinates.
Reimplemented in Garfield::ComponentAnalyticField, Garfield::ComponentChargedRing, Garfield::ComponentConstant, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem2d, Garfield::ComponentNeBem3dMap, Garfield::ComponentParallelPlate, Garfield::ComponentTcad2d, Garfield::ComponentTcad3d, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
virtual |
Get the coordinates of the elementary cell.
Reimplemented in Garfield::ComponentAnalyticField, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem2d, Garfield::ComponentTcad2d, Garfield::ComponentTcad3d, and Garfield::ComponentVoxel.
|
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.
|
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.
|
virtual |
Get the medium at a given location (x, y, z).
Reimplemented 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.
|
virtual |
Get the coordinates of a mesh node.
Reimplemented in Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentTcad2d, and Garfield::ComponentTcad3d.
|
inlinevirtual |
Return the number of mesh elements.
Reimplemented in Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentNeBem2d, Garfield::ComponentNeBem3d, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, and Garfield::ComponentTcadBase< 3 >.
Definition at line 148 of file Component.hh.
|
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.
|
pure virtual |
Calculate the voltage range [V].
Implemented in Garfield::ComponentAnalyticField, Garfield::ComponentChargedRing, Garfield::ComponentConstant, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem2d, Garfield::ComponentNeBem3d, Garfield::ComponentNeBem3dMap, Garfield::ComponentParallelPlate, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
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.
|
virtual |
Does the component have a non-zero magnetic field?
Reimplemented in Garfield::ComponentGrid, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
|
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.
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.
xc,yc | centre of the circle [cm] |
r | radius [cm] |
nI | number of intervals for the integration |
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).
x0,y0,z0 | coordinates of the starting point |
x1,y1,z1 | coordinates of the end point |
xp,yp,zp | normal vector |
nI | number of intervals for the integration |
isign | include both negative and positive contributions (0) or only contributions with a given polarity (+1,-1) |
|
private |
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.
x0,y0,z0 | coordinates of one of the corners [cm] |
dx1,dy1,dz1 | vector to one of the adjacent corners [cm] |
dx2,dy2,dz2 | vector to the other adjacent corner [cm] |
nU,nV | number of integration points in the two directions |
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.
xc,yc,zc | centre of the sphere [cm] |
r | radius of the sphere [cm] |
nI | number of integration intervals in phi and theta |
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.
|
virtual |
Determine whether a particle is inside the trap radius of a wire.
q0 | charge of the particle [in elementary charges]. |
x0,y0,z0 | position [cm] of the particle. |
xw,yw | coordinates of the wire (if applicable). |
rw | radius of the wire (if applicable). |
Reimplemented in Garfield::ComponentAnalyticField, and Garfield::ComponentNeBem2d.
|
inlinevirtual |
Does the component have a 3D field (map)?
Reimplemented in Garfield::ComponentFieldMap.
Definition at line 131 of file Component.hh.
|
inline |
Return axial periodicity flags.
Definition at line 306 of file Component.hh.
|
inline |
Return mirror periodicity flags.
Definition at line 284 of file Component.hh.
|
inline |
Return periodicity flags.
Definition at line 262 of file Component.hh.
|
inlinevirtual |
|
inline |
Return rotation symmetry flags.
Definition at line 328 of file Component.hh.
|
virtual |
Calculate the magnetic field at a given point.
x,y,z | coordinates [cm]. |
bx,by,bz | components of the magnetic field [Tesla]. |
status | status flag. |
Reimplemented in Garfield::ComponentGrid, Garfield::ComponentNeBem3dMap, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
protectedpure virtual |
Reset the component.
Implemented in Garfield::ComponentAnalyticField, Garfield::ComponentChargedRing, Garfield::ComponentConstant, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem2d, Garfield::ComponentNeBem3d, Garfield::ComponentNeBem3dMap, Garfield::ComponentParallelPlate, Garfield::ComponentTcad2d, Garfield::ComponentTcad3d, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
virtual |
Define the geometry.
void Garfield::Component::SetMagneticField | ( | const double | bx, |
const double | by, | ||
const double | bz ) |
Set a constant magnetic field.
|
inlinevirtual |
Reimplemented in Garfield::ComponentAnalyticField.
Definition at line 426 of file Component.hh.
|
protectedpure virtual |
Verify periodicities.
Implemented in Garfield::ComponentAnalyticField, Garfield::ComponentChargedRing, Garfield::ComponentConstant, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem2d, Garfield::ComponentNeBem3d, Garfield::ComponentNeBem3dMap, Garfield::ComponentParallelPlate, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
virtual |
Calculate the weighting field at a given point and for a given electrode.
x,y,z | coordinates [cm]. |
wx,wy,wz | components of the weighting field [1/cm]. |
label | name of the electrode |
Reimplemented 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.
|
virtual |
Calculate the weighting potential at a given point.
x,y,z | coordinates [cm]. |
label | name of the electrode. |
Reimplemented in Garfield::ComponentAnalyticField, Garfield::ComponentConstant, Garfield::ComponentCST, Garfield::ComponentFieldMap, Garfield::ComponentGrid, Garfield::ComponentNeBem3d, Garfield::ComponentNeBem3dMap, Garfield::ComponentParallelPlate, Garfield::ComponentTcadBase< N >, Garfield::ComponentTcadBase< 2 >, Garfield::ComponentTcadBase< 3 >, Garfield::ComponentUser, and Garfield::ComponentVoxel.
|
protected |
Axial periodicity in x, y, z.
Definition at line 451 of file Component.hh.
|
protected |
|
protected |
|
protected |
|
protected |
|
protected |
Mirror periodicity in x, y, z.
Definition at line 449 of file Component.hh.
|
protected |
Definition at line 460 of file Component.hh.
|
protected |
Simple periodicity in x, y, z.
Definition at line 447 of file Component.hh.
|
protected |
|
protected |
Rotation symmetry around x-axis, y-axis, z-axis.
Definition at line 453 of file Component.hh.
|
protected |
|
protected |
Triangle symmetric in the xy, xz, and yz plane.
Definition at line 455 of file Component.hh.
|
protected |
Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).
Definition at line 457 of file Component.hh.
|
protected |
Time steps at which the delayed weighting potentials/fields are stored.
Definition at line 462 of file Component.hh.