Garfield++ v1r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ComponentAnalyticField Class Reference

#include <ComponentAnalyticField.hh>

+ Inheritance diagram for Garfield::ComponentAnalyticField:

Public Member Functions

 ComponentAnalyticField ()
 
 ~ComponentAnalyticField ()
 
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)
 
void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status)
 
bool GetVoltageRange (double &pmin, double &pmax)
 
void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string label)
 
double WeightingPotential (const double x, const double y, const double z, const std::string label)
 
bool GetBoundingBox (double &x0, double &y0, double &z0, double &x1, double &y1, double &z1)
 
bool IsWireCrossed (double x0, double y0, double z0, double x1, double y1, double z1, double &xc, double &yc, double &zc)
 
bool IsInTrapRadius (double x0, double y0, double z0, double &xw, double &yx, double &rw)
 
void AddWire (const double x, const double y, const double diameter, const double voltage, const std::string label, const double length=100., const double tension=50., const double rho=19.3, const int ntrap=5)
 
void AddTube (const double radius, const double voltage, const int nEdges, const std::string label)
 
void AddPlaneX (const double x, const double voltage, const std::string label)
 
void AddPlaneY (const double y, const double voltage, const std::string label)
 
void AddStripOnPlaneX (const char direction, const double x, const double smin, const double smax, const std::string label, const double gap=-1.)
 
void AddStripOnPlaneY (const char direction, const double y, const double smin, const double smax, const std::string label, const double gap=-1.)
 
void SetPeriodicityX (const double s)
 
void SetPeriodicityY (const double s)
 
bool GetPeriodicityX (double &s)
 
bool GetPeriodicityY (double &s)
 
void AddCharge (const double x, const double y, const double z, const double q)
 
void ClearCharges ()
 
void PrintCharges ()
 
std::string GetCellType ()
 
void AddReadout (const std::string label)
 
void EnableChargeCheck ()
 
void DisableChargeCheck ()
 
int GetNumberOfWires ()
 
bool GetWire (const int i, double &x, double &y, double &diameter, double &voltage, std::string &label, double &length, double &charge, int &ntrap)
 
int GetNumberOfPlanesX ()
 
int GetNumberOfPlanesY ()
 
bool GetPlaneX (const int i, double &x, double &voltage, std::string &label)
 
bool GetPlaneY (const int i, double &y, double &voltage, std::string &label)
 
bool GetTube (double &r, double &voltage, int &nEdges, std::string &label)
 
- Public Member Functions inherited from Garfield::ComponentBase
 ComponentBase ()
 
virtual ~ComponentBase ()
 
virtual void SetGeometry (GeometryBase *geo)
 
virtual void Clear ()
 
virtual MediumGetMedium (const double &x, const double &y, const double &z)
 
virtual void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)=0
 
virtual void ElectricField (const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&m, int &status)=0
 
virtual bool GetVoltageRange (double &vmin, double &vmax)=0
 
virtual void WeightingField (const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string label)
 
virtual double WeightingPotential (const double x, const double y, const double z, const std::string label)
 
virtual void MagneticField (const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
 
void SetMagneticField (const double bx, const double by, const double bz)
 
virtual bool IsReady ()
 
virtual bool GetBoundingBox (double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)
 
virtual bool IsWireCrossed (const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, double &xc, double &yc, double &zc)
 
virtual bool IsInTrapRadius (double x0, double y0, double z0, double &xw, double &yw, double &rw)
 
void EnablePeriodicityX ()
 
void DisablePeriodicityX ()
 
void EnablePeriodicityY ()
 
void DisablePeriodicityY ()
 
void EnablePeriodicityZ ()
 
void DisablePeriodicityZ ()
 
void EnableMirrorPeriodicityX ()
 
void DisableMirrorPeriodicityX ()
 
void EnableMirrorPeriodicityY ()
 
void DisableMirrorPeriodicityY ()
 
void EnableMirrorPeriodicityZ ()
 
void DisableMirrorPeriodicityZ ()
 
void EnableAxialPeriodicityX ()
 
void DisableAxialPeriodicityX ()
 
void EnableAxialPeriodicityY ()
 
void DisableAxialPeriodicityY ()
 
void EnableAxialPeriodicityZ ()
 
void DisableAxialPeriodicityZ ()
 
void EnableRotationSymmetryX ()
 
void DisableRotationSymmetryX ()
 
void EnableRotationSymmetryY ()
 
void DisableRotationSymmetryY ()
 
void EnableRotationSymmetryZ ()
 
void DisableRotationSymmetryZ ()
 
void EnableDebugging ()
 
void DisableDebugging ()
 

Additional Inherited Members

virtual void Reset ()=0
 
virtual void UpdatePeriodicity ()=0
 
- Protected Attributes inherited from Garfield::ComponentBase
std::string m_className
 
GeometryBasetheGeometry
 
bool ready
 
bool xPeriodic
 
bool yPeriodic
 
bool zPeriodic
 
bool xMirrorPeriodic
 
bool yMirrorPeriodic
 
bool zMirrorPeriodic
 
bool xAxiallyPeriodic
 
bool yAxiallyPeriodic
 
bool zAxiallyPeriodic
 
bool xRotationSymmetry
 
bool yRotationSymmetry
 
bool zRotationSymmetry
 
double bx0
 
double by0
 
double bz0
 
bool debug
 

Detailed Description

Definition at line 12 of file ComponentAnalyticField.hh.

Constructor & Destructor Documentation

◆ ComponentAnalyticField()

Garfield::ComponentAnalyticField::ComponentAnalyticField ( )

Definition at line 11 of file ComponentAnalyticField.cc.

11 {
12
13 m_className = "ComponentAnalyticField";
14 chargeCheck = false;
15 CellInit();
16}

◆ ~ComponentAnalyticField()

Garfield::ComponentAnalyticField::~ComponentAnalyticField ( )
inline

Definition at line 18 of file ComponentAnalyticField.hh.

18{}

Member Function Documentation

◆ AddCharge()

void Garfield::ComponentAnalyticField::AddCharge ( const double  x,
const double  y,
const double  z,
const double  q 
)

Definition at line 682 of file ComponentAnalyticField.cc.

683 {
684
685 // Convert from fC to internal units (division by 4 pi epsilon0).
686 charge3d newCharge;
687 newCharge.x = x;
688 newCharge.y = y;
689 newCharge.z = z;
690 newCharge.e = q / FourPiEpsilon0;
691 ch3d.push_back(newCharge);
692 ++n3d;
693}

◆ AddPlaneX()

void Garfield::ComponentAnalyticField::AddPlaneX ( const double  x,
const double  voltage,
const std::string  label 
)

Definition at line 415 of file ComponentAnalyticField.cc.

416 {
417
418 if (ynplan[0] && ynplan[1]) {
419 std::cerr << m_className << "::AddPlaneX:\n";
420 std::cerr << " There are already two x planes defined.\n";
421 return;
422 }
423
424 if (ynplan[0]) {
425 ynplan[1] = true;
426 coplan[1] = x;
427 vtplan[1] = v;
428 planes[1].type = lab;
429 planes[1].ind = -1;
430 } else {
431 ynplan[0] = true;
432 coplan[0] = x;
433 vtplan[0] = v;
434 planes[0].type = lab;
435 planes[0].ind = -1;
436 }
437
438 // Force recalculation of the capacitance and signal matrices.
439 cellset = false;
440 sigset = false;
441}

◆ AddPlaneY()

void Garfield::ComponentAnalyticField::AddPlaneY ( const double  y,
const double  voltage,
const std::string  label 
)

Definition at line 443 of file ComponentAnalyticField.cc.

444 {
445
446 if (ynplan[2] && ynplan[3]) {
447 std::cerr << m_className << "::AddPlaneY:\n";
448 std::cerr << " There are already two y planes defined.\n";
449 return;
450 }
451
452 if (ynplan[2]) {
453 ynplan[3] = true;
454 coplan[3] = y;
455 vtplan[3] = v;
456 planes[3].type = lab;
457 planes[3].ind = -1;
458 } else {
459 ynplan[2] = true;
460 coplan[2] = y;
461 vtplan[2] = v;
462 planes[2].type = lab;
463 planes[2].ind = -1;
464 }
465
466 // Force recalculation of the capacitance and signal matrices.
467 cellset = false;
468 sigset = false;
469}

◆ AddReadout()

void Garfield::ComponentAnalyticField::AddReadout ( const std::string  label)

Definition at line 1938 of file ComponentAnalyticField.cc.

1938 {
1939
1940 // Check if this readout group already exists.
1941 for (int i = 0; i < nReadout; ++i) {
1942 if (readout[i] == label) {
1943 std::cout << m_className << "::AddReadout:\n";
1944 std::cout << " Readout group " << label << " already exists.\n";
1945 return;
1946 }
1947 }
1948
1949 readout.push_back(label);
1950 ++nReadout;
1951
1952 int nWiresFound = 0;
1953 for (int i = 0; i < nWires; ++i) {
1954 if (w[i].type == label) ++nWiresFound;
1955 }
1956
1957 int nPlanesFound = 0;
1958 int nStripsFound = 0;
1959 for (int i = 0; i < 5; ++i) {
1960 if (planes[i].type == label) ++nPlanesFound;
1961 for (int j = 0; j < planes[i].nStrips1; ++j) {
1962 if (planes[i].strips1[j].type == label) ++nStripsFound;
1963 }
1964 for (int j = 0; j < planes[i].nStrips2; ++j) {
1965 if (planes[i].strips2[j].type == label) ++nStripsFound;
1966 }
1967 }
1968
1969 if (nWiresFound == 0 && nPlanesFound == 0 && nStripsFound == 0) {
1970 std::cerr << m_className << "::AddReadout:\n";
1971 std::cerr << " At present there are no wires, planes or strips\n";
1972 std::cerr << " associated to readout group " << label << ".\n";
1973 } else {
1974 std::cout << m_className << "::AddReadout:\n";
1975 std::cout << " Readout group " << label << " comprises:\n";
1976 if (nWiresFound > 1) {
1977 std::cout << " " << nWiresFound << " wires\n";
1978 } else if (nWiresFound == 1) {
1979 std::cout << " 1 wire\n";
1980 }
1981 if (nPlanesFound > 1) {
1982 std::cout << " " << nPlanesFound << " planes\n";
1983 } else if (nPlanesFound == 1) {
1984 std::cout << " 1 plane\n";
1985 }
1986 if (nStripsFound > 1) {
1987 std::cout << " " << nStripsFound << " strips\n";
1988 } else if (nStripsFound == 1) {
1989 std::cout << " 1 strip\n";
1990 }
1991 }
1992
1993 sigset = false;
1994}

◆ AddStripOnPlaneX()

void Garfield::ComponentAnalyticField::AddStripOnPlaneX ( const char  direction,
const double  x,
const double  smin,
const double  smax,
const std::string  label,
const double  gap = -1. 
)

Definition at line 471 of file ComponentAnalyticField.cc.

475 {
476
477 if (!ynplan[0] && !ynplan[1]) {
478 std::cerr << m_className << "::AddStripOnPlaneX:\n";
479 std::cerr << " There are no planes at constant x defined.\n";
480 return;
481 }
482
483 if (direction != 'y' && direction != 'Y' && direction != 'z' &&
484 direction != 'Z') {
485 std::cerr << m_className << "::AddStripOnPlaneX:\n";
486 std::cerr << " Invalid direction (" << direction << ").\n";
487 std::cerr << " Only strips in y or z direction are possible.\n";
488 return;
489 }
490
491 if (fabs(smax - smin) < Small) {
492 std::cerr << m_className << "::AddStripOnPlaneX:\n";
493 std::cerr << " Strip width must be greater than zero.\n";
494 return;
495 }
496
497 strip newStrip;
498 newStrip.type = label;
499 newStrip.ind = -1;
500 newStrip.smin = std::min(smin, smax);
501 newStrip.smax = std::max(smin, smax);
502 if (gap > Small) {
503 newStrip.gap = gap;
504 } else {
505 newStrip.gap = -1.;
506 }
507
508 int iplane = 0;
509 if (ynplan[1]) {
510 const double d0 = fabs(coplan[0] - x);
511 const double d1 = fabs(coplan[1] - x);
512 if (d1 < d0) iplane = 1;
513 }
514
515 if (direction == 'y' || direction == 'Y') {
516 planes[iplane].nStrips1++;
517 planes[iplane].strips1.push_back(newStrip);
518 } else {
519 planes[iplane].nStrips2++;
520 planes[iplane].strips2.push_back(newStrip);
521 }
522}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:616

◆ AddStripOnPlaneY()

void Garfield::ComponentAnalyticField::AddStripOnPlaneY ( const char  direction,
const double  y,
const double  smin,
const double  smax,
const std::string  label,
const double  gap = -1. 
)

Definition at line 524 of file ComponentAnalyticField.cc.

528 {
529
530 if (!ynplan[2] && !ynplan[3]) {
531 std::cerr << m_className << "::AddStripOnPlaneY:\n";
532 std::cerr << " There are no planes at constant y defined.\n";
533 return;
534 }
535
536 if (direction != 'x' && direction != 'X' && direction != 'z' &&
537 direction != 'Z') {
538 std::cerr << m_className << "::AddStripOnPlaneY:\n";
539 std::cerr << " Invalid direction (" << direction << ").\n";
540 std::cerr << " Only strips in x or z direction are possible.\n";
541 return;
542 }
543
544 if (fabs(smax - smin) < Small) {
545 std::cerr << m_className << "::AddStripOnPlaneY:\n";
546 std::cerr << " Strip width must be greater than zero.\n";
547 return;
548 }
549
550 strip newStrip;
551 newStrip.type = label;
552 newStrip.ind = -1;
553 newStrip.smin = std::min(smin, smax);
554 newStrip.smax = std::max(smin, smax);
555 if (gap > Small) {
556 newStrip.gap = gap;
557 } else {
558 newStrip.gap = -1.;
559 }
560
561 int iplane = 2;
562 if (ynplan[3]) {
563 const double d2 = fabs(coplan[2] - y);
564 const double d3 = fabs(coplan[3] - y);
565 if (d3 < d2) iplane = 3;
566 }
567
568 if (direction == 'x' || direction == 'X') {
569 planes[iplane].nStrips1++;
570 planes[iplane].strips1.push_back(newStrip);
571 } else {
572 planes[iplane].nStrips2++;
573 planes[iplane].strips2.push_back(newStrip);
574 }
575}

◆ AddTube()

void Garfield::ComponentAnalyticField::AddTube ( const double  radius,
const double  voltage,
const int  nEdges,
const std::string  label 
)

Definition at line 374 of file ComponentAnalyticField.cc.

376 {
377
378 // Check if the provided parameters make sense.
379 if (radius <= 0.0) {
380 std::cerr << m_className << "::AddTube:\n";
381 std::cerr << " Unphysical tube dimension.\n";
382 return;
383 }
384
385 if (nEdges < 3 && nEdges != 0) {
386 std::cerr << m_className << "::AddTube:\n";
387 std::cerr << " Unphysical number of tube edges (" << nEdges << ")\n";
388 return;
389 }
390
391 // If there is already a tube defined, print a warning message.
392 if (tube) {
393 std::cout << m_className << "::AddTube:\n";
394 std::cout << " Warning: Existing tube settings will be overwritten.\n";
395 }
396
397 // Set the coordinate system.
398 tube = true;
399 polar = false;
400
401 // Set the tube parameters.
402 cotube = radius;
403 vttube = voltage;
404
405 ntube = nEdges;
406
407 planes[4].type = label;
408 planes[4].ind = -1;
409
410 // Force recalculation of the capacitance and signal matrices.
411 cellset = false;
412 sigset = false;
413}

Referenced by GarfieldPhysics::CreateGeometry().

◆ AddWire()

void Garfield::ComponentAnalyticField::AddWire ( const double  x,
const double  y,
const double  diameter,
const double  voltage,
const std::string  label,
const double  length = 100.,
const double  tension = 50.,
const double  rho = 19.3,
const int  ntrap = 5 
)

Definition at line 317 of file ComponentAnalyticField.cc.

322 {
323
324 // Check if the provided parameters make sense.
325 if (diameter <= 0.) {
326 std::cerr << m_className << "::AddWire:\n";
327 std::cerr << " Unphysical wire diameter.\n";
328 return;
329 }
330
331 if (tension <= 0.) {
332 std::cerr << m_className << "::AddWire:\n";
333 std::cerr << " Unphysical wire tension.\n";
334 return;
335 }
336
337 if (rho <= 0.0) {
338 std::cerr << m_className << "::AddWire:\n";
339 std::cerr << " Unphysical wire density.\n";
340 return;
341 }
342
343 if (length <= 0.0) {
344 std::cerr << m_className << "::AddWire:\n";
345 std::cerr << " Unphysical wire length.\n";
346 return;
347 }
348
349 if (ntrap <= 0) {
350 std::cerr << m_className << "::AddWire:\n";
351 std::cerr << " Number of trap radii must be > 0.\n";
352 return;
353 }
354 // Create a new wire
355 wire newWire;
356 newWire.x = x;
357 newWire.y = y;
358 newWire.d = diameter;
359 newWire.v = voltage;
360 newWire.u = length;
361 newWire.type = label;
362 newWire.e = 0.;
363 newWire.ind = -1;
364 newWire.nTrap = ntrap;
365 // Add the wire to the list
366 w.push_back(newWire);
367 ++nWires;
368
369 // Force recalculation of the capacitance and signal matrices.
370 cellset = false;
371 sigset = false;
372}

Referenced by GarfieldPhysics::CreateGeometry().

◆ ClearCharges()

void Garfield::ComponentAnalyticField::ClearCharges ( )

Definition at line 695 of file ComponentAnalyticField.cc.

695 {
696
697 n3d = 0;
698 ch3d.clear();
699 nTermBessel = 10;
700 nTermPoly = 100;
701}

◆ DisableChargeCheck()

void Garfield::ComponentAnalyticField::DisableChargeCheck ( )
inline

Definition at line 96 of file ComponentAnalyticField.hh.

96{ chargeCheck = false; }

◆ ElectricField() [1/2]

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

Implements Garfield::ComponentBase.

Definition at line 53 of file ComponentAnalyticField.cc.

56 {
57
58 // Initialize electric field and medium.
59 ex = ey = ez = v = 0.;
60 m = NULL;
61
62 // Make sure the charges have been calculated.
63 if (!cellset) {
64 if (!Prepare()) {
65 status = -11;
66 return;
67 }
68 }
69
70 // Request calculation of the potential.
71 const bool opt = true;
72
73 // Calculate the field.
74 status = Field(x, y, z, ex, ey, ez, v, opt);
75
76 // If the field is ok, get the medium.
77 if (status == 0) {
78 m = GetMedium(x, y, z);
79 if (m == NULL) {
80 status = -6;
81 } else if (!m->IsDriftable()) {
82 status = -5;
83 }
84 }
85}
virtual Medium * GetMedium(const double &x, const double &y, const double &z)

◆ ElectricField() [2/2]

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

Implements Garfield::ComponentBase.

Definition at line 18 of file ComponentAnalyticField.cc.

21 {
22
23 // Initialize electric field and medium.
24 ex = ey = ez = 0.;
25 m = NULL;
26
27 // Make sure the charges have been calculated.
28 if (!cellset) {
29 if (!Prepare()) {
30 status = -11;
31 return;
32 }
33 }
34
35 // Disable calculation of the potential.
36 const bool opt = false;
37 double v = 0.;
38
39 // Calculate the field.
40 status = Field(x, y, z, ex, ey, ez, v, opt);
41
42 // If the field is ok, get the medium.
43 if (status == 0) {
44 m = GetMedium(x, y, z);
45 if (m == NULL) {
46 status = -6;
47 } else if (!m->IsDriftable()) {
48 status = -5;
49 }
50 }
51}

◆ EnableChargeCheck()

void Garfield::ComponentAnalyticField::EnableChargeCheck ( )
inline

Definition at line 95 of file ComponentAnalyticField.hh.

95{ chargeCheck = true; }

◆ GetBoundingBox()

bool Garfield::ComponentAnalyticField::GetBoundingBox ( double &  x0,
double &  y0,
double &  z0,
double &  x1,
double &  y1,
double &  z1 
)
virtual

Reimplemented from Garfield::ComponentBase.

Definition at line 165 of file ComponentAnalyticField.cc.

167 {
168
169 // If a geometry is present, try to get the bounding box from there.
170 if (theGeometry != 0) {
171 if (theGeometry->GetBoundingBox(x0, y0, z0, x1, y1, z1)) return true;
172 }
173 // Otherwise, return the cell dimensions.
174 if (!cellset) return false;
175 x0 = xmin;
176 y0 = ymin;
177 z0 = zmin;
178 x1 = xmax;
179 y1 = ymax;
180 z1 = zmax;
181 return true;
182}
GeometryBase * theGeometry
virtual bool GetBoundingBox(double &xmin, double &ymin, double &zmin, double &xmax, double &ymax, double &zmax)=0

◆ GetCellType()

std::string Garfield::ComponentAnalyticField::GetCellType ( )
inline

Definition at line 74 of file ComponentAnalyticField.hh.

74{ return cellType; }

◆ GetNumberOfPlanesX()

int Garfield::ComponentAnalyticField::GetNumberOfPlanesX ( )

Definition at line 718 of file ComponentAnalyticField.cc.

718 {
719
720 if (ynplan[0] && ynplan[1]) {
721 return 2;
722 } else if (ynplan[0] || ynplan[1]) {
723 return 1;
724 }
725 return 0;
726}

◆ GetNumberOfPlanesY()

int Garfield::ComponentAnalyticField::GetNumberOfPlanesY ( )

Definition at line 728 of file ComponentAnalyticField.cc.

728 {
729
730 if (ynplan[2] && ynplan[3]) {
731 return 2;
732 } else if (ynplan[2] || ynplan[3]) {
733 return 1;
734 }
735 return 0;
736}

◆ GetNumberOfWires()

int Garfield::ComponentAnalyticField::GetNumberOfWires ( )
inline

Definition at line 98 of file ComponentAnalyticField.hh.

98{ return nWires; }

◆ GetPeriodicityX()

bool Garfield::ComponentAnalyticField::GetPeriodicityX ( double &  s)

Definition at line 603 of file ComponentAnalyticField.cc.

603 {
604
605 if (!xPeriodic) {
606 s = 0.;
607 return false;
608 }
609
610 s = sx;
611 return true;
612}

◆ GetPeriodicityY()

bool Garfield::ComponentAnalyticField::GetPeriodicityY ( double &  s)

Definition at line 614 of file ComponentAnalyticField.cc.

614 {
615
616 if (!yPeriodic) {
617 s = 0.;
618 return false;
619 }
620
621 s = sy;
622 return true;
623}

◆ GetPlaneX()

bool Garfield::ComponentAnalyticField::GetPlaneX ( const int  i,
double &  x,
double &  voltage,
std::string &  label 
)

Definition at line 760 of file ComponentAnalyticField.cc.

761 {
762
763 if (i < 0 || i >= 2 || (i == 1 && !ynplan[1])) {
764 std::cerr << m_className << "::GetPlaneX:\n";
765 std::cerr << " Plane index is out of range.\n";
766 return false;
767 }
768
769 x = coplan[i];
770 voltage = vtplan[i];
771 label = planes[i].type;
772 return true;
773}

◆ GetPlaneY()

bool Garfield::ComponentAnalyticField::GetPlaneY ( const int  i,
double &  y,
double &  voltage,
std::string &  label 
)

Definition at line 775 of file ComponentAnalyticField.cc.

776 {
777
778 if (i < 0 || i >= 2 || (i == 1 && !ynplan[3])) {
779 std::cerr << m_className << "::GetPlaneY:\n";
780 std::cerr << " Plane index is out of range.\n";
781 return false;
782 }
783
784 y = coplan[i + 2];
785 voltage = vtplan[i + 2];
786 label = planes[i + 2].type;
787 return true;
788}

◆ GetTube()

bool Garfield::ComponentAnalyticField::GetTube ( double &  r,
double &  voltage,
int &  nEdges,
std::string &  label 
)

Definition at line 790 of file ComponentAnalyticField.cc.

791 {
792
793 if (!tube) return false;
794 r = cotube;
795 voltage = vttube;
796 nEdges = ntube;
797 label = planes[4].type;
798 return true;
799}

◆ GetVoltageRange()

bool Garfield::ComponentAnalyticField::GetVoltageRange ( double &  pmin,
double &  pmax 
)
virtual

Implements Garfield::ComponentBase.

Definition at line 87 of file ComponentAnalyticField.cc.

87 {
88
89 // Make sure the cell is prepared.
90 if (!cellset) {
91 if (!Prepare()) {
92 std::cerr << m_className << "::GetVoltageRange:\n";
93 std::cerr << " Unable to return voltage range.\n";
94 std::cerr << " Cell could not be setup.\n";
95 return false;
96 }
97 }
98
99 pmin = vmin;
100 pmax = vmax;
101 return true;
102}

◆ GetWire()

bool Garfield::ComponentAnalyticField::GetWire ( const int  i,
double &  x,
double &  y,
double &  diameter,
double &  voltage,
std::string &  label,
double &  length,
double &  charge,
int &  ntrap 
)

Definition at line 738 of file ComponentAnalyticField.cc.

741 {
742
743 if (i < 0 || i >= nWires) {
744 std::cerr << m_className << "::GetWire:\n";
745 std::cerr << " Wire index is out of range.\n";
746 return false;
747 }
748
749 x = w[i].x;
750 y = w[i].y;
751 diameter = w[i].d;
752 voltage = w[i].v;
753 label = w[i].type;
754 length = w[i].u;
755 charge = w[i].e;
756 ntrap = w[i].nTrap;
757 return true;
758}

◆ IsInTrapRadius()

bool Garfield::ComponentAnalyticField::IsInTrapRadius ( double  x0,
double  y0,
double  z0,
double &  xw,
double &  yx,
double &  rw 
)
virtual

Reimplemented from Garfield::ComponentBase.

Definition at line 256 of file ComponentAnalyticField.cc.

258 {
259
260 // In case of periodicity, move the point into the basic cell.
261 double x0 = xin, y0 = yin;
262 int nX = 0, nY = 0, nPhi = 0;
263 if (perx) {
264 nX = int(round(xin / sx));
265 x0 -= sx * nX;
266 }
267 if (pery && tube) {
268 Cartesian2Polar(xin, yin, x0, y0);
269 nPhi = int(round((Pi * y0) / (sy * 180.)));
270 y0 -= 180. * sy * nPhi / Pi;
271 Polar2Cartesian(x0, y0, x0, y0);
272 } else if (pery) {
273 nY = int(round(yin / sy));
274 y0 -= sy * nY;
275 }
276
277 // Move the point to the correct side of the plane.
278 if (perx && ynplan[0] && x0 <= coplan[0]) x0 += sx;
279 if (perx && ynplan[1] && x0 >= coplan[1]) x0 -= sx;
280 if (pery && ynplan[2] && y0 <= coplan[2]) y0 += sy;
281 if (pery && ynplan[3] && y0 >= coplan[3]) y0 -= sy;
282
283 for (int i = 0; i < nWires; i++) {
284 const double xwc = w[i].x;
285 const double yxc = w[i].y;
286 const double r = sqrt(pow(xwc - x0, 2) + pow(yxc - y0, 2));
287 const double rTrap = 0.5 * w[i].d * w[i].nTrap;
288 if (r < rTrap) {
289 xw = w[i].x;
290 yw = w[i].y;
291 rw = w[i].d * 0.5;
292 if (perx && ynplan[0] && x0 <= coplan[0]) x0 -= sx;
293 if (perx && ynplan[1] && x0 >= coplan[1]) x0 += sx;
294 if (pery && ynplan[2] && y0 <= coplan[2]) y0 -= sy;
295 if (pery && ynplan[3] && y0 >= coplan[3]) y0 += sy;
296 if (pery && tube) {
297 double rhow, phiw;
298 Cartesian2Polar(xw, yw, rhow, phiw);
299 phiw += 180. * sy * nPhi / Pi;
300 Polar2Cartesian(rhow, phiw, xw, yw);
301 } else if (pery) {
302 y0 += sy * nY;
303 }
304 if (perx) xw += sx * nX;
305 if (debug) {
306 std::cout << m_className << "::IsInTrapRadius:\n";
307 std::cout << " (" << xin << ", " << yin << ", " << zin << ")"
308 << " within trap radius of wire " << i << ".\n";
309 }
310 return true;
311 }
312 }
313
314 return false;
315}
DoubleAc pow(const DoubleAc &f, double p)
Definition: DoubleAc.cpp:336
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313

◆ IsWireCrossed()

bool Garfield::ComponentAnalyticField::IsWireCrossed ( double  x0,
double  y0,
double  z0,
double  x1,
double  y1,
double  z1,
double &  xc,
double &  yc,
double &  zc 
)
virtual

Reimplemented from Garfield::ComponentBase.

Definition at line 184 of file ComponentAnalyticField.cc.

186 {
187
188 xc = x0;
189 yc = y0;
190 zc = z0;
191
192 if (nWires <= 0) return false;
193
194 const double dx = x1 - x0;
195 const double dy = y1 - y0;
196 const double d2 = dx * dx + dy * dy;
197 // Check that the step length is non-zero.
198 if (d2 < Small) return false;
199
200 // Check if a whole period has been crossed.
201 if ((perx && fabs(dx) >= sx) || (pery && fabs(dy) >= sy)) {
202 std::cerr << m_className << "::IsWireCrossed:\n";
203 std::cerr << " Particle crossed more than one period.\n";
204 return false;
205 }
206
207 // Both coordinates are assumed to be located inside
208 // the drift area and inside a drift medium.
209 // This should have been checked before this call.
210
211 const double xm = 0.5 * (x0 + x1);
212 const double ym = 0.5 * (y0 + y1);
213 double dMin2 = 0.;
214 for (int i = nWires; i--;) {
215 double xw = w[i].x, yw = w[i].y;
216 if (perx) {
217 xw += sx * int(round((xm - xw) / sx));
218 }
219 if (pery) {
220 yw += sy * int(round((ym - yw) / sy));
221 }
222 // Calculate the smallest distance between track and wire.
223 const double xIn0 = dx * (xw - x0) + dy * (yw - y0);
224 // Check if the minimum is located before (x0, y0).
225 if (xIn0 < 0.) continue;
226 const double xIn1 = -(dx * (xw - x1) + dy * (yw - y1));
227 // Check if the minimum is located behind (x1, y1).
228 if (xIn1 < 0.) continue;
229 // Minimum is located between (x0, y0) and (x1, y1).
230 const double dw02 = pow(xw - x0, 2) + pow(yw - y0, 2);
231 const double dw12 = pow(xw - x1, 2) + pow(yw - y1, 2);
232 if (xIn1 * xIn1 * dw02 > xIn0 * xIn0 * dw12) {
233 dMin2 = dw02 - xIn0 * xIn0 / d2;
234 } else {
235 dMin2 = dw12 - xIn1 * xIn1 / d2;
236 }
237 // Add in the times nTrap to account for the trap radius.
238 const double r2 = 0.25 * w[i].d * w[i].d;
239 if (dMin2 < r2) {
240 // Wire has been crossed.
241 // Find the point of intersection.
242 const double p = -xIn0 / d2;
243 const double q = (dw02 - r2) / d2;
244 const double t1 = -p + sqrt(p * p - q);
245 const double t2 = -p - sqrt(p * p - q);
246 const double t = std::min(t1, t2);
247 xc = x0 + t * dx;
248 yc = y0 + t * dy;
249 zc = z0 + t * (z1 - z0);
250 return true;
251 }
252 }
253 return false;
254}

◆ PrintCharges()

void Garfield::ComponentAnalyticField::PrintCharges ( )

Definition at line 703 of file ComponentAnalyticField.cc.

703 {
704
705 std::cout << m_className << "::PrintCharges:\n";
706 if (n3d <= 0) {
707 std::cout << " No charges present.\n";
708 return;
709 }
710 std::cout << " x [cm] y [cm] z [cm] charge [fC]\n";
711 for (int i = 0; i < n3d; ++i) {
712 std::cout << " " << std::setw(9) << ch3d[i].x << " " << std::setw(9)
713 << ch3d[i].y << " " << std::setw(9) << ch3d[i].z << " "
714 << std::setw(11) << ch3d[i].e * FourPiEpsilon0 << "\n";
715 }
716}

◆ SetPeriodicityX()

void Garfield::ComponentAnalyticField::SetPeriodicityX ( const double  s)

Definition at line 577 of file ComponentAnalyticField.cc.

577 {
578
579 if (s < Small) {
580 std::cerr << m_className << "::SetPeriodicityX:\n";
581 std::cerr << " Periodic length must be greater than zero.\n";
582 return;
583 }
584
585 xPeriodic = true;
586 sx = s;
587 UpdatePeriodicity();
588}

◆ SetPeriodicityY()

void Garfield::ComponentAnalyticField::SetPeriodicityY ( const double  s)

Definition at line 590 of file ComponentAnalyticField.cc.

590 {
591
592 if (s < Small) {
593 std::cerr << m_className << "::SetPeriodicityY:\n";
594 std::cerr << " Periodic length must be greater than zero.\n";
595 return;
596 }
597
598 yPeriodic = true;
599 sy = s;
600 UpdatePeriodicity();
601}

◆ WeightingField()

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

Reimplemented from Garfield::ComponentBase.

Definition at line 104 of file ComponentAnalyticField.cc.

107 {
108
109 wx = wy = wz = 0.;
110
111 if (nReadout <= 0) return;
112 if (!sigset) {
113 if (!PrepareSignals()) {
114 std::cerr << m_className << "::WeightingField::\n";
115 std::cerr << " Unable to calculate weighting fields.\n";
116 return;
117 }
118 }
119
120 if (label.empty()) return;
121 int index = -1;
122 for (int i = nReadout; i--;) {
123 if (readout[i] == label) {
124 index = i;
125 break;
126 }
127 }
128 if (index < 0) return;
129
130 double volt = 0.;
131 Wfield(x, y, z, wx, wy, wz, volt, index, false);
132}

◆ WeightingPotential()

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

Reimplemented from Garfield::ComponentBase.

Definition at line 134 of file ComponentAnalyticField.cc.

137 {
138
139 double volt = 0.;
140
141 if (nReadout <= 0) return volt;
142 if (!sigset) {
143 if (!PrepareSignals()) {
144 std::cerr << m_className << "::WeightingPotential::\n";
145 std::cerr << " Unable to calculate weighting fields.\n";
146 return volt;
147 }
148 }
149
150 if (label.empty()) return volt;
151 int index = -1;
152 for (int i = nReadout; i--;) {
153 if (readout[i] == label) {
154 index = i;
155 break;
156 }
157 }
158 if (index < 0) return volt;
159
160 double wx = 0., wy = 0., wz = 0.;
161 Wfield(x, y, z, wx, wy, wz, volt, index, true);
162 return volt;
163}

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