Garfield++ v2r0
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
7#include <TCanvas.h>
8#include <TF1.h>
9#include <TGraph.h>
10
12
13namespace Garfield {
14
15class Medium;
16
17/// Plot transport coefficients as function of electric and magnetic field.
18
20
21 public:
22 /// Constructor
23 ViewMedium();
24 /// Destructor
26
27 /// Set the canvas to be painted on.
28 void SetCanvas(TCanvas* c);
29 // Set the medium from which to retrieve the transport coefficients.
30 void SetMedium(Medium* m);
31
32 /// Set the limits of the electric field.
33 void SetElectricFieldRange(const double emin, const double emax);
34 /// Set the limits of the magnetic field.
35 void SetMagneticFieldRange(const double bmin, const double bmax);
36 /// Set the limits of the angle between electric and magnetic field.
37 void SetBAngleRange(const double amin, const double amax);
38 /// Set the range of the function (velocity etc.) to be plotted.
39 void SetFunctionRange(const double vmin, const double vmax);
40 /// Use the default function range.
41 void SetFunctionRange();
42
43 void PlotElectronVelocity(const char xaxis, const double e, const double b,
44 const double a = HalfPi);
45 void PlotHoleVelocity(const char xaxis, const double e, const double b,
46 const double a = HalfPi);
47 void PlotIonVelocity(const char xaxis, const double e, const double b,
48 const double a = HalfPi);
49 void PlotElectronDiffusion(const char xaxis, const double e, const double b,
50 const double a = HalfPi);
51 void PlotHoleDiffusion(const char xaxis, const double e, const double b,
52 const double a = HalfPi);
53 void PlotIonDiffusion(const char xaxis, const double e, const double b,
54 const double a = HalfPi);
55 void PlotElectronTownsend(const char xaxis, const double e, const double b,
56 const double a = HalfPi);
57 void PlotHoleTownsend(const char xaxis, const double e, const double b,
58 const double a = HalfPi);
59 void PlotElectronAttachment(const char xaxis, const double e, const double b,
60 const double a = HalfPi);
61 void PlotHoleAttachment(const char xaxis, const double e, const double b,
62 const double a = HalfPi);
63 void PlotElectronLorentzAngle(const char xaxis, const double e, const double b,
64 const double a = HalfPi);
66 double EvaluateFunction(double* pos, double* par);
67
68 enum Property {
87 };
88
89 private:
90 std::string m_className;
91
92 // Options
93 bool m_debug;
94
95 // Canvas
96 TCanvas* m_canvas;
97 bool m_hasExternalCanvas;
98
99 Medium* m_medium;
100
101 // Ranges for variable parameters
102 double m_eMin, m_eMax;
103 double m_bMin, m_bMax;
104 double m_aMin, m_aMax;
105 double m_vMin, m_vMax;
106
107 // Fixed parameters
108 double m_efield;
109 double m_bfield;
110 double m_angle;
111
112 // Tolerances for marker plotting
113 double m_etolerance;
114 double m_btolerance;
115 double m_atolerance;
116
117 // Labels
118 std::string m_labele;
119 std::string m_labelb;
120 std::string m_labela;
121 std::string m_labelv;
122 std::string m_labeld;
123
124 // Functions
125 std::vector<TF1> m_functions;
126 // Graphs
127 std::vector<TGraph> m_graphs;
128
129 void SetupCanvas();
130 void AddFunction(const double xmin, const double xmax, const double ymin,
131 const double ymax, const bool keep,
132 const std::string& xlabel, const std::string& ylabel,
133 const int type, const char xaxis,
134 const double e, const double b, const double a);
135 int GetColor(const unsigned int property) const;
136};
137}
138#endif
Abstract base class for media.
Definition: Medium.hh:11
Plot transport coefficients as function of electric and magnetic field.
Definition: ViewMedium.hh:19
double EvaluateFunction(double *pos, double *par)
Definition: ViewMedium.cc:619
void PlotElectronVelocity(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:119
void SetMedium(Medium *m)
Definition: ViewMedium.cc:63
void PlotIonDiffusion(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:218
void SetElectricFieldRange(const double emin, const double emax)
Set the limits of the electric field.
Definition: ViewMedium.cc:73
void PlotHoleTownsend(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:241
void PlotIonVelocity(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:181
void PlotElectronLorentzAngle(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:272
void PlotElectronTownsend(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:231
void PlotHoleVelocity(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:150
void PlotHoleAttachment(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:262
void PlotElectronCrossSections()
Definition: ViewMedium.cc:283
void SetCanvas(TCanvas *c)
Set the canvas to be painted on.
Definition: ViewMedium.cc:49
~ViewMedium()
Destructor.
Definition: ViewMedium.cc:44
void PlotElectronAttachment(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:251
void SetFunctionRange()
Use the default function range.
Definition: ViewMedium.cc:117
ViewMedium()
Constructor.
Definition: ViewMedium.cc:15
void SetMagneticFieldRange(const double bmin, const double bmax)
Set the limits of the magnetic field.
Definition: ViewMedium.cc:84
void PlotElectronDiffusion(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:191
void PlotHoleDiffusion(const char xaxis, const double e, const double b, const double a=HalfPi)
Definition: ViewMedium.cc:205
void SetBAngleRange(const double amin, const double amax)
Set the limits of the angle between electric and magnetic field.
Definition: ViewMedium.cc:95