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

#include <DriftLineRKF.hh>

Public Member Functions

 DriftLineRKF ()
 
 ~DriftLineRKF ()
 
void SetSensor (Sensor *s)
 
void EnablePlotting (ViewDrift *view)
 
void DisablePlotting ()
 
void SetIntegrationAccuracy (const double a)
 
void SetMaximumStepSize (const double ms)
 
void EnableStepSizeLimit ()
 
void DisableStepSizeLimit ()
 
void SetMaxSteps (const unsigned int m)
 
void EnableRejectKinks ()
 
void DisableRejectKinks ()
 
void SetElectronSignalScalingFactor (const double scale)
 
void SetHoleSignalScalingFactor (const double scale)
 
void SetIonSignalScalingFactor (const double scale)
 
bool DriftElectron (const double x0, const double y0, const double z0, const double t0)
 
bool DriftHole (const double x0, const double y0, const double z0, const double t0)
 
bool DriftIon (const double x0, const double y0, const double z0, const double t0)
 
void GetEndPoint (double &x, double &y, double &z, double &t, int &st) const
 
unsigned int GetNumberOfDriftLinePoints () const
 
void GetDriftLinePoint (const unsigned int i, double &x, double &y, double &z, double &t) const
 
double GetArrivalTimeSpread ()
 
double GetGain ()
 
double GetDriftTime () const
 
void EnableDebugging ()
 
void DisableDebugging ()
 
void EnableVerbose ()
 
void DisableVerbose ()
 

Detailed Description

Calculation of drift lines based on macroscopic transport coefficients using Runge-Kutta-Fehlberg integration.

Definition at line 17 of file DriftLineRKF.hh.

Constructor & Destructor Documentation

◆ DriftLineRKF()

Garfield::DriftLineRKF::DriftLineRKF ( )

Definition at line 11 of file DriftLineRKF.cc.

12 : m_sensor(NULL),
13 m_medium(NULL),
14 m_maxStepSize(1.e8),
15 m_accuracy(1.e-8),
16 m_maxSteps(1000),
17 m_maxStepsToWire(1000),
18 m_rejectKinks(true),
19 m_useStepSizeLimit(false),
20 m_usePlotting(false),
21 m_view(NULL),
22 m_status(0),
23 m_nPoints(0),
24 m_scaleElectronSignal(1.),
25 m_scaleHoleSignal(1.),
26 m_scaleIonSignal(1.),
27 m_debug(false),
28 m_verbose(false) {
29
30 const unsigned int nMaxPoints = m_maxSteps + m_maxStepsToWire + 10;
31 m_path.resize(nMaxPoints);
32 m_className = "DriftLineRKF";
33}

◆ ~DriftLineRKF()

Garfield::DriftLineRKF::~DriftLineRKF ( )
inline

Definition at line 21 of file DriftLineRKF.hh.

21{}

Member Function Documentation

◆ DisableDebugging()

void Garfield::DriftLineRKF::DisableDebugging ( )
inline

Definition at line 60 of file DriftLineRKF.hh.

60{ m_debug = false; }

◆ DisablePlotting()

void Garfield::DriftLineRKF::DisablePlotting ( )

Definition at line 78 of file DriftLineRKF.cc.

78 {
79
80 m_view = NULL;
81 m_usePlotting = false;
82}

◆ DisableRejectKinks()

void Garfield::DriftLineRKF::DisableRejectKinks ( )
inline

Definition at line 36 of file DriftLineRKF.hh.

36{ m_rejectKinks = false; }

◆ DisableStepSizeLimit()

void Garfield::DriftLineRKF::DisableStepSizeLimit ( )
inline

Definition at line 32 of file DriftLineRKF.hh.

32{ m_useStepSizeLimit = false; }

◆ DisableVerbose()

void Garfield::DriftLineRKF::DisableVerbose ( )
inline

Definition at line 63 of file DriftLineRKF.hh.

63{ m_verbose = false; }

◆ DriftElectron()

bool Garfield::DriftLineRKF::DriftElectron ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Definition at line 84 of file DriftLineRKF.cc.

85 {
86
87 m_particleType = ParticleTypeElectron;
88 if (!DriftLine(x0, y0, z0, t0)) return false;
89 GetGain();
90 ComputeSignal(-1., m_scaleElectronSignal);
91 return true;
92}

◆ DriftHole()

bool Garfield::DriftLineRKF::DriftHole ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Definition at line 94 of file DriftLineRKF.cc.

95 {
96
97 m_particleType = ParticleTypeHole;
98 if (!DriftLine(x0, y0, z0, t0)) return false;
99 ComputeSignal(1., m_scaleHoleSignal);
100 return true;
101}

◆ DriftIon()

bool Garfield::DriftLineRKF::DriftIon ( const double  x0,
const double  y0,
const double  z0,
const double  t0 
)

Definition at line 103 of file DriftLineRKF.cc.

104 {
105
106 m_particleType = ParticleTypeIon;
107 if (!DriftLine(x0, y0, z0, t0)) return false;
108 ComputeSignal(1., m_scaleIonSignal);
109 return true;
110}

◆ EnableDebugging()

void Garfield::DriftLineRKF::EnableDebugging ( )
inline

Definition at line 59 of file DriftLineRKF.hh.

59{ m_debug = true; }

◆ EnablePlotting()

void Garfield::DriftLineRKF::EnablePlotting ( ViewDrift view)

Definition at line 68 of file DriftLineRKF.cc.

68 {
69
70 if (!view) {
71 std::cerr << m_className << "::EnablePlotting:\n Null pointer.\n";
72 return;
73 }
74 m_usePlotting = true;
75 m_view = view;
76}

◆ EnableRejectKinks()

void Garfield::DriftLineRKF::EnableRejectKinks ( )
inline

Definition at line 35 of file DriftLineRKF.hh.

35{ m_rejectKinks = true; }

◆ EnableStepSizeLimit()

void Garfield::DriftLineRKF::EnableStepSizeLimit ( )
inline

Definition at line 31 of file DriftLineRKF.hh.

31{ m_useStepSizeLimit = true; }

◆ EnableVerbose()

void Garfield::DriftLineRKF::EnableVerbose ( )
inline

Definition at line 62 of file DriftLineRKF.hh.

62{ m_verbose = true; }

◆ GetArrivalTimeSpread()

double Garfield::DriftLineRKF::GetArrivalTimeSpread ( )

Definition at line 462 of file DriftLineRKF.cc.

462 {
463
464 if (m_nPoints < 2) return 0.;
465 double sum = 0.;
466 for (unsigned int i = 0; i < m_nPoints - 1; ++i) {
467 const unsigned int j = i + 1;
468 sum += IntegrateDiffusion(m_path[i].x, m_path[i].y, m_path[i].z,
469 m_path[j].x, m_path[j].y, m_path[j].z);
470 }
471 return sqrt(sum);
472}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

◆ GetDriftLinePoint()

void Garfield::DriftLineRKF::GetDriftLinePoint ( const unsigned int  i,
double &  x,
double &  y,
double &  z,
double &  t 
) const

Definition at line 870 of file DriftLineRKF.cc.

871 {
872
873 if (i >= m_nPoints) {
874 std::cerr << m_className << "::GetDriftLinePoint:\n";
875 std::cerr << " Index is outside the range.\n";
876 return;
877 }
878
879 // Return midpoint of drift line stage
880 x = m_path[i].x;
881 y = m_path[i].y;
882 z = m_path[i].z;
883 t = m_path[i].t;
884}

◆ GetDriftTime()

double Garfield::DriftLineRKF::GetDriftTime ( ) const
inline

Definition at line 55 of file DriftLineRKF.hh.

55 {
56 return m_nPoints > 0 ? m_path[m_nPoints - 1].t : 0.;
57 }

◆ GetEndPoint()

void Garfield::DriftLineRKF::GetEndPoint ( double &  x,
double &  y,
double &  z,
double &  t,
int &  st 
) const

Definition at line 855 of file DriftLineRKF.cc.

856 {
857
858 if (m_nPoints == 0) {
859 x = y = z = t = 0.;
860 stat = m_status;
861 return;
862 }
863 x = m_path[m_nPoints - 1].x;
864 y = m_path[m_nPoints - 1].y;
865 z = m_path[m_nPoints - 1].z;
866 t = m_path[m_nPoints - 1].t;
867 stat = m_status;
868}

◆ GetGain()

double Garfield::DriftLineRKF::GetGain ( )

Definition at line 474 of file DriftLineRKF.cc.

474 {
475
476 if (m_nPoints < 2) return 0.;
477 if (m_status == StatusCalculationAbandoned) return 0.;
478 // First get a rough estimate of the result.
479 double crude = 0.;
480 double alphaPrev = 0.;
481 for (unsigned int i = 0; i < m_nPoints; ++i) {
482 // Get the Townsend coefficient at this step.
483 const double x = m_path[i].x;
484 const double y = m_path[i].y;
485 const double z = m_path[i].z;
486 double ex, ey, ez;
487 double bx, by, bz;
488 int status;
489 m_sensor->MagneticField(x, y, z, bx, by, bz, status);
490 m_sensor->ElectricField(x, y, z, ex, ey, ez, m_medium, status);
491 if (status != 0) {
492 std::cerr << m_className << "::GetGain:\n"
493 << " Invalid drift line point " << i << ".\n";
494 return 0.;
495 }
496 double alpha = 0.;
497 if (!GetTownsend(ex, ey, ez, bx, by, bz, alpha)) {
498 std::cerr << m_className << "::GetGain:\n"
499 << " Unable to retrieve Townsend coefficient.\n";
500 return 0.;
501 }
502 if (i == 0) {
503 alphaPrev = alpha;
504 continue;
505 }
506 const double dx = x - m_path[i - 1].x;
507 const double dy = y - m_path[i - 1].y;
508 const double dz = z - m_path[i - 1].z;
509 const double d = sqrt(dx * dx + dy * dy + dz * dz);
510 crude += 0.5 * d * (alpha + alphaPrev);
511 alphaPrev = alpha;
512 }
513 // Stop if the rough estimate is negligibly small.
514 if (crude < Small) return 1.;
515
516 // Calculate the integration tolerance based on the rough estimate.
517 const double tol = 1.e-4 * crude;
518 double sum = 0.;
519 m_path[0].alphaint = 0.;
520 for (unsigned int i = 0; i < m_nPoints - 1; ++i) {
521 const unsigned int j = i + 1;
522 sum += IntegrateTownsend(m_path[i].x, m_path[i].y, m_path[i].z,
523 m_path[j].x, m_path[j].y, m_path[j].z, tol);
524 m_path[j].alphaint = sum;
525 }
526 return exp(sum);
527}
void MagneticField(const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
Get the magnetic field at (x, y, z).
Definition: Sensor.cc:101
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&medium, int &status)
Get the drift field and potential at (x, y, z).
Definition: Sensor.cc:48
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:377

Referenced by DriftElectron().

◆ GetNumberOfDriftLinePoints()

unsigned int Garfield::DriftLineRKF::GetNumberOfDriftLinePoints ( ) const
inline

Definition at line 50 of file DriftLineRKF.hh.

50{ return m_nPoints; }

◆ SetElectronSignalScalingFactor()

void Garfield::DriftLineRKF::SetElectronSignalScalingFactor ( const double  scale)
inline

Definition at line 38 of file DriftLineRKF.hh.

38{ m_scaleElectronSignal = scale; }

◆ SetHoleSignalScalingFactor()

void Garfield::DriftLineRKF::SetHoleSignalScalingFactor ( const double  scale)
inline

Definition at line 39 of file DriftLineRKF.hh.

39{ m_scaleHoleSignal = scale; }

◆ SetIntegrationAccuracy()

void Garfield::DriftLineRKF::SetIntegrationAccuracy ( const double  a)

Definition at line 44 of file DriftLineRKF.cc.

44 {
45
46 if (a > 0.) {
47 m_accuracy = a;
48 } else {
49 std::cerr << m_className << "::SetIntegrationAccuracy:\n"
50 << " Accuracy must be greater than zero.\n";
51 }
52}

◆ SetIonSignalScalingFactor()

void Garfield::DriftLineRKF::SetIonSignalScalingFactor ( const double  scale)
inline

Definition at line 40 of file DriftLineRKF.hh.

40{ m_scaleIonSignal = scale; }

◆ SetMaximumStepSize()

void Garfield::DriftLineRKF::SetMaximumStepSize ( const double  ms)

Definition at line 54 of file DriftLineRKF.cc.

54 {
55
56 if (ms > 0.) {
57 m_maxStepSize = ms;
58 if (!m_useStepSizeLimit) {
59 std::cout << m_className << "::SetMaximumStepSize:\n"
60 << " Don't forget to call EnableStepSizeLimit.\n";
61 }
62 } else {
63 std::cerr << m_className << "::SetMaximumStepSize:\n"
64 << " Step size must be greater than zero.\n";
65 }
66}

◆ SetMaxSteps()

void Garfield::DriftLineRKF::SetMaxSteps ( const unsigned int  m)
inline

Definition at line 34 of file DriftLineRKF.hh.

34{ m_maxSteps = m; }

◆ SetSensor()

void Garfield::DriftLineRKF::SetSensor ( Sensor s)

Definition at line 35 of file DriftLineRKF.cc.

35 {
36
37 if (!s) {
38 std::cerr << m_className << "::SetSensor:\n Null pointer.\n";
39 return;
40 }
41 m_sensor = s;
42}

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