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

Plot transport coefficients as function of electric and magnetic field. More...

#include <ViewMedium.hh>

+ Inheritance diagram for Garfield::ViewMedium:

Public Member Functions

 ViewMedium ()
 Constructor.
 
 ~ViewMedium ()=default
 Destructor.
 
void SetMedium (Medium *m)
 Set the medium from which to retrieve the transport coefficients.
 
void EnableAutoRangeX (const bool on=true)
 Try to choose the x-axis range based on the field grid.
 
void SetRangeE (const double emin, const double emax, const bool logscale)
 Set the limits of the electric field.
 
void SetRangeB (const double bmin, const double bmax, const bool logscale)
 Set the limits of the magnetic field.
 
void SetRangeA (const double amin, const double amax, const bool logscale)
 Set the limits of the angle between electric and magnetic field.
 
void EnableAutoRangeY (const bool on=true)
 Choose the y-axis range based on the function's minima/maxima.
 
void SetRangeY (const double ymin, const double ymax, const bool logscale=false)
 Set the range of the function (velocity etc.) to be plotted.
 
void SetElectricField (const double efield)
 Set the electric field to use when plotting as function of B or angle.
 
void SetMagneticField (const double bfield)
 Set the magnetic field to use when plotting as function of E or angle.
 
void SetAngle (const double angle)
 Set the angle to use when plotting as function of E or B.
 
void EnableExport (const std::string &txtfile)
 
void DisableExport ()
 
void PlotVelocity (const std::string &carriers, const char xaxis)
 
void PlotDiffusion (const std::string &carriers, const char xaxis)
 Plot the transverse and longitudinal diffusion coefficients.
 
void PlotTownsend (const std::string &carriers, const char xaxis)
 Plot the Townsend coefficient.
 
void PlotAttachment (const std::string &carriers, const char xaxis)
 Plot the attachment coefficient.
 
void PlotAlphaEta (const std::string &carriers, const char xaxis)
 Plot Townsend and attachment coefficients.
 
void PlotElectronVelocity (const char xaxis='e', const bool same=false)
 
void PlotHoleVelocity (const char xaxis='e', const bool same=false)
 Plot the drift velocity components of holes in the medium.
 
void PlotIonVelocity (const char xaxis='e', const bool same=false)
 Plot the ion drift velocity in the gas.
 
void PlotElectronDiffusion (const char xaxis='e', const bool same=false)
 Plot the diffusion coefficients in the medium.
 
void PlotHoleDiffusion (const char xaxis='e', const bool same=false)
 Plot the diffusion coefficients of holes in the medium.
 
void PlotIonDiffusion (const char xaxis='e', const bool same=false)
 Plot the diffusion coefficients of ions in the gas.
 
void PlotElectronTownsend (const char xaxis='e', const bool same=false)
 Plot the Townsend coefficient for electrons.
 
void PlotHoleTownsend (const char xaxis='e', const bool same=false)
 Plot the Townsend coefficient for holes.
 
void PlotElectronAttachment (const char xaxis='e', const bool same=false)
 Plot the attachment coefficient for electrons.
 
void PlotHoleAttachment (const char xaxis='e', const bool same=false)
 Plot the attachment coefficient for holes.
 
void PlotElectronLorentzAngle (const char xaxis='e', const bool same=false)
 Plot the angle between drift velocity and field.
 
void SetColours (const std::vector< short > &cols)
 Set the (ROOT) colours to be used in the plots.
 
void SetLabels (const std::vector< std::string > &labels)
 Set user-defined plot labels.
 
- Public Member Functions inherited from Garfield::ViewBase
 ViewBase ()=delete
 Default constructor.
 
 ViewBase (const std::string &name)
 Constructor.
 
virtual ~ViewBase ()=default
 Destructor.
 
void SetCanvas (TPad *pad)
 Set the canvas to be painted on.
 
void SetCanvas ()
 Unset an external canvas.
 
TPad * GetCanvas ()
 Retrieve the canvas.
 
void SetArea (const double xmin, const double ymin, const double xmax, const double ymax)
 
virtual void SetArea (const double xmin, const double ymin, const double zmin, const double xmax, const double ymax, const double zmax)
 Set a bounding box (if applicable).
 
void SetArea ()
 
virtual void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
 
virtual void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0, const double hx, const double hy, const double hz)
 Set the projection plane specifying a hint for the in-plane x axis.
 
void Rotate (const double angle)
 Rotate the viewing plane (angle in radian).
 
void SetPlaneXY ()
 Set the viewing plane to x-y.
 
void SetPlaneXZ ()
 Set the viewing plane to x-z.
 
void SetPlaneYZ ()
 Set the viewing plane to y-z.
 
void SetPlaneZX ()
 Set the viewing plane to z-x.
 
void SetPlaneZY ()
 Set the viewing plane to z-y.
 
void EnableDebugging (const bool on=true)
 Switch on/off debugging output.
 

Additional Inherited Members

- Static Public Member Functions inherited from Garfield::ViewBase
static std::string FindUnusedFunctionName (const std::string &s)
 Find an unused function name.
 
static std::string FindUnusedHistogramName (const std::string &s)
 Find an unused histogram name.
 
static std::string FindUnusedCanvasName (const std::string &s)
 Find an unused canvas name.
 
- Protected Member Functions inherited from Garfield::ViewBase
void UpdateProjectionMatrix ()
 
template<typename T>
void ToPlane (const T x, const T y, const T z, T &xp, T &yp) const
 
template<typename T>
bool InBox (const std::array< T, 3 > &x) const
 
void Clip (const std::array< float, 3 > &x0, const std::array< float, 3 > &x1, std::array< float, 3 > &xc) const
 
void DrawLine (const std::vector< std::array< float, 3 > > &xl, const short col, const short lw)
 
std::string LabelX ()
 
std::string LabelY ()
 
std::string PlaneDescription ()
 
bool PlotLimits (Sensor *sensor, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimits (Component *cmp, double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimitsFromUserBox (double &xmin, double &ymin, double &xmax, double &ymax) const
 
bool PlotLimits (std::array< double, 3 > &bbmin, std::array< double, 3 > &bbmax, double &xmin, double &ymin, double &xmax, double &ymax) const
 
- Static Protected Member Functions inherited from Garfield::ViewBase
static bool RangeSet (TVirtualPad *)
 
static void SetRange (TVirtualPad *pad, const double x0, const double y0, const double x1, const double y1)
 
- Protected Attributes inherited from Garfield::ViewBase
std::string m_className = "ViewBase"
 
bool m_debug = false
 
bool m_userPlotLimits = false
 
double m_xMinPlot = -1.
 
double m_xMaxPlot = 1.
 
double m_yMinPlot = -1.
 
double m_yMaxPlot = 1.
 
bool m_userBox = false
 
double m_xMinBox = -1.
 
double m_xMaxBox = 1.
 
double m_yMinBox = -1.
 
double m_yMaxBox = 1.
 
double m_zMinBox = -1.
 
double m_zMaxBox = 1.
 
std::array< std::array< double, 3 >, 3 > m_proj
 
std::array< double, 4 > m_plane {{0, 0, 1, 0}}
 
std::array< std::array< double, 3 >, 3 > m_prmat
 

Detailed Description

Plot transport coefficients as function of electric and magnetic field.

Definition at line 17 of file ViewMedium.hh.

Constructor & Destructor Documentation

◆ ViewMedium()

Garfield::ViewMedium::ViewMedium ( )

Constructor.

Definition at line 43 of file ViewMedium.cc.

43: ViewBase("ViewMedium") {}
ViewBase()=delete
Default constructor.

◆ ~ViewMedium()

Garfield::ViewMedium::~ViewMedium ( )
default

Destructor.

Member Function Documentation

◆ DisableExport()

void Garfield::ViewMedium::DisableExport ( )
inline

Definition at line 49 of file ViewMedium.hh.

49{ m_outfile = ""; }

◆ EnableAutoRangeX()

void Garfield::ViewMedium::EnableAutoRangeX ( const bool on = true)
inline

Try to choose the x-axis range based on the field grid.

Definition at line 28 of file ViewMedium.hh.

28{ m_autoRangeX = on; }

◆ EnableAutoRangeY()

void Garfield::ViewMedium::EnableAutoRangeY ( const bool on = true)
inline

Choose the y-axis range based on the function's minima/maxima.

Definition at line 36 of file ViewMedium.hh.

36{ m_autoRangeY = on; }

◆ EnableExport()

void Garfield::ViewMedium::EnableExport ( const std::string & txtfile)
inline

Definition at line 48 of file ViewMedium.hh.

48{ m_outfile = txtfile; }

◆ PlotAlphaEta()

void Garfield::ViewMedium::PlotAlphaEta ( const std::string & carriers,
const char xaxis )

Plot Townsend and attachment coefficients.

Definition at line 644 of file ViewMedium.cc.

644 {
645
646 const auto ax = GetAxis(xaxis);
647 std::vector<Charge> carriers;
648 if (opt.find('e') != std::string::npos) carriers.push_back(Charge::Electron);
649 if (opt.find('h') != std::string::npos) carriers.push_back(Charge::Hole);
650 if (opt.find('i') != std::string::npos) carriers.push_back(Charge::Ion);
651 bool same = false;
652 for (const auto c : carriers) {
653 Plot(ax, c, Parameter::Townsend, same);
654 Plot(ax, c, Parameter::Attachment, true);
655 same = true;
656 }
657}

Referenced by Garfield::Medium::PlotAlphaEta().

◆ PlotAttachment()

void Garfield::ViewMedium::PlotAttachment ( const std::string & carriers,
const char xaxis )

Plot the attachment coefficient.

Definition at line 630 of file ViewMedium.cc.

630 {
631
632 const auto ax = GetAxis(xaxis);
633 std::vector<Charge> carriers;
634 if (opt.find('e') != std::string::npos) carriers.push_back(Charge::Electron);
635 if (opt.find('h') != std::string::npos) carriers.push_back(Charge::Hole);
636 if (opt.find('i') != std::string::npos) carriers.push_back(Charge::Ion);
637 const size_t nC = carriers.size();
638 for (size_t i = 0; i < nC; ++i) {
639 const bool same = i > 0 ? true : false;
640 Plot(ax, carriers[i], Parameter::Attachment, same);
641 }
642}

Referenced by Garfield::Medium::PlotAttachment().

◆ PlotDiffusion()

void Garfield::ViewMedium::PlotDiffusion ( const std::string & carriers,
const char xaxis )

Plot the transverse and longitudinal diffusion coefficients.

Definition at line 602 of file ViewMedium.cc.

602 {
603
604 const auto ax = GetAxis(xaxis);
605 std::vector<Charge> carriers;
606 if (opt.find('e') != std::string::npos) carriers.push_back(Charge::Electron);
607 if (opt.find('h') != std::string::npos) carriers.push_back(Charge::Hole);
608 if (opt.find('i') != std::string::npos) carriers.push_back(Charge::Ion);
609 const size_t nC = carriers.size();
610 for (size_t i = 0; i < nC; ++i) {
611 const bool same = i > 0 ? true : false;
612 PlotDiffusion(ax, carriers[i], same);
613 }
614}
void PlotDiffusion(const std::string &carriers, const char xaxis)
Plot the transverse and longitudinal diffusion coefficients.

Referenced by Garfield::Medium::PlotDiffusion(), PlotDiffusion(), PlotElectronDiffusion(), PlotHoleDiffusion(), and PlotIonDiffusion().

◆ PlotElectronAttachment()

void Garfield::ViewMedium::PlotElectronAttachment ( const char xaxis = 'e',
const bool same = false )
inline

Plot the attachment coefficient for electrons.

Definition at line 112 of file ViewMedium.hh.

112 {
113 Plot(GetAxis(xaxis), Charge::Electron, Parameter::Attachment, same);
114 }

◆ PlotElectronDiffusion()

void Garfield::ViewMedium::PlotElectronDiffusion ( const char xaxis = 'e',
const bool same = false )
inline

Plot the diffusion coefficients in the medium.

Definition at line 92 of file ViewMedium.hh.

92 {
93 PlotDiffusion(GetAxis(xaxis), Charge::Electron, same);
94 }

◆ PlotElectronLorentzAngle()

void Garfield::ViewMedium::PlotElectronLorentzAngle ( const char xaxis = 'e',
const bool same = false )
inline

Plot the angle between drift velocity and field.

Definition at line 120 of file ViewMedium.hh.

120 {
121 PlotLorentzAngle(GetAxis(xaxis), Charge::Electron, same);
122 }

◆ PlotElectronTownsend()

void Garfield::ViewMedium::PlotElectronTownsend ( const char xaxis = 'e',
const bool same = false )
inline

Plot the Townsend coefficient for electrons.

Definition at line 104 of file ViewMedium.hh.

104 {
105 Plot(GetAxis(xaxis), Charge::Electron, Parameter::Townsend, same);
106 }

◆ PlotElectronVelocity()

void Garfield::ViewMedium::PlotElectronVelocity ( const char xaxis = 'e',
const bool same = false )
inline

Plot the drift velocity components of electrons in the medium.

Parameters
xaxisabscissa.
  • 'e': electric field,
  • 'b': magnetic field,
  • 'a': angle between E and B.
sameflag to keep existing plots (true) or not.

Definition at line 80 of file ViewMedium.hh.

80 {
81 PlotVelocity(GetAxis(xaxis), Charge::Electron, same);
82 }
void PlotVelocity(const std::string &carriers, const char xaxis)

◆ PlotHoleAttachment()

void Garfield::ViewMedium::PlotHoleAttachment ( const char xaxis = 'e',
const bool same = false )
inline

Plot the attachment coefficient for holes.

Definition at line 116 of file ViewMedium.hh.

116 {
117 Plot(GetAxis(xaxis), Charge::Hole, Parameter::Attachment, same);
118 }

◆ PlotHoleDiffusion()

void Garfield::ViewMedium::PlotHoleDiffusion ( const char xaxis = 'e',
const bool same = false )
inline

Plot the diffusion coefficients of holes in the medium.

Definition at line 96 of file ViewMedium.hh.

96 {
97 PlotDiffusion(GetAxis(xaxis), Charge::Hole, same);
98 }

◆ PlotHoleTownsend()

void Garfield::ViewMedium::PlotHoleTownsend ( const char xaxis = 'e',
const bool same = false )
inline

Plot the Townsend coefficient for holes.

Definition at line 108 of file ViewMedium.hh.

108 {
109 Plot(GetAxis(xaxis), Charge::Hole, Parameter::Townsend, same);
110 }

◆ PlotHoleVelocity()

void Garfield::ViewMedium::PlotHoleVelocity ( const char xaxis = 'e',
const bool same = false )
inline

Plot the drift velocity components of holes in the medium.

Definition at line 84 of file ViewMedium.hh.

84 {
85 PlotVelocity(GetAxis(xaxis), Charge::Hole, same);
86 }

◆ PlotIonDiffusion()

void Garfield::ViewMedium::PlotIonDiffusion ( const char xaxis = 'e',
const bool same = false )
inline

Plot the diffusion coefficients of ions in the gas.

Definition at line 100 of file ViewMedium.hh.

100 {
101 PlotDiffusion(GetAxis(xaxis), Charge::Ion, same);
102 }

◆ PlotIonVelocity()

void Garfield::ViewMedium::PlotIonVelocity ( const char xaxis = 'e',
const bool same = false )
inline

Plot the ion drift velocity in the gas.

Definition at line 88 of file ViewMedium.hh.

88 {
89 PlotVelocity(GetAxis(xaxis), Charge::Ion, same);
90 }

◆ PlotTownsend()

void Garfield::ViewMedium::PlotTownsend ( const std::string & carriers,
const char xaxis )

Plot the Townsend coefficient.

Definition at line 616 of file ViewMedium.cc.

616 {
617
618 const auto ax = GetAxis(xaxis);
619 std::vector<Charge> carriers;
620 if (opt.find('e') != std::string::npos) carriers.push_back(Charge::Electron);
621 if (opt.find('h') != std::string::npos) carriers.push_back(Charge::Hole);
622 if (opt.find('i') != std::string::npos) carriers.push_back(Charge::Ion);
623 const size_t nC = carriers.size();
624 for (size_t i = 0; i < nC; ++i) {
625 const bool same = i > 0 ? true : false;
626 Plot(ax, carriers[i], Parameter::Townsend, same);
627 }
628}

Referenced by Garfield::Medium::PlotTownsend().

◆ PlotVelocity()

void Garfield::ViewMedium::PlotVelocity ( const std::string & carriers,
const char xaxis )

Plot the drift velocity components.

Parameters
optionstring indicating the carriers for which to plot the drift velocity
  • "e": electrons,
  • "h": holes,
  • "i": ions Options can be concatenated (e. g. "ei", "eh").
xaxisabscissa.
  • 'e': electric field,
  • 'b': magnetic field,
  • 'a': angle between E and B.

Definition at line 588 of file ViewMedium.cc.

588 {
589
590 const auto ax = GetAxis(xaxis);
591 std::vector<Charge> carriers;
592 if (opt.find('e') != std::string::npos) carriers.push_back(Charge::Electron);
593 if (opt.find('h') != std::string::npos) carriers.push_back(Charge::Hole);
594 if (opt.find('i') != std::string::npos) carriers.push_back(Charge::Ion);
595 const size_t nC = carriers.size();
596 for (size_t i = 0; i < nC; ++i) {
597 const bool same = i > 0 ? true : false;
598 PlotVelocity(ax, carriers[i], same);
599 }
600}

Referenced by PlotElectronVelocity(), PlotHoleVelocity(), PlotIonVelocity(), Garfield::Medium::PlotVelocity(), and PlotVelocity().

◆ SetAngle()

void Garfield::ViewMedium::SetAngle ( const double angle)
inline

Set the angle to use when plotting as function of E or B.

Definition at line 46 of file ViewMedium.hh.

46{ m_angle = angle; }

◆ SetColours()

void Garfield::ViewMedium::SetColours ( const std::vector< short > & cols)
inline

Set the (ROOT) colours to be used in the plots.

Definition at line 125 of file ViewMedium.hh.

125{ m_colours = cols; }

◆ SetElectricField()

void Garfield::ViewMedium::SetElectricField ( const double efield)
inline

Set the electric field to use when plotting as function of B or angle.

Definition at line 42 of file ViewMedium.hh.

42{ m_efield = efield; }

◆ SetLabels()

void Garfield::ViewMedium::SetLabels ( const std::vector< std::string > & labels)
inline

Set user-defined plot labels.

Definition at line 127 of file ViewMedium.hh.

127{ m_labels = labels; }

◆ SetMagneticField()

void Garfield::ViewMedium::SetMagneticField ( const double bfield)
inline

Set the magnetic field to use when plotting as function of E or angle.

Definition at line 44 of file ViewMedium.hh.

44{ m_bfield = bfield; }

◆ SetMedium()

void Garfield::ViewMedium::SetMedium ( Medium * m)

Set the medium from which to retrieve the transport coefficients.

Definition at line 45 of file ViewMedium.cc.

45 {
46 if (!m) {
47 std::cerr << m_className << "::SetMedium: Null pointer.\n";
48 return;
49 }
50
51 m_medium = m;
52}
std::string m_className
Definition ViewBase.hh:82

Referenced by Garfield::Medium::PlotAlphaEta(), Garfield::Medium::PlotAttachment(), Garfield::Medium::PlotDiffusion(), Garfield::Medium::PlotTownsend(), and Garfield::Medium::PlotVelocity().

◆ SetRangeA()

void Garfield::ViewMedium::SetRangeA ( const double amin,
const double amax,
const bool logscale )

Set the limits of the angle between electric and magnetic field.

Definition at line 77 of file ViewMedium.cc.

77 {
78 if (amin >= amax || amin < 0.) {
79 std::cerr << m_className << "::SetRangeA: Incorrect range.\n";
80 return;
81 }
82
83 m_aMin = amin;
84 m_aMax = amax;
85 m_logA = logscale;
86}

◆ SetRangeB()

void Garfield::ViewMedium::SetRangeB ( const double bmin,
const double bmax,
const bool logscale )

Set the limits of the magnetic field.

Definition at line 66 of file ViewMedium.cc.

66 {
67 if (bmin >= bmax || bmin < 0.) {
68 std::cerr << m_className << "::SetRangeB: Incorrect range.\n";
69 return;
70 }
71
72 m_bMin = bmin;
73 m_bMax = bmax;
74 m_logB = logscale;
75}

◆ SetRangeE()

void Garfield::ViewMedium::SetRangeE ( const double emin,
const double emax,
const bool logscale )

Set the limits of the electric field.

Definition at line 54 of file ViewMedium.cc.

55 {
56 if (emin >= emax || emin < 0.) {
57 std::cerr << m_className << "::SetRangeE: Incorrect range.\n";
58 return;
59 }
60
61 m_eMin = emin;
62 m_eMax = emax;
63 m_logE = logscale;
64}

◆ SetRangeY()

void Garfield::ViewMedium::SetRangeY ( const double ymin,
const double ymax,
const bool logscale = false )

Set the range of the function (velocity etc.) to be plotted.

Definition at line 88 of file ViewMedium.cc.

89 {
90 if (ymin >= ymax || ymin < 0.) {
91 std::cerr << m_className << "::SetRangeY: Incorrect range.\n";
92 return;
93 }
94
95 m_yMin = ymin;
96 m_yMax = ymax;
97 m_logY = logscale;
98}

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