Belle II Software  release-08-01-10
SVDDQMClustersOnTrackModule.cc
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #include "svd/modules/svdDQM/SVDDQMClustersOnTrackModule.h"
10 
11 #include <hlt/softwaretrigger/core/FinalTriggerDecisionCalculator.h>
12 #include <framework/datastore/DataStore.h>
13 #include <framework/datastore/StoreObjPtr.h>
14 #include <framework/datastore/StoreArray.h>
15 #include <framework/dataobjects/EventMetaData.h>
16 #include <svd/dataobjects/SVDShaperDigit.h>
17 #include <svd/dataobjects/SVDRecoDigit.h>
18 #include <svd/dataobjects/SVDCluster.h>
19 #include <tracking/dataobjects/RecoTrack.h>
20 #include <vxd/geometry/GeoTools.h>
21 
22 #include "TDirectory.h"
23 
24 using namespace std;
25 using namespace Belle2;
26 using namespace SoftwareTrigger;
27 
28 //-----------------------------------------------------------------
29 // Register the Module
30 //-----------------------------------------------------------------
31 REG_MODULE(SVDDQMClustersOnTrack);
32 
33 
34 //-----------------------------------------------------------------
35 // Implementation
36 //-----------------------------------------------------------------
37 
38 SVDDQMClustersOnTrackModule::SVDDQMClustersOnTrackModule() : HistoModule()
39 {
40  //Set module properties
41  setDescription("SVD DQM module for clusters related to tracks.");
42 
43  setPropertyFlags(c_ParallelProcessingCertified); // specify this flag if you need parallel processing
44  addParam("skipHLTRejectedEvents", m_skipRejectedEvents, "If True, skip events rejected by HLT.", bool(true));
45  addParam("TriggerBin", m_tb, "select events for a specific trigger bin, if -1 then no selection is applied (default)", int(-1));
46  addParam("histogramDirectoryName", m_histogramDirectoryName, "Name of the directory where histograms will be placed.",
47  std::string("SVDClsTrk"));
48  addParam("desynchronizeSVDTime", m_desynchSVDTime,
49  "if True, svd time back in SVD time reference", bool(false));
50  addParam("EventInfo", m_svdEventInfoName, "SVDEventInfo StoreArray name.", std::string(""));
51  addParam("Clusters", m_svdClustersName, "SVDCluster StoreArray name.", std::string(""));
52  addParam("RecoDigits", m_svdRecoDigitsName, "SVDRecoDigits StoreArray name.", std::string(""));
53  addParam("ShaperDigits", m_svdShaperDigitsName, "SVDShaperDigits StoreArray name.", std::string(""));
54 
55  m_histoList = new TList();
56 
57  m_ladderMap = {
58  {{3, 1}, 0}, {{3, 2}, 1},
59  {{4, 1}, 2}, {{4, 2}, 3}, {{4, 3}, 4},
60  {{5, 1}, 5}, {{5, 2}, 6}, {{5, 3}, 7}, {{5, 4}, 8},
61  {{6, 1}, 9}, {{6, 2}, 10}, {{6, 3}, 11}, {{6, 4}, 12}, {{6, 5}, 13}
62  };
63 }
64 
65 
66 SVDDQMClustersOnTrackModule::~SVDDQMClustersOnTrackModule()
67 {
68 }
69 
70 //------------------------------------------------------------------
71 // Function to define histograms
72 //-----------------------------------------------------------------
73 
75 {
76  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
77  if (gTools->getNumberOfLayers() == 0) {
78  B2FATAL("Missing geometry for VXD, check steering file.");
79  }
80  if (gTools->getNumberOfSVDLayers() == 0) {
81  B2WARNING("Missing geometry for SVD, SVD-DQM is skipped.");
82  return;
83  }
84 
85  // Create a separate histogram directories and cd into it.
86  TDirectory* oldDir = gDirectory;
87  if (m_histogramDirectoryName != "") {
88  oldDir->mkdir(m_histogramDirectoryName.c_str());// do not use return value with ->cd(), its ZERO if dir already exists
89  oldDir->cd(m_histogramDirectoryName.c_str());
90  }
91 
92  int ChargeBins = 80;
93  float ChargeMax = 160;
94  int SNRBins = 50;
95  float SNRMax = 100;
96  int TimeBins = 300;
97  float TimeMin = -150;
98  float TimeMax = 150;
99 
100  int MaxBinBins = 6;
101  int MaxBinMax = 6;
102 
103  TString refFrame = "in FTSW reference";
104  if (m_desynchSVDTime)
105  refFrame = "in SVD reference";
106 
107 
108  //----------------------------------------------------------------
109  // Charge of clusters for L3/L456 sensors
110  //----------------------------------------------------------------
111  TString name = "SVDTRK_ClusterChargeU3";
112  TString title = "SVD U-Cluster-on-Track Charge for layer 3 sensors";
113  m_clsTrkChargeU3 = new TH1F(name.Data(), title.Data(), ChargeBins, 0, ChargeMax);
114  m_clsTrkChargeU3->GetXaxis()->SetTitle("cluster charge [ke-]");
115  m_clsTrkChargeU3->GetYaxis()->SetTitle("count");
117  name = "SVDTRK_ClusterChargeV3";
118  title = "SVD V-Cluster-on-Track Charge for layer 3 sensors";
119  m_clsTrkChargeV3 = new TH1F(name.Data(), title.Data(), ChargeBins, 0, ChargeMax);
120  m_clsTrkChargeV3->GetXaxis()->SetTitle("cluster charge [ke-]");
121  m_clsTrkChargeV3->GetYaxis()->SetTitle("count");
123 
124  name = "SVDTRK_ClusterChargeU456";
125  title = "SVD U-Cluster-on-Track Charge for layers 4,5,6 sensors";
126  m_clsTrkChargeU456 = new TH1F(name.Data(), title.Data(), ChargeBins, 0, ChargeMax);
127  m_clsTrkChargeU456->GetXaxis()->SetTitle("cluster charge [ke-]");
128  m_clsTrkChargeU456->GetYaxis()->SetTitle("count");
130 
131  name = "SVDTRK_ClusterChargeV456";
132  title = "SVD V-Cluster-on-Track Charge for layers 4,5,6 sensors";
133  m_clsTrkChargeV456 = new TH1F(name.Data(), title.Data(), ChargeBins, 0, ChargeMax);
134  m_clsTrkChargeV456->GetXaxis()->SetTitle("cluster charge [ke-]");
135  m_clsTrkChargeV456->GetYaxis()->SetTitle("count");
137 
138  m_clsTrkChargeL3 = new TH1F*[4];
139  m_clsTrkSNRL3 = new TH1F*[4];
140 
141  int ind = 0;
142  for (int ladder = 1; ladder <= 2; ++ladder) {
143  for (int sensor = 1; sensor <= 2; ++sensor) {
144 
145  name = Form("SVDTRK_ClusterCharge_L3.%d.%d", ladder, sensor);
146  title = Form("SVD Cluster-on-Track Charge for L3.%d.%d", ladder, sensor);
147  m_clsTrkChargeL3[ind] = new TH1F(name.Data(), title.Data(), ChargeBins, 0, ChargeMax);
148  m_clsTrkChargeL3[ind]->GetXaxis()->SetTitle("cluster charge [ke-]");
149  m_clsTrkChargeL3[ind]->GetYaxis()->SetTitle("count");
150  m_histoList->Add(m_clsTrkChargeL3[ind]);
151 
152  name = Form("SVDTRK_ClusterSNR_L3.%d.%d", ladder, sensor);
153  title = Form("SVD Cluster-on-Track SNR for L3.%d.%d", ladder, sensor);
154  m_clsTrkSNRL3[ind] = new TH1F(name.Data(), title.Data(), SNRBins, 0, SNRMax);
155  m_clsTrkSNRL3[ind]->GetXaxis()->SetTitle("cluster SNR");
156  m_clsTrkSNRL3[ind]->GetYaxis()->SetTitle("count");
157  m_histoList->Add(m_clsTrkSNRL3[ind]);
158  ind++;
159  }
160  }
161 
162  m_clsTrkCharge = new TH1F*[m_ladderMap.size()];
163  m_clsTrkSNR = new TH1F*[m_ladderMap.size()];
164 
165  for (const auto& it : m_ladderMap) {
166  std::pair<int, int> p = it.first;
167  int layer = p.first;
168  int sensor = p.second;
169  int idx = it.second;
170  name = Form("SVDTRK_ClusterCharge_L%d.x.%d", layer, sensor);
171  title = Form("SVD Cluster-on-Track Charge for L%d.x.%d", layer, sensor);
172  m_clsTrkCharge[idx] = new TH1F(name.Data(), title.Data(), ChargeBins, 0, ChargeMax);
173  m_clsTrkCharge[idx]->GetXaxis()->SetTitle("cluster charge [ke-]");
174  m_clsTrkCharge[idx]->GetYaxis()->SetTitle("count");
175  m_histoList->Add(m_clsTrkCharge[idx]);
176 
177  //printf("name %s layer %d sensor %d index %d\n", name.Data(), layer, sensor, idx);
178  name = Form("SVDTRK_ClusterSNR_L%d.x.%d", layer, sensor);
179  title = Form("SVD Cluster-on-Track SNR for L%d.x.%d", layer, sensor);
180  m_clsTrkSNR[idx] = new TH1F(name.Data(), title.Data(), SNRBins, 0, SNRMax);
181  m_clsTrkSNR[idx]->GetXaxis()->SetTitle("cluster SNR");
182  m_clsTrkSNR[idx]->GetYaxis()->SetTitle("count");
183  m_histoList->Add(m_clsTrkSNR[idx]);
184  }
185 
186  //----------------------------------------------------------------
187  // SNR of clusters for L3/L456 sensors
188  //----------------------------------------------------------------
189  name = "SVDTRK_ClusterSNRU3";
190  title = "SVD U-Cluster-on-Track SNR for layer 3 sensors";
191  m_clsTrkSNRU3 = new TH1F(name.Data(), title.Data(), SNRBins, 0, SNRMax);
192  m_clsTrkSNRU3->GetXaxis()->SetTitle("cluster SNR");
193  m_clsTrkSNRU3->GetYaxis()->SetTitle("count");
195  name = "SVDTRK_ClusterSNRV3";
196  title = "SVD V-Cluster-on-Track SNR for layer 3 sensors";
197  m_clsTrkSNRV3 = new TH1F(name.Data(), title.Data(), SNRBins, 0, SNRMax);
198  m_clsTrkSNRV3->GetXaxis()->SetTitle("cluster SNR");
199  m_clsTrkSNRV3->GetYaxis()->SetTitle("count");
201 
202  name = "SVDTRK_ClusterSNRU456";
203  title = "SVD U-Cluster-on-Track SNR for layers 4,5,6 sensors";
204  m_clsTrkSNRU456 = new TH1F(name.Data(), title.Data(), SNRBins, 0, SNRMax);
205  m_clsTrkSNRU456->GetXaxis()->SetTitle("cluster SNR");
206  m_clsTrkSNRU456->GetYaxis()->SetTitle("count");
208  name = "SVDTRK_ClusterSNRV456";
209  title = "SVD V-Cluster-on-Track SNR for layers 4,5,6 sensors";
210  m_clsTrkSNRV456 = new TH1F(name.Data(), title.Data(), SNRBins, 0, SNRMax);
211  m_clsTrkSNRV456->GetXaxis()->SetTitle("cluster SNR");
212  m_clsTrkSNRV456->GetYaxis()->SetTitle("count");
214 
215  //----------------------------------------------------------------
216  // Time of clusters for L3/L456 sensors
217  //----------------------------------------------------------------
218  name = "SVDTRK_ClusterTimeU3";
219  title = Form("SVD U-Cluster-on-Track Time %s for layer 3 sensors", refFrame.Data());
220  m_clsTrkTimeU3 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
221  m_clsTrkTimeU3->GetXaxis()->SetTitle("cluster time (ns)");
222  m_clsTrkTimeU3->GetYaxis()->SetTitle("count");
224  name = "SVDTRK_ClusterTimeV3";
225  title = Form("SVD V-Cluster-on-Track Time %s for layer 3 sensors", refFrame.Data());
226  m_clsTrkTimeV3 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
227  m_clsTrkTimeV3->GetXaxis()->SetTitle("cluster time (ns)");
228  m_clsTrkTimeV3->GetYaxis()->SetTitle("count");
230 
231  name = "SVDTRK_Cluster3TimeU3";
232  title = Form("SVD U-Cluster-on-Track Time %s for layer 3 sensors for 3 samples", refFrame.Data());
233  m_cls3TrkTimeU3 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
234  m_cls3TrkTimeU3->GetXaxis()->SetTitle("cluster time (ns)");
235  m_cls3TrkTimeU3->GetYaxis()->SetTitle("count");
237  name = "SVDTRK_Cluster3TimeV3";
238  title = Form("SVD V-Cluster-on-Track Time %s for layer 3 sensors for 3 samples", refFrame.Data());
239  m_cls3TrkTimeV3 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
240  m_cls3TrkTimeV3->GetXaxis()->SetTitle("cluster time (ns)");
241  m_cls3TrkTimeV3->GetYaxis()->SetTitle("count");
243 
244  name = "SVDTRK_Cluster6TimeU3";
245  title = Form("SVD U-Cluster-on-Track Time %s for layer 3 sensors for 6 samples", refFrame.Data());
246  m_cls6TrkTimeU3 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
247  m_cls6TrkTimeU3->GetXaxis()->SetTitle("cluster time (ns)");
248  m_cls6TrkTimeU3->GetYaxis()->SetTitle("count");
250  name = "SVDTRK_Cluster6TimeV3";
251  title = Form("SVD V-Cluster-on-Track Time %s for layer 3 sensors for 6 samples", refFrame.Data());
252  m_cls6TrkTimeV3 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
253  m_cls6TrkTimeV3->GetXaxis()->SetTitle("cluster time (ns)");
254  m_cls6TrkTimeV3->GetYaxis()->SetTitle("count");
256 
257 
258  name = "SVDTRK_ClusterTimeU456";
259  title = Form("SVD U-Cluster-on-Track Time %s for layers 4,5,6 sensors", refFrame.Data());
260  m_clsTrkTimeU456 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
261  m_clsTrkTimeU456->GetXaxis()->SetTitle("cluster time (ns)");
262  m_clsTrkTimeU456->GetYaxis()->SetTitle("count");
264  name = "SVDTRK_ClusterTimeV456";
265  title = Form("SVD V-Cluster-on-Track Time %s for layers 4,5,6 sensors", refFrame.Data());
266  m_clsTrkTimeV456 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
267  m_clsTrkTimeV456->GetXaxis()->SetTitle("cluster time (ns)");
268  m_clsTrkTimeV456->GetYaxis()->SetTitle("count");
270 
271  name = "SVDTRK_Cluster3TimeU456";
272  title = Form("SVD U-Cluster-on-Track Time %s for layers 4,5,6 sensors for 3 samples", refFrame.Data());
273  m_cls3TrkTimeU456 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
274  m_cls3TrkTimeU456->GetXaxis()->SetTitle("cluster time (ns)");
275  m_cls3TrkTimeU456->GetYaxis()->SetTitle("count");
277  name = "SVDTRK_Cluster3TimeV456";
278  title = Form("SVD V-Cluster-on-Track Time %s for layers 4,5,6 sensors for 3 samples", refFrame.Data());
279  m_cls3TrkTimeV456 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
280  m_cls3TrkTimeV456->GetXaxis()->SetTitle("cluster time (ns)");
281  m_cls3TrkTimeV456->GetYaxis()->SetTitle("count");
283 
284  name = "SVDTRK_Cluster6TimeU456";
285  title = Form("SVD U-Cluster-on-Track Time %s for layers 4,5,6 sensors for 6 samples", refFrame.Data());
286  m_cls6TrkTimeU456 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
287  m_cls6TrkTimeU456->GetXaxis()->SetTitle("cluster time (ns)");
288  m_cls6TrkTimeU456->GetYaxis()->SetTitle("count");
290  name = "SVDTRK_Cluster6TimeV456";
291  title = Form("SVD V-Cluster-on-Track Time %s for layers 4,5,6 sensors for 6 samples", refFrame.Data());
292  m_cls6TrkTimeV456 = new TH1F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax);
293  m_cls6TrkTimeV456->GetXaxis()->SetTitle("cluster time (ns)");
294  m_cls6TrkTimeV456->GetYaxis()->SetTitle("count");
296 
297  //----------------------------------------------------------------
298  // EventT0 vs Time of clusters for U and V sides
299  //----------------------------------------------------------------
300  name = "SVDTRK_ClusterTimeUvsEventT0";
301  title = Form("SVD U-Cluster-on-Track Time vs EventT0 %s for layer 3 sensors", refFrame.Data());
302  m_clsTrkTimeUEvtT0 = new TH2F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax, 100, -50, 50);
303  m_clsTrkTimeUEvtT0->GetXaxis()->SetTitle("clusters time (ns)");
304  m_clsTrkTimeUEvtT0->GetYaxis()->SetTitle("EventT0 (ns)");
306  name = "SVDTRK_ClusterTimeVvsEventT0";
307  title = Form("SVD V-Cluster-on-Track Time vs EventT0 %s for layer 3 sensors", refFrame.Data());
308  m_clsTrkTimeVEvtT0 = new TH2F(name.Data(), title.Data(), TimeBins, TimeMin, TimeMax, 100, -50, 50);
309  m_clsTrkTimeVEvtT0->GetXaxis()->SetTitle("cluster time (ns)");
310  m_clsTrkTimeVEvtT0->GetYaxis()->SetTitle("EventT0 (ns)");
312 
313  //----------------------------------------------------------------
314  // MaxBin of strips for all sensors (offline ZS)
315  //----------------------------------------------------------------
316  name = "SVDTRK_StripMaxBinUAll";
317  title = "SVD U-Strip-on-Track MaxBin for all sensors";
318  m_stripMaxBinUAll = new TH1F(name.Data(), title.Data(), MaxBinBins, 0, MaxBinMax);
319  m_stripMaxBinUAll->GetXaxis()->SetTitle("max bin");
320  m_stripMaxBinUAll->GetYaxis()->SetTitle("count");
322  name = "SVDTRK_StripMaxBinVAll";
323  title = "SVD V-Strip-on-Track MaxBin for all sensors";
324  m_stripMaxBinVAll = new TH1F(name.Data(), title.Data(), MaxBinBins, 0, MaxBinMax);
325  m_stripMaxBinVAll->GetXaxis()->SetTitle("max bin");
326  m_stripMaxBinVAll->GetYaxis()->SetTitle("count");
328 
329 
330  oldDir->cd();
331 }
332 
333 
335 {
336  // Register histograms (calls back defineHisto)
337  REG_HISTOGRAM
338 
339  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
340  if (gTools->getNumberOfSVDLayers() != 0) {
341 
343  m_eventT0.isOptional();
344  m_tracks.isOptional();
345  m_resultStoreObjectPointer.isOptional();
346 
347  }
348 }
349 
351 {
352  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
353  if (gTools->getNumberOfSVDLayers() == 0) return;
354 
355  // Add experiment and run number to the title of selected histograms (CR shifter plots)
356  StoreObjPtr<EventMetaData> evtMetaData;
357  m_expNumber = evtMetaData->getExperiment();
358  m_runNumber = evtMetaData->getRun();
359  TString runID = TString::Format(" ~ Exp%d Run%d", m_expNumber, m_runNumber);
360 
361  //reset histograms
362  TObject* obj;
363  TIter nextH(m_histoList);
364  while ((obj = nextH()))
365  if (obj->InheritsFrom("TH1")) {
366  ((TH1F*)obj)->Reset();
367 
368  TString tmp = (TString)obj->GetTitle();
369  Int_t pos = tmp.Last('~');
370  if (pos == -1) pos = tmp.Length() + 2;
371 
372  TString title = tmp(0, pos - 2);
373  ((TH1F*)obj)->SetTitle(title + runID);
374  }
375 }
376 
378 {
379 
380  if (!m_tracks.isValid()) {
381  B2WARNING("Missing Tracks StoreArray. Skipping SVDDQMClustersOnTrack");
382  return;
383  }
384 
385  if (!m_svdEventInfo.isValid())
386  m_tb = -1;
387  else {
388  if (m_tb != -1)
389  if (m_svdEventInfo->getModeByte().getTriggerBin() != m_tb)
390  return;
391  }
392 
393  int nSamples = 0;
394  if (m_svdEventInfo.isValid())
395  nSamples = m_svdEventInfo->getNSamples();
396  else
397  return;
398 
399  // get EventT0 if present and valid
400  double eventT0 = -1000;
401  if (m_eventT0.isOptional())
402  if (m_eventT0.isValid())
403  if (m_eventT0->hasEventT0())
404  eventT0 = m_eventT0->getEventT0();
405 
406  // if svd time in SVD time reference is shown, eventT0 is also synchronized with SVD reference frame, firstFrame = 0
407  if (m_desynchSVDTime && m_svdEventInfo.isValid())
408  eventT0 = eventT0 - m_svdEventInfo->getSVD2FTSWTimeShift(0);
409 
410  //check HLT decision and increase number of events only if the event has been accepted
411 
413  const bool eventAccepted = FinalTriggerDecisionCalculator::getFinalTriggerDecision(*m_resultStoreObjectPointer);
414  if (!eventAccepted) return;
415  }
416 
417  auto gTools = VXD::GeoCache::getInstance().getGeoTools();
418  if (gTools->getNumberOfSVDLayers() == 0) return;
419 
420 
421  for (const Track& track : m_tracks) {
422 
423  const TrackFitResult* tfr = track.getTrackFitResultWithClosestMass(Const::pion);
424  if (not tfr) continue;
425 
426  const RecoTrack* recoTrack = track.getRelated<RecoTrack>();
427  if (not recoTrack) continue;
428 
429  for (const SVDCluster& svdCluster : recoTrack->getRelationsWith<SVDCluster>(m_svdClustersName)) {
430 
431  int iLayer = svdCluster.getSensorID().getLayerNumber();
432  int iLadder = svdCluster.getSensorID().getLadderNumber();
433  int iSensor = svdCluster.getSensorID().getSensorNumber();
434 
435  if (iLayer == 3) {
436  int ind = -1;
437  if (iLadder == 1 && iSensor == 1) {
438  ind = 0;
439  } else if (iLadder == 1 && iSensor == 2) {
440  ind = 1;
441  } else if (iLadder == 2 && iSensor == 1) {
442  ind = 2;
443  } else if (iLadder == 2 && iSensor == 2) {
444  ind = 3;
445  }
446 
447  if (ind != -1) {
448  if (m_clsTrkChargeL3[ind] != nullptr) m_clsTrkChargeL3[ind]->Fill(svdCluster.getCharge() / 1000.0); // in kelectrons
449  if (m_clsTrkSNRL3[ind] != nullptr) m_clsTrkSNRL3[ind]->Fill(svdCluster.getSNR());
450  }
451 
452  }
453 
454  std::pair<int, int> p(iLayer, iSensor);
455  int idx = m_ladderMap[p];
456 
457  if (m_clsTrkCharge[idx] != nullptr) m_clsTrkCharge[idx]->Fill(svdCluster.getCharge() / 1000.0); // in kelectrons
458  if (m_clsTrkSNR[idx] != nullptr) m_clsTrkSNR[idx]->Fill(svdCluster.getSNR());
459 
460  float time = svdCluster.getClsTime();
461  if (m_desynchSVDTime && m_svdEventInfo.isValid())
462  time = time - m_svdEventInfo->getSVD2FTSWTimeShift(svdCluster.getFirstFrame());
463 
464  if (svdCluster.isUCluster()) {
465 
466  m_clsTrkTimeUEvtT0->Fill(time, eventT0);
467 
468  if (iLayer == 3) {
469  if (m_clsTrkChargeU3 != nullptr) m_clsTrkChargeU3->Fill(svdCluster.getCharge() / 1000.0); // in kelectrons
470  if (m_clsTrkSNRU3 != nullptr) m_clsTrkSNRU3->Fill(svdCluster.getSNR());
471  if (m_clsTrkTimeU3 != nullptr) m_clsTrkTimeU3->Fill(time);
472 
473  if (nSamples == 3) {
474  if (m_cls3TrkTimeU3 != nullptr) m_cls3TrkTimeU3->Fill(time);
475  } else {
476  if (m_cls6TrkTimeU3 != nullptr) m_cls6TrkTimeU3->Fill(time);
477  }
478  } else {
479  if (m_clsTrkChargeU456 != nullptr) m_clsTrkChargeU456->Fill(svdCluster.getCharge() / 1000.0); // in kelectrons
480  if (m_clsTrkSNRU456 != nullptr) m_clsTrkSNRU456->Fill(svdCluster.getSNR());
481  if (m_clsTrkTimeU456 != nullptr) m_clsTrkTimeU456->Fill(time);
482 
483  if (nSamples == 3) {
484  if (m_cls3TrkTimeU456 != nullptr) m_cls3TrkTimeU456->Fill(time);
485  } else {
486  if (m_cls6TrkTimeU456 != nullptr) m_cls6TrkTimeU456->Fill(time);
487  }
488  }
489 
490  for (const SVDRecoDigit& recoDigit : svdCluster.getRelationsTo<SVDRecoDigit>(m_svdRecoDigitsName)) {
491 
493  if (m_stripMaxBinUAll != nullptr and shaper != nullptr) m_stripMaxBinUAll->Fill(shaper->getMaxTimeBin());
494  }
495 
496 
497  } else {
498 
499  m_clsTrkTimeVEvtT0->Fill(time, eventT0);
500 
501  if (iLayer == 3) {
502  if (m_clsTrkChargeV3 != nullptr) m_clsTrkChargeV3->Fill(svdCluster.getCharge() / 1000.0); // in kelectrons
503  if (m_clsTrkSNRV3 != nullptr) m_clsTrkSNRV3->Fill(svdCluster.getSNR());
504  if (m_clsTrkTimeV3 != nullptr) m_clsTrkTimeV3->Fill(time);
505  if (nSamples == 3) {
506  if (m_cls3TrkTimeV3 != nullptr) m_cls3TrkTimeV3->Fill(time);
507  } else {
508  if (m_cls6TrkTimeV3 != nullptr) m_cls6TrkTimeV3->Fill(time);
509  }
510  } else {
511  if (m_clsTrkChargeV456 != nullptr) m_clsTrkChargeV456->Fill(svdCluster.getCharge() / 1000.0); // in kelectrons
512  if (m_clsTrkSNRV456 != nullptr) m_clsTrkSNRV456->Fill(svdCluster.getSNR());
513  if (m_clsTrkTimeV456 != nullptr) m_clsTrkTimeV456->Fill(time);
514  if (nSamples == 3) {
515  if (m_cls3TrkTimeV456 != nullptr) m_cls3TrkTimeV456->Fill(time);
516  } else {
517  if (m_cls6TrkTimeV456 != nullptr) m_cls6TrkTimeV456->Fill(time);
518  }
519  }
520 
521  for (const SVDRecoDigit& recoDigit : svdCluster.getRelationsTo<SVDRecoDigit>(m_svdRecoDigitsName)) {
522 
524  if (m_stripMaxBinVAll != nullptr and shaper != nullptr) m_stripMaxBinVAll->Fill(shaper->getMaxTimeBin());
525  }
526 
527  }
528 
529  }
530  }
531 }
532 
533 
535 {
536 
537  delete m_histoList;
538 
539 }
static const ChargedStable pion
charged pion particle
Definition: Const.h:652
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Definition: HistoModule.h:29
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
void setPropertyFlags(unsigned int propertyFlags)
Sets the flags for the module properties.
Definition: Module.cc:208
@ c_ParallelProcessingCertified
This module can be run in parallel processing mode safely (All I/O must be done through the data stor...
Definition: Module.h:80
This is the Reconstruction Event-Data Model Track.
Definition: RecoTrack.h:79
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
RelationVector< T > getRelationsWith(const std::string &name="", const std::string &namedRelation="") const
Get the relations between this object and another store array.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
The SVD Cluster class This class stores all information about reconstructed SVD clusters.
Definition: SVDCluster.h:29
TH1F * m_cls6TrkTimeV456
v Time of clusters related to tracks for layer 4,5,6 sensors for 6 samples
TH1F ** m_clsTrkSNRL3
SNR of clusters related to tracks per layer 3.
TH1F * m_cls6TrkTimeU3
u Time of clusters related to tracks for layer 3 sensors for 3 samples
std::string m_svdClustersName
SVDClusters data object name.
void initialize() override final
Module function initialize.
TH1F ** m_clsTrkChargeL3
charge of clusters related to tracks per layer 3
TH1F ** m_clsTrkCharge
charge of clusters related to tracks per ladder
TH1F * m_clsTrkSNRV456
v SNR of clusters related to tracks for layer 4,5,6 sensors
std::string m_svdShaperDigitsName
SVDShaperDigits data object name.
StoreObjPtr< EventT0 > m_eventT0
EventT0 data object.
bool m_skipRejectedEvents
if true skip events rejected by HLT (default)
TH1F * m_cls3TrkTimeU456
u Time of clusters related to tracks for layer 4,5,6 sensors for 3 samples
TH1F * m_clsTrkTimeU3
u Time of clusters related to tracks for layer 3 sensors
TH1F * m_cls3TrkTimeV456
v Time of clusters related to tracks for layer 4,5,6 sensors for 3 samples
StoreObjPtr< SVDEventInfo > m_svdEventInfo
SVDEventInfo data object.
TH1F * m_clsTrkChargeV3
v charge of clusters related to tracks for layer 3 sensors
void defineHisto() override final
Contains the Histogram definitions
TH1F * m_clsTrkSNRU456
u SNR of clusters related to tracks for layer 4,5,6 sensors
TH1F * m_clsTrkSNRV3
v SNR of clusters related to tracks for layer 3 sensors
TH1F * m_cls6TrkTimeU456
u Time of clusters related to tracks for layer 4,5,6 sensors for 6 samples
TH2F * m_clsTrkTimeUEvtT0
u Time of clusters related to tracks vs EventT0
void terminate() override final
Module function terminate.
TH1F * m_clsTrkChargeU456
u charge of clusters related to tracks for layer 4,5,6 sensors
TH1F * m_cls3TrkTimeV3
v Time of clusters related to tracks for layer 3 sensors for 3 sampes
void event() override final
Module function event.
TH1F * m_clsTrkTimeV3
v Time of clusters related to tracks for layer 3 sensors
std::string m_histogramDirectoryName
Name of the histogram directory in ROOT file.
StoreArray< Track > m_tracks
StoreArray of the Tracks.
int m_tb
choose one trigger bin, or none if the value is -1
TH1F * m_stripMaxBinUAll
u MaxBin of strips related to tracks for all sensors
TList * m_histoList
list of cumulative histograms
std::string m_svdRecoDigitsName
SVDRecoDigits data object name.
std::string m_svdEventInfoName
SVDEventInfo data object name.
TH1F * m_cls6TrkTimeV3
v Time of clusters related to tracks for layer 3 sensors for 3 sampes
TH1F * m_cls3TrkTimeU3
u Time of clusters related to tracks for layer 3 sensors for 3 samples
void beginRun() override final
Module function beginRun.
TH1F * m_clsTrkSNRU3
u SNR of clusters related to tracks for layer 3 sensors
bool m_desynchSVDTime
if TRUE: svdTime back in SVD time reference
TH1F * m_clsTrkChargeU3
u charge of clusters related to tracks for layer 3 sensors
TH1F * m_stripMaxBinVAll
v MaxBin of strips related to tracks for all sensors
StoreObjPtr< SoftwareTriggerResult > m_resultStoreObjectPointer
Store Object for reading the trigger decision.
TH1F * m_clsTrkTimeU456
u Time of clusters related to tracks for layer 4,5,6 sensors
TH1F ** m_clsTrkSNR
SNR of clusters related to tracks per ladder.
TH2F * m_clsTrkTimeVEvtT0
v Time of clusters related to tracks vs EventT0
TH1F * m_clsTrkChargeV456
v charge of clusters related to tracks for layer 4,5,6 sensors
TH1F * m_clsTrkTimeV456
v Time of clusters related to tracks for layer 4,5,6 sensors
std::map< std::pair< int, int >, int > m_ladderMap
map of ladder index
The SVD RecoDigit class.
Definition: SVDRecoDigit.h:43
The SVD ShaperDigit class.
int getMaxTimeBin() const
Get the max bin.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
Values of the result of a track fit with a given particle hypothesis.
Class that bundles various TrackFitResults.
Definition: Track.h:25
static GeoCache & getInstance()
Return a reference to the singleton instance.
Definition: GeoCache.cc:214
const GeoTools * getGeoTools()
Return a raw pointer to a GeoTools object.
Definition: GeoCache.h:147
void addParam(const std::string &name, T &paramVariable, const std::string &description, const T &defaultValue)
Adds a new parameter to the module.
Definition: Module.h:560
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:650
Abstract base class for different kinds of events.