Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
Garfield::ViewField Class Reference

Visualize the potential or electric field of a component or sensor. More...

#include <ViewField.hh>

+ Inheritance diagram for Garfield::ViewField:

Public Member Functions

 ViewField ()
 Constructor.
 
 ~ViewField ()
 Destructor.
 
void SetSensor (Sensor *s)
 Set the sensor from which to retrieve the field.
 
void SetComponent (ComponentBase *c)
 Set the component from which to retrieve the field.
 
void SetVoltageRange (const double vmin, const double vmax)
 Set the plot limits for the potential.
 
void SetElectricFieldRange (const double emin, const double emax)
 Set the plot limits for the electric field.
 
void SetWeightingFieldRange (const double wmin, const double wmax)
 Set the plot limits for the weighting field.
 
void SetArea (const double xmin, const double ymin, const double xmax, const double ymax)
 Set the viewing area (in local coordinates of the current viewing plane).
 
void SetArea ()
 Set the viewing area based on the bounding box of the sensor/component.
 
void SetPlane (const double fx, const double fy, const double fz, const double x0, const double y0, const double z0)
 
void SetDefaultProjection ()
 Set the default viewing plane ( $x$- $y$ at $z = 0$).
 
void Rotate (const double angle)
 Rotate the viewing plane (angle in radian).
 
void SetNumberOfContours (const unsigned int n)
 Set the number of contour levels (at most 50).
 
void SetNumberOfSamples1d (const unsigned int n)
 Set the number of points used for drawing 1D functions.
 
void SetNumberOfSamples2d (const unsigned int nx, const unsigned int ny)
 Set the number of points used for drawing 2D functions.
 
void PlotContour (const std::string &option="v")
 
void PlotSurface (const std::string &option="v")
 
void Plot (const std::string &option="v", const std::string &drawopt="arr")
 
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")
 
void PlotContourWeightingField (const std::string &label, const std::string &option)
 
void PlotSurfaceWeightingField (const std::string &label, const std::string &option)
 
void PlotWeightingField (const std::string &label, const std::string &option, const std::string &drawopt)
 
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")
 
void EnableAutoRange (const bool on=true)
 
void EnableAcknowledgeStatus (const double v0=0.)
 
void DisableAcknowledgeStatus ()
 Ignore the status flag returned by the sensor/component.
 
- Public Member Functions inherited from Garfield::ViewBase
 ViewBase ()=delete
 Default constructor.
 
 ViewBase (const std::string &name)
 Constructor.
 
virtual ~ViewBase ()
 Destructor.
 
void SetCanvas (TCanvas *c)
 Set the canvas to be painted on.
 
TCanvas * GetCanvas ()
 Retrieve the canvas.
 
void EnableDebugging (const bool on=true)
 Switch on/off debugging output.
 

Protected Member Functions

double Evaluate2D (double *pos, double *par)
 
double EvaluateProfile (double *pos, double *par)
 
- Protected Member Functions inherited from Garfield::ViewBase
std::string FindUnusedFunctionName (const std::string &s) const
 
std::string FindUnusedHistogramName (const std::string &s) const
 

Friends

class TF1
 
class TF2
 

Additional Inherited Members

- Protected Attributes inherited from Garfield::ViewBase
std::string m_className = "ViewBase"
 
bool m_debug = false
 
TCanvas * m_canvas = nullptr
 
bool m_hasExternalCanvas = false
 
double m_proj [3][3]
 

Detailed Description

Visualize the potential or electric field of a component or sensor.

Definition at line 16 of file ViewField.hh.

Constructor & Destructor Documentation

◆ ViewField()

Garfield::ViewField::ViewField ( )

Constructor.

Definition at line 47 of file ViewField.cc.

47 : ViewBase("ViewField") {
49}
ViewBase()=delete
Default constructor.
void SetDefaultProjection()
Set the default viewing plane ( - at ).
Definition: ViewField.cc:214

◆ ~ViewField()

Garfield::ViewField::~ViewField ( )

Destructor.

Definition at line 51 of file ViewField.cc.

51 {
52 if (m_f2d) delete m_f2d;
53 if (m_f2dW) delete m_f2dW;
54 if (m_fProfile) delete m_fProfile;
55 if (m_fProfileW) delete m_fProfileW;
56}

Member Function Documentation

◆ DisableAcknowledgeStatus()

void Garfield::ViewField::DisableAcknowledgeStatus ( )
inline

Ignore the status flag returned by the sensor/component.

Definition at line 131 of file ViewField.hh.

131{ m_useStatus = false; }

◆ EnableAcknowledgeStatus()

void Garfield::ViewField::EnableAcknowledgeStatus ( const double  v0 = 0.)
inline

Make use of the status flag returned by the sensor/component.

Parameters
v0Value to be used for regions with status != 0.

Definition at line 126 of file ViewField.hh.

126 {
127 m_useStatus = true;
128 m_vBkg = v0;
129 }

◆ EnableAutoRange()

void Garfield::ViewField::EnableAutoRange ( const bool  on = true)
inline

Definition at line 121 of file ViewField.hh.

121{ m_useAutoRange = on; }

◆ Evaluate2D()

double Garfield::ViewField::Evaluate2D ( double *  pos,
double *  par 
)
protected

Definition at line 236 of file ViewField.cc.

236 {
237 if (!m_sensor && !m_component) return 0.;
238
239 // Transform to global coordinates.
240 const double u = pos[0];
241 const double v = pos[1];
242 const double x = m_project[0][0] * u + m_project[1][0] * v + m_project[2][0];
243 const double y = m_project[0][1] * u + m_project[1][1] * v + m_project[2][1];
244 const double z = m_project[0][2] * u + m_project[1][2] * v + m_project[2][2];
245
246 // Determine which quantity to plot.
247 const PlotType plotType = static_cast<PlotType>(int(fabs(par[0]) / 10.));
248
249 // Compute the field.
250 double ex = 0., ey = 0., ez = 0., volt = 0.;
251 int status = 0;
252 if (par[0] > 0.) {
253 // "Drift" electric field.
254 Medium* medium = nullptr;
255 if (!m_sensor) {
256 m_component->ElectricField(x, y, z, ex, ey, ez, volt, medium, status);
257 } else {
258 m_sensor->ElectricField(x, y, z, ex, ey, ez, volt, medium, status);
259 }
260 } else {
261 // Weighting field.
262 if (!m_sensor) {
263 if (plotType == Potential) {
264 volt = m_component->WeightingPotential(x, y, z, m_electrode);
265 } else {
266 m_component->WeightingField(x, y, z, ex, ey, ez, m_electrode);
267 }
268 } else {
269 if (plotType == Potential) {
270 volt = m_sensor->WeightingPotential(x, y, z, m_electrode);
271 } else {
272 m_sensor->WeightingField(x, y, z, ex, ey, ez, m_electrode);
273 }
274 }
275 }
276
277 if (m_debug) {
278 std::cout << m_className << "::Evaluate2D:\n"
279 << " At (u, v) = (" << u << ", " << v << "), "
280 << " (x,y,z) = (" << x << "," << y << "," << z << ")\n"
281 << " E = " << ex << ", " << ey << ", " << ez
282 << "), V = " << volt << ", status = " << status << "\n";
283 }
284 if (m_useStatus && status != 0) return m_vBkg;
285
286 switch (plotType) {
287 case Potential:
288 return volt;
289 break;
290 case Magnitude:
291 return sqrt(ex * ex + ey * ey + ez * ez);
292 break;
293 case Ex:
294 return ex;
295 break;
296 case Ey:
297 return ey;
298 break;
299 case Ez:
300 return ez;
301 break;
302 default:
303 break;
304 }
305 return volt;
306}
virtual void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, Medium *&m, int &status)=0
virtual double WeightingPotential(const double x, const double y, const double z, const std::string &label)
virtual void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
void WeightingField(const double x, const double y, const double z, double &wx, double &wy, double &wz, const std::string &label)
Get the weighting field at (x, y, z).
Definition: Sensor.cc:138
void ElectricField(const double x, const double y, const double z, double &ex, double &ey, double &ez, double &v, Medium *&medium, int &status)
Get the drift field and potential at (x, y, z).
Definition: Sensor.cc:75
double WeightingPotential(const double x, const double y, const double z, const std::string &label)
Get the weighting potential at (x, y, z).
Definition: Sensor.cc:154
std::string m_className
Definition: ViewBase.hh:28
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615
DoubleAc sqrt(const DoubleAc &f)
Definition: DoubleAc.cpp:314

◆ EvaluateProfile()

double Garfield::ViewField::EvaluateProfile ( double *  pos,
double *  par 
)
protected

Definition at line 308 of file ViewField.cc.

308 {
309 if (!m_sensor && !m_component) return 0.;
310
311 // Get the start and end position.
312 const double x0 = par[0];
313 const double y0 = par[1];
314 const double z0 = par[2];
315 const double x1 = par[3];
316 const double y1 = par[4];
317 const double z1 = par[5];
318 // Compute the direction.
319 const double dx = x1 - x0;
320 const double dy = y1 - y0;
321 const double dz = z1 - z0;
322 // Get the position.
323 const double t = pos[0];
324 const double x = x0 + t * dx;
325 const double y = y0 + t * dy;
326 const double z = z0 + t * dz;
327 // Determine which quantity to plot.
328 const PlotType plotType = static_cast<PlotType>(int(fabs(par[6]) / 10.));
329
330 // Compute the field.
331 double ex = 0., ey = 0., ez = 0., volt = 0.;
332 int status = 0;
333 if (par[6] > 0.) {
334 Medium* medium = nullptr;
335 // "Drift" electric field.
336 if (!m_sensor) {
337 m_component->ElectricField(x, y, z, ex, ey, ez, volt, medium, status);
338 } else {
339 m_sensor->ElectricField(x, y, z, ex, ey, ez, volt, medium, status);
340 }
341 } else {
342 // Weighting field.
343 if (!m_sensor) {
344 if (plotType == Potential) {
345 volt = m_component->WeightingPotential(x, y, z, m_electrode);
346 } else {
347 m_component->WeightingField(x, y, z, ex, ey, ez, m_electrode);
348 }
349 } else {
350 if (plotType == Potential) {
351 volt = m_sensor->WeightingPotential(x, y, z, m_electrode);
352 } else {
353 m_sensor->WeightingField(x, y, z, ex, ey, ez, m_electrode);
354 }
355 }
356 }
357 if (m_useStatus && status != 0) volt = m_vBkg;
358
359 switch (plotType) {
360 case Potential:
361 return volt;
362 break;
363 case Magnitude:
364 return sqrt(ex * ex + ey * ey + ez * ez);
365 break;
366 case Ex:
367 return ex;
368 break;
369 case Ey:
370 return ey;
371 break;
372 case Ez:
373 return ez;
374 break;
375 default:
376 break;
377 }
378 return volt;
379}

◆ Plot()

void Garfield::ViewField::Plot ( const std::string &  option = "v",
const std::string &  drawopt = "arr" 
)

Make a 2D plot of the electric potential or field.

Parameters
optionquantity to be plotted (see PlotContour)
drawoptoption string passed to TF2::Draw

Definition at line 163 of file ViewField.cc.

163 {
164 SetupCanvas();
165 if (!SetupFunction(option, m_f2d, false, false)) return;
166 m_f2d->Draw(drawopt.c_str());
167 gPad->SetRightMargin(0.15);
168 gPad->Update();
169}

Referenced by PlotSurface().

◆ PlotContour()

void Garfield::ViewField::PlotContour ( const std::string &  option = "v")

Make a contour plot of the electric potential or field.

Parameters
optionquantity to be plotted
  • potential: "v", "voltage", "p", "potential"
  • magnitude of the electric field: "e", "field"
  • x-component of the electric field: "ex"
  • y-component of the electric field: "ey"
  • z-component of the electric field: "ez"

Definition at line 155 of file ViewField.cc.

155 {
156 SetupCanvas();
157 if (!SetupFunction(option, m_f2d, true, false)) return;
158 m_f2d->Draw("CONT4Z");
159 gPad->SetRightMargin(0.15);
160 gPad->Update();
161}

Referenced by main().

◆ PlotContourWeightingField()

void Garfield::ViewField::PlotContourWeightingField ( const std::string &  label,
const std::string &  option 
)

Make a contour plot of the weighting potential or field.

Parameters
labelidentifier of the electrode
optionquantity to be plotted (see PlotContour)

Definition at line 191 of file ViewField.cc.

192 {
193 m_electrode = label;
194 SetupCanvas();
195 if (!SetupFunction(option, m_f2dW, true, true)) return;
196 m_f2dW->Draw("CONT4Z");
197 gPad->SetRightMargin(0.15);
198 gPad->Update();
199}

◆ PlotProfile()

void Garfield::ViewField::PlotProfile ( const double  x0,
const double  y0,
const double  z0,
const double  x1,
const double  y1,
const double  z1,
const std::string &  option = "v" 
)

Make a 1D plot of the electric potential or field along a line.

Parameters
x0,y0,z0starting point
x1,y1,z1end point
optionquantity to be plotted (see PlotContour)

Definition at line 171 of file ViewField.cc.

173 {
174 SetupCanvas();
175 if (!SetupProfile(x0, y0, z0, x1, y1, z1, option, m_fProfile, false)) return;
176 m_fProfile->Draw();
177 gPad->Update();
178}

◆ PlotProfileWeightingField()

void Garfield::ViewField::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" 
)

Make a 1D plot of the weighting potential or field along a line.

Parameters
labelidentifier of the electrode
x0,y0,z0starting point
x1,y1,z1end point
optionquantity to be plotted (see PlotContour)

Definition at line 201 of file ViewField.cc.

205 {
206 m_electrode = label;
207 SetupCanvas();
208 if (!SetupProfile(x0, y0, z0, x1, y1, z1, option, m_fProfileW, true)) return;
209 m_fProfileW->Draw();
210 gPad->Update();
211}

◆ PlotSurface()

void Garfield::ViewField::PlotSurface ( const std::string &  option = "v")
inline

Make a surface plot ("SURF4") of the electric potential or field.

Parameters
optionquantity to be plotted (see PlotContour)

Definition at line 71 of file ViewField.hh.

71{ Plot(option, "SURF4"); }
void Plot(const std::string &option="v", const std::string &drawopt="arr")
Definition: ViewField.cc:163

◆ PlotSurfaceWeightingField()

void Garfield::ViewField::PlotSurfaceWeightingField ( const std::string &  label,
const std::string &  option 
)
inline

Make a surface plot ("SURF4") of the weighting potential or field.

Parameters
labelidentifier of the electrode
optionquantity to be plotted (see PlotContour)

Definition at line 97 of file ViewField.hh.

98 {
99 PlotWeightingField(label, option, "SURF4");
100 }
void PlotWeightingField(const std::string &label, const std::string &option, const std::string &drawopt)
Definition: ViewField.cc:180

◆ PlotWeightingField()

void Garfield::ViewField::PlotWeightingField ( const std::string &  label,
const std::string &  option,
const std::string &  drawopt 
)

Make a 2D plot of the weighting potential or field.

Parameters
labelidentifier of the electrode
optionquantity to be plotted (see PlotContour)
drawoptoption string passed to TF2::Draw

Definition at line 180 of file ViewField.cc.

182 {
183 m_electrode = label;
184 SetupCanvas();
185 if (!SetupFunction(option, m_f2dW, false, true)) return;
186 m_f2dW->Draw(drawopt.c_str());
187 gPad->SetRightMargin(0.15);
188 gPad->Update();
189}

Referenced by PlotSurfaceWeightingField().

◆ Rotate()

void Garfield::ViewField::Rotate ( const double  angle)

Rotate the viewing plane (angle in radian).

Definition at line 653 of file ViewField.cc.

653 {
654 // Rotate the axes
655 double auxu[3], auxv[3];
656 const double ctheta = cos(theta);
657 const double stheta = sin(theta);
658 for (int i = 0; i < 3; ++i) {
659 auxu[i] = ctheta * m_project[0][i] - stheta * m_project[1][i];
660 auxv[i] = stheta * m_project[0][i] + ctheta * m_project[1][i];
661 }
662 for (int i = 0; i < 3; ++i) {
663 m_project[0][i] = auxu[i];
664 m_project[1][i] = auxv[i];
665 }
666
667 // Make labels to be placed along the axes
668 Labels();
669}
DoubleAc cos(const DoubleAc &f)
Definition: DoubleAc.cpp:432
DoubleAc sin(const DoubleAc &f)
Definition: DoubleAc.cpp:384

◆ SetArea() [1/2]

void Garfield::ViewField::SetArea ( )
inline

Set the viewing area based on the bounding box of the sensor/component.

Definition at line 39 of file ViewField.hh.

39{ m_hasUserArea = false; }

◆ SetArea() [2/2]

void Garfield::ViewField::SetArea ( const double  xmin,
const double  ymin,
const double  xmax,
const double  ymax 
)

Set the viewing area (in local coordinates of the current viewing plane).

Definition at line 78 of file ViewField.cc.

79 {
80 // Check range, assign if non-null.
81 if (xmin == xmax || ymin == ymax) {
82 std::cerr << m_className << "::SetArea: Null area range is not permitted.\n"
83 << " " << xmin << " < x < " << xmax << "\n"
84 << " " << ymin << " < y < " << ymax << "\n";
85 return;
86 }
87 m_xmin = std::min(xmin, xmax);
88 m_ymin = std::min(ymin, ymax);
89 m_xmax = std::max(xmin, xmax);
90 m_ymax = std::max(ymin, ymax);
91 m_hasUserArea = true;
92}

Referenced by main().

◆ SetComponent()

void Garfield::ViewField::SetComponent ( ComponentBase c)

Set the component from which to retrieve the field.

Definition at line 68 of file ViewField.cc.

68 {
69 if (!c) {
70 std::cerr << m_className << "::SetComponent: Null pointer.\n";
71 return;
72 }
73
74 m_component = c;
75 m_sensor = nullptr;
76}

◆ SetDefaultProjection()

void Garfield::ViewField::SetDefaultProjection ( )

Set the default viewing plane ( $x$- $y$ at $z = 0$).

Definition at line 214 of file ViewField.cc.

214 {
215 // Default projection: x-y at z=0
216 m_project[0][0] = 1;
217 m_project[1][0] = 0;
218 m_project[2][0] = 0;
219 m_project[0][1] = 0;
220 m_project[1][1] = 1;
221 m_project[2][1] = 0;
222 m_project[0][2] = 0;
223 m_project[1][2] = 0;
224 m_project[2][2] = 0;
225
226 // Plane description
227 m_plane[0] = 0;
228 m_plane[1] = 0;
229 m_plane[2] = 1;
230 m_plane[3] = 0;
231
232 // Prepare axis labels.
233 Labels();
234}

Referenced by ViewField().

◆ SetElectricFieldRange()

void Garfield::ViewField::SetElectricFieldRange ( const double  emin,
const double  emax 
)

Set the plot limits for the electric field.

Definition at line 100 of file ViewField.cc.

100 {
101 m_emin = std::min(emin, emax);
102 m_emax = std::max(emin, emax);
103 m_useAutoRange = false;
104}

◆ SetNumberOfContours()

void Garfield::ViewField::SetNumberOfContours ( const unsigned int  n)

Set the number of contour levels (at most 50).

Definition at line 112 of file ViewField.cc.

112 {
113 if (n <= m_nMaxContours) {
114 m_nContours = n;
115 } else {
116 std::cerr << m_className << "::SetNumberOfContours:\n"
117 << " Max. number of contours is " << m_nMaxContours << ".\n";
118 }
119}

Referenced by main().

◆ SetNumberOfSamples1d()

void Garfield::ViewField::SetNumberOfSamples1d ( const unsigned int  n)

Set the number of points used for drawing 1D functions.

Definition at line 121 of file ViewField.cc.

121 {
122 constexpr unsigned int nmin = 10;
123 constexpr unsigned int nmax = 100000;
124 if (n < nmin || n > nmax) {
125 std::cerr << m_className << "::SetNumberOfSamples1d:\n"
126 << " Number of points (" << n << ") out of range.\n"
127 << " " << nmin << " <= n <= " << nmax << "\n";
128 return;
129 }
130
131 m_nSamples1d = n;
132}

◆ SetNumberOfSamples2d()

void Garfield::ViewField::SetNumberOfSamples2d ( const unsigned int  nx,
const unsigned int  ny 
)

Set the number of points used for drawing 2D functions.

Definition at line 134 of file ViewField.cc.

135 {
136 constexpr unsigned int nmin = 10;
137 constexpr unsigned int nmax = 10000;
138 if (nx < nmin || nx > nmax) {
139 std::cerr << m_className << "::SetNumberOfSamples2d:\n"
140 << " Number of x-points (" << nx << ") out of range.\n"
141 << " " << nmin << " <= nx <= " << nmax << "\n";
142 } else {
143 m_nSamples2dX = nx;
144 }
145
146 if (ny < nmin || ny > nmax) {
147 std::cerr << m_className << "::SetNumberOfSamples2d:\n"
148 << " Number of y-points (" << ny << ") out of range.\n"
149 << " " << nmin << " <= ny <= " << nmax << "\n";
150 } else {
151 m_nSamples2dY = ny;
152 }
153}

Referenced by main().

◆ SetPlane()

void Garfield::ViewField::SetPlane ( const double  fx,
const double  fy,
const double  fz,
const double  x0,
const double  y0,
const double  z0 
)

Set the projection (viewing plane).

Parameters
fx,fy,fznormal vector
x0,y0,z0in-plane point

Definition at line 612 of file ViewField.cc.

613 {
614 // Calculate two in-plane vectors for the normal vector
615 const double fnorm = sqrt(fx * fx + fy * fy + fz * fz);
616 if (fnorm > 0 && fx * fx + fz * fz > 0) {
617 const double fxz = sqrt(fx * fx + fz * fz);
618 m_project[0][0] = fz / fxz;
619 m_project[0][1] = 0;
620 m_project[0][2] = -fx / fxz;
621 m_project[1][0] = -fx * fy / (fxz * fnorm);
622 m_project[1][1] = (fx * fx + fz * fz) / (fxz * fnorm);
623 m_project[1][2] = -fy * fz / (fxz * fnorm);
624 m_project[2][0] = x0;
625 m_project[2][1] = y0;
626 m_project[2][2] = z0;
627 } else if (fnorm > 0 && fy * fy + fz * fz > 0) {
628 const double fyz = sqrt(fy * fy + fz * fz);
629 m_project[0][0] = (fy * fy + fz * fz) / (fyz * fnorm);
630 m_project[0][1] = -fx * fz / (fyz * fnorm);
631 m_project[0][2] = -fy * fz / (fyz * fnorm);
632 m_project[1][0] = 0;
633 m_project[1][1] = fz / fyz;
634 m_project[1][2] = -fy / fyz;
635 m_project[2][0] = x0;
636 m_project[2][1] = y0;
637 m_project[2][2] = z0;
638 } else {
639 std::cout << m_className << "::SetPlane:\n"
640 << " Normal vector has zero norm. No new projection set.\n";
641 }
642
643 // Store the plane description
644 m_plane[0] = fx;
645 m_plane[1] = fy;
646 m_plane[2] = fz;
647 m_plane[3] = fx * x0 + fy * y0 + fz * z0;
648
649 // Make labels to be placed along the axes
650 Labels();
651}

Referenced by main().

◆ SetSensor()

void Garfield::ViewField::SetSensor ( Sensor s)

Set the sensor from which to retrieve the field.

Definition at line 58 of file ViewField.cc.

58 {
59 if (!s) {
60 std::cerr << m_className << "::SetSensor: Null pointer.\n";
61 return;
62 }
63
64 m_sensor = s;
65 m_component = nullptr;
66}

Referenced by main().

◆ SetVoltageRange()

void Garfield::ViewField::SetVoltageRange ( const double  vmin,
const double  vmax 
)

Set the plot limits for the potential.

Definition at line 94 of file ViewField.cc.

94 {
95 m_vmin = std::min(vmin, vmax);
96 m_vmax = std::max(vmin, vmax);
97 m_useAutoRange = false;
98}

◆ SetWeightingFieldRange()

void Garfield::ViewField::SetWeightingFieldRange ( const double  wmin,
const double  wmax 
)

Set the plot limits for the weighting field.

Definition at line 106 of file ViewField.cc.

106 {
107 m_wmin = std::min(wmin, wmax);
108 m_wmax = std::max(wmin, wmax);
109 m_useAutoRange = false;
110}

Friends And Related Function Documentation

◆ TF1

friend class TF1
friend

Definition at line 133 of file ViewField.hh.

◆ TF2

friend class TF2
friend

Definition at line 134 of file ViewField.hh.


The documentation for this class was generated from the following files: