Garfield++ 5.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 <Rtypes.h>
5#include <TH1D.h>
6
7#include <array>
8#include <memory>
9#include <string>
10
11#include "ViewBase.hh"
12
13namespace Garfield {
14
15class Sensor;
16
17/// Plot the signal computed by a sensor as a ROOT histogram.
18
19class ViewSignal : public ViewBase {
20 public:
21 /// Constructor
22 ViewSignal();
23 /// Destructor
24 ~ViewSignal() = default;
25
26 /// Set the sensor from which to retrieve the signal.
27 void SetSensor(Sensor* s);
28
29 /** Plot the signal.
30 * \param label Identifier (weighting field) of the signal to be plotted.
31 * \param optTotal String containing information about the total signals
32 * you want to plot. The syntax is the first letter of
33 * the charge carrier signal component you want to plot:
34 * "t" for total, "e" for electron and "i" for ion/hole.
35 * "tei" enables all three components.
36 * The total signal is always plotted.
37 * \param optPrompt String containing information about the
38 * prompt signal components you want to plot.
39 * The syntax is identical to the one described above.
40 * \param optDelayed String containing information about the delayed
41 * signal components you want to plot.
42 * The syntax is identical to the one described above.
43 * \param same Flag to keep existing plots on the canvas or not.
44 */
45
46 void PlotSignal(const std::string& label,
47 const std::string& optTotal = "t",
48 const std::string& optPrompt = "",
49 const std::string& optDelayed = "",
50 const bool same = false);
51
52 /** Retrieve the histogram for the total, prompt and delayed induced charge or
53 signal.
54 * \param h histogram to be returned
55 ('t': total, 'e': electron-induced, 'i': ion/hole-induced).
56 **/
57
58 TH1D* GetHistogram(const char h = 't') {
59 return h == 'e' ? m_hSignalElectrons.get()
60 : h == 'i' ? m_hSignalIons.get() : m_hSignal.get();
61 }
62
63 /// Set the x-axis limits explicitly.
64 void SetRangeX(const double xmin, const double xmax);
65 /// Remove the user-defined x-axis limits.
66 void UnsetRangeX() { m_userRangeX = false; }
67
68 /// Set the y-axis limits explicitly.
69 void SetRangeY(const double ymin, const double ymax);
70 /// Remove the user-defined y-axis limits.
71 void UnsetRangeY() { m_userRangeY = false; }
72
73 /// Override the default y-axis label.
74 void SetLabelY(const std::string& label) { m_labelY = label; }
75
76 /// Draw a legend on the plot or not.
77 void EnableLegend(const bool on = true) { m_legend = on; }
78
79 /// Set the (ROOT) colour with which to draw the total signal.
80 void SetColourTotal(const short col) { m_colTotal = col; }
81 /// Set the (ROOT) colour with which to draw the electron component.
82 void SetColourElectrons(const short col) { m_colElectrons = col; }
83 /// Set the (ROOT) colour with which to draw the hole/ion component.
84 void SetColourIons(const short col) { m_colIons = col; }
85 /// Set the (ROOT) colour with which to draw the hole/ion component.
86 void SetColourHoles(const short col) { m_colIons = col; }
87
88 /// Set the (ROOT) colours with which to draw the delayed signal(s).
89 void SetColourDelayed(const short colTotal,
90 const short colElectrons = kYellow - 7,
91 const short colIons = kRed - 9) {
92 m_colDelayed = {colTotal, colElectrons, colIons};
93 }
94
95 private:
96 // Sensor.
97 Sensor* m_sensor = nullptr;
98
99 // Axis range.
100 double m_xmin = 0.;
101 double m_xmax = 0.;
102 bool m_userRangeX = false;
103 double m_ymin = 0.;
104 double m_ymax = 0.;
105 bool m_userRangeY = false;
106
107 // Axis label.
108 std::string m_labelY = "";
109
110 // Histograms.
111 std::unique_ptr<TH1D> m_hSignal;
112 std::unique_ptr<TH1D> m_hSignalElectrons;
113 std::unique_ptr<TH1D> m_hSignalIons;
114 std::unique_ptr<TH1D> m_hPromptSignal;
115 std::unique_ptr<TH1D> m_hPromptElectrons;
116 std::unique_ptr<TH1D> m_hPromptIons;
117 std::unique_ptr<TH1D> m_hDelayedSignal;
118 std::unique_ptr<TH1D> m_hDelayedElectrons;
119 std::unique_ptr<TH1D> m_hDelayedIons;
120
121 bool m_legend = false;
122
123 // Colours.
124 short m_colTotal = kBlue + 3;
125 short m_colElectrons = kOrange - 3;
126 short m_colIons = kRed + 1;
127 std::array<short, 6> m_colDelayed{
128 {kCyan + 2, kYellow - 7, kRed - 9, kGreen + 1, kYellow - 4, kRed - 9}};
129
130 std::array<short, 3> m_colPrompt{{kAzure + 10, kRed - 4, kMagenta + 2}};
131
132 void DrawHistogram(TH1D* h, const std::string& opt, const std::string& xlabel,
133 const std::string& ylabel);
134
135};
136} // namespace Garfield
137#endif
ViewBase()=delete
Default constructor.
void UnsetRangeY()
Remove the user-defined y-axis limits.
Definition ViewSignal.hh:71
ViewSignal()
Constructor.
Definition ViewSignal.cc:16
TH1D * GetHistogram(const char h='t')
Definition ViewSignal.hh:58
void SetColourTotal(const short col)
Set the (ROOT) colour with which to draw the total signal.
Definition ViewSignal.hh:80
~ViewSignal()=default
Destructor.
void SetColourDelayed(const short colTotal, const short colElectrons=kYellow - 7, const short colIons=kRed - 9)
Set the (ROOT) colours with which to draw the delayed signal(s).
Definition ViewSignal.hh:89
void EnableLegend(const bool on=true)
Draw a legend on the plot or not.
Definition ViewSignal.hh:77
void SetColourIons(const short col)
Set the (ROOT) colour with which to draw the hole/ion component.
Definition ViewSignal.hh:84
void SetLabelY(const std::string &label)
Override the default y-axis label.
Definition ViewSignal.hh:74
void SetSensor(Sensor *s)
Set the sensor from which to retrieve the signal.
Definition ViewSignal.cc:18
void SetColourHoles(const short col)
Set the (ROOT) colour with which to draw the hole/ion component.
Definition ViewSignal.hh:86
void UnsetRangeX()
Remove the user-defined x-axis limits.
Definition ViewSignal.hh:66
void SetRangeY(const double ymin, const double ymax)
Set the y-axis limits explicitly.
Definition ViewSignal.cc:36
void SetRangeX(const double xmin, const double xmax)
Set the x-axis limits explicitly.
Definition ViewSignal.cc:26
void SetColourElectrons(const short col)
Set the (ROOT) colour with which to draw the electron component.
Definition ViewSignal.hh:82
void PlotSignal(const std::string &label, const std::string &optTotal="t", const std::string &optPrompt="", const std::string &optDelayed="", const bool same=false)
Definition ViewSignal.cc:61