Garfield++ 5.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 void EnableExport(const std::string& txtfile) { m_outfile = txtfile; }
49 void DisableExport() { m_outfile = ""; }
50
51 /** Plot the drift velocity components.
52 * \param option string indicating the carriers for which to plot
53 * the drift velocity
54 * - "e": electrons,
55 * - "h": holes,
56 * - "i": ions
57 * Options can be concatenated (e. g. "ei", "eh").
58 * \param xaxis abscissa.
59 * - 'e': electric field,
60 * - 'b': magnetic field,
61 * - 'a': angle between E and B.
62 */
63 void PlotVelocity(const std::string& carriers, const char xaxis);
64 /// Plot the transverse and longitudinal diffusion coefficients.
65 void PlotDiffusion(const std::string& carriers, const char xaxis);
66 /// Plot the Townsend coefficient.
67 void PlotTownsend(const std::string& carriers, const char xaxis);
68 /// Plot the attachment coefficient.
69 void PlotAttachment(const std::string& carriers, const char xaxis);
70 /// Plot Townsend and attachment coefficients.
71 void PlotAlphaEta(const std::string& carriers, const char xaxis);
72
73 /** Plot the drift velocity components of electrons in the medium.
74 * \param xaxis abscissa.
75 * - 'e': electric field,
76 * - 'b': magnetic field,
77 * - 'a': angle between E and B.
78 * \param same flag to keep existing plots (true) or not.
79 */
80 void PlotElectronVelocity(const char xaxis = 'e', const bool same = false) {
81 PlotVelocity(GetAxis(xaxis), Charge::Electron, same);
82 }
83 /// Plot the drift velocity components of holes in the medium.
84 void PlotHoleVelocity(const char xaxis = 'e', const bool same = false) {
85 PlotVelocity(GetAxis(xaxis), Charge::Hole, same);
86 }
87 /// Plot the ion drift velocity in the gas.
88 void PlotIonVelocity(const char xaxis = 'e', const bool same = false) {
89 PlotVelocity(GetAxis(xaxis), Charge::Ion, same);
90 }
91 /// Plot the diffusion coefficients in the medium.
92 void PlotElectronDiffusion(const char xaxis = 'e', const bool same = false) {
93 PlotDiffusion(GetAxis(xaxis), Charge::Electron, same);
94 }
95 /// Plot the diffusion coefficients of holes in the medium.
96 void PlotHoleDiffusion(const char xaxis = 'e', const bool same = false) {
97 PlotDiffusion(GetAxis(xaxis), Charge::Hole, same);
98 }
99 /// Plot the diffusion coefficients of ions in the gas.
100 void PlotIonDiffusion(const char xaxis = 'e', const bool same = false) {
101 PlotDiffusion(GetAxis(xaxis), Charge::Ion, same);
102 }
103 /// Plot the Townsend coefficient for electrons.
104 void PlotElectronTownsend(const char xaxis = 'e', const bool same = false) {
105 Plot(GetAxis(xaxis), Charge::Electron, Parameter::Townsend, same);
106 }
107 /// Plot the Townsend coefficient for holes.
108 void PlotHoleTownsend(const char xaxis = 'e', const bool same = false) {
109 Plot(GetAxis(xaxis), Charge::Hole, Parameter::Townsend, same);
110 }
111 /// Plot the attachment coefficient for electrons.
112 void PlotElectronAttachment(const char xaxis = 'e', const bool same = false) {
113 Plot(GetAxis(xaxis), Charge::Electron, Parameter::Attachment, same);
114 }
115 /// Plot the attachment coefficient for holes.
116 void PlotHoleAttachment(const char xaxis = 'e', const bool same = false) {
117 Plot(GetAxis(xaxis), Charge::Hole, Parameter::Attachment, same);
118 }
119 /// Plot the angle between drift velocity and field.
120 void PlotElectronLorentzAngle(const char xaxis = 'e', const bool same = false) {
121 PlotLorentzAngle(GetAxis(xaxis), Charge::Electron, same);
122 }
123
124 /// Set the (ROOT) colours to be used in the plots.
125 void SetColours(const std::vector<short>& cols) { m_colours = cols; }
126 /// Set user-defined plot labels.
127 void SetLabels(const std::vector<std::string>& labels) { m_labels = labels; }
128
129 private:
130
131 enum class Parameter {
132 VelocityE,
133 VelocityB,
134 VelocityExB,
135 TransverseDiffusion,
136 LongitudinalDiffusion,
137 Townsend,
138 Attachment,
139 LorentzAngle
140 };
141
142 enum class Charge {
143 Electron,
144 Hole,
145 Ion
146 };
147
148 enum class Axis {
149 E,
150 B,
151 Angle,
152 None
153 };
154
155 Medium* m_medium = nullptr;
156
157 // X axis
158 double m_eMin = 100., m_eMax = 100000.;
159 double m_bMin = 0., m_bMax = 2.;
160 double m_aMin = 0., m_aMax = Pi;
161 bool m_logE = true;
162 bool m_logB = false;
163 bool m_logA = false;
164 bool m_logX = true;
165 bool m_autoRangeX = true;
166 Axis m_xaxis = Axis::None;
167
168 // Y axis
169 double m_yMin = 0., m_yMax = 1.;
170 bool m_logY = false;
171 bool m_autoRangeY = true;
172
173 // E-field to use when plotting as function of B-field or angle.
174 double m_efield = 1000.;
175 // B-field to use when plotting as function of E-field or angle.
176 double m_bfield = 0.;
177 // Angle to use when plotting as function of E-field or B-field.
178 double m_angle = HalfPi;
179
180 std::vector<double> m_xPlot;
181 std::vector<std::vector<double> > m_yPlot;
182 std::vector<Parameter> m_par;
183 std::vector<Charge> m_q;
184
185 std::vector<std::vector<double> > m_xGraph;
186 std::vector<std::vector<double> > m_yGraph;
187
188 std::vector<short> m_colours;
189 std::vector<std::string> m_labels;
190
191 std::string m_outfile;
192
193 void PlotVelocity(const Axis xaxis, const Charge particle,
194 const bool same);
195 void PlotDiffusion(const Axis xaxis, const Charge particle,
196 const bool same);
197 void Plot(const Axis xaxis, const Charge particle,
198 const Parameter par, const bool same);
199 void PlotLorentzAngle(const Axis xaxis, const Charge particle,
200 const bool same);
201
202 void ResetY();
203 void ResetX(const Axis xaxis);
204 void Draw();
205 void Export();
206
207 Axis GetAxis(const char xaxis) const;
208 bool GetGrid(std::array<std::vector<double>, 3>& grid,
209 int& ie, int& ib, int& ia, const Axis xaxis) const;
210};
211}
212#endif
Abstract base class for media.
Definition Medium.hh:16
ViewBase()=delete
Default constructor.
void PlotHoleTownsend(const char xaxis='e', const bool same=false)
Plot the Townsend coefficient for holes.
void PlotElectronVelocity(const char xaxis='e', const bool same=false)
Definition ViewMedium.hh:80
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:88
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:66
void PlotElectronTownsend(const char xaxis='e', const bool same=false)
Plot the Townsend coefficient for electrons.
void PlotAttachment(const std::string &carriers, const char xaxis)
Plot the attachment coefficient.
void PlotElectronAttachment(const char xaxis='e', const bool same=false)
Plot the attachment coefficient for electrons.
void SetLabels(const std::vector< std::string > &labels)
Set user-defined plot labels.
void SetMedium(Medium *m)
Set the medium from which to retrieve the transport coefficients.
Definition ViewMedium.cc:45
void SetMagneticField(const double bfield)
Set the magnetic field to use when plotting as function of E or angle.
Definition ViewMedium.hh:44
void PlotAlphaEta(const std::string &carriers, const char xaxis)
Plot Townsend and attachment coefficients.
void PlotElectronLorentzAngle(const char xaxis='e', const bool same=false)
Plot the angle between drift velocity and field.
void PlotIonDiffusion(const char xaxis='e', const bool same=false)
Plot the diffusion coefficients of ions in the gas.
void PlotHoleAttachment(const char xaxis='e', const bool same=false)
Plot the attachment coefficient for holes.
void PlotTownsend(const std::string &carriers, const char xaxis)
Plot the Townsend coefficient.
~ViewMedium()=default
Destructor.
void EnableExport(const std::string &txtfile)
Definition ViewMedium.hh:48
void PlotHoleVelocity(const char xaxis='e', const bool same=false)
Plot the drift velocity components of holes in the medium.
Definition ViewMedium.hh:84
void SetColours(const std::vector< short > &cols)
Set the (ROOT) colours to be used in the plots.
void PlotDiffusion(const std::string &carriers, const char xaxis)
Plot the transverse and longitudinal diffusion coefficients.
void PlotIonVelocity(const char xaxis='e', const bool same=false)
Plot the ion drift velocity in the gas.
Definition ViewMedium.hh:88
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:54
void SetElectricField(const double efield)
Set the electric field to use when plotting as function of B or angle.
Definition ViewMedium.hh:42
void PlotElectronDiffusion(const char xaxis='e', const bool same=false)
Plot the diffusion coefficients in the medium.
Definition ViewMedium.hh:92
void PlotVelocity(const std::string &carriers, const char xaxis)
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:77
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:43
void PlotHoleDiffusion(const char xaxis='e', const bool same=false)
Plot the diffusion coefficients of holes in the medium.
Definition ViewMedium.hh:96