BmnRoot
Loading...
Searching...
No Matches
BmnEventSelector.cxx
Go to the documentation of this file.
1#include "BmnEventSelector.h"
2
3#include "BmnEventHeader.h"
4#include "BmnStripDigit.h"
5#include "BmnTof1Digit.h"
6#include "BmnTrigDigit.h"
7#include "BmnTrigWaveDigit.h"
8#include "FairRootManager.h"
9
10#include <TClonesArray.h>
11#include <TF1.h>
12#include <TFile.h>
13#include <TH2F.h>
14#include <THnSparse.h>
15#include <TProfile.h>
16#include <numeric>
17#include <string>
18#include <vector>
19
20// ---- Default constructor -------------------------------------------
22 : FairTask("BmnEventSelector")
23{
24 LOG(debug) << "Default Constructor of BmnEventSelector";
25}
26
27// ---- Destructor ----------------------------------------------------
29{
30 LOG(debug) << "Destructor of BmnEventSelector";
31}
32
33// ---- Initialisation ----------------------------------------------
35{
36 LOG(debug) << "SetParContainers of BmnEventSelector";
37 // Load all necessary parameter containers from the runtime data base
38 /*
39 FairRunAna* ana = FairRunAna::Instance();
40 FairRuntimeDb* rtdb = ana->GetRuntimeDb();
41
42 <BmnEventSelectorDataMember> = (<ClassPointer>*)
43 (rtdb->getContainer("<ContainerName>"));
44 */
45}
46
47template<typename T>
48T BmnEventSelector::ReadBranch(const char* name)
49{
50 auto branch = (T)ioman->GetObject(name);
51 if (!branch)
52 LOG(fatal) << Form("BmnEventSelector: %s branch not found!!!\n", name);
53 return branch;
54}
55
56bool BmnEventSelector::FitRunDistr(std::string histName, std::vector<double>& fitPar, float fitMin, float fitMax)
57{
58 cout << endl;
59 LOG(info) << __func__ << ": " << histName;
60 fitPar.clear();
61 auto h2 = static_cast<TH2*>(fInputFile->Get(histName.c_str()));
62 if (!h2) {
63 fNoCalibrationData = true;
64 LOG(error) << "No input histogram: " << histName;
65 return false;
66 }
67 auto runBin = h2->GetXaxis()->FindBin(fRunId);
68 auto h1 = h2->ProjectionY(Form("%s_fit", h2->GetName()), runBin, runBin);
69 if (!h1 || h1->GetEntries() == 0) {
70 fNoCalibrationData = true;
71 LOG(error) << "No input histogram " << histName << " for run number " << fRunId;
72 }
73 h1->Fit("gaus", "", "", fitMin, fitMax);
74 auto f = (TF1*)h1->GetListOfFunctions()->FindObject("gaus");
75 if (!f) {
76 fNoCalibrationData = true;
77 LOG(error) << "ERROR fitting distribution for run number " << fRunId;
78 return false;
79 }
80 if (fFillCalibrationHists)
81 fOutputList->Add(h1);
82 fitPar = {f->GetParameter(0), f->GetParameter(1), f->GetParameter(2)};
83 return true;
84}
85
86// ---- Init ----------------------------------------------------------
88{
89 LOG(debug) << "Initilization of BmnEventSelector";
90
91 if (fFillCalibrationHists) {
92 fOutputFile = new TFile(fOutputFileName.c_str(), "recreate");
93 fOutputList = new TList;
95 }
96
97 fInputFile = new TFile(fInputFileName.c_str());
98 if (fInputFile->IsOpen()) {
99 vector<double> fitPar;
100 if (FitRunDistr("h2runId_time_SingleBC1S", fitPar, 2000, 2500))
101 fSingleBC1hitTimeMean = fitPar.at(1);
102 LOG(info) << "fSingleBC1hitTimeMean=" << fSingleBC1hitTimeMean;
103
104 if (FitRunDistr("h2runId_shapeShifted_SingleBC1S", fitPar, 0, 100)) {
105 fBC1SpeakTimeMin = fitPar.at(1) - 2 * fitPar.at(2);
106 fBC1SpeakTimeMax = fitPar.at(1) + 2 * fitPar.at(2);
107 }
108 LOG(info) << "fBC1SpeakTimeMin=" << fBC1SpeakTimeMin << "\tfBC1SpeakTimeMax=" << fBC1SpeakTimeMax;
109
110 if (FitRunDistr("h2runId_time_SingleVCS_BC1S", fitPar, -10, 10)) {
111 fSingleVCS_BC1S_hitTimeMin = fitPar.at(1) - 10;
112 fSingleVCS_BC1S_hitTimeMax = fitPar.at(1) + 10;
113 }
114 LOG(info) << "fSingleVCS_BC1S_hitTimeMin=" << fSingleVCS_BC1S_hitTimeMin
115 << "\tfSingleVCS_BC1S_hitTimeMax=" << fSingleVCS_BC1S_hitTimeMax;
116
117 if (FitRunDistr("h2runId_shapeShiftedBC1S_SingleVCS", fitPar, 0, 80)) {
118 fVCSpeakTimeMin = fitPar.at(1) - 2 * fitPar.at(2);
119 fVCSpeakTimeMax = fitPar.at(1) + 2 * fitPar.at(2);
120 }
121 LOG(info) << "fVCSpeakTimeMin=" << fVCSpeakTimeMin << "\tfVCSpeakTimeMax=" << fVCSpeakTimeMax;
122
123 if (FitRunDistr("h2runId_time_SingleBC2AS_BC1S", fitPar, -20, 0)) {
124 fSingleBC2AS_BC1S_hitTimeMin = fitPar.at(1) - 2;
125 fSingleBC2AS_BC1S_hitTimeMax = fitPar.at(1) + 2;
126 }
127 LOG(info) << "fSingleBC2AS_BC1S_hitTimeMin=" << fSingleBC2AS_BC1S_hitTimeMin
128 << "\tfSingleBC2AS_BC1S_hitTimeMax=" << fSingleBC2AS_BC1S_hitTimeMax;
129
130 if (FitRunDistr("h2runId_shapeShiftedBC1S_SingleBC2AS", fitPar, 0, 80)) {
131 fBC2ASpeakTimeMin = fitPar.at(1) - 2 * fitPar.at(2);
132 fBC2ASpeakTimeMax = fitPar.at(1) + 2 * fitPar.at(2);
133 }
134 LOG(info) << "fBC2ASpeakTimeMin=" << fBC2ASpeakTimeMin << "\tfBC2ASpeakTimeMax=" << fBC2ASpeakTimeMax;
135
136 if (FitRunDistr("h2runId_peak_BC2AS", fitPar, 0., 1.e4))
137 fSingleBC2ASpeakMin = fitPar.at(1) - 3 * fitPar.at(2);
138 LOG(info) << "\tfSingleBC2ASpeakMin=" << fSingleBC2ASpeakMin;
139
140 if (FitRunDistr("h2runId_time_SingleFD_BC1S", fitPar, 100, 200)) {
141 fSingleFD_BC1S_hitTimeMin = fitPar.at(1) - 10;
142 fSingleFD_BC1S_hitTimeMax = fitPar.at(1) + 20;
143 }
144 LOG(info) << "fSingleFD_BC1S_hitTimeMin=" << fSingleFD_BC1S_hitTimeMin
145 << "\tfSingleFD_BC1S_hitTimeMax=" << fSingleFD_BC1S_hitTimeMax;
146
147 if (FitRunDistr("h2runId_shapeShiftedBC1S_SingleFD", fitPar, 150, 250)) {
148 fFDpeakTimeMin = fitPar.at(1) - 2 * fitPar.at(2);
149 fFDpeakTimeMax = fitPar.at(1) + 2 * fitPar.at(2);
150 }
151 LOG(info) << "fFDpeakTimeMin=" << fFDpeakTimeMin << "\tfFDpeakTimeMax=" << fFDpeakTimeMax;
152
153 if (FitRunDistr("h2runId_peak_FD_nc", fitPar, 4.e3, 7.e3))
154 fSingleFDpeakMin = fitPar.at(1) - 3 * fitPar.at(2);
155 LOG(info) << "\tfSingleFDpeakMin=" << fSingleFDpeakMin;
156
157 if (FitRunDistr("h2runId_shapeShiftedBC1S_BD", fitPar, -150, -50)) {
158 fBDtimeMin = fitPar.at(1) - 3 * fitPar.at(2);
159 fBDtimeMax = fitPar.at(1) + 3 * fitPar.at(2);
160 }
161 LOG(info) << "fBDtimeMin=" << fBDtimeMin << "\tfBDtimeMax=" << fBDtimeMax;
162
163 if (FitRunDistr("h2runId_shapeShiftedBC1S_TOF400digit", fitPar)) {
164 fTOF400digitTimeMin = fitPar.at(1) - 3 * fitPar.at(2);
165 fTOF400digitTimeMax = fitPar.at(1) + 3 * fitPar.at(2);
166 }
167 LOG(info) << "fTOF400digitTimeMin=" << fTOF400digitTimeMin << "\tfTOF400digitTimeMax=" << fTOF400digitTimeMax;
168
169 if (FitRunDistr("h2runId_shapeShiftedBC1S_TOF700digit", fitPar)) {
170 fTOF700digitTimeMin = fitPar.at(1) - 3 * fitPar.at(2);
171 fTOF700digitTimeMax = fitPar.at(1) + 3 * fitPar.at(2);
172 }
173 LOG(info) << "fTOF700digitTimeMin=" << fTOF700digitTimeMin << "\tfTOF700digitTimeMax=" << fTOF700digitTimeMax;
174
175 } else {
176 fNoCalibrationData = true;
177 LOG(error) << "No input histogram file!!!";
178 }
179 // Get a handle from the IO manager
180 ioman = FairRootManager::Instance();
181 if (!ioman)
182 LOG(fatal) << "Init" << "FairRootManager is not instantiated";
183
184 // Get a pointer to the previous already existing data level
185 fEventHeader = ReadBranch<BmnEventHeader*>("BmnEventHeader.");
186 fArrTQDC_BC1S = ReadBranch<TClonesArray*>("TQDC_BC1S");
187 fArrTQDC_VCS = ReadBranch<TClonesArray*>("TQDC_VCS");
188 fArrTQDC_BC2AS = ReadBranch<TClonesArray*>("TQDC_BC2AS");
189 fArrTQDC_FD = ReadBranch<TClonesArray*>("TQDC_FD");
190 fArrBDdigits = ReadBranch<TClonesArray*>("BD");
191 fArrTOF400digits = ReadBranch<TClonesArray*>("TOF400");
192 fArrTOF700digits = ReadBranch<TClonesArray*>(
193 "TOF701"); // new tof700 digits are stored in TOF701 branch (09.2025). Will be fixed later
194 fArrSILICONdigits = ReadBranch<TClonesArray*>("SILICON");
195 fArrGEMdigits = ReadBranch<TClonesArray*>("GEM");
196
197 ioman->Register("BmnBC1hitInfo.", "BC1hitInfo", &fBC1hitInfo, fWriteBC1hitInfo);
198
199 return kSUCCESS;
200}
201
203{
204 h2runId_eventClass = AddHist(new TH2F("h2runId_eventClass", "Event class;runId;Event class", 1, fRunId, fRunId + 1,
206 for (int i = 0; i < BmnEventClass::kNclasses; i++)
207 h2runId_eventClass->GetYaxis()->SetBinLabel(i + 1, BmnEventClass::names.at(i).c_str());
208 auto h2runId_time =
209 new TH2F("h2runId_time", "TDC time;runId;TDC time (ns);N entries", 1, fRunId, fRunId + 1, 900, 0, 3600);
210 h2runId_time_SingleBC1S = AddHistClone(h2runId_time, "_SingleBC1S");
211 h2runId_time_BC1S_k0 = AddHistClone(h2runId_time, "_BC1S_k0");
212 h2runId_time_BC1S_k1 = AddHistClone(h2runId_time, "_BC1S_k1");
213 h2runId_time_BC1S_kV0 = AddHistClone(h2runId_time, "_BC1S_kV0");
214 h2runId_time_BC1S_kV1 = AddHistClone(h2runId_time, "_BC1S_kV1");
215 h2runId_time_SingleVCS_BC1S = AddHistClone(h2runId_time, "_SingleVCS_BC1S");
216 h2runId_time_SingleVCS_BC1S->GetYaxis()->SetLimits(-50, 50);
217 h2runId_time_SingleBC2AS_BC1S = AddHistClone(h2runId_time, "_SingleBC2AS_BC1S");
218 h2runId_time_SingleBC2AS_BC1S->GetYaxis()->SetLimits(-20, 0);
219 h2runId_time_SingleFD_BC1S = AddHistClone(h2runId_time, "_SingleFD_BC1S");
220 h2runId_time_SingleFD_BC1S->GetYaxis()->SetLimits(-50, 200);
221 auto h2runId_shape =
222 new TH2F("h2runId_shape", "TQDC shape;runId;Time (ns);TQDC signal", 1, fRunId, fRunId + 1, 500, -400, 3600);
223 h2runId_shapeShifted_SingleBC1S = AddHistClone(h2runId_shape, "Shifted_SingleBC1S");
224 h2runId_shapeShiftedBC1S_SingleVCS = AddHistClone(h2runId_shape, "ShiftedBC1S_SingleVCS");
225 h2runId_shapeShiftedBC1S_SingleBC2AS = AddHistClone(h2runId_shape, "ShiftedBC1S_SingleBC2AS");
226 h2runId_shapeShiftedBC1S_SingleFD = AddHistClone(h2runId_shape, "ShiftedBC1S_SingleFD");
227 h2runId_shapeShiftedBC1S_BD = AddHistClone(h2runId_shape, "ShiftedBC1S_BD");
228 h2runId_shapeShiftedBC1S_BD->GetYaxis()->SetLimits(-450, 0);
229 h2runId_NdigitsBD =
230 AddHist(new TH2F("h2runId_NdigitsBD", "NdigitsBD;runId;NdigitsBD", 1, fRunId, fRunId + 1, 50, 0, 50));
231 h2runId_SumAmpBD =
232 AddHist(new TH2F("h2runId_SumAmpBD", "SumAmpBD;runId;SumAmpBD", 1, fRunId, fRunId + 1, 200, 0, 800));
233 h2runId_NdigitsBD_SingleBC1S = AddHistClone(h2runId_NdigitsBD, "_SingleBC1S");
234 h2runId_SumAmpBD_SingleBC1S = AddHistClone(h2runId_SumAmpBD, "_SingleBC1S");
235 h2runId_NdigitsBD_nc = AddHistClone(h2runId_NdigitsBD, "_nc");
236 h2runId_SumAmpBD_nc = AddHistClone(h2runId_SumAmpBD, "_nc");
237 auto h2runId_NdigitsTOF =
238 new TH2F("h2runId_NdigitsTOF", "NdigitsTOF;runId;NdigitsTOF", 1, fRunId, fRunId + 1, 400, 0, 400);
239 h2runId_NdigitsTOF400 = AddHistClone(h2runId_NdigitsTOF, "400");
240 h2runId_NdigitsTOF400_SingleBC1S = AddHistClone(h2runId_NdigitsTOF, "400_SingleBC1S");
241 h2runId_NdigitsTOF400_nc = AddHistClone(h2runId_NdigitsTOF, "400_nc");
242 h2runId_NdigitsTOF700 = AddHistClone(h2runId_NdigitsTOF, "700");
243 h2runId_NdigitsTOF700_SingleBC1S = AddHistClone(h2runId_NdigitsTOF, "700_SingleBC1S");
244 h2runId_NdigitsTOF700_nc = AddHistClone(h2runId_NdigitsTOF, "700_nc");
245 auto h2runId_peak = new TH2F("h2runId_peak", "Peak;runId;Peak;N entries", 1, fRunId, fRunId + 1, 300, 0, 3e4);
246 h2runId_peak_SingleBC1S = AddHistClone(h2runId_peak, "_SingleBC1S");
247 h2runId_peak_BC1S = AddHistClone(h2runId_peak, "_BC1S");
248 h2runId_peak_SingleVCS = AddHistClone(h2runId_peak, "_SingleVCS");
249 h2runId_peak_VCS = AddHistClone(h2runId_peak, "_VCS");
250 h2runId_peak_SingleBC2AS = AddHistClone(h2runId_peak, "_SingleBC2AS");
251 h2runId_peak_BC2AS = AddHistClone(h2runId_peak, "_BC2AS");
252 h2runId_peak_SingleFD = AddHistClone(h2runId_peak, "_SingleFD");
253 h2runId_peak_FD = AddHistClone(h2runId_peak, "_FD");
254 h2runId_peak_FD_nc = AddHistClone(h2runId_peak, "_FD_nc");
255
256 h2runId_shapeShiftedBC1S_TOF400digit =
257 AddHist(new TH2F("h2runId_shapeShiftedBC1S_TOF400digit", "TOF400 digit times;runId;Time (ns)", 1, fRunId,
258 fRunId + 1, 1000, -700, -600));
259 h2runId_shapeShiftedBC1S_TOF700digit =
260 AddHist(new TH2F("h2runId_shapeShiftedBC1S_TOF700digit", "TOF700 digit times;runId;Time (ns)", 1, fRunId,
261 fRunId + 1, 1000, -2300, -1800));
262
263 if (fFillPerformanceHists) {
264 vector<string> trigNames{"BT", "CCT1", "CCT2", "MBT"};
265 vector<int> bins{
266 fRunMax - fRunMin, int(trigNames.size()), BmnEventClass::kNclasses, 185, 10000,
267 };
268 vector<double> xmin{
269 double(fRunMin), 0, 0, -2200, 0,
270 };
271 vector<double> xmax{
272 double(fRunMax), double(trigNames.size()), double(BmnEventClass::kNclasses), 1500, 10000,
273 };
274 auto h5runId_trigger_eventClass_dt_n =
275 new THnSparseI("h5runId_trigger_eventClass_dt_n", "Number of;runId;trigger;Event class;#Deltat (ns);N",
276 bins.size(), bins.data(), xmin.data(), xmax.data());
277 for (uint i = 0; i < trigNames.size(); i++)
278 h5runId_trigger_eventClass_dt_n->GetAxis(1)->SetBinLabel(i + 1, trigNames.at(i).c_str());
279 for (int i = 0; i < BmnEventClass::kNclasses; i++)
280 h5runId_trigger_eventClass_dt_n->GetAxis(2)->SetBinLabel(i + 1, BmnEventClass::names.at(i).c_str());
281
282 int n = 400;
283 vector<int> rebin = {1, 1, 1, 1, int(xmax.at(4) / n)};
284 auto h5runId_trigger_eventClass_dt_nTOF =
285 static_cast<THnSparse*>(h5runId_trigger_eventClass_dt_n->Clone())->Rebin(rebin.data());
286 h5runId_trigger_eventClass_dt_nTOF->SetName("h5runId_trigger_eventClass_dt_nTOF");
287 h5runId_trigger_eventClass_dt_nTOF->SetTitle("Number of TOF");
288 h5runId_trigger_eventClass_dt_nTOF->GetAxis(4)->SetLimits(0, n);
289 h5runId_trigger_eventClass_dt_nTOF400digits = AddHistClone(h5runId_trigger_eventClass_dt_nTOF, "400digits");
290 h5runId_trigger_eventClass_dt_nTOF700digits = AddHistClone(h5runId_trigger_eventClass_dt_nTOF, "700digits");
291
292 n = 5000;
293 rebin = {1, 1, 1, 1, int(xmax.at(4) / n)};
294 auto h5runId_trigger_eventClass_dt_nSILICON =
295 static_cast<THnSparse*>(h5runId_trigger_eventClass_dt_n->Clone())->Rebin(rebin.data());
296 h5runId_trigger_eventClass_dt_nSILICON->SetName("h5runId_trigger_eventClass_dt_nSILICON");
297 h5runId_trigger_eventClass_dt_nSILICON->SetTitle("Number of SILICON");
298 h5runId_trigger_eventClass_dt_nSILICON->GetAxis(4)->SetLimits(0, n);
299 h5runId_trigger_eventClass_dt_nSILICONdigits = AddHistClone(h5runId_trigger_eventClass_dt_nSILICON, "digits");
300
301 h5runId_trigger_eventClass_dt_nGEMdigits = AddHistClone(h5runId_trigger_eventClass_dt_n, "GEMdigits");
302 }
303}
304
305// ---- ReInit -------------------------------------------------------
307{
308 LOG(debug) << "Initilization of BmnEventSelector";
309 return kSUCCESS;
310}
311
312int BmnEventSelector::GetClosestHitIndex(std::vector<double> hitTimes, float refTime)
313{
314 if (hitTimes.size() == 0)
315 return -1;
316 int index = 0;
317 float min = fabs(hitTimes.at(0) - refTime);
318 for (uint i = 1; i < hitTimes.size(); i++) {
319 auto dt = fabs(hitTimes.at(i) - refTime);
320 if (dt < min) {
321 min = dt;
322 index = i;
323 }
324 }
325 return index;
326}
327
328int BmnEventSelector::GetPeak(BmnTrigWaveDigit* digit, float timeMin, float timeMax)
329{
330 if (!digit)
331 return 0;
332 int binMin = timeMin / fNsPerChannel;
333 uint binMax = timeMax / fNsPerChannel;
334 if (binMin < 0 || binMax >= digit->GetNSamples())
335 return -999;
336 return digit->GetPeak(binMin, binMax);
337}
338
339template<typename T>
340std::vector<float> BmnEventSelector::GetSumAmp(TClonesArray* array,
341 std::vector<double> hitTimes,
342 float relTimeMin,
343 float relTimeMax)
344{
345 std::vector<float> sum(hitTimes.size(), 0.);
346 for (uint i = 0; i < hitTimes.size(); i++) {
347 for (auto obj : *array) {
348 auto digit = static_cast<T>(obj);
349 auto relTime = digit->GetTime() - hitTimes.at(i);
350 if (relTimeMin < relTime && relTime < relTimeMax)
351 sum.at(i) += digit->GetAmp();
352 }
353 }
354 return sum;
355}
356
357template<typename T>
358std::vector<int> BmnEventSelector::GetNdigits(TClonesArray* array,
359 std::vector<double> hitTimes,
360 float relTimeMin,
361 float relTimeMax)
362{
363 std::vector<int> n(hitTimes.size(), 0);
364 for (uint i = 0; i < hitTimes.size(); i++) {
365 for (auto obj : *array) {
366 auto relTime = static_cast<T>(obj)->GetTime() - hitTimes.at(i);
367 if (relTimeMin < relTime && relTime < relTimeMax)
368 n.at(i)++;
369 }
370 }
371 return n;
372}
373
374std::vector<int> BmnEventSelector::GetPeaks(BmnTrigWaveDigit* digit,
375 std::vector<double> hitTimes,
376 float relTimeMin,
377 float relTimeMax,
378 float relHitTimeMin,
379 float relHitTimeMax)
380{
381 std::vector<int> ret(hitTimes.size(), -999);
382 if (!digit)
383 return ret;
384 auto hitTimesDigit = digit->TdcVector();
385 auto tMax = fNsPerChannel * digit->GetNSamples();
386 for (uint i = 0; i < hitTimes.size(); i++) {
387 auto hitTime = hitTimes.at(i);
388 auto hitTimeMin = hitTime + relHitTimeMin;
389 auto hitTimeMax = hitTime + relHitTimeMax;
390 if (hitTimeMin < 0 || hitTimeMax > tMax)
391 continue;
392 auto hitFound = true; // FIX ME!!! -> false
393 for (uint j = 0; j < hitTimesDigit.size(); j++) {
394 auto hitTimeDigit = hitTimesDigit.at(j);
395 if (hitTimeMin <= hitTimeDigit && hitTimeDigit <= hitTimeMax) {
396 hitFound = true;
397 break;
398 }
399 }
400 auto timeMin = hitTime + relTimeMin;
401 auto timeMax = hitTime + relTimeMax;
402 if (hitFound)
403 ret.at(i) = GetPeak(digit, timeMin, timeMax);
404 else
405 ret.at(i) = 0;
406 }
407 return ret;
408}
409
411{
412 fTQDC_BC1S = (BmnTrigWaveDigit*)fArrTQDC_BC1S->At(0);
413 fTQDC_VCS = (BmnTrigWaveDigit*)fArrTQDC_VCS->At(0);
414 fTQDC_BC2AS = (BmnTrigWaveDigit*)fArrTQDC_BC2AS->At(0);
415 fTQDC_FD = (BmnTrigWaveDigit*)fArrTQDC_FD->At(0);
416
417 if (!fTQDC_BC1S || !fTQDC_VCS || !fTQDC_FD)
418 return;
419
420 fTimesBC1S = fTQDC_BC1S->TdcVector();
421 fNhitsBC1S = fTimesBC1S.size();
422 if (fNhitsBC1S == 0) {
423 if (fFillCalibrationHists)
424 h2runId_eventClass->Fill(fRunId, BmnEventClass::kNull);
425 return;
426 }
427 fTimesVCS = fTQDC_VCS->TdcVector();
428 fNhitsVCS = fTimesVCS.size();
429 fTimesBC2AS = fTQDC_BC2AS->TdcVector();
430 fNhitsBC2AS = fTimesBC2AS.size();
431 fTimesFD = fTQDC_FD->TdcVector();
432 fNhitsFD = fTimesFD.size();
433
434 fBC1Speaks = GetPeaks(fTQDC_BC1S, fTimesBC1S, fBC1SpeakTimeMin, fBC1SpeakTimeMax);
435 fVCSpeaks = GetPeaks(fTQDC_VCS, fTimesBC1S, fVCSpeakTimeMin, fVCSpeakTimeMax, fSingleVCS_BC1S_hitTimeMin,
436 fSingleVCS_BC1S_hitTimeMax);
437 fBC2ASpeaks = GetPeaks(fTQDC_BC2AS, fTimesBC1S, fBC2ASpeakTimeMin, fBC2ASpeakTimeMax, fSingleBC2AS_BC1S_hitTimeMin,
438 fSingleBC2AS_BC1S_hitTimeMax);
439 fFDpeaks = GetPeaks(fTQDC_FD, fTimesBC1S, fFDpeakTimeMin, fFDpeakTimeMax, fSingleFD_BC1S_hitTimeMin,
440 fSingleFD_BC1S_hitTimeMax);
441 fNdigitsBD = GetNdigits<BmnTrigDigit*>(fArrBDdigits, fTimesBC1S, fBDtimeMin, fBDtimeMax);
442 fSumAmpBD = GetSumAmp<BmnTrigDigit*>(fArrBDdigits, fTimesBC1S, fBDtimeMin, fBDtimeMax);
443 fNdigitsTOF400 = GetNdigits<BmnTof1Digit*>(fArrTOF400digits, fTimesBC1S, fTOF400digitTimeMin, fTOF400digitTimeMax);
444 fNdigitsTOF700 = GetNdigits<BmnTof1Digit*>(fArrTOF700digits, fTimesBC1S, fTOF700digitTimeMin, fTOF700digitTimeMax);
445}
446
448{
449 using namespace BmnEventClass;
450 BmnEventClass::id status;
451
452 for (int i = 0; i < fNhitsBC1S; i++) {
453 status = k0;
454 if (fNoCalibrationData || fBC1Speaks.at(i) < -998 || fVCSpeaks.at(i) < -998 || fBC2ASpeaks.at(i) < -998
455 || fFDpeaks.at(i) < -998)
456 {
457 fBC1hitClasses.push_back(kUndef);
458 continue;
459 }
460
461 if (fFDpeaks.at(i) < fSingleFDpeakMin)
462 status = k1;
463 if (fVCSpeaks.at(i) > fSingleVCSpeakMin || fBC2ASpeaks.at(i) < fSingleBC2ASpeakMin) {
464 if (status == k0)
465 status = kV0;
466 if (status == k1)
467 status = kV1;
468 }
469 fBC1hitClasses.push_back(status);
470 }
471
472 fCentralHitIndexBC1S = GetClosestHitIndex(fTimesBC1S, fSingleBC1hitTimeMean);
473 if (fCentralHitIndexBC1S < 0)
474 return kNull;
475
476 if (fBC1hitClasses.at(fCentralHitIndexBC1S) == k0
477 || fBC1hitClasses.at(fCentralHitIndexBC1S) == kV0) // search backwards
478 {
479 const float DTmax = 300; // ns around mean for single BC1S hit
480 for (int i = fCentralHitIndexBC1S - 1; i >= 0 && fSingleBC1hitTimeMean - fTimesBC1S.at(i) < DTmax; i--) {
481 if (fBC1hitClasses.at(i) == k1 || fBC1hitClasses.at(i) == kV1) {
482 fCentralHitIndexBC1S = i;
483 break;
484 }
485 }
486 }
487 fCentralHitTimeBC1S = fTimesBC1S.at(fCentralHitIndexBC1S);
488 fRelTimesBC1S = fTimesBC1S;
489 for (auto& t : fRelTimesBC1S)
490 t -= fCentralHitTimeBC1S;
491
492 if (fNhitsBC1S == 1)
493 return status;
494 if (fNhitsBC1S == 2) {
495 int c = fCentralHitIndexBC1S;
496 int nc = abs(c - 1);
497
498 if (fBC1hitClasses.at(c) == k0 && fBC1hitClasses.at(nc) == k0)
499 return k00;
500 if (fBC1hitClasses.at(c) == k0 && fBC1hitClasses.at(nc) == k1)
501 return k01;
502 if (fBC1hitClasses.at(c) == k0 && fBC1hitClasses.at(nc) == kV0)
503 return k0V0;
504 if (fBC1hitClasses.at(c) == k0 && fBC1hitClasses.at(nc) == kV1)
505 return k0V1;
506 if (fBC1hitClasses.at(c) == k1 && fBC1hitClasses.at(nc) == k0)
507 return k10;
508 if (fBC1hitClasses.at(c) == k1 && fBC1hitClasses.at(nc) == k1)
509 return k11;
510 if (fBC1hitClasses.at(c) == k1 && fBC1hitClasses.at(nc) == kV0)
511 return k1V0;
512 if (fBC1hitClasses.at(c) == k1 && fBC1hitClasses.at(nc) == kV1)
513 return k1V1;
514 }
515 if (fNhitsBC1S > 2)
516 return kMult;
517 return kUndef;
518}
519
520// ---- FillHistograms
521// ----------------------------------------------------------
522void BmnEventSelector::FillHistograms()
523{
524 using namespace BmnTriggerBits;
525 bool isBT_AR = fTrigAR & BIT(kBT);
526 bool isMBT_AR = fTrigAR & BIT(kMBT);
527 bool isCCT1_AR = fTrigAR & BIT(kCCT1);
528 bool isCCT2_AR = fTrigAR & BIT(kCCT2);
529 vector<bool> trigAR = {isBT_AR, isCCT1_AR, isCCT2_AR, isMBT_AR};
530
531 fNdigitsSILICON = 0;
532 for (auto digit : *fArrSILICONdigits)
533 if (static_cast<BmnStripDigit*>(digit)->IsGoodDigit())
534 fNdigitsSILICON++;
535
536 fNdigitsGEM = 0;
537 for (auto digit : *fArrGEMdigits)
538 if (static_cast<BmnStripDigit*>(digit)->IsGoodDigit())
539 fNdigitsGEM++;
540
541 h2runId_eventClass->Fill(fRunId, fEventClass);
542 for (int i = 0; i < fNhitsBC1S; i++) {
543 if (fBC1hitClasses.at(i) == BmnEventClass::k0)
544 h2runId_time_BC1S_k0->Fill(fRunId, fTimesBC1S.at(i));
545 else if (fBC1hitClasses.at(i) == BmnEventClass::k1)
546 h2runId_time_BC1S_k1->Fill(fRunId, fTimesBC1S.at(i));
547 else if (fBC1hitClasses.at(i) == BmnEventClass::kV0)
548 h2runId_time_BC1S_kV0->Fill(fRunId, fTimesBC1S.at(i));
549 else if (fBC1hitClasses.at(i) == BmnEventClass::kV1)
550 h2runId_time_BC1S_kV1->Fill(fRunId, fTimesBC1S.at(i));
551 h2runId_peak_BC1S->Fill(fRunId, fBC1Speaks.at(i));
552 h2runId_peak_VCS->Fill(fRunId, fVCSpeaks.at(i));
553 h2runId_peak_BC2AS->Fill(fRunId, fBC2ASpeaks.at(i));
554 h2runId_peak_FD->Fill(fRunId, fFDpeaks.at(i));
555 h2runId_NdigitsBD->Fill(fRunId, fNdigitsBD.at(i));
556 h2runId_SumAmpBD->Fill(fRunId, fSumAmpBD.at(i));
557 h2runId_NdigitsTOF400->Fill(fRunId, fNdigitsTOF400.at(i));
558 h2runId_NdigitsTOF700->Fill(fRunId, fNdigitsTOF700.at(i));
559 if (i != fCentralHitIndexBC1S) {
560 h2runId_peak_FD_nc->Fill(fRunId, fFDpeaks.at(i));
561 h2runId_NdigitsBD_nc->Fill(fRunId, fNdigitsBD.at(i));
562 h2runId_SumAmpBD_nc->Fill(fRunId, fSumAmpBD.at(i));
563 h2runId_NdigitsTOF400_nc->Fill(fRunId, fNdigitsTOF400.at(i));
564 h2runId_NdigitsTOF700_nc->Fill(fRunId, fNdigitsTOF700.at(i));
565 }
566 }
567
568 if (fNhitsBC1S == 1) {
569 h2runId_time_SingleBC1S->Fill(fRunId, fTimesBC1S.at(0));
570 auto pulseBC1S = fTQDC_BC1S->GetShortValue();
571 for (uint i = 0; i < fTQDC_BC1S->GetNSamples(); i++) {
572 auto channelTime = i * fNsPerChannel;
573 h2runId_shapeShifted_SingleBC1S->Fill(fRunId, channelTime - fTimesBC1S.at(0), pulseBC1S[i]);
574 }
575 h2runId_peak_SingleBC1S->Fill(fRunId, fBC1Speaks.at(0));
576
577 for (auto digit : *fArrBDdigits) {
578 auto time = static_cast<BmnTrigDigit*>(digit)->GetTime() - fTimesBC1S.at(0);
579 h2runId_shapeShiftedBC1S_BD->Fill(fRunId, time);
580 }
581 h2runId_NdigitsBD_SingleBC1S->Fill(fRunId, fNdigitsBD.at(0));
582 h2runId_SumAmpBD_SingleBC1S->Fill(fRunId, fSumAmpBD.at(0));
583 h2runId_NdigitsTOF400_SingleBC1S->Fill(fRunId, fNdigitsTOF400.at(0));
584 h2runId_NdigitsTOF700_SingleBC1S->Fill(fRunId, fNdigitsTOF700.at(0));
585
586 for (auto digit : *fArrTOF400digits) {
587 auto time = static_cast<BmnTof1Digit*>(digit)->GetTime() - fTimesBC1S.at(0);
588 h2runId_shapeShiftedBC1S_TOF400digit->Fill(fRunId, time);
589 }
590
591 for (auto digit : *fArrTOF700digits) {
592 auto time = static_cast<BmnTof1Digit*>(digit)->GetTime() - fTimesBC1S.at(0);
593 h2runId_shapeShiftedBC1S_TOF700digit->Fill(fRunId, time);
594 }
595
596 if (fNhitsFD == 1) {
597 h2runId_time_SingleFD_BC1S->Fill(fRunId, fTimesFD.at(0) - fTimesBC1S.at(0));
598 h2runId_peak_SingleFD->Fill(fRunId, fFDpeaks.at(0));
599 auto pulseFD = fTQDC_FD->GetShortValue();
600 for (uint i = 0; i < fTQDC_FD->GetNSamples(); i++) {
601 auto channelTime = i * fNsPerChannel;
602 h2runId_shapeShiftedBC1S_SingleFD->Fill(fRunId, channelTime - fTimesBC1S.at(0), pulseFD[i]);
603 }
604 }
605 if (fNhitsBC2AS == 1) {
606 h2runId_time_SingleBC2AS_BC1S->Fill(fRunId, fTimesBC2AS.at(0) - fTimesBC1S.at(0));
607 h2runId_peak_SingleBC2AS->Fill(fRunId, fBC2ASpeaks.at(0));
608 auto pulseBC2AS = fTQDC_BC2AS->GetShortValue();
609 for (uint i = 0; i < fTQDC_BC2AS->GetNSamples(); i++) {
610 auto channelTime = i * fNsPerChannel;
611 h2runId_shapeShiftedBC1S_SingleBC2AS->Fill(fRunId, channelTime - fTimesBC1S.at(0), pulseBC2AS[i]);
612 }
613 }
614 }
615
616 if (fNhitsBC1S == 2 && fNhitsVCS == 1) {
617 h2runId_time_SingleVCS_BC1S->Fill(fRunId, fTimesVCS.at(0) - fTimesBC1S.at(0));
618 h2runId_peak_SingleVCS->Fill(fRunId, fVCSpeaks.at(0));
619 auto closestBC1Sindex = GetClosestHitIndex(fTimesBC1S, fTimesVCS.at(0));
620 auto closestBC1Stime = fTimesBC1S.at(closestBC1Sindex);
621 auto pulseVCS = fTQDC_VCS->GetShortValue();
622 for (uint i = 0; i < fTQDC_VCS->GetNSamples(); i++) {
623 auto channelTime = i * fNsPerChannel;
624 h2runId_shapeShiftedBC1S_SingleVCS->Fill(fRunId, channelTime - closestBC1Stime, pulseVCS[i]);
625 }
626 }
627
628 if (fFillPerformanceHists) {
629 float dt = 0;
630 if (fNhitsBC1S > 1) {
631 float dtLeft = 1e10, dtRight = 1e10;
632 if (fCentralHitIndexBC1S > 0)
633 dtLeft = fRelTimesBC1S.at(fCentralHitIndexBC1S - 1);
634 if (fCentralHitIndexBC1S < fNhitsBC1S - 1)
635 dtRight = fRelTimesBC1S.at(fCentralHitIndexBC1S + 1);
636 if (fabs(dtLeft) < dtRight)
637 dt = dtLeft;
638 else
639 dt = dtRight;
640 }
641
642 vector<double> x = {fRunId + 0.5, 0., fEventClass + 0.5, dt, 0.};
643 for (uint i = 0; i < trigAR.size(); i++) {
644 if (trigAR.at(i)) {
645 x.at(1) = i + 0.5;
646 x.at(4) = fNdigitsTOF400.at(fCentralHitIndexBC1S);
647 h5runId_trigger_eventClass_dt_nTOF400digits->Fill(x.data());
648 x.at(4) = fNdigitsTOF700.at(fCentralHitIndexBC1S);
649 h5runId_trigger_eventClass_dt_nTOF700digits->Fill(x.data());
650 x.at(4) = fNdigitsSILICON;
651 h5runId_trigger_eventClass_dt_nSILICONdigits->Fill(x.data());
652 x.at(4) = fNdigitsGEM;
653 h5runId_trigger_eventClass_dt_nGEMdigits->Fill(x.data());
654 }
655 }
656 }
657}
658
659// ---- Exec ----------------------------------------------------------
660void BmnEventSelector::Exec(Option_t* /*option*/)
661{
662 // cout << "Event " << ++iEvent << endl;
663 fTrigAR = fEventHeader->GetInputSignalsAR();
664 fEventClass = BmnEventClass::kUndef;
665 fNhitsBC1S = 0;
666 fCentralHitIndexBC1S = -1;
667 fTimesBC1S.clear();
668 fRelTimesBC1S.clear();
669 fBC1Speaks.clear();
670 fBC1hitClasses.clear();
671 fNdigitsBD.clear();
672 fSumAmpBD.clear();
673 fBC2ASpeaks.clear();
674 fVCSpeaks.clear();
675 fFDpeaks.clear();
676 fNdigitsTOF400.clear();
677 fNdigitsTOF700.clear();
678 fNdigitsSILICON = -999;
679 fNdigitsGEM = -999;
680 Calculate();
681 fEventClass = Classify();
682
683 fBC1hitInfo.fEventClass = fEventClass;
684 fBC1hitInfo.fCentralHitIndexBC1S = fCentralHitIndexBC1S;
685 fBC1hitInfo.fTimesBC1S = fTimesBC1S;
686 fBC1hitInfo.fRelTimesBC1S = fRelTimesBC1S;
687 fBC1hitInfo.fBC1hitClasses = fBC1hitClasses;
688 fBC1hitInfo.fBC1Speaks = fBC1Speaks;
689 fBC1hitInfo.fVCSpeaks = fVCSpeaks;
690 fBC1hitInfo.fBC2ASpeaks = fBC2ASpeaks;
691 fBC1hitInfo.fFDpeaks = fFDpeaks;
692 fBC1hitInfo.fNdigitsBD = fNdigitsBD;
693 fBC1hitInfo.fSumAmpBD = fSumAmpBD;
694 fBC1hitInfo.fNdigitsTOF400 = fNdigitsTOF400;
695 fBC1hitInfo.fNdigitsTOF700 = fNdigitsTOF700;
696
697 if (fFillCalibrationHists && fEventClass != BmnEventClass::kNull)
698 FillHistograms();
699}
700
701// ---- Finish --------------------------------------------------------
703{
704 LOG(debug) << "Finish of BmnEventSelector";
705 if (fFillCalibrationHists) {
706 fOutputFile->cd();
707 for (auto h : *fOutputList) {
708 string className = h->ClassName();
709 if (className.find("THn") == 0) {
710 vector<int> axes = {0, 1, 2, 4};
711 static_cast<THnSparse*>(h)->Projection(axes.size(), axes.data())->Write();
712 axes = {1, 2, 3, 4};
713 static_cast<THnSparse*>(h)->Projection(axes.size(), axes.data())->Write();
714 } else
715 h->Write();
716 }
717 for (auto h : *fOutputList) {
718 string name = h->GetName();
719 string className = h->ClassName();
720 if (className.find("TH2") == 0 && name.find("h2runId") == 0) {
721 // static_cast<TH2*>(h)->ProfileX()->Write();
722 static_cast<TH2*>(h)->ProjectionY()->Write();
723 }
724 }
725 fOutputFile->Close();
726 }
727}
friend F32vec4 fabs(const F32vec4 &a)
Definition P4_F32vec4.h:52
int i
Definition P4_F32vec4.h:22
friend F32vec4 min(const F32vec4 &a, const F32vec4 &b)
Definition P4_F32vec4.h:30
float f
Definition P4_F32vec4.h:21
UInt_t GetInputSignalsAR()
virtual InitStatus ReInit()
int GetClosestHitIndex(std::vector< double > hitTimes, float refTime)
virtual void Finish()
virtual void Exec(Option_t *opt)
virtual void SetParContainers()
virtual InitStatus Init()
BmnEventClass::id Classify()
Short_t * GetShortValue() const
int GetPeak(UInt_t start=0, UInt_t stop=1e9) const
UInt_t GetNSamples() const
Double_t GetTime() const
vector< Double_t > & TdcVector()
const std::map< int, std::string > names
@ array
array (ordered collection of values)