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