Garfield++ v1r0
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 SetIntegrationAccuracy (const double a)
 
void SetMaximumStepSize (const double ms)
 
void EnablePlotting (ViewDrift *view)
 
void DisablePlotting ()
 
void SetMaxSteps (const unsigned int &m)
 
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

Definition at line 14 of file DriftLineRKF.hh.

Constructor & Destructor Documentation

◆ DriftLineRKF()

Garfield::DriftLineRKF::DriftLineRKF ( )

Definition at line 10 of file DriftLineRKF.cc.

11 : m_sensor(0),
12 m_medium(0),
13 m_maxStepSize(1.e8),
14 m_accuracy(1.e-8),
15 m_maxSteps(1000),
16 m_usePlotting(false),
17 m_view(0),
18 m_scaleElectronSignal(1.),
19 m_scaleHoleSignal(1.),
20 m_scaleIonSignal(1.),
21 m_debug(false),
22 m_verbose(false) {
23
24 m_className = "DriftLineRKF";
25 m_path.clear();
26}

◆ ~DriftLineRKF()

Garfield::DriftLineRKF::~DriftLineRKF ( )
inline

Definition at line 18 of file DriftLineRKF.hh.

18{}

Member Function Documentation

◆ DisableDebugging()

void Garfield::DriftLineRKF::DisableDebugging ( )
inline

Definition at line 49 of file DriftLineRKF.hh.

49{ m_debug = false; }

◆ DisablePlotting()

void Garfield::DriftLineRKF::DisablePlotting ( )

Definition at line 69 of file DriftLineRKF.cc.

69 {
70
71 m_view = 0;
72 m_usePlotting = false;
73}

◆ DisableVerbose()

void Garfield::DriftLineRKF::DisableVerbose ( )
inline

Definition at line 52 of file DriftLineRKF.hh.

52{ m_verbose = false; }

◆ DriftElectron()

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

Definition at line 75 of file DriftLineRKF.cc.

76 {
77
78 m_particleType = ParticleTypeElectron;
79 if (!DriftLine(x0, y0, z0, t0)) return false;
80 GetGain();
81 ComputeSignal();
82 return true;
83}

◆ DriftHole()

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

Definition at line 85 of file DriftLineRKF.cc.

86 {
87
88 m_particleType = ParticleTypeHole;
89 if (!DriftLine(x0, y0, z0, t0)) return false;
90 ComputeSignal();
91 return true;
92}

◆ DriftIon()

bool Garfield::DriftLineRKF::DriftIon ( 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 = ParticleTypeIon;
98 if (!DriftLine(x0, y0, z0, t0)) return false;
99 ComputeSignal();
100 return true;
101}

◆ EnableDebugging()

void Garfield::DriftLineRKF::EnableDebugging ( )
inline

Definition at line 48 of file DriftLineRKF.hh.

48{ m_debug = true; }

◆ EnablePlotting()

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

Definition at line 58 of file DriftLineRKF.cc.

58 {
59
60 if (!view) {
61 std::cerr << m_className << "::EnablePlotting:\n";
62 std::cerr << " Viewer pointer is null.\n";
63 return;
64 }
65 m_usePlotting = true;
66 m_view = view;
67}

◆ EnableVerbose()

void Garfield::DriftLineRKF::EnableVerbose ( )
inline

Definition at line 51 of file DriftLineRKF.hh.

51{ m_verbose = true; }

◆ GetArrivalTimeSpread()

double Garfield::DriftLineRKF::GetArrivalTimeSpread ( )

Definition at line 401 of file DriftLineRKF.cc.

401 {
402
403 const unsigned int nSteps = m_path.size();
404 double sum = 0.;
405 for (unsigned int i = 0; i < nSteps; ++i) {
406 sum += IntegrateDiffusion(m_path[i].xi, m_path[i].yi, m_path[i].zi,
407 m_path[i].xf, m_path[i].yf, m_path[i].zf);
408 }
409 return sqrt(sum);
410}
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:313

◆ GetDriftLinePoint()

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

Definition at line 816 of file DriftLineRKF.cc.

817 {
818
819 if (i >= m_path.size()) {
820 std::cerr << m_className << "::GetDriftLinePoint:\n";
821 std::cerr << " Index is outside the range.\n";
822 return;
823 }
824
825 // Return midpoint of drift line stage
826 x = 0.5*(m_path.at(i).xi+m_path.at(i).xf);
827 y = 0.5*(m_path.at(i).yi+m_path.at(i).yf);
828 z = 0.5*(m_path.at(i).zi+m_path.at(i).zf);
829 t = 0.5*(m_path.at(i).ti+m_path.at(i).tf);
830}

◆ GetDriftTime()

double Garfield::DriftLineRKF::GetDriftTime ( ) const

Definition at line 466 of file DriftLineRKF.cc.

466 {
467
468 if (m_path.empty()) return 0.;
469 return m_path.back().tf;
470}

◆ GetEndPoint()

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

Definition at line 806 of file DriftLineRKF.cc.

807 {
808
809 x = m_path.back().xi;
810 y = m_path.back().yi;
811 z = m_path.back().zi;
812 t = m_path.back().ti;
813 stat = m_path.back().status;
814}

◆ GetGain()

double Garfield::DriftLineRKF::GetGain ( )

Definition at line 412 of file DriftLineRKF.cc.

412 {
413
414 if (m_path.empty()) return 0.;
415 const unsigned int nSteps = m_path.size();
416 // First get a rough estimate of the result.
417 double crude = 0.;
418 double alpha0 = 0.;
419 double alpha1 = 0.;
420 for (unsigned int i = 0; i < nSteps; ++i) {
421 // Get the Townsend coefficient at this step.
422 const double x = m_path[i].xi;
423 const double y = m_path[i].yi;
424 const double z = m_path[i].zi;
425 double ex, ey, ez;
426 double bx, by, bz;
427 int status;
428 m_sensor->MagneticField(x, y, z, bx, by, bz, status);
429 m_sensor->ElectricField(x, y, z, ex, ey, ez, m_medium, status);
430 if (status != 0) {
431 std::cerr << m_className << "::GetGain:\n";
432 std::cerr << " Initial drift line point.\n";
433 return 0.;
434 }
435 if (0 == i) {
436 if (!GetTownsend(ex, ey, ez, bx, by, bz, alpha0)) {
437 std::cerr << m_className << "::GetGain:\n";
438 std::cerr << " Unable to retrieve Townsend coefficient.\n";
439 return 0.;
440 }
441 } else {
442 if (!GetTownsend(ex, ey, ez, bx, by, bz, alpha1)) {
443 std::cerr << m_className << "::GetGain:\n";
444 std::cerr << " Unable to retrieve Townsend coefficient.\n";
445 return 0.;
446 }
447 const double dx = x - m_path[i - 1].xi;
448 const double dy = y - m_path[i - 1].yi;
449 const double dz = z - m_path[i - 1].zi;
450 const double d = sqrt(dx * dx + dy * dy + dz * dz);
451 crude += 0.5 * d * (alpha1 + alpha0);
452 alpha0 = alpha1;
453 }
454 }
455 // Calculate the integration tolerance based on the rough estimate.
456 const double tol = 1.e-4 * crude;
457 double sum = 0.;
458 for (unsigned int i = 0; i < nSteps; ++i) {
459 sum += IntegrateTownsend(m_path[i].xi, m_path[i].yi, m_path[i].zi,
460 m_path[i].xf, m_path[i].yf, m_path[i].zf, tol);
461 m_path[i].alphaint = sum;
462 }
463 return exp(sum);
464}
DoubleAc exp(const DoubleAc &f)
Definition: DoubleAc.cpp:376
void MagneticField(const double x, const double y, const double z, double &bx, double &by, double &bz, int &status)
Definition: Sensor.cc:95
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&medium, int &status)
Definition: Sensor.cc:44

Referenced by DriftElectron().

◆ GetNumberOfDriftLinePoints()

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

Definition at line 41 of file DriftLineRKF.hh.

41{ return m_path.size(); }

◆ SetElectronSignalScalingFactor()

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

Definition at line 29 of file DriftLineRKF.hh.

29{ m_scaleElectronSignal = scale; }

◆ SetHoleSignalScalingFactor()

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

Definition at line 30 of file DriftLineRKF.hh.

30{ m_scaleHoleSignal = scale; }

◆ SetIntegrationAccuracy()

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

Definition at line 38 of file DriftLineRKF.cc.

38 {
39
40 if (a > 0.) {
41 m_accuracy = a;
42 } else {
43 std::cerr << m_className << "::SetIntegrationAccuracy:\n";
44 std::cerr << " Accuracy must be greater than zero.\n";
45 }
46}

◆ SetIonSignalScalingFactor()

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

Definition at line 31 of file DriftLineRKF.hh.

31{ m_scaleIonSignal = scale; }

◆ SetMaximumStepSize()

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

Definition at line 48 of file DriftLineRKF.cc.

48 {
49
50 if (ms > 0.) {
51 m_maxStepSize = ms;
52 } else {
53 std::cerr << m_className << "::SetMaximumStepSize:\n";
54 std::cerr << " Step size must be greater than zero.\n";
55 }
56}

◆ SetMaxSteps()

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

Definition at line 27 of file DriftLineRKF.hh.

27{ m_maxSteps = m; }

◆ SetSensor()

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

Definition at line 28 of file DriftLineRKF.cc.

28 {
29
30 if (!s) {
31 std::cerr << m_className << "::SetSensor:\n";
32 std::cerr << " Sensor pointer is null.\n";
33 return;
34 }
35 m_sensor = s;
36}

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