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.cc
Go to the documentation of this file.
2
3#include <TAxis.h>
4#include <TGraph.h>
5
6#include <iostream>
7
10#include "Garfield/Sensor.hh"
11#include "TLegend.h"
12#include "TPaveLabel.h"
13
14namespace Garfield {
15
16ViewSignal::ViewSignal() : ViewBase("ViewSignal") {}
17
19 if (!s) {
20 std::cerr << m_className << "::SetSensor: Null pointer.\n";
21 return;
22 }
23 m_sensor = s;
24}
25
26void ViewSignal::SetRangeX(const double xmin, const double xmax) {
27 if (fabs(xmax - xmin) < Small) {
28 std::cerr << m_className << "::SetRangeX: Invalid range.\n";
29 return;
30 }
31 m_xmin = std::min(xmin, xmax);
32 m_xmax = std::max(xmin, xmax);
33 m_userRangeX = true;
34}
35
36void ViewSignal::SetRangeY(const double ymin, const double ymax) {
37 if (fabs(ymax - ymin) < Small) {
38 std::cerr << m_className << "::SetRangeY: Invalid range.\n";
39 return;
40 }
41 m_ymin = std::min(ymin, ymax);
42 m_ymax = std::max(ymin, ymax);
43 m_userRangeY = true;
44}
45
46void ViewSignal::DrawHistogram(TH1D* h, const std::string& opt,
47 const std::string& xlabel,
48 const std::string& ylabel) {
49 if (!h) return;
50 h->GetXaxis()->SetTitle(xlabel.c_str());
51 h->GetYaxis()->SetTitle(ylabel.c_str());
52 h->SetStats(0);
53 auto hCopy = h->DrawCopy(opt.c_str());
54 if (m_userRangeX) hCopy->SetAxisRange(m_xmin, m_xmax, "X");
55 if (m_userRangeY) {
56 hCopy->SetMinimum(m_ymin);
57 hCopy->SetMaximum(m_ymax);
58 }
59}
60
61void ViewSignal::PlotSignal(const std::string& label,
62 const std::string& settingTotal,
63 const std::string& settingPrompt,
64 const std::string& settingDelayed,
65 const bool same) {
66 const bool totalT = true;
67 const bool totalP =
68 settingPrompt.find("t") != std::string::npos ? true : false;
69 const bool totalD =
70 settingDelayed.find("t") != std::string::npos ? true : false;
71 const bool electronT =
72 settingTotal.find("e") != std::string::npos ? true : false;
73 const bool electronP =
74 settingPrompt.find("e") != std::string::npos ? true : false;
75 const bool electronD =
76 settingDelayed.find("e") != std::string::npos ? true : false;
77 const bool ionT = settingTotal.find("i") != std::string::npos ? true : false;
78 const bool ionP = settingPrompt.find("i") != std::string::npos ? true : false;
79 const bool ionD =
80 settingDelayed.find("i") != std::string::npos ? true : false;
81
82 constexpr int lineThickness = 5;
83 constexpr double tol = 1e-50;
84
85 if (!m_sensor) {
86 std::cerr << m_className << "::PlotSignal: Sensor is not defined.\n";
87 return;
88 }
89
90 auto canvas = GetCanvas();
91 canvas->cd();
92 canvas->SetTitle("Signal");
93
94 unsigned int nBins = 100;
95 double t0 = 0., dt = 1.;
96 m_sensor->GetTimeWindow(t0, dt, nBins);
97 const double t1 = t0 + nBins * dt;
98
99 std::string xlabel = "time [ns]";
100 std::string ylabel = m_labelY;
101 if (ylabel.empty()) {
102 ylabel = m_sensor->IsIntegrated(label) ? "signal [fC]" : "signal [fC / ns]";
103 }
104 unsigned int nPlots = same ? 1 : 0;
105 if (!RangeSet(gPad)) nPlots = 0;
106
107 if (totalT) {
108 const auto hname = FindUnusedHistogramName("hSignal_");
109 m_hSignal.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
110 m_hSignal->SetDirectory(nullptr);
111 m_hSignal->SetLineColor(m_colTotal);
112 for (unsigned int i = 0; i < nBins; ++i) {
113 const double sig = m_sensor->GetSignal(label, i, 0);
114 if (!std::isnan(sig) && std::abs(sig) < tol) {
115 m_hSignal->SetBinContent(i + 1, 0.);
116 } else {
117 m_hSignal->SetBinContent(i + 1, sig);
118 }
119 }
120
121 m_hSignal->SetLineWidth(lineThickness);
122
123 const std::string opt = nPlots > 0 ? "same" : "";
124 ++nPlots;
125
126 // Get and plot threshold crossings.
127 const auto nCrossings = m_sensor->GetNumberOfThresholdCrossings();
128 if (nCrossings > 0) {
129 TGraph gCrossings;
130 gCrossings.SetMarkerStyle(20);
131 gCrossings.SetMarkerColor(m_colTotal);
132 std::vector<double> xp;
133 std::vector<double> yp;
134 double time = 0., level = 0.;
135 bool rise = true;
136 for (unsigned int i = 0; i < nCrossings; ++i) {
137 if (m_sensor->GetThresholdCrossing(i, time, level, rise)) {
138 xp.push_back(time);
139 yp.push_back(level);
140 }
141 }
142 gCrossings.DrawGraph(xp.size(), xp.data(), yp.data(), "psame");
143 }
144 DrawHistogram(m_hSignal.get(), opt, xlabel, ylabel);
145 }
146
147 if (totalD) {
148 const auto hname = FindUnusedHistogramName("hDelayedSignal_");
149 m_hDelayedSignal.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
150 m_hDelayedSignal->SetDirectory(nullptr);
151 m_hDelayedSignal->SetLineColor(m_colDelayed[3]);
152 m_hDelayedSignal->SetLineStyle(7);
153 m_hDelayedSignal->SetStats(0);
154 for (unsigned int i = 0; i < nBins; ++i) {
155 const double sig = m_sensor->GetSignal(label, i, 2);
156 if (!std::isnan(sig) && std::abs(sig) < tol) {
157 m_hDelayedSignal->SetBinContent(i + 1, 0.);
158 } else {
159 m_hDelayedSignal->SetBinContent(i + 1, sig);
160 }
161 }
162 m_hDelayedSignal->SetLineWidth(lineThickness);
163 m_hDelayedSignal->DrawCopy("same");
164 }
165
166 if (totalP) {
167 const auto hname = FindUnusedHistogramName("hPromptSignal_");
168 m_hPromptSignal.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
169 m_hPromptSignal->SetDirectory(nullptr);
170 m_hPromptSignal->SetLineColor(m_colPrompt[0]);
171 m_hPromptSignal->SetLineStyle(2);
172 m_hPromptSignal->SetStats(0);
173 for (unsigned int i = 0; i < nBins; ++i) {
174 const double sig = m_sensor->GetSignal(label, i, 1);
175 if (!std::isnan(sig) && std::abs(sig) < tol) {
176 m_hPromptSignal->SetBinContent(i + 1, 0.);
177 } else {
178 m_hPromptSignal->SetBinContent(i + 1, sig);
179 }
180 }
181
182 m_hPromptSignal->SetLineWidth(lineThickness);
183 m_hPromptSignal->DrawCopy("same");
184 }
185
186 if (electronT) {
187 const auto hname = FindUnusedHistogramName("hSignalElectrons_");
188 m_hSignalElectrons.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
189 m_hSignalElectrons->SetDirectory(nullptr);
190 m_hSignalElectrons->SetLineColor(m_colElectrons);
191 for (unsigned int i = 0; i < nBins; ++i) {
192 const double sig = m_sensor->GetElectronSignal(label, i);
193 m_hSignalElectrons->SetBinContent(i + 1, sig);
194 }
195 m_hSignalElectrons->SetLineWidth(lineThickness);
196 m_hSignalElectrons->DrawCopy("same");
197 }
198
199 if (electronD) {
200 const auto hname = FindUnusedHistogramName("m_hDelayedElectrons_");
201 m_hDelayedElectrons.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
202 m_hDelayedElectrons->SetDirectory(nullptr);
203 m_hDelayedElectrons->SetLineColor(m_colDelayed[4]);
204 m_hDelayedElectrons->SetLineStyle(7);
205 m_hDelayedElectrons->SetStats(0);
206 for (unsigned int i = 0; i < nBins; ++i) {
207 const double sig = m_sensor->GetDelayedElectronSignal(label, i);
208 if (!std::isnan(sig) && std::abs(sig) < tol) {
209 m_hDelayedElectrons->SetBinContent(i + 1, 0.);
210 } else {
211 m_hDelayedElectrons->SetBinContent(i + 1, sig);
212 }
213 }
214 m_hDelayedElectrons->SetLineWidth(lineThickness);
215 m_hDelayedElectrons->DrawCopy("same");
216 }
217
218 if (electronP) {
219 const auto hname = FindUnusedHistogramName("m_hPromptElectrons_");
220 m_hPromptElectrons.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
221 m_hPromptElectrons->SetDirectory(nullptr);
222 m_hPromptElectrons->SetLineColor(m_colPrompt[1]);
223 m_hPromptElectrons->SetLineStyle(2);
224 m_hPromptElectrons->SetStats(0);
225 for (unsigned int i = 0; i < nBins; ++i) {
226 const double sig = m_sensor->GetElectronSignal(label, i) -
227 m_sensor->GetDelayedElectronSignal(label, i);
228 if (!std::isnan(sig) && std::abs(sig) < tol) {
229 m_hPromptElectrons->SetBinContent(i + 1, 0.);
230 } else {
231 m_hPromptElectrons->SetBinContent(i + 1, sig);
232 }
233 }
234
235 m_hPromptElectrons->SetLineWidth(lineThickness);
236 m_hPromptElectrons->DrawCopy("same");
237 }
238
239 if (ionT) {
240 const auto hname = FindUnusedHistogramName("m_hSignalIons_");
241 m_hSignalIons.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
242 m_hSignalIons->SetDirectory(nullptr);
243 m_hSignalIons->SetLineColor(m_colIons);
244 for (unsigned int i = 0; i < nBins; ++i) {
245 const double sig = m_sensor->GetIonSignal(label, i);
246 m_hSignalIons->SetBinContent(i + 1, sig);
247 }
248 m_hSignalIons->SetLineWidth(lineThickness);
249 m_hSignalIons->DrawCopy("same");
250 }
251
252 if (ionD) {
253 const auto hname = FindUnusedHistogramName("m_hDelayedIons_");
254 m_hDelayedIons.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
255 m_hDelayedIons->SetDirectory(nullptr);
256 m_hDelayedIons->SetLineColor(m_colDelayed[5]);
257 m_hDelayedIons->SetLineStyle(7);
258 m_hDelayedIons->SetStats(0);
259 for (unsigned int i = 0; i < nBins; ++i) {
260 const double sig = m_sensor->GetDelayedIonSignal(label, i);
261 if (!std::isnan(sig) && std::abs(sig) < tol) {
262 m_hDelayedIons->SetBinContent(i + 1, 0.);
263 } else {
264 m_hDelayedIons->SetBinContent(i + 1, sig);
265 }
266 }
267 m_hDelayedIons->SetLineWidth(lineThickness);
268 m_hDelayedIons->DrawCopy("same");
269 }
270
271 if (ionP) {
272 std::cerr << m_className << "::ionP.\n";
273 const auto hname = FindUnusedHistogramName("m_hPromptIons_");
274 m_hPromptIons.reset(new TH1D(hname.c_str(), "", nBins, t0, t1));
275 m_hPromptIons->SetDirectory(nullptr);
276 m_hPromptIons->SetLineColor(m_colPrompt[2]);
277 m_hPromptIons->SetLineStyle(2);
278 m_hPromptIons->SetStats(0);
279 for (unsigned int i = 0; i < nBins; ++i) {
280 const double sig = m_sensor->GetIonSignal(label, i) -
281 m_sensor->GetDelayedIonSignal(label, i);
282 if (!std::isnan(sig) && std::abs(sig) < tol) {
283 m_hPromptIons->SetBinContent(i + 1, 0.);
284 } else {
285 m_hPromptIons->SetBinContent(i + 1, sig);
286 }
287 }
288
289 m_hPromptIons->SetLineWidth(lineThickness);
290 m_hPromptIons->DrawCopy("same");
291 }
292
293 if (m_legend) {
294 TLegend* leg = new TLegend(0.7, 0.7, 0.9, 0.9);
295 leg->SetHeader("Induced current components");
296 if (totalT) leg->AddEntry(m_hSignal.get(), "Total induced signal", "l");
297 if (totalP) {
298 leg->AddEntry(m_hPromptSignal.get(), "Prompt induced signal", "l");
299 }
300 if (totalD) {
301 leg->AddEntry(m_hDelayedSignal.get(), "Delayed induced signal", "l");
302 }
303 if (electronT) {
304 leg->AddEntry(m_hSignalElectrons.get(), "Electron induced signal", "l");
305 }
306 if (electronP) {
307 leg->AddEntry(m_hPromptElectrons.get(), "Electron prompt induced signal",
308 "l");
309 }
310 if (electronD) {
311 leg->AddEntry(m_hDelayedElectrons.get(),
312 "Electron delayed induced signal", "l");
313 }
314 if (ionT) leg->AddEntry(m_hSignalIons.get(), "Ion induced signal", "l");
315 if (ionP) {
316 leg->AddEntry(m_hPromptIons.get(), "Ion/hole prompt induced signal", "l");
317 }
318 if (ionD) {
319 leg->AddEntry(m_hDelayedIons.get(), "Ion/hole delayed induced signal",
320 "l");
321 }
322 leg->Draw();
323 }
324 gPad->Update();
325}
326
327} // namespace Garfield
static std::string FindUnusedHistogramName(const std::string &s)
Find an unused histogram name.
Definition ViewBase.cc:327
std::string m_className
Definition ViewBase.hh:82
TPad * GetCanvas()
Retrieve the canvas.
Definition ViewBase.cc:74
static bool RangeSet(TVirtualPad *)
Definition ViewBase.cc:83
ViewBase()=delete
Default constructor.
ViewSignal()
Constructor.
Definition ViewSignal.cc:16
void SetSensor(Sensor *s)
Set the sensor from which to retrieve the signal.
Definition ViewSignal.cc:18
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 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