Garfield 0.3
Toolkit for the detailed simulation of particle detectors based on ionization 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
7
8namespace Garfield {
9
10class Sensor;
11class Component;
12
14
15class ViewField : public ViewBase {
16 public:
20 ViewField(Sensor* sensor);
22 ViewField(Component* component);
24 ~ViewField() = default;
25
27 void SetSensor(Sensor* s);
29 void SetComponent(Component* c);
30
32 void SetVoltageRange(const double vmin, const double vmax);
34 void SetElectricFieldRange(const double emin, const double emax);
36 void SetWeightingFieldRange(const double wmin, const double wmax);
38 void SetMagneticFieldRange(const double bmin, const double bmax);
39
41 void SetNumberOfContours(const unsigned int n);
43 void SetNumberOfSamples1d(const unsigned int n);
45 void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny);
46
60 void PlotContour(const std::string& option = "v");
66 void Plot(const std::string& option = "v",
67 const std::string& drawopt = "arr");
74 void PlotProfile(const double x0, const double y0, const double z0,
75 const double x1, const double y1, const double z1,
76 const std::string& option = "v",
77 const bool normalised = true);
78
83 void PlotContourWeightingField(const std::string& label,
84 const std::string& option);
91 void PlotWeightingField(const std::string& label, const std::string& option,
92 const std::string& drawopt, const double t = 0.);
93
101 void PlotProfileWeightingField(const std::string& label, const double x0,
102 const double y0, const double z0,
103 const double x1, const double y1,
104 const double z1,
105 const std::string& option = "v",
106 const bool normalised = true);
109 void EnableAutoRange(const bool on = true,
110 const bool samplePotential = true) {
111 m_useAutoRange = on;
112 m_samplePotential = samplePotential;
113 }
114
119 void AcknowledgeStatus(const bool on, const double v0 = 0.) {
120 m_useStatus = on;
121 m_vBkg = v0;
122 }
123
125 void PlotFieldLines(const std::vector<double>& x0,
126 const std::vector<double>& y0,
127 const std::vector<double>& z0, const bool electron = true,
128 const bool axis = true, const short col = kOrange - 3);
130 bool EqualFluxIntervals(const double x0, const double y0, const double z0,
131 const double x1, const double y1, const double z1,
132 std::vector<double>& xf, std::vector<double>& yf,
133 std::vector<double>& zf,
134 const unsigned int nPoints = 20) const;
136 bool FixedFluxIntervals(const double x0, const double y0, const double z0,
137 const double x1, const double y1, const double z1,
138 std::vector<double>& xf, std::vector<double>& yf,
139 std::vector<double>& zf,
140 const double interval = 10.) const;
141
142 private:
155
156 bool m_useAutoRange = true;
157 bool m_samplePotential = true;
158 bool m_useStatus = false;
159 double m_vBkg = 0.;
160
161 // Sensor
162 Sensor* m_sensor = nullptr;
163 Component* m_component = nullptr;
164
165 // Function range
166 double m_vmin = 0., m_vmax = 100.;
167 double m_emin = 0., m_emax = 10000.;
168 double m_wmin = 0., m_wmax = 100.;
169 double m_bmin = 0., m_bmax = 10.;
170
171 // Number of contours
172 unsigned int m_nContours = 20;
173 // Number of points used to draw the functions
174 unsigned int m_nSamples1d = 1000;
175 unsigned int m_nSamples2dX = 200;
176 unsigned int m_nSamples2dY = 200;
177
179 void Draw2d(const std::string& option, const bool contour, const bool wfield,
180 const std::string& electrode, const std::string& drawopt,
181 const double t = 0.);
182 void DrawProfile(const double x0, const double y0, const double z0,
183 const double x1, const double y1, const double z1,
184 const std::string& option, const bool wfield,
185 const std::string& electrode, const bool normalised);
186 Parameter GetPar(const std::string& option, std::string& title,
187 bool& bfield) const;
188 double Efield(const double x, const double y, const double z,
189 const Parameter par) const;
190 double Wfield(const double x, const double y, const double z,
191 const Parameter par, const std::string& electrode,
192 const double t = 0.) const;
193 double Bfield(const double x, const double y, const double z,
194 const Parameter par) const;
195};
196} // namespace Garfield
197#endif
ViewBase()=delete
Default constructor.
void Draw2d(const std::string &option, const bool contour, const bool wfield, const std::string &electrode, const std::string &drawopt, const double t=0.)
ViewField()
Default constructor.
void SetMagneticFieldRange(const double bmin, const double bmax)
Set the plot limits for the magnetic field.
unsigned int m_nContours
Definition ViewField.hh:172
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)
Make a 1D plot of the weighting potential or field along a line.
void SetElectricFieldRange(const double emin, const double emax)
Set the plot limits for the electric field.
void SetNumberOfSamples2d(const unsigned int nx, const unsigned int ny)
Set the number of points used for drawing 2D functions.
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.
void SetComponent(Component *c)
Set the component for which to plot the field.
ViewField(Component *component)
Constructor from component.
ViewField(Sensor *sensor)
Constructor from sensor.
void SetNumberOfSamples1d(const unsigned int n)
Set the number of points used for drawing 1D functions.
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.
unsigned int m_nSamples1d
Definition ViewField.hh:174
Parameter GetPar(const std::string &option, std::string &title, bool &bfield) const
void PlotWeightingField(const std::string &label, const std::string &option, const std::string &drawopt, const double t=0.)
Make a 2D plot of the weighting potential or field.
double Bfield(const double x, const double y, const double z, const Parameter par) const
void DrawProfile(const double x0, const double y0, const double z0, const double x1, const double y1, const double z1, const std::string &option, const bool wfield, const std::string &electrode, const bool normalised)
void SetVoltageRange(const double vmin, const double vmax)
Set the plot limits for the potential.
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)
Make a 1D plot of the potential or field along a line.
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.
void PlotContour(const std::string &option="v")
Make a contour plot of the electric potential, electric field, or magnetic field.
void Plot(const std::string &option="v", const std::string &drawopt="arr")
Make a 2D plot of the electric potential, electric field or magnetic field.
double Wfield(const double x, const double y, const double z, const Parameter par, const std::string &electrode, const double t=0.) const
unsigned int m_nSamples2dX
Definition ViewField.hh:175
void SetNumberOfContours(const unsigned int n)
Set the number of contour levels.
double Efield(const double x, const double y, const double z, const Parameter par) const
void SetWeightingFieldRange(const double wmin, const double wmax)
Set the plot limits for the weighting field.
unsigned int m_nSamples2dY
Definition ViewField.hh:176
void SetSensor(Sensor *s)
Set the sensor for which to plot the field.
void EnableAutoRange(const bool on=true, const bool samplePotential=true)
Determine the range of the potential/field automatically (true) or set it explicitly (false).
Definition ViewField.hh:109
~ViewField()=default
Destructor.
void PlotContourWeightingField(const std::string &label, const std::string &option)
Make a contour plot of the weighting potential or field.
void AcknowledgeStatus(const bool on, const double v0=0.)
Make use (or not) of the status flag returned by the sensor/component.
Definition ViewField.hh:119
Component * m_component
Definition ViewField.hh:163