Garfield++ v2r0
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>

Public Member Functions

 ViewSignal ()
 Constructor.
 
 ~ViewSignal ()
 Destructor.
 
void SetSensor (Sensor *s)
 Set the sensor from which to retrieve the signal.
 
void SetCanvas (TCanvas *c)
 Set the pad on which to draw the histogram.
 
void PlotSignal (const std::string &label, const bool total=true, const bool electron=false, const bool ion=false)
 
TH1D * GetHistogram (const char h='t')
 
void EnableDebugging (const bool on=true)
 Enable/disable debugging output.
 

Detailed Description

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

Definition at line 16 of file ViewSignal.hh.

Constructor & Destructor Documentation

◆ ViewSignal()

Garfield::ViewSignal::ViewSignal ( )

Constructor.

Definition at line 12 of file ViewSignal.cc.

13 : m_className("ViewSignal"),
14 m_debug(false),
15 m_sensor(NULL),
16 m_canvas(NULL),
17 m_hasExternalCanvas(false),
18 m_hSignal(NULL), m_hSignalElectrons(NULL), m_hSignalIons(NULL),
19 m_gCrossings(NULL) {
20
22}
PlottingEngineRoot plottingEngine

◆ ~ViewSignal()

Garfield::ViewSignal::~ViewSignal ( )

Destructor.

Definition at line 24 of file ViewSignal.cc.

24 {
25
26 if (!m_hasExternalCanvas && m_canvas) delete m_canvas;
27 if (m_hSignal) delete m_hSignal;
28 if (m_hSignalElectrons) delete m_hSignalElectrons;
29 if (m_hSignalIons) delete m_hSignalIons;
30 if (m_gCrossings) delete m_gCrossings;
31}

Member Function Documentation

◆ EnableDebugging()

void Garfield::ViewSignal::EnableDebugging ( const bool  on = true)
inline

Enable/disable debugging output.

Definition at line 46 of file ViewSignal.hh.

46{ m_debug = on; }

◆ 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, 'h': ion-induced).

Definition at line 41 of file ViewSignal.hh.

41 {
42 return h == 'e' ? m_hSignalElectrons : 'i' ? m_hSignalIons : m_hSignal;
43 }

◆ PlotSignal()

void Garfield::ViewSignal::PlotSignal ( const std::string &  label,
const bool  total = true,
const bool  electron = false,
const bool  ion = 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.

Definition at line 53 of file ViewSignal.cc.

54 {
55
56 if (!m_sensor) {
57 std::cerr << m_className << "::PlotSignal: Sensor is not defined.\n";
58 return;
59 }
60
61 // Setup the canvas
62 if (!m_canvas) {
63 m_canvas = new TCanvas();
64 m_canvas->SetTitle("Signal");
65 if (m_hasExternalCanvas) m_hasExternalCanvas = false;
66 }
67 m_canvas->cd();
68
69 unsigned int nBins = 100;
70 double t0 = 0., dt = 1.;
71 m_sensor->GetTimeWindow(t0, dt, nBins);
72
73 if (total) {
74 if (m_hSignal) {
75 delete m_hSignal;
76 m_hSignal = NULL;
77 }
78 const std::string hname = FindHistogramName("hSignal_");
79 m_hSignal = new TH1D(hname.c_str(), label.c_str(), nBins, t0, t0 + nBins * dt);
80 m_hSignal->SetLineColor(plottingEngine.GetRootColorLine1());
81 m_hSignal->GetXaxis()->SetTitle("time [ns]");
82 m_hSignal->GetYaxis()->SetTitle("signal [fC / ns]");
83
84 for (unsigned int i = 0; i < nBins; ++i) {
85 const double sig = m_sensor->GetSignal(label, i);
86 m_hSignal->SetBinContent(i + 1, sig);
87 }
88
89 if (m_gCrossings) {
90 delete m_gCrossings;
91 m_gCrossings = NULL;
92 }
93
94 // Get threshold crossings.
95 const int nCrossings = m_sensor->GetNumberOfThresholdCrossings();
96 if (nCrossings > 0) {
97 m_gCrossings = new TGraph(nCrossings);
98 m_gCrossings->SetMarkerStyle(20);
99 m_gCrossings->SetMarkerColor(plottingEngine.GetRootColorLine1());
100 double time = 0., level = 0.;
101 bool rise = true;
102 for (int i = nCrossings; i--;) {
103 if (m_sensor->GetThresholdCrossing(i, time, level, rise)) {
104 m_gCrossings->SetPoint(i, time, level);
105 }
106 }
107 }
108
109 m_hSignal->Draw("");
110 if (nCrossings > 0) m_gCrossings->Draw("psame");
111 m_canvas->Update();
112 }
113
114 // Plot the electron and ion signals if requested.
115 if (electron) {
116 if (m_hSignalElectrons) {
117 delete m_hSignalElectrons;
118 m_hSignalElectrons = NULL;
119 }
120 const std::string hname = FindHistogramName("hSignalElectrons_");
121 m_hSignalElectrons = new TH1D(hname.c_str(),
122 label.c_str(), nBins, t0, t0 + nBins * dt);
123 m_hSignalElectrons->SetLineColor(plottingEngine.GetRootColorElectron());
124 m_hSignalElectrons->GetXaxis()->SetTitle("time [ns]");
125 m_hSignalElectrons->GetYaxis()->SetTitle("signal [fC / ns]");
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 m_canvas->Update();
132 }
133 if (ion) {
134 if (m_hSignalIons) {
135 delete m_hSignalIons;
136 m_hSignalIons = NULL;
137 }
138 const std::string hname = FindHistogramName("hSignalIons_");
139 m_hSignalIons = new TH1D(hname.c_str(),
140 label.c_str(), nBins, t0, t0 + nBins * dt);
141 m_hSignalIons->SetLineColor(plottingEngine.GetRootColorIon());
142 m_hSignalIons->GetXaxis()->SetTitle("time [ns]");
143 m_hSignalIons->GetYaxis()->SetTitle("signal [fC / ns]");
144 for (unsigned int i = 0; i < nBins; ++i) {
145 const double sig = m_sensor->GetIonSignal(label, i);
146 m_hSignalIons->SetBinContent(i + 1, sig);
147 }
148 m_hSignalIons->Draw("same");
149 m_canvas->Update();
150 }
151}
double GetElectronSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:563
double GetSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:601
bool GetThresholdCrossing(const unsigned int i, double &time, double &level, bool &rise) const
Definition: Sensor.cc:960
unsigned int GetNumberOfThresholdCrossings() const
Definition: Sensor.hh:109
void GetTimeWindow(double &tstart, double &tstep, unsigned int &nsteps)
Definition: Sensor.hh:89
double GetIonSignal(const std::string &label, const unsigned int bin)
Definition: Sensor.cc:583

◆ SetCanvas()

void Garfield::ViewSignal::SetCanvas ( TCanvas *  c)

Set the pad on which to draw the histogram.

Definition at line 42 of file ViewSignal.cc.

42 {
43
44 if (!c) return;
45 if (!m_hasExternalCanvas && m_canvas) {
46 delete m_canvas;
47 m_canvas = NULL;
48 }
49 m_canvas = c;
50 m_hasExternalCanvas = true;
51}

◆ SetSensor()

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

Set the sensor from which to retrieve the signal.

Definition at line 33 of file ViewSignal.cc.

33 {
34
35 if (!s) {
36 std::cerr << m_className << "::SetSensor: Null pointer.\n";
37 return;
38 }
39 m_sensor = s;
40}

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