Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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 <array>
5#include <string>
6#include <vector>
7
10
11namespace Garfield {
12
13class Medium;
14
16
17class ViewMedium : public ViewBase {
18 public:
20 ViewMedium() : ViewMedium(nullptr) {}
22 ViewMedium(Medium* medium);
24 ~ViewMedium() = default;
25
27 void SetMedium(Medium* m);
28
30 void EnableAutoRangeX(const bool on = true) { m_autoRangeX = on; }
32 void SetRangeE(const double emin, const double emax, const bool logscale);
34 void SetRangeB(const double bmin, const double bmax, const bool logscale);
36 void SetRangeA(const double amin, const double amax, const bool logscale);
38 void SetRangeEN(const double emin, const double emax, const bool logscale);
40 void SetRangeEP(const double emin, const double emax, const bool logscale);
42 void EnableAutoRangeY(const bool on = true) { m_autoRangeY = on; }
44 void SetRangeY(const double ymin, const double ymax,
45 const bool logscale = false);
46
48 void SetElectricField(const double efield) { m_efield = efield; }
50 void SetMagneticField(const double bfield) { m_bfield = bfield; }
52 void SetAngle(const double angle) { m_angle = angle; }
53
54 void EnableExport(const std::string& txtfile) { m_outfile = txtfile; }
55 void DisableExport() { m_outfile = ""; }
56
71 void PlotVelocity(const std::string& carriers, const char xaxis);
73 void PlotDiffusion(const std::string& carriers, const char xaxis);
75 void PlotTownsend(const std::string& carriers, const char xaxis);
77 void PlotAttachment(const std::string& carriers, const char xaxis);
79 void PlotAlphaEta(const std::string& carriers, const char xaxis);
80
90 void PlotElectronVelocity(const char xaxis = 'e', const bool same = false) {
92 }
93
94 void PlotElectronVelocityFluxBulk(const char xaxis = 'e',
95 const bool same = false) {
97 }
98
99 void PlotHoleVelocity(const char xaxis = 'e', const bool same = false) {
100 PlotVelocity(GetAxis(xaxis), Charge::Hole, same);
101 }
102
103 void PlotIonVelocity(const char xaxis = 'e', const bool same = false) {
104 PlotVelocity(GetAxis(xaxis), Charge::Ion, same);
105 }
106
107 void PlotElectronDiffusion(const char xaxis = 'e', const bool same = false) {
109 }
110
111 void PlotHoleDiffusion(const char xaxis = 'e', const bool same = false) {
112 PlotDiffusion(GetAxis(xaxis), Charge::Hole, same);
113 }
114
115 void PlotIonDiffusion(const char xaxis = 'e', const bool same = false) {
116 PlotDiffusion(GetAxis(xaxis), Charge::Ion, same);
117 }
118
119 void PlotElectronTownsend(const char xaxis = 'e', const bool same = false) {
121 }
122
123 void PlotElectronReducedTownsendN(const char xaxis = 'r',
124 const bool same = false) {
126 }
127
128 void PlotElectronReducedTownsendP(const char xaxis = 'p',
129 const bool same = false) {
131 }
132
133 void PlotElectronTOFIonization(const char xaxis = 'e',
134 const bool same = false) {
136 }
137
138 void PlotHoleTownsend(const char xaxis = 'e', const bool same = false) {
140 }
141
142 void PlotElectronAttachment(const char xaxis = 'e', const bool same = false) {
144 }
145
146 void PlotElectronTOFAttachment(const char xaxis = 'e',
147 const bool same = false) {
149 }
150
151 void PlotHoleAttachment(const char xaxis = 'e', const bool same = false) {
153 }
154
156 void PlotElectronLorentzAngle(const char xaxis = 'e',
157 const bool same = false) {
159 }
160
162 void SetColours(const std::vector<short>& cols) { m_colours = cols; }
164 void SetLabels(const std::vector<std::string>& labels) { m_labels = labels; }
165
166 private:
183
184 enum class Charge { Electron, Hole, Ion };
185
186 enum class Axis { E, B, Angle, EoverN, EoverP, None };
187
188 Medium* m_medium = nullptr;
189
190 // X axis
191 double m_eMin = 100., m_eMax = 100000.;
192 double m_bMin = 0., m_bMax = 2.;
193 double m_aMin = 0., m_aMax = Pi;
194 double m_enMin = 0.25, m_enMax = 400.;
195 double m_epMin = 0.13, m_epMax = 132.;
196 bool m_logE = true;
197 bool m_logB = false;
198 bool m_logA = false;
199 bool m_logEN = true;
200 bool m_logEP = true;
201 bool m_logX = true;
202 bool m_autoRangeX = true;
204
205 // Y axis
206 double m_yMin = 0., m_yMax = 1.;
207 bool m_logY = false;
208 bool m_autoRangeY = true;
209
210 // E-field to use when plotting as function of B-field or angle.
211 double m_efield = 1000.;
212 // B-field to use when plotting as function of E-field or angle.
213 double m_bfield = 0.;
214 // Angle to use when plotting as function of E-field or B-field.
215 double m_angle = HalfPi;
216
217 std::vector<double> m_xPlot;
218 std::vector<std::vector<double> > m_yPlot;
219 std::vector<Parameter> m_par;
220 std::vector<Charge> m_q;
221
222 std::vector<std::vector<double> > m_xGraph;
223 std::vector<std::vector<double> > m_yGraph;
224
225 std::vector<short> m_colours;
226 std::vector<std::string> m_labels;
227
228 std::string m_outfile;
229
230 void PlotVelocity(const Axis xaxis, const Charge particle, const bool same);
231 void PlotVelocityFluxBulk(const Axis xaxis, const Charge particle,
232 const bool same);
233 void PlotDiffusion(const Axis xaxis, const Charge particle, const bool same);
234 void Plot(const Axis xaxis, const Charge particle, const Parameter par,
235 const bool same);
236 void PlotLorentzAngle(const Axis xaxis, const Charge particle,
237 const bool same);
238
239 void ResetY();
240 void ResetX(const Axis xaxis);
241 void Draw();
242 void Export();
243
244 Axis GetAxis(const char xaxis) const;
245 bool GetGrid(std::array<std::vector<double>, 3>& grid, int& ie, int& ib,
246 int& ia, const Axis xaxis) const;
247
248 double ConvertToEN(const double e0Vcm);
249 double UnConvertFromEN(const double e0Td);
250 double ConvertToEP(const double e0Vcm);
251 double UnConvertFromEP(const double e0VcmTorr);
252};
253} // namespace Garfield
254#endif
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)
Plot the drift velocity components of electrons in the medium.
Definition ViewMedium.hh:90
void SetRangeY(const double ymin, const double ymax, const bool logscale=false)
Set the range of the function (velocity etc.) to be plotted.
void SetAngle(const double angle)
Set the angle to use when plotting as function of E or B.
Definition ViewMedium.hh:52
void SetRangeB(const double bmin, const double bmax, const bool logscale)
Set the limits of the magnetic field.
void PlotElectronTownsend(const char xaxis='e', const bool same=false)
Plot the Townsend coefficient for electrons.
void PlotElectronReducedTownsendP(const char xaxis='p', 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 PlotElectronVelocityFluxBulk(const char xaxis='e', const bool same=false)
Plot Flux and Bulk drift velocity.
Definition ViewMedium.hh:94
void SetRangeEN(const double emin, const double emax, const bool logscale)
Set the limits of the reduced electric field (E/N).
std::vector< std::vector< double > > m_xGraph
double ConvertToEN(const double e0Vcm)
void SetMedium(Medium *m)
Set the medium from which to retrieve the transport coefficients.
std::vector< std::string > m_labels
void SetMagneticField(const double bfield)
Set the magnetic field to use when plotting as function of E or angle.
Definition ViewMedium.hh:50
void PlotLorentzAngle(const Axis xaxis, const Charge particle, const bool same)
double UnConvertFromEP(const double e0VcmTorr)
void PlotAlphaEta(const std::string &carriers, const char xaxis)
Plot Townsend and attachment coefficients.
void PlotElectronTOFIonization(const char xaxis='e', const bool same=false)
Plot the TOF ionization rate.
void PlotElectronLorentzAngle(const char xaxis='e', const bool same=false)
Plot the angle between drift velocity and field.
void PlotElectronReducedTownsendN(const char xaxis='r', const bool same=false)
Plot the Townsend coefficient for electrons.
void PlotIonDiffusion(const char xaxis='e', const bool same=false)
Plot the diffusion coefficients of ions in the gas.
std::vector< Parameter > m_par
ViewMedium(Medium *medium)
Constructor.
void PlotHoleAttachment(const char xaxis='e', const bool same=false)
Plot the attachment coefficient for holes.
std::string m_outfile
Axis GetAxis(const char xaxis) const
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:54
void PlotHoleVelocity(const char xaxis='e', const bool same=false)
Plot the drift velocity components of holes in the medium.
Definition ViewMedium.hh:99
void SetColours(const std::vector< short > &cols)
Set the (ROOT) colours to be used in the plots.
std::vector< Charge > m_q
void PlotDiffusion(const std::string &carriers, const char xaxis)
Plot the transverse and longitudinal diffusion coefficients.
double UnConvertFromEN(const double e0Td)
void ResetX(const Axis xaxis)
void PlotIonVelocity(const char xaxis='e', const bool same=false)
Plot the ion drift velocity in the gas.
void SetRangeEP(const double emin, const double emax, const bool logscale)
Set the limits of the reduced electric field (E/P).
void EnableAutoRangeY(const bool on=true)
Choose the y-axis range based on the function's minima/maxima.
Definition ViewMedium.hh:42
void SetRangeE(const double emin, const double emax, const bool logscale)
Set the limits of the electric field.
bool GetGrid(std::array< std::vector< double >, 3 > &grid, int &ie, int &ib, int &ia, const Axis xaxis) const
void PlotVelocityFluxBulk(const Axis xaxis, const Charge particle, const bool same)
void SetElectricField(const double efield)
Set the electric field to use when plotting as function of B or angle.
Definition ViewMedium.hh:48
double ConvertToEP(const double e0Vcm)
void PlotDiffusion(const Axis xaxis, const Charge particle, const bool same)
void PlotElectronDiffusion(const char xaxis='e', const bool same=false)
Plot the diffusion coefficients in the medium.
std::vector< double > m_xPlot
void Plot(const Axis xaxis, const Charge particle, const Parameter par, const bool same)
void PlotVelocity(const std::string &carriers, const char xaxis)
Plot the drift velocity components.
void SetRangeA(const double amin, const double amax, const bool logscale)
Set the limits of the angle between electric and magnetic field.
void EnableAutoRangeX(const bool on=true)
Try to choose the x-axis range based on the field grid.
Definition ViewMedium.hh:30
void PlotVelocity(const Axis xaxis, const Charge particle, const bool same)
std::vector< std::vector< double > > m_yGraph
ViewMedium()
Default constructor.
Definition ViewMedium.hh:20
void PlotElectronTOFAttachment(const char xaxis='e', const bool same=false)
Plot the TOF attachment rate.
std::vector< std::vector< double > > m_yPlot
void PlotHoleDiffusion(const char xaxis='e', const bool same=false)
Plot the diffusion coefficients of holes in the medium.
std::vector< short > m_colours