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::ViewSignal Class Reference

Plot the signal computed by a sensor as a ROOT histogram. More...

#include <ViewSignal.hh>

+ Inheritance diagram for Garfield::ViewSignal:

Public Member Functions

 ViewSignal ()
 Constructor.
 
 ~ViewSignal ()=default
 Destructor.
 
void SetSensor (Sensor *s)
 Set the sensor from which to retrieve the signal.
 
void PlotSignal (const std::string &label, const bool total=true, const bool electron=false, const bool ion=false, const bool delayed=false)
 
TH1D * GetHistogram (const char h='t')
 
void SetRangeX (const double xmin, const double xmax)
 Set the x-axis limits explicitly.
 
void UnsetRangeX ()
 Remove the user-defined x-axis limits.
 
void SetRangeY (const double ymin, const double ymax)
 Set the y-axis limits explicitly.
 
void UnsetRangeY ()
 Remove the user-defined y-axis limits.
 
void SetLabelY (const std::string &label)
 Override the default y-axis label.
 
- 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.
 

Additional Inherited Members

- Protected Member Functions inherited from Garfield::ViewBase
std::string FindUnusedFunctionName (const std::string &s) const
 
std::string FindUnusedHistogramName (const std::string &s) const
 
- 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

Plot the signal computed by a sensor as a ROOT histogram.

Definition at line 18 of file ViewSignal.hh.

Constructor & Destructor Documentation

◆ ViewSignal()

Garfield::ViewSignal::ViewSignal ( )

Constructor.

Definition at line 12 of file ViewSignal.cc.

12: ViewBase("ViewSignal") {}
ViewBase()=delete
Default constructor.

◆ ~ViewSignal()

Garfield::ViewSignal::~ViewSignal ( )
default

Destructor.

Member Function Documentation

◆ GetHistogram()

TH1D * Garfield::ViewSignal::GetHistogram ( const char  h = 't')
inline

Retrieve the histogram for the induced signal.

Parameters
hhistogram to be returned ('t': total, 'e': electron-induced, 'i': ion/hole-induced).

Definition at line 43 of file ViewSignal.hh.

43 {
44 return h == 'e' ? m_hSignalElectrons.get()
45 : h == 'i' ? m_hSignalIons.get() : m_hSignal.get();
46 }

◆ PlotSignal()

void Garfield::ViewSignal::PlotSignal ( const std::string &  label,
const bool  total = true,
const bool  electron = false,
const bool  ion = false,
const bool  delayed = false 
)

Plot the signal.

Parameters
labelIdentifier (weighting field) of the signal to be plotted.
totalFlag whether to plot the total induced signal.
electronFlag whether to plot the electron-induced signal.
ionFlag whether to plot the ion/hole-induced signal.
delayedFlag whether to plot the delayed signal component.

Definition at line 44 of file ViewSignal.cc.

46 {
47 if (!m_sensor) {
48 std::cerr << m_className << "::PlotSignal: Sensor is not defined.\n";
49 return;
50 }
51
52 // Setup the canvas
53 if (!m_canvas) {
54 m_canvas = new TCanvas();
55 m_canvas->SetTitle("Signal");
57 }
58 m_canvas->cd();
59
60 unsigned int nBins = 100;
61 double t0 = 0., dt = 1.;
62 m_sensor->GetTimeWindow(t0, dt, nBins);
63 const double t1 = t0 + nBins * dt;
64
65 const auto title = label.c_str();
66 std::string xlabel = "time [ns]";
67 std::string ylabel = m_labelY;
68 if (ylabel.empty()) {
69 ylabel = m_sensor->IsIntegrated() ? "signal [fC]" : "signal [fC / ns]";
70 }
71 if (total) {
72 auto hname = FindUnusedHistogramName("hSignal_").c_str();
73 m_hSignal.reset(new TH1D(hname, title, nBins, t0, t1));
74 m_hSignal->SetLineColor(plottingEngine.GetRootColorLine1());
75 m_hSignal->GetXaxis()->SetTitle(xlabel.c_str()),
76 m_hSignal->GetYaxis()->SetTitle(ylabel.c_str());
77 m_hSignal->SetStats(0);
78 for (unsigned int i = 0; i < nBins; ++i) {
79 const double sig = m_sensor->GetSignal(label, i);
80 m_hSignal->SetBinContent(i + 1, sig);
81 }
82 m_hSignal->Draw("");
83 if (m_userRangeX) m_hSignal->SetAxisRange(m_xmin, m_xmax, "X");
84 if (m_userRangeY) m_hSignal->SetAxisRange(m_ymin, m_ymax, "Y");
85
86 // Get and plot threshold crossings.
87 const auto nCrossings = m_sensor->GetNumberOfThresholdCrossings();
88 if (nCrossings > 0) {
89 m_gCrossings.reset(new TGraph(nCrossings));
90 m_gCrossings->SetMarkerStyle(20);
91 m_gCrossings->SetMarkerColor(plottingEngine.GetRootColorLine1());
92 double time = 0., level = 0.;
93 bool rise = true;
94 for (unsigned int i = 0; i < nCrossings; ++i) {
95 if (m_sensor->GetThresholdCrossing(i, time, level, rise)) {
96 m_gCrossings->SetPoint(i, time, level);
97 }
98 }
99 m_gCrossings->Draw("psame");
100 }
101
102 if (delayed) {
103 hname = FindUnusedHistogramName("hDelayedSignal_").c_str();
104 m_hDelayedSignal.reset(new TH1D(hname, title, nBins, t0, t1));
105 m_hDelayedSignal->SetLineColor(kCyan + 2);
106 m_hDelayedSignal->SetLineStyle(2);
107 m_hDelayedSignal->SetStats(0);
108 for (unsigned int i = 0; i < nBins; ++i) {
109 const double sig = m_sensor->GetDelayedElectronSignal(label, i) +
110 m_sensor->GetDelayedIonSignal(label, i);
111 m_hDelayedSignal->SetBinContent(i + 1, sig);
112 }
113 m_hDelayedSignal->Draw("same");
114 }
115 m_canvas->Update();
116 }
117
118 // Plot the electron and ion signals if requested.
119 if (electron) {
120 auto hname = FindUnusedHistogramName("hSignalElectrons_").c_str();
121 m_hSignalElectrons.reset(new TH1D(hname, title, nBins, t0, t1));
122 m_hSignalElectrons->SetLineColor(plottingEngine.GetRootColorElectron());
123 m_hSignalElectrons->GetXaxis()->SetTitle(xlabel.c_str());
124 m_hSignalElectrons->GetYaxis()->SetTitle(ylabel.c_str());
125 m_hSignalElectrons->SetStats(0);
126 for (unsigned int i = 0; i < nBins; ++i) {
127 const double sig = m_sensor->GetElectronSignal(label, i);
128 m_hSignalElectrons->SetBinContent(i + 1, sig);
129 }
130 m_hSignalElectrons->Draw("same");
131 if (!total) {
132 if (m_userRangeX) m_hSignalElectrons->SetAxisRange(m_xmin, m_xmax, "X");
133 if (m_userRangeY) m_hSignalElectrons->SetAxisRange(m_ymin, m_ymax, "Y");
134 }
135 if (delayed) {
136 hname = FindUnusedHistogramName("hDelayedSignalElectrons_").c_str();
137 m_hDelayedSignalElectrons.reset(new TH1D(hname, title, nBins, t0, t1));
138 m_hDelayedSignalElectrons->SetLineColor(kYellow - 7);
139 m_hDelayedSignalElectrons->SetLineStyle(2);
140 m_hDelayedSignalElectrons->SetStats(0);
141 for (unsigned int i = 0; i < nBins; ++i) {
142 const double sig = m_sensor->GetDelayedElectronSignal(label, i);
143 m_hDelayedSignalElectrons->SetBinContent(i + 1, sig);
144 }
145 m_hDelayedSignalElectrons->Draw("same");
146 }
147 m_canvas->Update();
148 }
149 if (ion) {
150 auto hname = FindUnusedHistogramName("hSignalIons_").c_str();
151 m_hSignalIons.reset(new TH1D(hname, title, nBins, t0, t1));
152 m_hSignalIons->SetLineColor(plottingEngine.GetRootColorIon());
153 m_hSignalIons->GetXaxis()->SetTitle(xlabel.c_str());
154 m_hSignalIons->GetYaxis()->SetTitle(ylabel.c_str());
155 m_hSignalIons->SetStats(0);
156 for (unsigned int i = 0; i < nBins; ++i) {
157 const double sig = m_sensor->GetIonSignal(label, i);
158 m_hSignalIons->SetBinContent(i + 1, sig);
159 }
160 m_hSignalIons->Draw("same");
161 if (!(total || electron)) {
162 if (m_userRangeX) m_hSignalIons->SetAxisRange(m_xmin, m_xmax, "X");
163 if (m_userRangeY) m_hSignalIons->SetAxisRange(m_ymin, m_ymax, "Y");
164 }
165 if (delayed) {
166 hname = FindUnusedHistogramName("hDelayedSignalIons_").c_str();
167 m_hDelayedSignalIons.reset(new TH1D(hname, title, nBins, t0, t1));
168 m_hDelayedSignalIons->SetLineColor(kRed - 9);
169 m_hDelayedSignalIons->SetLineStyle(2);
170 m_hDelayedSignalIons->SetStats(0);
171 for (unsigned int i = 0; i < nBins; ++i) {
172 const double sig = m_sensor->GetDelayedIonSignal(label, i);
173 m_hDelayedSignalIons->SetBinContent(i + 1, sig);
174 }
175 m_hDelayedSignalIons->Draw("same");
176 }
177 m_canvas->Update();
178 }
179}
bool IsIntegrated() const
Return whether the signal has been integrated/convoluted.
Definition: Sensor.hh:143
double GetElectronSignal(const std::string &label, const unsigned int bin)
Retrieve the electron signal for a given electrode and time bin.
Definition: Sensor.cc:733
double GetDelayedIonSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed ion/hole signal for a given electrode and time bin.
Definition: Sensor.cc:765
void GetTimeWindow(double &tstart, double &tstep, unsigned int &nsteps) const
Retrieve the time window and binning.
Definition: Sensor.hh:89
double GetSignal(const std::string &label, const unsigned int bin)
Retrieve the total signal for a given electrode and time bin.
Definition: Sensor.cc:788
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
Definition: Sensor.cc:1216
unsigned int GetNumberOfThresholdCrossings() const
Definition: Sensor.hh:174
double GetDelayedElectronSignal(const std::string &label, const unsigned int bin)
Retrieve the delayed electron signal for a given electrode and time bin.
Definition: Sensor.cc:754
double GetIonSignal(const std::string &label, const unsigned int bin)
Retrieve the ion or hole signal for a given electrode and time bin.
Definition: Sensor.cc:744
std::string m_className
Definition: ViewBase.hh:28
std::string FindUnusedHistogramName(const std::string &s) const
Definition: ViewBase.cc:40
bool m_hasExternalCanvas
Definition: ViewBase.hh:35
TCanvas * m_canvas
Definition: ViewBase.hh:34
PlottingEngineRoot plottingEngine

Referenced by main().

◆ SetLabelY()

void Garfield::ViewSignal::SetLabelY ( const std::string &  label)
inline

Override the default y-axis label.

Definition at line 59 of file ViewSignal.hh.

59{ m_labelY = label; }

◆ SetRangeX()

void Garfield::ViewSignal::SetRangeX ( const double  xmin,
const double  xmax 
)

Set the x-axis limits explicitly.

Definition at line 22 of file ViewSignal.cc.

22 {
23
24 if (fabs(xmax - xmin) < Small) {
25 std::cerr << m_className << "::SetRangeX: Invalid range.\n";
26 return;
27 }
28 m_xmin = std::min(xmin, xmax);
29 m_xmax = std::max(xmin, xmax);
30 m_userRangeX = true;
31}
DoubleAc fabs(const DoubleAc &f)
Definition: DoubleAc.h:615

◆ SetRangeY()

void Garfield::ViewSignal::SetRangeY ( const double  ymin,
const double  ymax 
)

Set the y-axis limits explicitly.

Definition at line 33 of file ViewSignal.cc.

33 {
34
35 if (fabs(ymax - ymin) < Small) {
36 std::cerr << m_className << "::SetRangeY: Invalid range.\n";
37 return;
38 }
39 m_ymin = std::min(ymin, ymax);
40 m_ymax = std::max(ymin, ymax);
41 m_userRangeY = true;
42}

◆ SetSensor()

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

Set the sensor from which to retrieve the signal.

Definition at line 14 of file ViewSignal.cc.

14 {
15 if (!s) {
16 std::cerr << m_className << "::SetSensor: Null pointer.\n";
17 return;
18 }
19 m_sensor = s;
20}

Referenced by main().

◆ UnsetRangeX()

void Garfield::ViewSignal::UnsetRangeX ( )
inline

Remove the user-defined x-axis limits.

Definition at line 51 of file ViewSignal.hh.

51{ m_userRangeX = false; }

◆ UnsetRangeY()

void Garfield::ViewSignal::UnsetRangeY ( )
inline

Remove the user-defined y-axis limits.

Definition at line 56 of file ViewSignal.hh.

56{ m_userRangeY = false; }

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