Garfield++ v2r0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewField.hh
Go to the documentation of this file.
1#ifndef G_VIEW_FIELD
2#define G_VIEW_FIELD
3
4#include <TCanvas.h>
5#include <TF2.h>
6#include <TF1.h>
7
8namespace Garfield {
9
10class Sensor;
11class ComponentBase;
12
13/// Visualize the potential or electric field of a component or sensor.
14
15class ViewField {
16
17 public:
18 /// Constructor
19 ViewField();
20 /// Destructor
21 ~ViewField();
22
23 /// Set the sensor from which to retrieve the field.
24 void SetSensor(Sensor* s);
25 /// Set the component from which to retrieve the field.
27 /// Set the canvas to be painted on.
28 void SetCanvas(TCanvas* c);
29
30 /// Set the plot limits for the potential.
31 void SetVoltageRange(const double vmin, const double vmax);
32 /// Set the plot limits for the electric field.
33 void SetElectricFieldRange(const double emin, const double emax);
34 /// Set the plot limits for the weighting field.
35 void SetWeightingFieldRange(const double wmin, const double wmax);
36
37 /// Set the viewing area (in local coordinates of the current viewing plane).
38 void SetArea(const double xmin, const double ymin,
39 const double xmax, const double ymax);
40 /// Set the viewing area based on the bounding box of the sensor/component.
41 void SetArea() { m_hasUserArea = false; }
42
43 /** Set the projection (viewing plane).
44 * \param fx,fy,fz normal vector
45 * \param x0,y0,z0 in-plane point
46 */
47 void SetPlane(const double fx, const double fy, const double fz,
48 const double x0, const double y0, const double z0);
49 /// Set the default viewing plane (\f$x\f$-\f$y\f$ at \f$z = 0\f$).
51 /// Rotate the viewing plane (angle in radian).
52 void Rotate(const double angle);
53
54 /// Set the number of contour levels (at most 50).
55 void SetNumberOfContours(const unsigned int n);
56 /// Set the number of points used for drawing 1D functions.
57 void SetNumberOfSamples1d(const unsigned int n);
58 /// Set the number of points used for drawing 2D functions.
59 void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny);
60
61 /** Make a contour plot of the electric potential or field.
62 * \param option quantity to be plotted
63 * - potential: "v", "voltage", "p", "potential"
64 * - magnitude of the electric field: "e", "field"
65 * - x-component of the electric field: "ex"
66 * - y-component of the electric field: "ey"
67 * - z-component of the electric field: "ez"
68 **/
69 void PlotContour(const std::string& option = "v");
70 /** Make a surface plot ("SURF4") of the electric potential or field.
71 * \param option quantity to be plotted (see PlotContour)
72 **/
73 void PlotSurface(const std::string& option = "v") { Plot(option, "SURF4"); }
74 /** Make a 2D plot of the electric potential or field.
75 * \param option quantity to be plotted (see PlotContour)
76 * \param drawopt option string passed to TF2::Draw
77 **/
78 void Plot(const std::string& option = "v",
79 const std::string& drawopt = "arr");
80 /** Make a 1D plot of the electric potential or field along a line.
81 * \param x0,y0,z0 starting point
82 * \param x1,y1,z1 end point
83 * \param quantity to be plotted (see PlotContour)
84 **/
85 void PlotProfile(const double x0, const double y0, const double z0,
86 const double x1, const double y1, const double z1,
87 const std::string& option = "v");
88
89 /** Make a contour plot of the weighting potential or field.
90 * \param label identifier of the electrode
91 * \param option quantity to be plotted (see PlotContour)
92 **/
93 void PlotContourWeightingField(const std::string& label,
94 const std::string& option);
95 /** Make a surface plot ("SURF4") of the weighting potential or field.
96 * \param label identifier of the electrode
97 * \param option quantity to be plotted (see PlotContour)
98 **/
99 void PlotSurfaceWeightingField(const std::string& label,
100 const std::string& option) {
101 PlotWeightingField(label, option, "SURF4");
102 }
103 /** Make a 2D plot of the weighting potential or field.
104 * \param label identifier of the electrode
105 * \param option quantity to be plotted (see PlotContour)
106 * \param drawopt option string passed to TF2::Draw
107 **/
108 void PlotWeightingField(const std::string& label,
109 const std::string& option,
110 const std::string& drawopt);
111
112 /** Make a 1D plot of the weighting potential or field along a line.
113 * \param label identifier of the electrode
114 * \param x0,y0,z0 starting point
115 * \param x1,y1,z1 end point
116 * \param quantity to be plotted (see PlotContour)
117 **/
118 void PlotProfileWeightingField(const std::string& label,
119 const double x0, const double y0, const double z0,
120 const double x1, const double y1, const double z1,
121 const std::string& option = "v");
122
123 void EnableAutoRange(const bool on = true) { m_useAutoRange = on; }
124
125 /** Make use of the status flag returned by the sensor/component.
126 * \param v0 Value to be used for regions with status != 0.
127 */
128 void EnableAcknowledgeStatus(const double v0 = 0.) {
129 m_useStatus = true;
130 m_vBkg = v0;
131 }
132 /// Ignore the status flag returned by the sensor/component.
133 void DisableAcknowledgeStatus() { m_useStatus = false; }
134
135 /// Switch on/off debugging output.
136 void EnableDebugging(const bool on = true) { m_debug = on; }
137
138 friend class TF1;
139 friend class TF2;
140
141 protected:
142 // Functions called by TF1/TF2.
143 double Evaluate2D(double* pos, double* par);
144 double EvaluateProfile(double* pos, double* par);
145
146 private:
147 enum PlotType {
148 Potential = 0,
149 Magnitude,
150 Ex,
151 Ey,
152 Ez,
153 Unknown
154 };
155
156 std::string m_className;
157
158 static const unsigned int m_nMaxContours = 50;
159
160 // Options
161 bool m_debug;
162
163 bool m_useAutoRange;
164 bool m_useStatus;
165 double m_vBkg;
166
167 // Sensor
168 Sensor* m_sensor;
169 ComponentBase* m_component;
170
171 // Projection for viewing
172 double m_project[3][3];
173 double m_plane[4];
174 char m_xLabel[50], m_yLabel[50], m_description[50];
175
176 // Plot area
177 bool m_hasUserArea;
178 double m_xmin, m_ymin;
179 double m_xmax, m_ymax;
180
181 // Function range
182 double m_vmin, m_vmax;
183 double m_emin, m_emax;
184 double m_wmin, m_wmax;
185
186 // Number of contours
187 unsigned int m_nContours;
188 // Number of points used to draw the functions
189 unsigned int m_nSamples1d;
190 unsigned int m_nSamples2dX, m_nSamples2dY;
191 // Weighting field label
192 std::string m_electrode;
193
194 // Canvas
195 TCanvas* m_canvas;
196 bool m_hasExternalCanvas;
197
198 // Potential function
199 TF2* m_f2d;
200 TF2* m_f2dW;
201 TF1* m_fProfile;
202 TF1* m_fProfileW;
203
204 void Labels();
205 void CreateFunction();
206 bool SetupFunction(const std::string& option, TF2*& f, const bool contour,
207 const bool wfield = false);
208 bool SetupProfile(const double x0, const double y0, const double z0,
209 const double x1, const double y1, const double z1,
210 const std::string& option, TF1*& f, const bool wfield);
211 void SetupCanvas();
212 PlotType GetPlotType(const std::string& option, std::string& title) const;
213 std::string FindUnusedFunctionName(const std::string& s);
214};
215}
216#endif
Abstract base class for components.
Visualize the potential or electric field of a component or sensor.
Definition: ViewField.hh:15
~ViewField()
Destructor.
Definition: ViewField.cc:87
ViewField()
Constructor.
Definition: ViewField.cc:52
void SetPlane(const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
Definition: ViewField.cc:694
void SetElectricFieldRange(const double emin, const double emax)
Set the plot limits for the electric field.
Definition: ViewField.cc:153
double Evaluate2D(double *pos, double *par)
Definition: ViewField.cc:314
void PlotWeightingField(const std::string &label, const std::string &option, const std::string &drawopt)
Definition: ViewField.cc:242
void PlotProfileWeightingField(const std::string &label, const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const std::string &option="v")
Definition: ViewField.cc:265
friend class TF2
Definition: ViewField.hh:139
void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny)
Set the number of points used for drawing 2D functions.
Definition: ViewField.cc:191
void Rotate(const double angle)
Rotate the viewing plane (angle in radian).
Definition: ViewField.cc:736
void EnableAcknowledgeStatus(const double v0=0.)
Definition: ViewField.hh:128
friend class TF1
Definition: ViewField.hh:138
void SetArea()
Set the viewing area based on the bounding box of the sensor/component.
Definition: ViewField.hh:41
void DisableAcknowledgeStatus()
Ignore the status flag returned by the sensor/component.
Definition: ViewField.hh:133
void SetComponent(ComponentBase *c)
Set the component from which to retrieve the field.
Definition: ViewField.cc:107
void SetNumberOfSamples1d(const unsigned int n)
Set the number of points used for drawing 1D functions.
Definition: ViewField.cc:177
void SetVoltageRange(const double vmin, const double vmax)
Set the plot limits for the potential.
Definition: ViewField.cc:146
void PlotContour(const std::string &option="v")
Definition: ViewField.cc:213
void PlotSurface(const std::string &option="v")
Definition: ViewField.hh:73
void Plot(const std::string &option="v", const std::string &drawopt="arr")
Definition: ViewField.cc:222
void SetNumberOfContours(const unsigned int n)
Set the number of contour levels (at most 50).
Definition: ViewField.cc:167
void EnableAutoRange(const bool on=true)
Definition: ViewField.hh:123
void EnableDebugging(const bool on=true)
Switch on/off debugging output.
Definition: ViewField.hh:136
void SetWeightingFieldRange(const double wmin, const double wmax)
Set the plot limits for the weighting field.
Definition: ViewField.cc:160
void PlotSurfaceWeightingField(const std::string &label, const std::string &option)
Definition: ViewField.hh:99
void SetSensor(Sensor *s)
Set the sensor from which to retrieve the field.
Definition: ViewField.cc:96
void SetDefaultProjection()
Set the default viewing plane ( - at ).
Definition: ViewField.cc:291
double EvaluateProfile(double *pos, double *par)
Definition: ViewField.cc:387
void PlotContourWeightingField(const std::string &label, const std::string &option)
Definition: ViewField.cc:254
void PlotProfile(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const std::string &option="v")
Definition: ViewField.cc:231
void SetCanvas(TCanvas *c)
Set the canvas to be painted on.
Definition: ViewField.cc:118