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// Include this header if we're compiling with the GPU or this is the first time
2// without
3#if defined(__GPUCOMPILE__) || !defined(G_COMPONENT_H)
4
5#if !defined(__GPUCOMPILE__) && !defined(G_COMPONENT_H)
6#define G_COMPONENT_H
7#endif
8
10
11#ifdef __GPUCOMPILE__
12#else
13
14#include <array>
15#include <string>
16#include <vector>
17
18#endif
19
20namespace Garfield {
21
22// setup class names depending on if this is compiling the GPU static version or
23// not
24#if !defined(__GPUCOMPILE__)
25class ComponentGPU;
26class Geometry;
27class Medium;
28#endif
29
31class GARFIELD_CLASS_NAME(Component) {
32 public:
33#ifdef __GPUCOMPILE__
34 GARFIELD_CLASS_NAME(Component)() = default;
35#else
37 GARFIELD_CLASS_NAME(Component)() = delete;
39 GARFIELD_CLASS_NAME(Component)(const std::string& name);
40#endif
42 virtual ~GARFIELD_CLASS_NAME(Component)() {};
43
44#ifndef __GPUCOMPILE__
46 virtual void SetGeometry(Geometry* geo);
48 virtual void Clear();
49
51 virtual Medium* GetMedium(const double x, const double y, const double z);
52
70 virtual void ElectricField(const double x, const double y, const double z,
71 double& ex, double& ey, double& ez, Medium*& m,
72 int& status) = 0;
74 virtual void ElectricField(const double x, const double y, const double z,
75 double& ex, double& ey, double& ez, double& v,
76 Medium*& m, int& status) = 0;
77#else
78 __device__ void ElectricField(const cuda_t xin, const cuda_t yin,
79 const cuda_t zin, cuda_t& ex, cuda_t& ey,
80 cuda_t& ez, MediumGPU*& m, int& status);
81#endif
82
83#ifndef __GPUCOMPILE__
85 std::array<double, 3> ElectricField(const double x, const double y,
86 const double z);
88 virtual double ElectricPotential(const double x, const double y,
89 const double z);
91 virtual bool GetVoltageRange(double& vmin, double& vmax) = 0;
92
98 virtual void WeightingField(const double x, const double y, const double z,
99 double& wx, double& wy, double& wz,
100 const std::string& label);
106 virtual double WeightingPotential(const double x, const double y,
107 const double z, const std::string& label);
108
112 virtual const std::vector<double>& DelayedSignalTimes(
113 const std::string& /*label*/) {
114 return m_wdtimes;
115 }
123 virtual void DelayedWeightingField(const double x, const double y,
124 const double z, const double t, double& wx,
125 double& wy, double& wz,
126 const std::string& label);
127
134 virtual double DelayedWeightingPotential(const double x, const double y,
135 const double z, const double t,
136 const std::string& label);
137
141 virtual void DelayedWeightingPotentials(const double x, const double y,
142 const double z,
143 const std::string& label,
144 std::vector<double>& dwp);
151 virtual void MagneticField(const double x, const double y, const double z,
152 double& bx, double& by, double& bz, int& status);
154 void SetMagneticField(const double bx, const double by, const double bz);
155
157 virtual bool IsReady() { return m_ready; }
159 virtual bool Is3d() { return true; }
160
162 virtual bool GetBoundingBox(double& xmin, double& ymin, double& zmin,
163 double& xmax, double& ymax, double& zmax);
164
166 virtual bool GetElementaryCell(double& xmin, double& ymin, double& zmin,
167 double& xmax, double& ymax, double& zmax);
168 // Get the x-length of the elementary cell.
169 double CellSizeX();
170 // Get the y-length of the elementary cell.
171 double CellSizeY();
172 // Get the z-length of the elementary cell.
173 double CellSizeZ();
174
176 virtual size_t GetNumberOfElements() const { return 0; }
178 virtual bool GetElementNodes(const size_t /*i*/,
179 std::vector<size_t>& /*nodes*/) const {
180 return false;
181 }
184 virtual bool GetElementRegion(const size_t /*i*/, size_t& /*mat*/,
185 bool& /*drift*/) const {
186 return false;
187 }
189 virtual size_t GetNumberOfNodes() const { return 0; }
191 virtual bool GetNode(const size_t i, double& x, double& y, double& z) const;
192
200 double IntegrateFluxCircle(const double xc, const double yc, const double r,
201 const unsigned int nI = 50);
209 double IntegrateFluxSphere(const double xc, const double yc, const double zc,
210 const double r, const unsigned int nI = 20);
211
220 double IntegrateFluxParallelogram(
221 const double x0, const double y0, const double z0, const double dx1,
222 const double dy1, const double dz1, const double dx2, const double dy2,
223 const double dz2, const unsigned int nU = 20, const unsigned int nV = 20);
224
227 double IntegrateWeightingFluxParallelogram(
228 const std::string& label, const double x0, const double y0,
229 const double z0, const double dx1, const double dy1, const double dz1,
230 const double dx2, const double dy2, const double dz2,
231 const unsigned int nU = 20, const unsigned int nV = 20);
232
242 double IntegrateFluxLine(const double x0, const double y0, const double z0,
243 const double x1, const double y1, const double z1,
244 const double xp, const double yp, const double zp,
245 const unsigned int nI, const int isign = 0);
246
256 virtual bool CrossedWire(const double x0, const double y0, const double z0,
257 const double x1, const double y1, const double z1,
258 double& xc, double& yc, double& zc,
259 const bool centre, double& rc);
266 virtual bool InTrapRadius(const double q0, const double x0, const double y0,
267 const double z0, double& xw, double& yw,
268 double& rw);
271 virtual bool CrossedPlane(const double x0, const double y0, const double z0,
272 const double x1, const double y1, const double z1,
273 double& xc, double& yc, double& zc);
274
276 void EnablePeriodicityX(const bool on = true) {
277 m_periodic[0] = on;
278 UpdatePeriodicity();
279 }
281 void EnablePeriodicityY(const bool on = true) {
282 m_periodic[1] = on;
283 UpdatePeriodicity();
284 }
286 void EnablePeriodicityZ(const bool on = true) {
287 m_periodic[2] = on;
288 UpdatePeriodicity();
289 }
291 void IsPeriodic(bool& perx, bool& pery, bool& perz) {
292 perx = m_periodic[0];
293 pery = m_periodic[1];
294 perz = m_periodic[2];
295 }
296
298 void EnableMirrorPeriodicityX(const bool on = true) {
299 m_mirrorPeriodic[0] = on;
300 UpdatePeriodicity();
301 }
303 void EnableMirrorPeriodicityY(const bool on = true) {
304 m_mirrorPeriodic[1] = on;
305 UpdatePeriodicity();
306 }
308 void EnableMirrorPeriodicityZ(const bool on = true) {
309 m_mirrorPeriodic[2] = on;
310 UpdatePeriodicity();
311 }
313 void IsMirrorPeriodic(bool& perx, bool& pery, bool& perz) {
314 perx = m_mirrorPeriodic[0];
315 pery = m_mirrorPeriodic[1];
316 perz = m_mirrorPeriodic[2];
317 }
318
320 void EnableAxialPeriodicityX(const bool on = true) {
321 m_axiallyPeriodic[0] = on;
322 UpdatePeriodicity();
323 }
325 void EnableAxialPeriodicityY(const bool on = true) {
326 m_axiallyPeriodic[1] = on;
327 UpdatePeriodicity();
328 }
330 void EnableAxialPeriodicityZ(const bool on = true) {
331 m_axiallyPeriodic[2] = on;
332 UpdatePeriodicity();
333 }
335 void IsAxiallyPeriodic(bool& perx, bool& pery, bool& perz) {
336 perx = m_axiallyPeriodic[0];
337 pery = m_axiallyPeriodic[1];
338 perz = m_axiallyPeriodic[2];
339 }
340
342 void EnableRotationSymmetryX(const bool on = true) {
343 m_rotationSymmetric[0] = on;
344 UpdatePeriodicity();
345 }
347 void EnableRotationSymmetryY(const bool on = true) {
348 m_rotationSymmetric[1] = on;
349 UpdatePeriodicity();
350 }
352 void EnableRotationSymmetryZ(const bool on = true) {
353 m_rotationSymmetric[2] = on;
354 UpdatePeriodicity();
355 }
357 void IsRotationSymmetric(bool& rotx, bool& roty, bool& rotz) {
358 rotx = m_rotationSymmetric[0];
359 roty = m_rotationSymmetric[1];
360 rotz = m_rotationSymmetric[2];
361 }
362
364 void EnableTriangleSymmetricXY(const bool on = true,
365 const bool oct = 2) {
366 m_triangleSymmetric[0] = on;
367 m_triangleSymmetricOct = oct;
368 m_mirrorPeriodic[0] = on;
369 m_mirrorPeriodic[1] = on;
370 }
372 void EnableTriangleSymmetricXZ(const bool on = true,
373 const bool oct = 2) {
374 m_triangleSymmetric[1] = on;
375 m_triangleSymmetricOct = oct;
376 m_mirrorPeriodic[0] = on;
377 m_mirrorPeriodic[2] = on;
378 }
380 void EnableTriangleSymmetricYZ(const bool on = true,
381 const bool oct = 2) {
382 m_triangleSymmetric[2] = on;
383 m_triangleSymmetricOct = oct;
384 m_mirrorPeriodic[1] = on;
385 m_mirrorPeriodic[2] = on;
386 }
387
389 void EnableDebugging(const bool on = true) { m_debug = on; }
391 void DisableDebugging() { m_debug = false; }
392
394 virtual bool HasMagneticField() const;
395
397 virtual bool HasTownsendMap() const { return false; }
399 virtual bool HasAttachmentMap() const { return false; }
401 virtual bool HasMobilityMap() const { return false; }
403 virtual bool HasVelocityMap() const { return false; }
404
406 virtual bool ElectronAttachment(const double /*x*/, const double /*y*/,
407 const double /*z*/, double& eta) {
408 eta = 0;
409 return false;
410 }
412 virtual bool HoleAttachment(const double /*x*/, const double /*y*/,
413 const double /*z*/, double& eta) {
414 eta = 0;
415 return false;
416 }
417
418 // Get the electron mobility coefficient.
419 virtual bool ElectronMobility(const double /*x*/, const double /*y*/,
420 const double /*z*/, double& mu) {
421 mu = 0;
422 return false;
423 }
425 virtual bool HoleMobility(const double /*x*/, const double /*y*/,
426 const double /*z*/, double& mu) {
427 mu = 0;
428 return false;
429 }
430
432 virtual bool ElectronTownsend(const double /*x*/, const double /*y*/,
433 const double /*z*/, double& alpha) {
434 alpha = 0;
435 return false;
436 }
438 virtual bool HoleTownsend(const double /*x*/, const double /*y*/,
439 const double /*z*/, double& alpha) {
440 alpha = 0;
441 return false;
442 }
444 virtual bool ElectronVelocity(const double /*x*/, const double /*y*/,
445 const double /*z*/, double& vx, double& vy,
446 double& vz) {
447 vx = vy = vz = 0;
448 return false;
449 }
451 virtual bool HoleVelocity(const double /*x*/, const double /*y*/,
452 const double /*z*/, double& vx, double& vy,
453 double& vz) {
454 vx = vy = vz = 0;
455 return false;
456 }
457
458 virtual double StepSizeHint() { return -1.; }
459
461 virtual double CreateGPUTransferObject(ComponentGPU*& comp_gpu);
462
463 protected:
465 std::string m_className = "Component";
466
468 Geometry* m_geometry = nullptr;
469
471 std::array<double, 3> m_b0 = {{0., 0., 0.}};
472#endif
474 bool m_ready = false;
475
476#ifdef __GPUCOMPILE__
478 bool m_periodic[3] = {false, false, false};
480 bool m_mirrorPeriodic[3] = {false, false, false};
482 bool m_axiallyPeriodic[3] = {false, false, false};
484 bool m_rotationSymmetric[3] = {false, false, false};
486 bool m_triangleSymmetric[3] = {false, false, false};
488 int m_triangleSymmetricOct = 0;
490 const int m_triangleOctRules[4] = {1, 4, 5, 8};
491 bool m_outsideCone = false;
492#else
494 bool m_debug = false;
495
497 std::array<bool, 3> m_periodic = {{false, false, false}};
499 std::array<bool, 3> m_mirrorPeriodic = {{false, false, false}};
501 std::array<bool, 3> m_axiallyPeriodic = {{false, false, false}};
503 std::array<bool, 3> m_rotationSymmetric = {{false, false, false}};
505 std::array<bool, 3> m_triangleSymmetric = {{false, false, false}};
507 int m_triangleSymmetricOct = 0;
509 const std::array<int, 4> m_triangleOctRules = {1, 4, 5, 8};
510 bool m_outsideCone = false;
511#endif
512
513#ifndef __GPUCOMPILE__
515 std::vector<double> m_wdtimes;
516
518 virtual void Reset() = 0;
520 virtual void UpdatePeriodicity() = 0;
521
522 private:
523 double IntegrateFluxParallelogram(const double x0, const double y0,
524 const double z0, const double dx1,
525 const double dy1, const double dz1,
526 const double dx2, const double dy2,
527 const double dz2, const unsigned int nU,
528 const unsigned int nV, const bool wfield,
529 const std::string& label);
530#else
531
532// include parts from derived class due to big performance hit from using
533// virtual methods
534// TODO GPU TN: It isn't clear to me that we actually need (at least) the Ansys
535// include - all the code is protected by an ifndef __GPUCOMPILE__, whereas it
536// will be defined when these includes are made
537#include "ComponentAnsys123.hh"
538#include "ComponentFieldMap.hh"
539
540 friend class ComponentAnsys123;
541 friend class ComponentComsol;
542 friend class ComponentElmer;
543 friend class ComponentFieldMap;
544 friend class Component;
545
546 // enum to mimic polymorphism
547 enum class ComponentType {
548 Component = 0,
553 };
554
555 ComponentType m_ComponentType{ComponentType::Component};
556
557#endif
558};
559} // namespace Garfield
560
561#endif
#define GARFIELD_CLASS_NAME(name)
#define __device__
Definition Vector.hh:12
Component for importing and interpolating three-dimensional ANSYS field maps.
Component for importing and interpolating Comsol field maps.
Component for importing field maps computed by Elmer.
Base class for components based on finite-element field maps.
Abstract base class for geometry classes.
Definition Geometry.hh:13
class GARFIELD_CLASS_NAME(Component)
Abstract base class for components.
Definition Component.hh:31