Garfield++ 3.0
A toolkit for the detailed simulation of particle detectors based on ionisation measurement in gases and semiconductors
Loading...
Searching...
No Matches
ViewSignal.hh
Go to the documentation of this file.
1#ifndef G_VIEW_SIGNAL
2#define G_VIEW_SIGNAL
3
4#include <memory>
5#include <string>
6
7#include <TGraph.h>
8#include <TH1D.h>
9
10#include "ViewBase.hh"
11
12namespace Garfield {
13
14class Sensor;
15
16/// Plot the signal computed by a sensor as a ROOT histogram.
17
18class ViewSignal : public ViewBase {
19 public:
20 /// Constructor
21 ViewSignal();
22 /// Destructor
23 ~ViewSignal() = default;
24
25 /// Set the sensor from which to retrieve the signal.
26 void SetSensor(Sensor* s);
27
28 /** Plot the signal.
29 * \param label Identifier (weighting field) of the signal to be plotted.
30 * \param total Flag whether to plot the total induced signal.
31 * \param electron Flag whether to plot the electron-induced signal.
32 * \param ion Flag whether to plot the ion/hole-induced signal.
33 * \param delayed Flag whether to plot the delayed signal component.
34 */
35 void PlotSignal(const std::string& label, const bool total = true,
36 const bool electron = false, const bool ion = false,
37 const bool delayed = false);
38
39 /** Retrieve the histogram for the induced signal.
40 * \param h histogram to be returned
41 ('t': total, 'e': electron-induced, 'i': ion/hole-induced).
42 **/
43 TH1D* GetHistogram(const char h = 't') {
44 return h == 'e' ? m_hSignalElectrons.get()
45 : h == 'i' ? m_hSignalIons.get() : m_hSignal.get();
46 }
47
48 /// Set the x-axis limits explicitly.
49 void SetRangeX(const double xmin, const double xmax);
50 /// Remove the user-defined x-axis limits.
51 void UnsetRangeX() { m_userRangeX = false; }
52
53 /// Set the y-axis limits explicitly.
54 void SetRangeY(const double ymin, const double ymax);
55 /// Remove the user-defined y-axis limits.
56 void UnsetRangeY() { m_userRangeY = false; }
57
58 /// Override the default y-axis label.
59 void SetLabelY(const std::string& label) { m_labelY = label; }
60
61 private:
62 // Sensor.
63 Sensor* m_sensor = nullptr;
64
65 // Axis range.
66 double m_xmin = 0.;
67 double m_xmax = 0.;
68 bool m_userRangeX = false;
69 double m_ymin = 0.;
70 double m_ymax = 0.;
71 bool m_userRangeY = false;
72
73 // Axis label.
74 std::string m_labelY = "";
75
76 // Histograms.
77 std::unique_ptr<TH1D> m_hSignal;
78 std::unique_ptr<TH1D> m_hSignalElectrons;
79 std::unique_ptr<TH1D> m_hSignalIons;
80 std::unique_ptr<TH1D> m_hDelayedSignal;
81 std::unique_ptr<TH1D> m_hDelayedSignalElectrons;
82 std::unique_ptr<TH1D> m_hDelayedSignalIons;
83
84 // Threshold crossings.
85 std::unique_ptr<TGraph> m_gCrossings;
86};
87}
88#endif
Base class for visualization classes.
Definition: ViewBase.hh:10
Plot the signal computed by a sensor as a ROOT histogram.
Definition: ViewSignal.hh:18
void UnsetRangeY()
Remove the user-defined y-axis limits.
Definition: ViewSignal.hh:56
ViewSignal()
Constructor.
Definition: ViewSignal.cc:12
TH1D * GetHistogram(const char h='t')
Definition: ViewSignal.hh:43
~ViewSignal()=default
Destructor.
void SetLabelY(const std::string &label)
Override the default y-axis label.
Definition: ViewSignal.hh:59
void SetSensor(Sensor *s)
Set the sensor from which to retrieve the signal.
Definition: ViewSignal.cc:14
void PlotSignal(const std::string &label, const bool total=true, const bool electron=false, const bool ion=false, const bool delayed=false)
Definition: ViewSignal.cc:44
void UnsetRangeX()
Remove the user-defined x-axis limits.
Definition: ViewSignal.hh:51
void SetRangeY(const double ymin, const double ymax)
Set the y-axis limits explicitly.
Definition: ViewSignal.cc:33
void SetRangeX(const double xmin, const double xmax)
Set the x-axis limits explicitly.
Definition: ViewSignal.cc:22