Garfield++ 4.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewMedium.hh
Go to the documentation of this file.
1#ifndef G_VIEW_MEDIUM
2#define G_VIEW_MEDIUM
3
4#include <string>
5#include <vector>
6#include <array>
7
8#include "ViewBase.hh"
10
11namespace Garfield {
12
13class Medium;
14
15/// Plot transport coefficients as function of electric and magnetic field.
16
17class ViewMedium : public ViewBase {
18 public:
19 /// Constructor.
20 ViewMedium();
21 /// Destructor.
22 ~ViewMedium() = default;
23
24 /// Set the medium from which to retrieve the transport coefficients.
25 void SetMedium(Medium* m);
26
27 /// Try to choose the x-axis range based on the field grid.
28 void EnableAutoRangeX(const bool on = true) { m_autoRangeX = on; }
29 /// Set the limits of the electric field.
30 void SetRangeE(const double emin, const double emax, const bool logscale);
31 /// Set the limits of the magnetic field.
32 void SetRangeB(const double bmin, const double bmax, const bool logscale);
33 /// Set the limits of the angle between electric and magnetic field.
34 void SetRangeA(const double amin, const double amax, const bool logscale);
35 /// Choose the y-axis range based on the function's minima/maxima.
36 void EnableAutoRangeY(const bool on = true) { m_autoRangeY = on; }
37 /// Set the range of the function (velocity etc.) to be plotted.
38 void SetRangeY(const double ymin, const double ymax,
39 const bool logscale = false);
40
41 /// Set the electric field to use when plotting as function of B or angle.
42 void SetElectricField(const double efield) { m_efield = efield; }
43 /// Set the magnetic field to use when plotting as function of E or angle.
44 void SetMagneticField(const double bfield) { m_bfield = bfield; }
45 /// Set the angle to use when plotting as function of E or B.
46 void SetAngle(const double angle) { m_angle = angle; }
47
48 /** Plot the drift velocity components of electrons in the medium.
49 * \param xaxis abscissa.
50 * - 'e': electric field,
51 * - 'b': magnetic field,
52 * - 'a': angle between E and B.
53 * \param same flag to keep existing plots (true) or not.
54 */
55 void PlotElectronVelocity(const char xaxis, const bool same = false) {
56 PlotVelocity(GetAxis(xaxis), Charge::Electron, same);
57 }
58 /// Plot the drift velocity components of holes in the medium.
59 void PlotHoleVelocity(const char xaxis, const bool same = false) {
60 PlotVelocity(GetAxis(xaxis), Charge::Hole, same);
61 }
62 /// Plot the ion drift velocity in the gas.
63 void PlotIonVelocity(const char xaxis, const bool same = false) {
64 PlotVelocity(GetAxis(xaxis), Charge::Ion, same);
65 }
66 /// Plot the diffusion coefficients in the medium.
67 void PlotElectronDiffusion(const char xaxis, const bool same = false) {
68 PlotDiffusion(GetAxis(xaxis), Charge::Electron, same);
69 }
70 /// Plot the diffusion coefficients of holes in the medium.
71 void PlotHoleDiffusion(const char xaxis, const bool same = false) {
72 PlotDiffusion(GetAxis(xaxis), Charge::Hole, same);
73 }
74 /// Plot the diffusion coefficients of ions in the gas.
75 void PlotIonDiffusion(const char xaxis, const bool same = false) {
76 PlotDiffusion(GetAxis(xaxis), Charge::Ion, same);
77 }
78 /// Plot the Townsend coefficient for electrons.
79 void PlotElectronTownsend(const char xaxis, const bool same = false) {
80 Plot(GetAxis(xaxis), Charge::Electron, Parameter::Townsend, same);
81 }
82 /// Plot the Townsend coefficient for holes.
83 void PlotHoleTownsend(const char xaxis, const bool same = false) {
84 Plot(GetAxis(xaxis), Charge::Hole, Parameter::Townsend, same);
85 }
86 /// Plot the attachment coefficient for electrons.
87 void PlotElectronAttachment(const char xaxis, const bool same = false) {
88 Plot(GetAxis(xaxis), Charge::Electron, Parameter::Attachment, same);
89 }
90 /// Plot the attachment coefficient for holes.
91 void PlotHoleAttachment(const char xaxis, const bool same = false) {
92 Plot(GetAxis(xaxis), Charge::Hole, Parameter::Attachment, same);
93 }
94 /// Plot the angle between drift velocity and field.
95 void PlotElectronLorentzAngle(const char xaxis, const bool same = false) {
96 PlotLorentzAngle(GetAxis(xaxis), Charge::Electron, same);
97 }
99
100 /// Set the (ROOT) colours to be used in the plots.
101 void SetColours(const std::vector<short>& cols) { m_colours = cols; }
102 /// Set user-defined plot labels.
103 void SetLabels(const std::vector<std::string>& labels) { m_labels = labels; }
104
105 private:
106
107 enum class Parameter {
108 VelocityE,
109 VelocityB,
110 VelocityExB,
111 TransverseDiffusion,
112 LongitudinalDiffusion,
113 Townsend,
114 Attachment,
115 LorentzAngle
116 };
117
118 enum class Charge {
119 Electron,
120 Hole,
121 Ion
122 };
123
124 enum class Axis {
125 E,
126 B,
127 Angle,
128 None
129 };
130
131 Medium* m_medium = nullptr;
132
133 // X axis
134 double m_eMin = 100., m_eMax = 100000.;
135 double m_bMin = 0., m_bMax = 2.;
136 double m_aMin = 0., m_aMax = Pi;
137 bool m_logE = true;
138 bool m_logB = false;
139 bool m_logA = false;
140 bool m_logX = true;
141 bool m_autoRangeX = true;
142 Axis m_xaxis = Axis::None;
143
144 // Y axis
145 double m_yMin = 0., m_yMax = 1.;
146 bool m_logY = false;
147 bool m_autoRangeY = true;
148
149 // E-field to use when plotting as function of B-field or angle.
150 double m_efield = 1000.;
151 // B-field to use when plotting as function of E-field or angle.
152 double m_bfield = 0.;
153 // Angle to use when plotting as function of E-field or B-field.
154 double m_angle = HalfPi;
155
156 std::vector<double> m_xPlot;
157 std::vector<std::vector<double> > m_yPlot;
158 std::vector<Parameter> m_par;
159 std::vector<Charge> m_q;
160
161 std::vector<std::vector<double> > m_xGraph;
162 std::vector<std::vector<double> > m_yGraph;
163
164 std::vector<short> m_colours;
165 std::vector<std::string> m_labels;
166
167 void PlotVelocity(const Axis xaxis, const Charge particle,
168 const bool same);
169 void PlotDiffusion(const Axis xaxis, const Charge particle,
170 const bool same);
171 void Plot(const Axis xaxis, const Charge particle,
172 const Parameter par, const bool same);
173 void PlotLorentzAngle(const Axis xaxis, const Charge particle,
174 const bool same);
175
176 void ResetY();
177 void ResetX(const Axis xaxis);
178 void Draw();
179
180 Axis GetAxis(const char xaxis) const;
181 bool GetGrid(std::array<std::vector<double>, 3>& grid,
182 int& ie, int& ib, int& ia, const Axis xaxis) const;
183};
184}
185#endif
Abstract base class for media.
Definition: Medium.hh:13
Base class for visualization classes.
Definition: ViewBase.hh:18
Plot transport coefficients as function of electric and magnetic field.
Definition: ViewMedium.hh:17
void SetRangeY(const double ymin, const double ymax, const bool logscale=false)
Set the range of the function (velocity etc.) to be plotted.
Definition: ViewMedium.cc:87
void SetAngle(const double angle)
Set the angle to use when plotting as function of E or B.
Definition: ViewMedium.hh:46
void SetRangeB(const double bmin, const double bmax, const bool logscale)
Set the limits of the magnetic field.
Definition: ViewMedium.cc:65
void SetLabels(const std::vector< std::string > &labels)
Set user-defined plot labels.
Definition: ViewMedium.hh:103
void PlotIonVelocity(const char xaxis, const bool same=false)
Plot the ion drift velocity in the gas.
Definition: ViewMedium.hh:63
void PlotElectronTownsend(const char xaxis, const bool same=false)
Plot the Townsend coefficient for electrons.
Definition: ViewMedium.hh:79
void SetMedium(Medium *m)
Set the medium from which to retrieve the transport coefficients.
Definition: ViewMedium.cc:44
void PlotHoleVelocity(const char xaxis, const bool same=false)
Plot the drift velocity components of holes in the medium.
Definition: ViewMedium.hh:59
void SetMagneticField(const double bfield)
Set the magnetic field to use when plotting as function of E or angle.
Definition: ViewMedium.hh:44
void PlotHoleDiffusion(const char xaxis, const bool same=false)
Plot the diffusion coefficients of holes in the medium.
Definition: ViewMedium.hh:71
~ViewMedium()=default
Destructor.
void PlotElectronAttachment(const char xaxis, const bool same=false)
Plot the attachment coefficient for electrons.
Definition: ViewMedium.hh:87
void SetColours(const std::vector< short > &cols)
Set the (ROOT) colours to be used in the plots.
Definition: ViewMedium.hh:101
void PlotIonDiffusion(const char xaxis, const bool same=false)
Plot the diffusion coefficients of ions in the gas.
Definition: ViewMedium.hh:75
void PlotElectronCrossSections()
Definition: ViewMedium.cc:99
void PlotElectronLorentzAngle(const char xaxis, const bool same=false)
Plot the angle between drift velocity and field.
Definition: ViewMedium.hh:95
void PlotElectronVelocity(const char xaxis, const bool same=false)
Definition: ViewMedium.hh:55
void EnableAutoRangeY(const bool on=true)
Choose the y-axis range based on the function's minima/maxima.
Definition: ViewMedium.hh:36
void SetRangeE(const double emin, const double emax, const bool logscale)
Set the limits of the electric field.
Definition: ViewMedium.cc:53
void SetElectricField(const double efield)
Set the electric field to use when plotting as function of B or angle.
Definition: ViewMedium.hh:42
void PlotHoleAttachment(const char xaxis, const bool same=false)
Plot the attachment coefficient for holes.
Definition: ViewMedium.hh:91
void SetRangeA(const double amin, const double amax, const bool logscale)
Set the limits of the angle between electric and magnetic field.
Definition: ViewMedium.cc:76
void EnableAutoRangeX(const bool on=true)
Try to choose the x-axis range based on the field grid.
Definition: ViewMedium.hh:28
ViewMedium()
Constructor.
Definition: ViewMedium.cc:42
void PlotElectronDiffusion(const char xaxis, const bool same=false)
Plot the diffusion coefficients in the medium.
Definition: ViewMedium.hh:67
void PlotHoleTownsend(const char xaxis, const bool same=false)
Plot the Townsend coefficient for holes.
Definition: ViewMedium.hh:83