CbmRoot
Loading...
Searching...
No Matches
CbmMcbm2024CheckEventsDt.cxx
Go to the documentation of this file.
1/* Copyright (C) 2024 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Pierre-Alain Loizeau [committer] */
4
5
7
8#include "CbmDefs.h"
9#include "CbmEvent.h"
10
11#include <FairRootManager.h>
12#include <Logger.h>
13
14#include <TCanvas.h>
15#include <TFile.h>
16#include <TH1.h>
17#include <TH2.h>
18#include <TROOT.h>
19#include <TStopwatch.h>
20
21#include <algorithm>
22#include <cassert>
23#include <iomanip>
24#include <vector>
25
26
27using namespace std;
28
29
30// ----- Constructor -----------------------------------------------------
32// ---------------------------------------------------------------------------
33
34
35// ----- Destructor ------------------------------------------------------
37// ---------------------------------------------------------------------------
38
39
40// ----- Execution -------------------------------------------------------
42{
43 // --- First TS: find boundaries of plots
45 if (0 == fNumTs) {
46 TDirectory* oldDir = gDirectory;
47 TFile* oldFile = gFile;
48 gROOT->cd();
49
50 for (std::string sDet : fvDets) {
51 fHistDt[sDet] = new TH1I(Form("histDt%s", sDet.data()),
52 Form("Time difference to event seed, %s; dt [ns]", sDet.data()), 1000, -500., 500.);
53 fHistDtEvo[sDet] = new TH2I(Form("histDtEvo%s", sDet.data()),
54 Form("Time difference to event seed vs TS, %s; TS []; dt [ns]", sDet.data()), //
55 100, 0., 1000., 1000, -500., 500.);
56
57 fHistDtToBmon[sDet] =
58 new TH1I(Form("histDtToBmon%s", sDet.data()),
59 Form("Time difference to first/single BMon, %s; dt [ns]", sDet.data()), 1000, -500., 500.);
60 fHistDtToBmonEvo[sDet] =
61 new TH2I(Form("histDtToBmonEvo%s", sDet.data()),
62 Form("Time difference to first/single BMon vs TS, %s; TS []; dt [ns]", sDet.data()), //
63 100, 0., 1000., 1000, -500., 500.);
64
65 fHistMul[sDet] =
66 new TH1I(Form("histMul%s", sDet.data()), Form("Nb %s digis per event; N []", sDet.data()), 100, -0.5, 99.5);
67 fHistDtMul[sDet] =
68 new TH2I(Form("histDtMul%s", sDet.data()), Form("Dt Bmon vs Nb %s digis per event; dt [ns]; N []", sDet.data()),
69 1000, -500., 500., 100, -0.5, 99.5);
70 }
71
72
73 fCanvDt = new TCanvas("canvDt", "Time differences to event seed");
74 fCanvDt->Divide(3, 3);
75
76 uint32_t uPadIdx = 1;
77 for (std::string sDet : fvDets) {
78 fCanvDt->cd(uPadIdx);
79 gPad->SetLogy();
80 gPad->SetGridx();
81 gPad->SetGridy();
82 fHistDt[sDet]->Draw("hist");
83 uPadIdx++;
84 }
85
86 fCanvDtEvo = new TCanvas("canvDtEvo", "Time differences to event seed vs TS");
87 fCanvDtEvo->Divide(3, 3);
88
89 uPadIdx = 1;
90 for (std::string sDet : fvDets) {
91 fCanvDtEvo->cd(uPadIdx);
92 gPad->SetLogz();
93 gPad->SetGridx();
94 gPad->SetGridy();
95 fHistDtEvo[sDet]->Draw("colz");
96 uPadIdx++;
97 }
98
99
100 fCanvDtToBmon = new TCanvas("canvDtToBmon", "Time differences to event seed");
101 fCanvDtToBmon->Divide(3, 3);
102
103 uPadIdx = 1;
104 for (std::string sDet : fvDets) {
105 fCanvDtToBmon->cd(uPadIdx);
106 gPad->SetLogy();
107 gPad->SetGridx();
108 gPad->SetGridy();
109 fHistDtToBmon[sDet]->Draw("hist");
110 uPadIdx++;
111 }
112
113 fCanvDtToBmonEvo = new TCanvas("canvDtToBmonEvo", "Time differences to event seed vs TS");
114 fCanvDtToBmonEvo->Divide(3, 3);
115
116 uPadIdx = 1;
117 for (std::string sDet : fvDets) {
118 fCanvDtToBmonEvo->cd(uPadIdx);
119 gPad->SetLogz();
120 gPad->SetGridx();
121 gPad->SetGridy();
122 fHistDtToBmonEvo[sDet]->Draw("colz");
123 uPadIdx++;
124 }
125
126 for (std::string sDet : fvDets) {
127 fCanvMul[sDet] = new TCanvas(Form("canvMul%s", sDet.data()), Form("Multiplicity %s", sDet.data()));
128 fCanvMul[sDet]->Divide(2);
129
130 fCanvMul[sDet]->cd(1);
131 gPad->SetLogy();
132 gPad->SetGridx();
133 gPad->SetGridy();
134 fHistMul[sDet]->Draw("hist");
135
136 fCanvMul[sDet]->cd(2);
137 gPad->SetLogz();
138 gPad->SetGridx();
139 gPad->SetGridy();
140 fHistDtMul[sDet]->Draw("colz");
141 }
142
143 gFile = oldFile;
144 gDirectory = oldDir;
145 }
146
147 size_t numEventsInTs = 0;
148 std::string sDet = "Bmon";
149 for (auto& event : *fEvents) {
150 bool bSingleBmon = false;
151 double_t dSingleBmonTime = 0.;
152 if (0 < event.fData.fBmon.Size()) {
153 if (1 == event.fData.fBmon.Size()) {
154 bSingleBmon = true;
155 }
156 dSingleBmonTime = event.fData.fBmon.fDigis[0].GetTime();
157 }
158
159 sDet = "Bmon";
160 fHistMul[sDet]->Fill(event.fData.fBmon.Size());
161 for (auto& digi : event.fData.fBmon.fDigis) {
162 double_t dDt = digi.GetTime() - event.fTime;
163 fHistDt[sDet]->Fill(dDt);
164 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
165 fHistDtMul[sDet]->Fill(dDt, event.fData.fBmon.Size());
166
167 // Special case: make internal Dt to first BMon instead of Dt to single Bmon
168 dDt = digi.GetTime() - dSingleBmonTime;
169 fHistDtToBmon[sDet]->Fill(dDt);
170 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
171 }
172
173 sDet = "Sts";
174 fHistMul[sDet]->Fill(event.fData.fSts.Size());
175 for (auto& digi : event.fData.fSts.fDigis) {
176 double_t dDt = digi.GetTime() - event.fTime;
177 fHistDt[sDet]->Fill(dDt);
178 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
179 fHistDtMul[sDet]->Fill(dDt, event.fData.fSts.Size());
180 if (bSingleBmon) {
181 dDt = digi.GetTime() - dSingleBmonTime;
182 fHistDtToBmon[sDet]->Fill(dDt);
183 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
184 }
185 }
186
187 sDet = "Much";
188 fHistMul[sDet]->Fill(event.fData.fMuch.Size());
189 for (auto& digi : event.fData.fMuch.fDigis) {
190 double_t dDt = digi.GetTime() - event.fTime;
191 fHistDt[sDet]->Fill(dDt);
192 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
193 fHistDtMul[sDet]->Fill(dDt, event.fData.fMuch.Size());
194 if (bSingleBmon) {
195 dDt = digi.GetTime() - dSingleBmonTime;
196 fHistDtToBmon[sDet]->Fill(dDt);
197 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
198 }
199 }
200
201 sDet = "Trd1d";
202 fHistMul[sDet]->Fill(event.fData.fTrd.Size());
203 for (auto& digi : event.fData.fTrd.fDigis) {
204 double_t dDt = digi.GetTime() - event.fTime;
205 fHistDt[sDet]->Fill(dDt);
206 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
207 fHistDtMul[sDet]->Fill(dDt, event.fData.fTrd.Size());
208 if (bSingleBmon) {
209 dDt = digi.GetTime() - dSingleBmonTime;
210 fHistDtToBmon[sDet]->Fill(dDt);
211 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
212 }
213 }
214
215 sDet = "Trd2d";
216 fHistMul[sDet]->Fill(event.fData.fTrd2d.Size());
217 for (auto& digi : event.fData.fTrd2d.fDigis) {
218 double_t dDt = digi.GetTime() - event.fTime;
219 fHistDt[sDet]->Fill(dDt);
220 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
221 fHistDtMul[sDet]->Fill(dDt, event.fData.fTrd2d.Size());
222 if (bSingleBmon) {
223 dDt = digi.GetTime() - dSingleBmonTime;
224 fHistDtToBmon[sDet]->Fill(dDt);
225 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
226 }
227 }
228
229 sDet = "Tof";
230 fHistMul[sDet]->Fill(event.fData.fTof.Size());
231 for (auto& digi : event.fData.fTof.fDigis) {
232 double_t dDt = digi.GetTime() - event.fTime;
233 fHistDt[sDet]->Fill(dDt);
234 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
235 fHistDtMul[sDet]->Fill(dDt, event.fData.fTof.Size());
236 if (bSingleBmon) {
237 dDt = digi.GetTime() - dSingleBmonTime;
238 fHistDtToBmon[sDet]->Fill(dDt);
239 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
240 }
241 }
242
243 sDet = "Rich";
244 fHistMul[sDet]->Fill(event.fData.fRich.Size());
245 for (auto& digi : event.fData.fRich.fDigis) {
246 double_t dDt = digi.GetTime() - event.fTime;
247 fHistDt[sDet]->Fill(dDt);
248 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
249 fHistDtMul[sDet]->Fill(dDt, event.fData.fRich.Size());
250 if (bSingleBmon) {
251 dDt = digi.GetTime() - dSingleBmonTime;
252 fHistDtToBmon[sDet]->Fill(dDt);
253 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
254 }
255 }
256
257 sDet = "Fsd";
258 fHistMul[sDet]->Fill(event.fData.fFsd.Size());
259 for (auto& digi : event.fData.fFsd.fDigis) {
260 double_t dDt = digi.GetTime() - event.fTime;
261 fHistDt[sDet]->Fill(dDt);
262 fHistDtEvo[sDet]->Fill(fNumTs, dDt);
263 fHistDtMul[sDet]->Fill(dDt, event.fData.fFsd.Size());
264 if (bSingleBmon) {
265 dDt = digi.GetTime() - dSingleBmonTime;
266 fHistDtToBmon[sDet]->Fill(dDt);
267 fHistDtToBmonEvo[sDet]->Fill(fNumTs, dDt);
268 }
269 }
270 numEventsInTs++;
271 }
272
273 // --- Run statistics
274 fNumTs++;
275 fNumEvents += fEvents->size();
276}
277// ----------------------------------------------------------------------------
278
279
280// ----- End-of-timeslice action ------------------------------------------
282{
283 LOG(info) << "=====================================";
284 LOG(info) << GetName() << ": Run summary";
285 LOG(info) << "Timeslices : " << fNumTs;
286 LOG(info) << "Events : " << fNumEvents;
287 LOG(info) << "=====================================";
288}
289// ----------------------------------------------------------------------------
290
291
292// ----- Initialisation ---------------------------------------------------
294{
295 // --- Get FairRootManager instance
296 FairRootManager* ioman = FairRootManager::Instance();
297 assert(ioman);
298
299 LOG(info) << "==================================================";
300 LOG(info) << GetName() << ": Initialising...";
301
302 // --- Input data
303 fEvents = ioman->InitObjectAs<const std::vector<CbmDigiEvent>*>("DigiEvent");
304 if (!fEvents) {
305 LOG(error) << GetName() << ": No input branch DigiEvent!";
306 return kFATAL;
307 }
308 LOG(info) << "--- Found branch DigiEvent at " << fEvents;
309
310 LOG(info) << "==================================================";
311 return kSUCCESS;
312}
313// ----------------------------------------------------------------------------
314
ClassImp(CbmConverterManager)
Create and fills plots of time differences to trigger in DigiEvents.
virtual void Finish()
Finish timeslice.
std::map< std::string, TCanvas * > fCanvMul
virtual ~CbmMcbm2024CheckEventsDt()
Destructor.
std::map< std::string, TH2 * > fHistDtEvo
size_t fNumEvents
Number of events.
std::map< std::string, TH2 * > fHistDtToBmonEvo
virtual void Exec(Option_t *opt)
Task execution.
std::map< std::string, TH1 * > fHistDt
virtual InitStatus Init()
Task initialisation.
std::map< std::string, TH2 * > fHistDtMul
std::vector< std::string > fvDets
std::map< std::string, TH1 * > fHistMul
const std::vector< CbmDigiEvent > * fEvents
size_t fNumTs
Input data (events)
std::map< std::string, TH1 * > fHistDtToBmon
Hash for CbmL1LinkKey.