Belle II Software  release-08-01-10
BKLMDigitAnalyzerModule.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 /* Own header. */
10 #include <klm/bklm/modules/bklmDigitAnalyzer/BKLMDigitAnalyzerModule.h>
11 
12 /* Basf2 headers. */
13 #include <framework/datastore/StoreObjPtr.h>
14 #include <framework/dataobjects/EventMetaData.h>
15 
16 /* ROOT headers. */
17 #include <TGaxis.h>
18 
19 using namespace Belle2;
20 
21 
22 REG_MODULE(BKLMDigitAnalyzer);
23 
24 
26  m_runNumber(0),
27  m_outputRootFile(nullptr),
28  m_histoList(nullptr)
29 {
30  for (int i = 0; i < 2; ++i) {
31  m_histoLayerVsSector[i] = nullptr;
32  m_histoLayerVsSectorPerPlane[i][0] = nullptr;
33  m_histoLayerVsSectorPerPlane[i][1] = nullptr;
34  for (int s = 0; s < 8; ++s) {
35  m_histoLayer[i][s] = nullptr;
36  for (int p = 0; p < 2; ++p) {
37  m_histoChannel[i][s][p] = nullptr;
38  m_histoStrip[i][s][p] = nullptr;
39  m_histoTdc[i][s][p] = nullptr;
40  m_histoCTimeDiff[i][s][p] = nullptr;
41  }
42  }
43  }
44 
45  // Set module properties
46  setDescription("Module useful to quickly analyze BKLM unpacked data.");
47 
48  // Parameter definitions
49  addParam("outputRootName", m_outputRootName, "Name of output .root file (without .root!)", std::string("bklmHitmap"));
50 
51 }
52 
54 {
55 }
56 
58 {
59  m_digit.isRequired();
60  m_digitRaw.isRequired();
61  m_digitOutOfRange.isRequired("KLMDigitsOutOfRange");
62  m_digitEventInfo.isRequired();
63 }
64 
66 {
67  time_t rawTime;
68  time(&rawTime);
69  struct tm* tm = gmtime(&rawTime);
70 
71  StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
72  m_runNumber = eventMetaData->getRun();
73 
74  TString runNumberTString(toString(m_runNumber).c_str());
75 
76  TString outputRootNameTString(m_outputRootName);
77  outputRootNameTString += "_run" + runNumberTString + ".root";
78 
79  const TString label[2] = {"BF", "BB"};
80 
81  m_outputRootFile = new TFile(outputRootNameTString, "RECREATE");
82  B2INFO("BKLMDigitAnalyzer:: the output file '" << outputRootNameTString.Data() << "' will be created for run " << m_runNumber);
83 
84  int exp = -1;
85  int run = -1;
86  int year = -1;
87  int month = -1;
88  int day = -1;
89  int hour = -1;
90  int min = -1;
91  int sec = -1;
92  m_extraInfo = new TTree("extraInfo", "Extra informations");
93  m_extraInfo->Branch("exp", &exp, "exp/I");
94  m_extraInfo->Branch("run", &run, "run/I");
95  m_extraInfo->Branch("year", &year, "year/I");
96  m_extraInfo->Branch("month", &month, "month/I");
97  m_extraInfo->Branch("day", &day, "day/I");
98  m_extraInfo->Branch("hour", &hour, "hour/I");
99  m_extraInfo->Branch("min", &min, "min/I");
100  m_extraInfo->Branch("sec", &sec, "sec/I");
101  exp = eventMetaData->getExperiment();
102  run = m_runNumber;
103  year = tm->tm_year + 1900;
104  month = tm->tm_mon + 1;
105  day = tm->tm_mday;
106  hour = tm->tm_hour;
107  min = tm->tm_min;
108  sec = tm->tm_sec;
109  m_extraInfo->Fill();
110 
111  m_histoList = new TList;
112 
113  for (int fb = 0; fb < 2; fb++) {
114 
115  m_histoLayerVsSector[fb] = createTH2("SectLay" + label[fb], label[fb] + " Hitmap -- run" + runNumberTString, 31, -0.5, 15.,
116  "Layer (0-based)", 17, -0.5, 8., "BF Sector", 0, m_histoList);
117 
118  for (int isRPCorPhi = 0; isRPCorPhi < 2; isRPCorPhi++) {
119 
120  if (isRPCorPhi == 0)
121  m_histoLayerVsSectorPerPlane[fb][isRPCorPhi] = createTH2("SectLayPlaneZ" + label[fb],
122  label[fb] + " Hitmap of plane z -- run" + runNumberTString, 31, -0.5, 15., "Layer (0-based)", 17, -0.5, 8., "BF Sector", 0,
123  m_histoList);
124  else
125  m_histoLayerVsSectorPerPlane[fb][isRPCorPhi] = createTH2("SectLayPlanePhi" + label[fb],
126  label[fb] + " Hitmap of plane phi -- run" + runNumberTString, 31, -0.5, 15., "Layer (0-based)", 17, -0.5, 8., "BF Sector", 0,
127  m_histoList);
128 
129  for (int iSector = 0; iSector < 8; iSector++) {
130 
131  TString iSectorTString(toString(iSector).c_str());
132  TString nameSector = label[fb] + iSectorTString;
133 
134  if (isRPCorPhi == 0) {
135 
136  // Create the histogram belowo only one time
137  m_histoLayer[fb][iSector] = createTH1("Layer" + nameSector, nameSector + " Layer -- run" + runNumberTString, 31, -0.5, 15.,
138  "Layer (0-based)", "Counts", 1, m_histoList);
139 
140  m_histoChannel[fb][iSector][isRPCorPhi] = createTH2("PlaneZ" + nameSector,
141  nameSector + " Plane z, electronic channels -- run" + runNumberTString, 31,
142  -0.5, 15., "Layer (0-based)", 130, -0.5, 64.5, "Channel", 1, m_histoList);
143 
144  m_histoStrip[fb][iSector][isRPCorPhi] = createTH2("PlaneZStrip" + nameSector,
145  nameSector + " Plane z, strips -- run" + runNumberTString,
146  31, -0.5, 15., "Layer (0-based)", 130, -0.5, 64.5, "Strip", 1, m_histoList);
147 
148  m_histoTdc[fb][iSector][isRPCorPhi] = createTH1("SciTdc" + nameSector,
149  nameSector + " TDC (Scintillators) -- run" + runNumberTString, 60, 0, 30, "TDC", "Counts", 1, m_histoList);
150 
151  m_histoCTimeDiff[fb][iSector][isRPCorPhi] = createTH1("SciCTimeDiff" + nameSector,
152  nameSector + " CTime diff. (Scintillators) -- run" + runNumberTString, 200, -200, -100, "Sci_CTime - Trg_CTime [ns]", "Counts", 1,
153  m_histoList);
154 
155  } else {
156 
157  m_histoChannel[fb][iSector][isRPCorPhi] = createTH2("PlanePhi" + nameSector,
158  nameSector + " Plane phi, electronic channels -- run" + runNumberTString,
159  31, -0.5, 15., "Layer (0-based)", 130, -0.5, 64.5, "Channel", 1, m_histoList);
160 
161  m_histoStrip[fb][iSector][isRPCorPhi] = createTH2("PlanePhiStrip" + nameSector,
162  nameSector + " Plane phi, strips -- run" + runNumberTString,
163  31, -0.5, 15., "Layer (0-based)", 130, -0.5, 64.5, "Strip", 1, m_histoList);
164 
165  m_histoTdc[fb][iSector][isRPCorPhi] = createTH1("RPCTdc" + nameSector, nameSector + " TDC (RPCs) -- run" + runNumberTString, 230,
166  -100, 2200, "TDC", "Counts", 1, m_histoList);
167 
168  m_histoCTimeDiff[fb][iSector][isRPCorPhi] = createTH1("RPCCTimeDiff" + nameSector,
169  nameSector + " CTime diff. (RPCs) -- run" + runNumberTString, 350, -3000, 500, "RPC_CTime - Trg_CTime [ns]", "Counts", 1,
170  m_histoList);
171 
172  }
173  }
174  }
175  }
176 }
177 
179 {
180  for (int i = 0; i < m_digitEventInfo.getEntries(); i++) {
181 
182  KLMDigitEventInfo* digitEventInfo = m_digitEventInfo[i];
183  // Some warnings (they should never appear, but it's better to be sure)
184  if ((digitEventInfo->getRPCHits() + digitEventInfo->getSciHits()) != (int)digitEventInfo->getRelationsFrom<KLMDigit>().size())
185  B2WARNING("BKLMDigitAnalyzer:: the total number of BKLMDigit differs from the sum of RPC and scintillator hits stored in BKLMEventDigitDebug!");
186  if (digitEventInfo->getOutOfRangeHits() != (int)digitEventInfo->getRelationsTo<KLMDigit>("BKLMDigitsOutOfRange").size())
187  B2WARNING("BKLMDigitAnalyzer:: the total number of BKLMDigit differs from the number of outOfRange-flagged hits stored in BKLMEventDigitDebug!");
188 
189  for (const KLMDigit& digit : digitEventInfo->getRelationsFrom<KLMDigit>()) {
190  if (digit.getSubdetector() != KLMElementNumbers::c_BKLM)
191  continue;
192 
193  KLMDigitRaw* digitRaw = digit.getRelatedTo<KLMDigitRaw>();
194 
195  m_histoLayerVsSector[1 - digit.getSection()]->Fill(digit.getLayer() - 1, digit.getSector() - 1);
196 
197  m_histoLayerVsSectorPerPlane[1 - digit.getSection()][digit.isPhiReadout()]->Fill(digit.getLayer() - 1, digit.getSector() - 1);
198 
199  m_histoLayer[1 - digit.getSection()][digit.getSector() - 1]->Fill(digit.getLayer() - 1);
200 
201  m_histoChannel[1 - digit.getSection()][digit.getSector() - 1][digit.isPhiReadout()]->Fill(digit.getLayer() - 1,
202  digitRaw->getChannel() - 1);
203 
204  m_histoStrip[1 - digit.getSection()][digit.getSector() - 1][digit.isPhiReadout()]->Fill(digit.getLayer() - 1,
205  digit.getStrip() - 1);
206 
207  // getTime() retruns the TDC
208  m_histoTdc[1 - digit.getSection()][digit.getSector() - 1][digit.inRPC()]->Fill(digit.getTime());
209 
210  m_histoCTimeDiff[1 - digit.getSection()][digit.getSector() - 1][digit.inRPC()]->Fill(digit.getCTime() -
211  digitEventInfo->getIntTriggerCTime());
212 
213  }
214  }
215 }
216 
218 {
219  StoreObjPtr<EventMetaData> eventMetaData("EventMetaData", DataStore::c_Event);
220 
221  // Save the .root file
222  if (m_outputRootFile != nullptr) {
223  m_outputRootFile->cd();
224 
225  m_extraInfo->Write();
226 
227  TIter nextHisto(m_histoList);
228  TObject* obj;
229  while ((obj = nextHisto()))
230  obj->Write("", TObject::kWriteDelete);
231  m_outputRootFile->Close();
232  }
233 }
234 
236 {
237 }
238 
239 TH1F* BKLMDigitAnalyzerModule::createTH1(const char* name, const char* title, Int_t nbinsX, Double_t minX, Double_t maxX,
240  const char* titleX, const char* titleY, bool drawStat, TList* histoList)
241 {
242  TH1F* hist = new TH1F(name, title, nbinsX, minX, maxX);
243  hist->GetXaxis()->SetTitle(titleX);
244  hist->GetYaxis()->SetTitle(titleY);
245  TGaxis::SetMaxDigits(3);
246  hist->SetLineWidth(1);
247  hist->SetLineColor(kRed);
248  if (!drawStat)
249  hist->SetStats(0);
250  if (histoList)
251  histoList->Add(hist);
252  return hist;
253 }
254 
255 TH2F* BKLMDigitAnalyzerModule::createTH2(const char* name, const char* title, Int_t nbinsX, Double_t minX, Double_t maxX,
256  const char* titleX, Int_t nbinsY, Double_t minY, Double_t maxY, const char* titleY, bool drawStat, TList* histoList)
257 {
258  TH2F* hist = new TH2F(name, title, nbinsX, minX, maxX, nbinsY, minY, maxY);
259  hist->GetXaxis()->SetTitle(titleX);
260  hist->GetYaxis()->SetTitle(titleY);
261  TGaxis::SetMaxDigits(3);
262  if (!drawStat)
263  hist->SetStats(0);
264  if (histoList)
265  histoList->Add(hist);
266  return hist;
267 }
int m_runNumber
Run number of the current data set.
TH2F * m_histoLayerVsSector[2]
Pointer to occupancy 2D histograms (numerator, denominator) of layer vs sector.
void initialize() override
Initializer.
void event() override
This method is called for each event.
void endRun() override
This method is called if the current run ends.
TFile * m_outputRootFile
Pointer to ROOT output file.
void terminate() override
This method is called at the end of the event processing.
StoreArray< KLMDigitEventInfo > m_digitEventInfo
Output data array of analyzed KLMDigit information.
TH2F * createTH2(const char *name, const char *title, Int_t nBinsX, Double_t minX, Double_t maxX, const char *titleX, Int_t nBinsY, Double_t minY, Double_t maxY, const char *titleY, bool drawStat, TList *histoList=nullptr)
Create a ROOT 2D histogram.
TH1F * m_histoLayer[2][8]
Pointer to occupancy 1D histogram of each layer for each sector.
StoreArray< KLMDigit > m_digitOutOfRange
Input BKLMDigits whose time is out of range.
TH1F * createTH1(const char *name, const char *title, Int_t nBinsX, Double_t minX, Double_t maxX, const char *titleX, const char *titleY, bool drawStat, TList *histoList=nullptr)
Create a ROOT 1D histogram.
void beginRun() override
Called when entering a new run.
StoreArray< KLMDigit > m_digit
Input KLMDigits.
~BKLMDigitAnalyzerModule() override
Destructor.
TH2F * m_histoChannel[2][8][2]
Pointer to occupancy 2D histogram of each channel.
TList * m_histoList
Pointer to ROOT list of histograms.
TH2F * m_histoLayerVsSectorPerPlane[2][2]
Pointer to occupancy 2D histogram of layer vs sector for each view (=plane)
TH1F * m_histoTdc[2][8][2]
Pointer to TDC 1D histogram of each channel for each sector and view.
TH2F * m_histoStrip[2][8][2]
Pointer to occupancy 2D histogram of each strip.
std::string toString(T val)
Convert a number of type T into a string.
std::string m_outputRootName
Output filename.
TH1F * m_histoCTimeDiff[2][8][2]
Pointer to TDC-difference 1D histogram of each channel for each sector and view.
TTree * m_extraInfo
Pointer to ROOT tree with extra info.
StoreArray< KLMDigitRaw > m_digitRaw
Input raw BLMDigits.
@ c_Event
Different object in each event, all objects/arrays are invalidated after event() function has been ca...
Definition: DataStore.h:59
Class to store debugging informations from the unpacker (event based).
int getRPCHits() const
Returns the number of RPC hits in the event.
int getIntTriggerCTime() const
Returns trigger CTIME as int.
int getOutOfRangeHits() const
Returns the number of OutOfRange-flagged hits in the event.
int getSciHits() const
Returns the number of scintillator hits in the event.
Class to store the raw words from the unpacker, digit-by-digit.
Definition: KLMDigitRaw.h:29
uint16_t getChannel() const
Get the channel number.
Definition: KLMDigitRaw.h:111
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:29
Base class for Modules.
Definition: Module.h:72
void setDescription(const std::string &description)
Sets the description of the module.
Definition: Module.cc:214
RelationVector< FROM > getRelationsFrom(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from another store array to this object.
RelationVector< TO > getRelationsTo(const std::string &name="", const std::string &namedRelation="") const
Get the relations that point from this object to another store array.
TO * getRelatedTo(const std::string &name="", const std::string &namedRelation="") const
Get the object to which this object has a relation.
Type-safe access to single objects in the data store.
Definition: StoreObjPtr.h:96
REG_MODULE(arichBtest)
Register the Module.
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
Abstract base class for different kinds of events.