Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization measurement in gases and semiconductors
Loading...
Searching...
No Matches
Component.hh
Go to the documentation of this file.
1#ifndef G_COMPONENT_H
2#define G_COMPONENT_H
3
4#include <array>
5#include <cstdint>
6#include <string>
7#include <vector>
8
9namespace Garfield {
10
11class ComponentGPU;
12class Geometry;
13class Medium;
14
16class Component {
17 public:
19 Component() = delete;
21 Component(const std::string& name);
23 virtual ~Component() = default;
24
26 virtual void SetGeometry(Geometry* geo);
28 virtual void Clear();
29
31 virtual Medium* GetMedium(const double x, const double y, const double z);
32
50 virtual void ElectricField(const double x, const double y, const double z,
51 double& ex, double& ey, double& ez, Medium*& m,
52 int& status) = 0;
54 virtual void ElectricField(const double x, const double y, const double z,
55 double& ex, double& ey, double& ez, double& v,
56 Medium*& m, int& status) = 0;
57
59 std::array<double, 3> ElectricField(const double x, const double y,
60 const double z);
62 virtual double ElectricPotential(const double x, const double y,
63 const double z);
65 virtual bool GetVoltageRange(double& vmin, double& vmax) = 0;
66
72 virtual void WeightingField(const double x, const double y, const double z,
73 double& wx, double& wy, double& wz,
74 const std::string& label);
80 virtual double WeightingPotential(const double x, const double y,
81 const double z, const std::string& label);
82
85 virtual const std::vector<double>& DelayedSignalTimes(
86 const std::string& /*label*/) {
87 return m_wdtimes;
88 }
89
96 virtual void DelayedWeightingField(const double x, const double y,
97 const double z, const double t, double& wx,
98 double& wy, double& wz,
99 const std::string& label);
100
107 virtual double DelayedWeightingPotential(const double x, const double y,
108 const double z, const double t,
109 const std::string& label);
110
113 virtual void DelayedWeightingPotentials(const double x, const double y,
114 const double z,
115 const std::string& label,
116 std::vector<double>& dwp);
123 virtual void MagneticField(const double x, const double y, const double z,
124 double& bx, double& by, double& bz, int& status);
126 void SetMagneticField(const double bx, const double by, const double bz);
127
129 virtual bool IsReady() { return m_ready; }
131 virtual bool Is3d() { return true; }
132
134 virtual bool GetBoundingBox(double& xmin, double& ymin, double& zmin,
135 double& xmax, double& ymax, double& zmax);
136
138 virtual bool GetElementaryCell(double& xmin, double& ymin, double& zmin,
139 double& xmax, double& ymax, double& zmax);
140 // Get the x-length of the elementary cell.
141 double CellSizeX();
142 // Get the y-length of the elementary cell.
143 double CellSizeY();
144 // Get the z-length of the elementary cell.
145 double CellSizeZ();
146
148 virtual size_t GetNumberOfElements() const { return 0; }
150 virtual bool GetElementNodes(const size_t /*i*/,
151 std::vector<size_t>& /*nodes*/) const {
152 return false;
153 }
154
156 virtual bool GetElementRegion(const size_t /*i*/, size_t& /*mat*/,
157 bool& /*drift*/) const {
158 return false;
159 }
160
161 virtual size_t GetNumberOfNodes() const { return 0; }
163 virtual bool GetNode(const size_t i, double& x, double& y, double& z) const;
164
172 double IntegrateFluxCircle(const double xc, const double yc, const double r,
173 const unsigned int nI = 50);
181 double IntegrateFluxSphere(const double xc, const double yc, const double zc,
182 const double r, const unsigned int nI = 20);
183
193 const double x0, const double y0, const double z0, const double dx1,
194 const double dy1, const double dz1, const double dx2, const double dy2,
195 const double dz2, const unsigned int nU = 20, const unsigned int nV = 20);
196
200 const std::string& label, const double x0, const double y0,
201 const double z0, const double dx1, const double dy1, const double dz1,
202 const double dx2, const double dy2, const double dz2,
203 const unsigned int nU = 20, const unsigned int nV = 20);
204
214 double IntegrateFluxLine(const double x0, const double y0, const double z0,
215 const double x1, const double y1, const double z1,
216 const double xp, const double yp, const double zp,
217 const unsigned int nI, const int isign = 0);
218
228 virtual bool CrossedWire(const double x0, const double y0, const double z0,
229 const double x1, const double y1, const double z1,
230 double& xc, double& yc, double& zc,
231 const bool centre, double& rc);
238 virtual bool InTrapRadius(const double q0, const double x0, const double y0,
239 const double z0, double& xw, double& yw,
240 double& rw);
242 virtual bool CrossedPlane(const double x0, const double y0, const double z0,
243 const double x1, const double y1, const double z1,
244 double& xc, double& yc, double& zc);
245
247 void EnablePeriodicityX(const bool on = true) {
248 m_periodic[0] = on;
250 }
251
252 void EnablePeriodicityY(const bool on = true) {
253 m_periodic[1] = on;
255 }
256
257 void EnablePeriodicityZ(const bool on = true) {
258 m_periodic[2] = on;
260 }
261
262 void IsPeriodic(bool& perx, bool& pery, bool& perz) {
263 perx = m_periodic[0];
264 pery = m_periodic[1];
265 perz = m_periodic[2];
266 }
267
269 void EnableMirrorPeriodicityX(const bool on = true) {
270 m_mirrorPeriodic[0] = on;
272 }
273
274 void EnableMirrorPeriodicityY(const bool on = true) {
275 m_mirrorPeriodic[1] = on;
277 }
278
279 void EnableMirrorPeriodicityZ(const bool on = true) {
280 m_mirrorPeriodic[2] = on;
282 }
283
284 void IsMirrorPeriodic(bool& perx, bool& pery, bool& perz) {
285 perx = m_mirrorPeriodic[0];
286 pery = m_mirrorPeriodic[1];
287 perz = m_mirrorPeriodic[2];
288 }
289
291 void EnableAxialPeriodicityX(const bool on = true) {
292 m_axiallyPeriodic[0] = on;
294 }
295
296 void EnableAxialPeriodicityY(const bool on = true) {
297 m_axiallyPeriodic[1] = on;
299 }
300
301 void EnableAxialPeriodicityZ(const bool on = true) {
302 m_axiallyPeriodic[2] = on;
304 }
305
306 void IsAxiallyPeriodic(bool& perx, bool& pery, bool& perz) {
307 perx = m_axiallyPeriodic[0];
308 pery = m_axiallyPeriodic[1];
309 perz = m_axiallyPeriodic[2];
310 }
311
313 void EnableRotationSymmetryX(const bool on = true) {
314 m_rotationSymmetric[0] = on;
316 }
317
318 void EnableRotationSymmetryY(const bool on = true) {
319 m_rotationSymmetric[1] = on;
321 }
322
323 void EnableRotationSymmetryZ(const bool on = true) {
324 m_rotationSymmetric[2] = on;
326 }
327
328 void IsRotationSymmetric(bool& rotx, bool& roty, bool& rotz) {
329 rotx = m_rotationSymmetric[0];
330 roty = m_rotationSymmetric[1];
331 rotz = m_rotationSymmetric[2];
332 }
333
335 void EnableTriangleSymmetricXY(const bool on = true, const bool oct = 2) {
336 m_triangleSymmetric[0] = on;
338 m_mirrorPeriodic[0] = on;
339 m_mirrorPeriodic[1] = on;
340 }
341
342 void EnableTriangleSymmetricXZ(const bool on = true, const bool oct = 2) {
343 m_triangleSymmetric[1] = on;
345 m_mirrorPeriodic[0] = on;
346 m_mirrorPeriodic[2] = on;
347 }
348
349 void EnableTriangleSymmetricYZ(const bool on = true, const bool oct = 2) {
350 m_triangleSymmetric[2] = on;
352 m_mirrorPeriodic[1] = on;
353 m_mirrorPeriodic[2] = on;
354 }
355
357 void EnableDebugging(const bool on = true) { m_debug = on; }
359 void DisableDebugging() { m_debug = false; }
360
362 virtual bool HasMagneticField() const;
363
365 virtual bool HasTownsendMap() const { return false; }
367 virtual bool HasAttachmentMap() const { return false; }
369 virtual bool HasMobilityMap() const { return false; }
371 virtual bool HasVelocityMap() const { return false; }
372
374 virtual bool ElectronAttachment(const double /*x*/, const double /*y*/,
375 const double /*z*/, double& eta) {
376 eta = 0;
377 return false;
378 }
379
380 virtual bool HoleAttachment(const double /*x*/, const double /*y*/,
381 const double /*z*/, double& eta) {
382 eta = 0;
383 return false;
384 }
385
386 // Get the electron mobility coefficient.
387 virtual bool ElectronMobility(const double /*x*/, const double /*y*/,
388 const double /*z*/, double& mu) {
389 mu = 0;
390 return false;
391 }
392
393 virtual bool HoleMobility(const double /*x*/, const double /*y*/,
394 const double /*z*/, double& mu) {
395 mu = 0;
396 return false;
397 }
398
400 virtual bool ElectronTownsend(const double /*x*/, const double /*y*/,
401 const double /*z*/, double& alpha) {
402 alpha = 0;
403 return false;
404 }
405
406 virtual bool HoleTownsend(const double /*x*/, const double /*y*/,
407 const double /*z*/, double& alpha) {
408 alpha = 0;
409 return false;
410 }
411
412 virtual bool ElectronVelocity(const double /*x*/, const double /*y*/,
413 const double /*z*/, double& vx, double& vy,
414 double& vz) {
415 vx = vy = vz = 0;
416 return false;
417 }
418
419 virtual bool HoleVelocity(const double /*x*/, const double /*y*/,
420 const double /*z*/, double& vx, double& vy,
421 double& vz) {
422 vx = vy = vz = 0;
423 return false;
424 }
425
426 virtual double StepSizeHint() { return -1.; }
427
429 virtual double CreateGPUTransferObject(ComponentGPU*& comp_gpu);
430
431 protected:
433 std::string m_className{"Component"};
434
437
439 std::array<double, 3> m_b0 = {{0., 0., 0.}};
441 bool m_ready{false};
442
444 bool m_debug{false};
445
447 std::array<bool, 3> m_periodic = {{false, false, false}};
449 std::array<bool, 3> m_mirrorPeriodic = {{false, false, false}};
451 std::array<bool, 3> m_axiallyPeriodic = {{false, false, false}};
453 std::array<bool, 3> m_rotationSymmetric = {{false, false, false}};
455 std::array<bool, 3> m_triangleSymmetric = {{false, false, false}};
459 const std::array<int, 4> m_triangleOctRules = {1, 4, 5, 8};
460 bool m_outsideCone = false;
462 std::vector<double> m_wdtimes;
463
465 virtual void Reset() = 0;
467 virtual void UpdatePeriodicity() = 0;
468
469 private:
470 double IntegrateFluxParallelogram(const double x0, const double y0,
471 const double z0, const double dx1,
472 const double dy1, const double dz1,
473 const double dx2, const double dy2,
474 const double dz2, const unsigned int nU,
475 const unsigned int nV, const bool wfield,
476 const std::string& label);
477};
478
479} // namespace Garfield
480#endif
virtual double CreateGPUTransferObject(ComponentGPU *&comp_gpu)
Create and initialise GPU Transfer class.
virtual size_t GetNumberOfNodes() const
Return the number of mesh nodes.
Definition Component.hh:161
Component(const std::string &name)
Constructor.
virtual void UpdatePeriodicity()=0
Verify periodicities.
void IsAxiallyPeriodic(bool &perx, bool &pery, bool &perz)
Return axial periodicity flags.
Definition Component.hh:306
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 bool GetNode(const size_t i, double &x, double &y, double &z) const
Get the coordinates of a mesh node.
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)
virtual bool HasMagneticField() const
Does the component have a non-zero magnetic field?
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.
std::array< bool, 3 > m_rotationSymmetric
Rotation symmetry around x-axis, y-axis, z-axis.
Definition Component.hh:453
void EnableMirrorPeriodicityX(const bool on=true)
Enable mirror periodicity in the direction.
Definition Component.hh:269
virtual bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Get the bounding box coordinates.
std::array< bool, 3 > m_mirrorPeriodic
Mirror periodicity in x, y, z.
Definition Component.hh:449
virtual bool ElectronMobility(const double, const double, const double, double &mu)
Definition Component.hh:387
const std::array< int, 4 > m_triangleOctRules
Octants where |x| >= |y|.
Definition Component.hh:459
int m_triangleSymmetricOct
Triangle symmetric octant of imported map (0 < phi < Pi/4 --> octant 1).
Definition Component.hh:457
virtual void Reset()=0
Reset the component.
void DisableDebugging()
Switch off debugging messages.
Definition Component.hh:359
void EnableTriangleSymmetricXZ(const bool on=true, const bool oct=2)
Enable triangular periodicity in the plane.
Definition Component.hh:342
void EnablePeriodicityX(const bool on=true)
Enable simple periodicity in the direction.
Definition Component.hh:247
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.
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.
void IsRotationSymmetric(bool &rotx, bool &roty, bool &rotz)
Return rotation symmetry flags.
Definition Component.hh:328
virtual bool GetElementNodes(const size_t, std::vector< size_t > &) const
Get the indices of the nodes constituting a given element.
Definition Component.hh:150
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).
virtual bool ElectronVelocity(const double, const double, const double, double &vx, double &vy, double &vz)
Get the electron drift velocity.
Definition Component.hh:412
void EnablePeriodicityY(const bool on=true)
Enable simple periodicity in the direction.
Definition Component.hh:252
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,...
bool m_debug
Switch on/off debugging messages.
Definition Component.hh:444
virtual bool Is3d()
Does the component have a 3D field (map)?
Definition Component.hh:131
virtual void Clear()
Reset.
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.
virtual bool ElectronAttachment(const double, const double, const double, double &eta)
Get the electron attachment coefficient.
Definition Component.hh:374
void EnablePeriodicityZ(const bool on=true)
Enable simple periodicity in the direction.
Definition Component.hh:257
virtual bool GetVoltageRange(double &vmin, double &vmax)=0
Calculate the voltage range [V].
Component()=delete
Default constructor.
void EnableRotationSymmetryY(const bool on=true)
Enable rotation symmetry around the axis.
Definition Component.hh:318
std::array< bool, 3 > m_periodic
Simple periodicity in x, y, z.
Definition Component.hh:447
void EnableAxialPeriodicityY(const bool on=true)
Enable axial periodicity in the direction.
Definition Component.hh:296
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.
void EnableTriangleSymmetricYZ(const bool on=true, const bool oct=2)
Enable triangular periodicity in the plane.
Definition Component.hh:349
virtual bool ElectronTownsend(const double, const double, const double, double &alpha)
Get the electron Townsend coefficient.
Definition Component.hh:400
virtual bool HasMobilityMap() const
Does the component have maps of the low-field mobility?
Definition Component.hh:369
virtual bool HasVelocityMap() const
Does the component have velocity maps?
Definition Component.hh:371
void EnableRotationSymmetryZ(const bool on=true)
Enable rotation symmetry around the axis.
Definition Component.hh:323
std::vector< double > m_wdtimes
Time steps at which the delayed weighting potentials/fields are stored.
Definition Component.hh:462
void EnableAxialPeriodicityX(const bool on=true)
Enable axial periodicity in the direction.
Definition Component.hh:291
void EnableMirrorPeriodicityZ(const bool on=true)
Enable mirror periodicity in the direction.
Definition Component.hh:279
void IsMirrorPeriodic(bool &perx, bool &pery, bool &perz)
Return mirror periodicity flags.
Definition Component.hh:284
void EnableDebugging(const bool on=true)
Switch on debugging messages.
Definition Component.hh:357
void EnableMirrorPeriodicityY(const bool on=true)
Enable mirror periodicity in the direction.
Definition Component.hh:274
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 bool HasTownsendMap() const
Does the component have maps of the Townsend coefficient?
Definition Component.hh:365
virtual bool HoleAttachment(const double, const double, const double, double &eta)
Get the hole attachment coefficient.
Definition Component.hh:380
virtual bool GetElementaryCell(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
Get the coordinates of the elementary cell.
virtual bool HoleMobility(const double, const double, const double, double &mu)
Get the hole Mobility coefficient.
Definition Component.hh:393
std::string m_className
Class name.
Definition Component.hh:433
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.
void EnableTriangleSymmetricXY(const bool on=true, const bool oct=2)
Enable triangular periodicity in the plane.
Definition Component.hh:335
Geometry * m_geometry
Pointer to the geometry.
Definition Component.hh:436
virtual bool HoleTownsend(const double, const double, const double, double &alpha)
Get the hole Townsend coefficient.
Definition Component.hh:406
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.
void IsPeriodic(bool &perx, bool &pery, bool &perz)
Return periodicity flags.
Definition Component.hh:262
bool m_ready
Ready for use?
Definition Component.hh:441
void SetMagneticField(const double bx, const double by, const double bz)
Set a constant magnetic field.
virtual bool HoleVelocity(const double, const double, const double, double &vx, double &vy, double &vz)
Get the hole drift velocity.
Definition Component.hh:419
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.
Definition Component.hh:156
std::array< double, 3 > m_b0
Constant magnetic field.
Definition Component.hh:439
std::array< bool, 3 > m_axiallyPeriodic
Axial periodicity in x, y, z.
Definition Component.hh:451
virtual bool IsReady()
Ready for use?
Definition Component.hh:129
std::array< bool, 3 > m_triangleSymmetric
Triangle symmetric in the xy, xz, and yz plane.
Definition Component.hh:455
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 StepSizeHint()
Definition Component.hh:426
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.
virtual bool HasAttachmentMap() const
Does the component have maps of the attachment coefficient?
Definition Component.hh:367
void EnableAxialPeriodicityZ(const bool on=true)
Enable axial periodicity in the direction.
Definition Component.hh:301
virtual Medium * GetMedium(const double x, const double y, const double z)
Get the medium at a given location (x, y, z).
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.
virtual ~Component()=default
Destructor.
virtual void SetGeometry(Geometry *geo)
Define the geometry.
virtual const std::vector< double > & DelayedSignalTimes(const std::string &)
Return the time steps at which the delayed weighting potential/field are stored/evaluated.
Definition Component.hh:85
void EnableRotationSymmetryX(const bool on=true)
Enable rotation symmetry around the axis.
Definition Component.hh:313
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,...
virtual size_t GetNumberOfElements() const
Return the number of mesh elements.
Definition Component.hh:148
virtual double ElectricPotential(const double x, const double y, const double z)
Calculate the (drift) electrostatic potential [V] at (x, y, z).
Abstract base class for geometry classes.
Definition Geometry.hh:13
Abstract base class for components.
Definition Medium.hh:17