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