Garfield++ 5.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 <Rtypes.h>
5
6#include "ViewBase.hh"
7
8namespace Garfield {
9
10class Sensor;
11class Component;
12
13/// Visualize the potential or electric field of a component or sensor.
14
15class ViewField : public ViewBase {
16 public:
17 /// Constructor.
18 ViewField();
19 /// Destructor.
20 ~ViewField() = default;
21
22 /// Set the sensor for which to plot the field.
23 void SetSensor(Sensor* s);
24 /// Set the component for which to plot the field.
25 void SetComponent(Component* c);
26
27 /// Set the plot limits for the potential.
28 void SetVoltageRange(const double vmin, const double vmax);
29 /// Set the plot limits for the electric field.
30 void SetElectricFieldRange(const double emin, const double emax);
31 /// Set the plot limits for the weighting field.
32 void SetWeightingFieldRange(const double wmin, const double wmax);
33 /// Set the plot limits for the magnetic field.
34 void SetMagneticFieldRange(const double bmin, const double bmax);
35
36 /// Set the number of contour levels.
37 void SetNumberOfContours(const unsigned int n);
38 /// Set the number of points used for drawing 1D functions.
39 void SetNumberOfSamples1d(const unsigned int n);
40 /// Set the number of points used for drawing 2D functions.
41 void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny);
42
43 /** Make a contour plot of the electric potential, electric field,
44 * or magnetic field.
45 * \param option quantity to be plotted
46 * - potential: "v", "voltage", "p", "potential"
47 * - magnitude of the electric field: "emag", "field"
48 * - x-component of the electric field: "ex"
49 * - y-component of the electric field: "ey"
50 * - z-component of the electric field: "ez"
51 * - magnitude of the magnetic field: "bmag"
52 * - x-component of the magnetic field: "bx"
53 * - y-component of the magnetic field: "by"
54 * - z-component of the magnetic field: "bz"
55 **/
56 void PlotContour(const std::string& option = "v");
57 /** Make a 2D plot of the electric potential, electric field or
58 * magnetic field.
59 * \param option quantity to be plotted (see PlotContour)
60 * \param drawopt option string passed to TF2::Draw
61 **/
62 void Plot(const std::string& option = "v",
63 const std::string& drawopt = "arr");
64 /** Make a 1D plot of the potential or field along a line.
65 * \param x0,y0,z0 starting point
66 * \param x1,y1,z1 end point
67 * \param option quantity to be plotted (see PlotContour)
68 * \param normalised flag whether to use normalised x-axis coordinates
69 **/
70 void PlotProfile(const double x0, const double y0, const double z0,
71 const double x1, const double y1, const double z1,
72 const std::string& option = "v",
73 const bool normalised = true);
74
75 /** Make a contour plot of the weighting potential or field.
76 * \param label identifier of the electrode
77 * \param option quantity to be plotted (see PlotContour)
78 **/
79 void PlotContourWeightingField(const std::string& label,
80 const std::string& option);
81 /** Make a 2D plot of the weighting potential or field.
82 * \param label identifier of the electrode
83 * \param option quantity to be plotted (see PlotContour)
84 * \param drawopt option string passed to TF2::Draw
85 *\param t time slice of dynamic weighting potential [ns].
86 **/
87 void PlotWeightingField(const std::string& label, const std::string& option,
88 const std::string& drawopt, const double t = 0.);
89
90 /** Make a 1D plot of the weighting potential or field along a line.
91 * \param label identifier of the electrode
92 * \param x0,y0,z0 starting point
93 * \param x1,y1,z1 end point
94 * \param option quantity to be plotted (see PlotContour)
95 * \param normalised flag whether to use normalised x-axis coordinates
96 **/
97 void PlotProfileWeightingField(const std::string& label, const double x0,
98 const double y0, const double z0,
99 const double x1, const double y1,
100 const double z1,
101 const std::string& option = "v",
102 const bool normalised = true);
103 /// Determine the range of the potential/field automatically (true)
104 /// or set it explicitly (false).
105 void EnableAutoRange(const bool on = true,
106 const bool samplePotential = true) {
107 m_useAutoRange = on;
108 m_samplePotential = samplePotential;
109 }
110
111 /** Make use (or not) of the status flag returned by the sensor/component.
112 * \param on Take status flag into account (true) or ignore it (false).
113 * \param v0 Value to be used for regions with status != 0.
114 */
115 void AcknowledgeStatus(const bool on, const double v0 = 0.) {
116 m_useStatus = on;
117 m_vBkg = v0;
118 }
119
120 /// Draw electric field lines from a set of starting points.
121 void PlotFieldLines(const std::vector<double>& x0,
122 const std::vector<double>& y0,
123 const std::vector<double>& z0,
124 const bool electron = true, const bool axis = true,
125 const short col = kOrange - 3);
126 /// Generates point along a line, spaced by equal flux intervals.
127 bool EqualFluxIntervals(const double x0, const double y0, const double z0,
128 const double x1, const double y1, const double z1,
129 std::vector<double>& xf, std::vector<double>& yf,
130 std::vector<double>& zf,
131 const unsigned int nPoints = 20) const;
132 /// Generate points along a line, spaced by a given flux interval.
133 bool FixedFluxIntervals(const double x0, const double y0, const double z0,
134 const double x1, const double y1, const double z1,
135 std::vector<double>& xf, std::vector<double>& yf,
136 std::vector<double>& zf,
137 const double interval = 10.) const;
138
139 private:
140 enum class Parameter {
141 Potential = 0,
142 Emag,
143 Ex,
144 Ey,
145 Ez,
146 Bmag,
147 Bx,
148 By,
149 Bz,
150 Unknown };
151
152 bool m_useAutoRange = true;
153 bool m_samplePotential = true;
154 bool m_useStatus = false;
155 double m_vBkg = 0.;
156
157 // Sensor
158 Sensor* m_sensor = nullptr;
159 Component* m_component = nullptr;
160
161 // Function range
162 double m_vmin = 0., m_vmax = 100.;
163 double m_emin = 0., m_emax = 10000.;
164 double m_wmin = 0., m_wmax = 100.;
165 double m_bmin = 0., m_bmax = 10.;
166
167 // Number of contours
168 unsigned int m_nContours = 20;
169 // Number of points used to draw the functions
170 unsigned int m_nSamples1d = 1000;
171 unsigned int m_nSamples2dX = 200;
172 unsigned int m_nSamples2dY = 200;
173
174 bool SetPlotLimits();
175 void Draw2d(const std::string& option, const bool contour,
176 const bool wfield, const std::string& electrode,
177 const std::string& drawopt, const double t = 0.);
178 void DrawProfile(const double x0, const double y0, const double z0,
179 const double x1, const double y1, const double z1,
180 const std::string& option, const bool wfield,
181 const std::string& electrode, const bool normalised);
182 Parameter GetPar(const std::string& option, std::string& title,
183 bool& bfield) const;
184 double Efield(const double x, const double y, const double z,
185 const Parameter par) const;
186 double Wfield(const double x, const double y, const double z,
187 const Parameter par, const std::string& electrode, const double t = 0.) const;
188 double Bfield(const double x, const double y, const double z,
189 const Parameter par) const;
190
191};
192}
193#endif
Abstract base class for components.
Definition Component.hh:13
ViewBase()=delete
Default constructor.
ViewField()
Constructor.
Definition ViewField.cc:71
void SetMagneticFieldRange(const double bmin, const double bmax)
Set the plot limits for the magnetic field.
Definition ViewField.cc:109
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", const bool normalised=true)
Definition ViewField.cc:161
void SetElectricFieldRange(const double emin, const double emax)
Set the plot limits for the electric field.
Definition ViewField.cc:97
void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny)
Set the number of points used for drawing 2D functions.
Definition ViewField.cc:123
bool FixedFluxIntervals(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, std::vector< double > &xf, std::vector< double > &yf, std::vector< double > &zf, const double interval=10.) const
Generate points along a line, spaced by a given flux interval.
Definition ViewField.cc:832
void SetComponent(Component *c)
Set the component for which to plot the field.
Definition ViewField.cc:82
void SetNumberOfSamples1d(const unsigned int n)
Set the number of points used for drawing 1D functions.
Definition ViewField.cc:119
void PlotFieldLines(const std::vector< double > &x0, const std::vector< double > &y0, const std::vector< double > &z0, const bool electron=true, const bool axis=true, const short col=kOrange - 3)
Draw electric field lines from a set of starting points.
Definition ViewField.cc:666
void PlotWeightingField(const std::string &label, const std::string &option, const std::string &drawopt, const double t=0.)
Definition ViewField.cc:150
void SetVoltageRange(const double vmin, const double vmax)
Set the plot limits for the potential.
Definition ViewField.cc:91
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", const bool normalised=true)
Definition ViewField.cc:144
bool EqualFluxIntervals(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, std::vector< double > &xf, std::vector< double > &yf, std::vector< double > &zf, const unsigned int nPoints=20) const
Generates point along a line, spaced by equal flux intervals.
Definition ViewField.cc:736
void PlotContour(const std::string &option="v")
Definition ViewField.cc:129
void Plot(const std::string &option="v", const std::string &drawopt="arr")
Definition ViewField.cc:133
void SetNumberOfContours(const unsigned int n)
Set the number of contour levels.
Definition ViewField.cc:115
void SetWeightingFieldRange(const double wmin, const double wmax)
Set the plot limits for the weighting field.
Definition ViewField.cc:103
void SetSensor(Sensor *s)
Set the sensor for which to plot the field.
Definition ViewField.cc:73
void EnableAutoRange(const bool on=true, const bool samplePotential=true)
Definition ViewField.hh:105
~ViewField()=default
Destructor.
void PlotContourWeightingField(const std::string &label, const std::string &option)
Definition ViewField.cc:156
void AcknowledgeStatus(const bool on, const double v0=0.)
Definition ViewField.hh:115