9 #include <framework/core/HistoModule.h>
10 #include <top/modules/TOPTBCComparator/TOPTBCComparatorModule.h>
11 #include <top/geometry/TOPGeometryPar.h>
12 #include <top/geometry/FrontEndMapper.h>
13 #include "TDirectory.h"
14 #include "TSystemDirectory.h"
15 #include "TSystemFile.h"
18 #include <boost/format.hpp>
34 setDescription(
"TOP Time Base Correction monitor module ");
36 addParam(
"inputDirectorList", m_inputDirectoryList,
37 "List of the directories (one per IOV) in which the files with the calibration constants of the SCODS are stored. The format of each line must be: folderpath/ label ",
39 addParam(
"compareToPreviousSet", m_compareToPreviousSet,
40 "If true, each CalSet is compared to the previous one in the order determined by the inputDirectorList file. If false, the reference CalSet is the first of the list",
42 addParam(
"outputFile", m_outputFile,
43 "output file containing the monitoring plots",
string(
"TBCComparisonResults.root"));
44 addParam(
"minCalPulses", m_minCalPulses,
45 "minimum number of calibration pulses need to flag a sample as non-empty",
short(200));
46 addParam(
"numSamples", m_numSamples,
47 "number of samples that have been calibrated",
short(256));
52 void TOPTBCComparatorModule::defineHisto()
55 ifstream inputDirectoryListFile(m_inputDirectoryList.c_str());
58 if (!inputDirectoryListFile) {
59 B2ERROR(
"Unable to open the input file with the list of CalSets to analyze");
63 B2INFO(
"Initializing histograms");
66 std::string inputString;
70 while (std::getline(inputDirectoryListFile, inputString)) {
73 parseInputDirectoryLine(inputString);
75 B2INFO(
"Initializing histograms for Calibration set located at " << m_calSetDirectory <<
" with label " << m_calSetLabel);
88 name = str(format(
"TOPAverageDeltaT_CalSet_%1%") % m_calSetLabel);
89 title = str(format(
"Average value of #Delta T VS global channel number. CalSet %1%") % m_calSetLabel);
90 m_topAverageDeltaT.push_back(
new TH1F(name.c_str(), title.c_str(), 512 * 16, 0., 512 * 16.));
93 name = str(format(
"TOPSigmaDeltaT_CalSet_%1%") % m_calSetLabel);
94 title = str(format(
"Sigma value of #Delta T VS global channel number. CalSet %1%") % m_calSetLabel);
95 m_topSigmaDeltaT.push_back(
new TH1F(name.c_str(), title.c_str(), 512 * 16, 0., 512 * 16.));
98 name = str(format(
"TOPSampleOccupancy_CalSet_%1%") % m_calSetLabel);
99 title = str(format(
"Average number of calpulses per sample VS global channel number. CalSet %1%") % m_calSetLabel);
100 m_topSampleOccupancy.push_back(
new TH1F(name.c_str(), title.c_str(), 512 * 16, 0., 512 * 16.));
104 name = str(format(
"TOPAverageDeltaTComparison_CalSet_%1%") % m_calSetLabel);
105 title = str(format(
"Ratio of the average #Delta T in CalSet %1% over the previous one, over the whole detector") % m_calSetLabel);
106 m_topAverageDeltaTComparison.push_back(
new TH1F(name.c_str(), title.c_str(), 512 * 16, 0., 512 * 16.));
109 name = str(format(
"TOPSigmaDeltaTComparison_CalSet_%1%") % m_calSetLabel);
110 title = str(format(
"Ratio of the st. dev. of #Delta T in CalSet %1% over the previous one, , over the whole detector") %
112 m_topSigmaDeltaTComparison.push_back(
new TH1F(name.c_str(), title.c_str(), 512 * 16, 0., 512 * 16.));
115 name = str(format(
"TOPSampleOccupancyComparison_CalSet_%1%") % m_calSetLabel);
117 format(
"Ratio of the average number of calpulses per sample in CalSet %1% over the previous one, over the whole detector") %
119 m_topSampleOccupancyComparison.push_back(
new TH1F(name.c_str(), title.c_str(), 512 * 16, 0., 512 * 16.));
126 for (
int iSlot = 0; iSlot < 16; iSlot ++) {
129 name = str(format(
"SlotAverageDeltaT_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
130 title = str(format(
"Average value of #Delta T VS Channel number. Slot %1%, CalSet %2%") % (iSlot) % m_calSetLabel);
131 m_slotAverageDeltaT[iSlot].push_back(
new TH1F(name.c_str(), title.c_str(), 512, 0., 512.));
133 name = str(format(
"SlotSigmaDeltaT_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
134 title = str(format(
"Standard deviation of #Delta T VS Channel number. Slot %1%, CalSet %2%") % (iSlot) % m_calSetLabel);
135 m_slotSigmaDeltaT[iSlot].push_back(
new TH1F(name.c_str(), title.c_str(), 512, 0., 512.));
137 name = str(format(
"SlotAverageDeltaTMap_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
138 title = str(format(
"Map of the average value of #Delta T on the 256 samples. Slot %1%, CalSet %2%") % (iSlot) % m_calSetLabel);
139 m_slotAverageDeltaTMap[iSlot].push_back(
new TH2F(name.c_str(), title.c_str(), 64, 0., 64., 8, 0., 8.));
141 name = str(format(
"SlotSigmaDeltaTMap_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
142 title = str(format(
"Map of the RMS of #Delta T on the 256 samples . Slot %1%, CalSet %2%") % (iSlot) % m_calSetLabel);
143 m_slotSigmaDeltaTMap[iSlot].push_back(
new TH2F(name.c_str(), title.c_str(), 64, 0., 64., 8, 0., 8.));
147 name = str(format(
"SlotSampleOccupancy_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
148 title = str(format(
"Average occupancy per sample. Slot %1%, CalSet %2%") % (iSlot) % m_calSetLabel);
149 m_slotSampleOccupancy[iSlot].push_back(
new TH1F(name.c_str(), title.c_str(), 512, 0., 512.));
151 name = str(format(
"SlotEmptySamples_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
152 title = str(format(
"Number samples with less than %3% calpulses per each slot channel. Slot %1%, CalSet %2%") %
153 (iSlot) % m_calSetLabel % m_minCalPulses);
154 m_slotEmptySamples[iSlot].push_back(
new TH1F(name.c_str(), title.c_str(), 512, 0., 512.));
156 name = str(format(
"SlotSampleOccupancyMap_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
157 title = str(format(
"Map of the average occupancy per sample. Slot %1%, CalSet %2%") %
158 (iSlot) % m_calSetLabel);
159 m_slotSampleOccupancyMap[iSlot].push_back(
new TH2F(name.c_str(), title.c_str(), 64, 0., 64., 8, 0., 8.));
161 name = str(format(
"SlotEmptySamplesMap_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
162 title = str(format(
"Map of the number samples with less than %3% calpulses per each. Slot %1%, CalSet %2%") %
163 (iSlot) % m_calSetLabel % m_minCalPulses);
164 m_slotEmptySamplesMap[iSlot].push_back(
new TH2F(name.c_str(), title.c_str(), 64, 0., 64., 8, 0., 8.));
169 name = str(format(
"SlotAverageDeltaTComparison_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
170 title = str(format(
"Ratio of the average #Delta T in CalSet %2% over the previous one. Slot %1%") % (iSlot) % m_calSetLabel);
171 m_slotAverageDeltaTComparison[iSlot].push_back(
new TH1F(name.c_str(), title.c_str(), 512, 0., 512.));
173 name = str(format(
"SlotRMSDeltaTComparison_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
174 title = str(format(
"Ratio of the RMS of #Delta T in CalSet %2% over the previous one. Slot %1%") % (iSlot) % m_calSetLabel);
175 m_slotSigmaDeltaTComparison[iSlot].push_back(
new TH1F(name.c_str(), title.c_str(), 512, 0., 512.));
177 name = str(format(
"SlotAverageDeltaTMapComparison_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
178 title = str(format(
"Map of the ratio of the average #Delta T in CalSet %2% over the previous one. Slot %1%") %
179 (iSlot) % m_calSetLabel);
180 m_slotAverageDeltaTMapComparison[iSlot].push_back(
new TH2F(name.c_str(), title.c_str(), 64, 0., 64., 8, 0., 8.));
182 name = str(format(
"SlotRMSDeltaTMapComparison_Slot%1%_CalSet_%2%") % (iSlot) % m_calSetLabel);
183 title = str(format(
"Map of the ratio of the RMS of #Delta T in CalSet %2% over the previous one. Slot %1%") %
184 (iSlot) % m_calSetLabel);
185 m_slotSigmaDeltaTMapComparison[iSlot].push_back(
new TH2F(name.c_str(), title.c_str(), 64, 0., 64., 8, 0., 8.));
192 B2INFO(
"Initialization of the histograms for " << m_totCalSets <<
" CalSets done.");
198 int TOPTBCComparatorModule::analyzeCalFile()
204 if (m_slotID < 0 || m_boardstackID < 0 || m_scrodID < 0) {
205 B2WARNING(
"Negative slot, BS or scrod ID found while calling analyzeCalFile(). Looks like they have not been initialized, or that a function re-initialized them");
210 for (
short iChannel = 0; iChannel < 128; iChannel++) {
213 if (!m_calSetFile->Get(str(format(
"timeDiff_ch%1%") % (iChannel)).c_str())) {
223 short hardwareChannel = iChannel + 128 * m_boardstackID;
224 auto& chMapper = TOP::TOPGeometryPar::Instance()->getChannelMapper();
225 short pixelID = chMapper.getPixelID(hardwareChannel);
226 short colNum = (pixelID - 1) % 64 + 1;
227 short rowNum = (pixelID - 1) / 64 + 1 ;
228 short globalChannel = hardwareChannel + 512 * (m_slotID -
237 if (!m_calSetFile->Get(str(format(
"timeDiffcal_ch%1%") % (iChannel)).c_str())) {
238 B2WARNING(
"Error opening " << str(format(
"timeDiffcal_ch%1%") % (iChannel)));
240 TH2F* h_timeDiffcal = (TH2F*)m_calSetFile->Get(str(format(
"timeDiffcal_ch%1%") % (iChannel)).c_str());
241 TH1D* h_projection = h_timeDiffcal->ProjectionY(
"h_projection", 1, m_numSamples);
243 m_slotAverageDeltaT[m_slotID - 1][m_calSetID]->SetBinContent(hardwareChannel + 1, h_projection->GetMean());
244 m_slotAverageDeltaT[m_slotID - 1][m_calSetID]->SetBinError(hardwareChannel + 1,
245 h_projection->GetMeanError());
246 m_slotSigmaDeltaT[m_slotID - 1][m_calSetID]->SetBinContent(hardwareChannel + 1,
247 h_projection->GetRMS());
248 m_slotSigmaDeltaT[m_slotID - 1][m_calSetID]->SetBinError(hardwareChannel + 1,
249 h_projection->GetRMSError());
251 m_slotAverageDeltaTMap[m_slotID - 1][m_calSetID]->SetBinContent(colNum, rowNum, h_projection->GetMean());
252 m_slotSigmaDeltaTMap[m_slotID - 1][m_calSetID]->SetBinContent(colNum, rowNum,
253 h_projection->GetRMS());
255 m_topAverageDeltaT[m_calSetID]->SetBinContent(globalChannel + 1, h_projection->GetMean());
256 m_topSigmaDeltaT[m_calSetID]->SetBinContent(globalChannel + 1, h_projection->GetRMS());
264 if (!m_calSetFile->Get(str(format(
"sampleOccup_ch%1%") % (iChannel)).c_str())) {
265 B2WARNING(
"Error opening " << str(format(
"sampleOccup_ch%1%") % (iChannel)));
267 TH1F* h_sampleOccup = (TH1F*)m_calSetFile->Get(str(format(
"sampleOccup_ch%1%") % (iChannel)).c_str());
271 for (
int iSample = 1; iSample < m_numSamples + 1 ; iSample++) {
272 if (h_sampleOccup->GetBinContent(iSample) < m_minCalPulses) nEmpty++;
275 m_slotSampleOccupancy[m_slotID - 1][m_calSetID]->SetBinContent(hardwareChannel + 1, h_sampleOccup->Integral() / m_numSamples);
276 m_slotEmptySamples[m_slotID - 1][m_calSetID]->SetBinContent(hardwareChannel + 1, nEmpty);
278 m_slotSampleOccupancyMap[m_slotID - 1][m_calSetID]->SetBinContent(colNum, rowNum, h_sampleOccup->Integral() / m_numSamples);
279 m_slotEmptySamplesMap[m_slotID - 1][m_calSetID]->SetBinContent(colNum, rowNum, nEmpty);
281 m_topSampleOccupancy[m_calSetID]->SetBinContent(globalChannel + 1, h_sampleOccup->Integral() / m_numSamples);
292 int TOPTBCComparatorModule::makeComparisons()
296 B2INFO(
"Making comparisons for " << m_totCalSets <<
" sets.");
299 for (
int iSet = 1; iSet < m_totCalSets; iSet++) {
300 if (m_compareToPreviousSet) refSet = iSet - 1;
302 m_topAverageDeltaTComparison[iSet] = calculateHistoRatio(m_topAverageDeltaTComparison[iSet], m_topAverageDeltaT[iSet],
303 m_topAverageDeltaT[refSet]);
304 m_topSigmaDeltaTComparison[iSet] = calculateHistoRatio(m_topSigmaDeltaTComparison[iSet], m_topSigmaDeltaT[iSet],
305 m_topSigmaDeltaT[refSet]);
306 m_topSampleOccupancyComparison[iSet] = calculateHistoRatio(m_topSampleOccupancyComparison[iSet], m_topSampleOccupancy[iSet],
307 m_topSampleOccupancy[refSet]);
310 for (
int iSlot = 0; iSlot < 16; iSlot++) {
311 m_slotAverageDeltaTComparison[iSlot][iSet] = calculateHistoRatio(m_slotAverageDeltaTComparison[iSlot][iSet],
312 m_slotAverageDeltaT[iSlot][iSet], m_slotAverageDeltaT[iSlot][refSet]);
313 m_slotAverageDeltaTMapComparison[iSlot][iSet] = calculateHistoRatio(m_slotAverageDeltaTMapComparison[iSlot][iSet],
314 m_slotAverageDeltaTMap[iSlot][iSet], m_slotAverageDeltaTMap[iSlot][refSet]);
315 m_slotSigmaDeltaTComparison[iSlot][iSet] = calculateHistoRatio(m_slotSigmaDeltaTComparison[iSlot][iSet],
316 m_slotSigmaDeltaT[iSlot][iSet], m_slotSigmaDeltaT[iSlot][refSet]);
317 m_slotSigmaDeltaTMapComparison[iSlot][iSet] = calculateHistoRatio(m_slotSigmaDeltaTMapComparison[iSlot][iSet],
318 m_slotSigmaDeltaTMap[iSlot][iSet], m_slotSigmaDeltaTMap[iSlot][refSet]);
321 B2INFO(
"Comparisons done");
327 void TOPTBCComparatorModule::initialize()
332 void TOPTBCComparatorModule::beginRun()
338 void TOPTBCComparatorModule::event()
347 void TOPTBCComparatorModule::endRun()
350 ifstream inputDirectoryListFile(m_inputDirectoryList.c_str());
353 if (!inputDirectoryListFile) {
354 B2ERROR(
"Unable to open the input file with the list of CalSets to analyze");
358 std::string inputString;
360 while (std::getline(inputDirectoryListFile, inputString)) {
363 parseInputDirectoryLine(inputString);
365 B2INFO(
"Processing the calibration set located in " << m_calSetDirectory <<
" and labelled " << m_calSetLabel);
366 TSystemDirectory calSetDir(m_calSetDirectory.c_str(), (m_calSetDirectory).c_str());
369 TList* calSetDirContent = calSetDir.GetListOfFiles();
370 if (calSetDirContent) {
372 TIter next(calSetDirContent);
373 while ((entry = (TSystemFile*)next())) {
376 std::string entryName = entry->GetName();
379 if (!entry->IsDirectory() && TString(entryName.c_str()).EndsWith(
".root")) {
381 m_calSetFile =
new TFile((m_calSetDirectory + entryName).c_str(),
" ");
383 int parserStatus = parseSlotAndScrodIDfromFileName(entryName);
385 if (parserStatus == 1)
388 m_calSetFile->Close();
395 if (entry->IsDirectory() && entryName !=
"." && entryName !=
"..") {
398 TSystemDirectory calSetSubDir((m_calSetDirectory + entryName).c_str(), (m_calSetDirectory + entryName).c_str());
401 TList* calSetSubDirContent = calSetSubDir.GetListOfFiles();
402 if (calSetSubDirContent) {
404 TIter nextSub(calSetSubDirContent);
405 while ((file = (TSystemFile*)nextSub())) {
407 std::string fileName = file->GetName();
409 if (!file->IsDirectory() && TString(fileName.c_str()).EndsWith(
".root")) {
411 m_calSetFile =
new TFile((m_calSetDirectory + entryName + fileName).c_str(),
" ");
413 int parserStatus = parseSlotAndScrodIDfromFileName(fileName);
415 if (parserStatus == 1)
418 m_calSetFile->Close();
419 }
else if (file->IsDirectory() && fileName !=
"." && fileName !=
"..") {
420 B2WARNING(
"Additional subdirectory " << fileName <<
" found in " << m_calSetDirectory + entryName <<
421 ", where only .root and .log files are expected ");
430 B2WARNING(
"Error in creating the TList form the directory " << m_calSetDirectory);
435 B2INFO(
"Analisys concluded.");
444 void TOPTBCComparatorModule::terminate()
447 B2INFO(
"Creating output file " << m_outputFile);
448 TFile outfile(m_outputFile.c_str(),
"recreate");
450 B2INFO(
"Writing histograms ");
453 ifstream inputDirectoryListFile(m_inputDirectoryList.c_str());
454 std::string inputString;
457 while (std::getline(inputDirectoryListFile, inputString)) {
458 parseInputDirectoryLine(inputString);
459 TDirectory* dirSet = outfile.mkdir(m_calSetLabel.c_str());
462 m_topAverageDeltaT[iSet]->Write();
463 m_topSigmaDeltaT[iSet]->Write();
464 m_topSampleOccupancy[iSet]->Write();
465 m_topAverageDeltaTComparison[iSet]->Write();
466 m_topSigmaDeltaTComparison[iSet]->Write();
467 m_topSampleOccupancyComparison[iSet]->Write();
469 for (
int iSlot = 0; iSlot < 16; iSlot ++) {
470 TDirectory* dirSlot = dirSet->mkdir(str(format(
"slot%1%") % (iSlot + 1)).c_str());
473 m_slotAverageDeltaT[iSlot][iSet]->Write();
474 m_slotSigmaDeltaT[iSlot][iSet]->Write();
475 m_slotAverageDeltaTMap[iSlot][iSet]->Write();
476 m_slotSigmaDeltaTMap[iSlot][iSet]->Write();
477 m_slotSampleOccupancy[iSlot][iSet]->Write();
478 m_slotEmptySamples[iSlot][iSet]->Write();
479 m_slotSampleOccupancyMap[iSlot][iSet]->Write();
480 m_slotEmptySamplesMap[iSlot][iSet]->Write();
482 m_slotAverageDeltaTComparison[iSlot][iSet]->Write();
483 m_slotSigmaDeltaTComparison[iSlot][iSet]->Write();
484 m_slotAverageDeltaTMapComparison[iSlot][iSet]->Write();
485 m_slotSigmaDeltaTMapComparison[iSlot][iSet]->Write();
500 int TOPTBCComparatorModule::parseInputDirectoryLine(std::string inputString)
503 m_calSetDirectory.clear();
504 m_calSetLabel.clear();
507 bool isDirectoryNameChar =
true;
508 bool isAtBeginning =
true;
510 for (std::string::size_type i = 0; i < inputString.size(); ++i) {
511 char c = inputString[i];
513 if (c ==
' ' && isAtBeginning)
continue;
514 if (c ==
' ' && !isAtBeginning) {
515 isDirectoryNameChar =
false;
516 isAtBeginning =
false;
519 if (c !=
' ' && isDirectoryNameChar) {
520 m_calSetDirectory += c;
521 isAtBeginning =
false;
524 if (c !=
' ' && !isDirectoryNameChar) {
526 isAtBeginning =
false;
529 B2WARNING(
"Uncaught exception in parsing the input string. Ending the parsing.");
538 int TOPTBCComparatorModule::parseSlotAndScrodIDfromFileName(std::string inputString)
548 std::string stringSlot =
"";
549 std::string stringBS =
"";
550 std::string stringScrod =
"";
560 if (!(inputString[0] ==
't' && inputString[1] ==
'b' && inputString[2] ==
'c')) {
561 B2WARNING(inputString <<
" is not a valid TBC file. Skipping it.");
564 stringSlot += inputString[7];
565 stringSlot += inputString[8];
566 stringBS += inputString[10];
567 stringScrod += inputString[17];
568 stringScrod += inputString[18];
569 if (inputString[19] !=
'.')
570 stringScrod += inputString[19];
573 m_slotID = std::stoi(stringSlot);
574 m_scrodID = std::stoi(stringScrod);
575 m_boardstackID = std::stoi(stringBS);
582 TH1F* TOPTBCComparatorModule::calculateHistoRatio(TH1F* hRatio, TH1F* hNum, TH1F* hDen)
585 const char* name = hRatio->GetName();
586 const char* title = hRatio->GetTitle();
588 const char* xAxisTitle = hRatio->GetXaxis()->GetTitle();
589 const char* yAxisTitle = hRatio->GetYaxis()->GetTitle();
592 hRatio = (TH1F*)hNum->Clone();
596 hRatio->Divide(hDen);
599 hRatio->SetName(name);
600 hRatio->SetTitle(title);
601 hRatio->GetXaxis()->SetTitle(xAxisTitle);
602 hRatio->GetYaxis()->SetTitle(yAxisTitle);
606 TH2F* TOPTBCComparatorModule::calculateHistoRatio(TH2F* hRatio, TH2F* hNum, TH2F* hDen)
609 const char* name = hRatio->GetName();
610 const char* title = hRatio->GetTitle();
612 const char* xAxisTitle = hRatio->GetXaxis()->GetTitle();
613 const char* yAxisTitle = hRatio->GetYaxis()->GetTitle();
616 hRatio = (TH2F*)hNum->Clone();
620 hRatio->Divide(hDen);
623 hRatio->SetName(name);
624 hRatio->SetTitle(title);
625 hRatio->GetXaxis()->SetTitle(xAxisTitle);
626 hRatio->GetYaxis()->SetTitle(yAxisTitle);
HistoModule.h is supposed to be used instead of Module.h for the modules with histogram definitions t...
Module for the comparison of different sets of time base correction (TBC) constants and to produce mo...
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Abstract base class for different kinds of events.