Belle II Software development
KLMChannelStatusCollectorModule.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/modules/KLMChannelStatusCollector/KLMChannelStatusCollectorModule.h>
11
12/* KLM headers. */
13#include <klm/dataobjects/KLMChannelIndex.h>
14#include <klm/dataobjects/KLMChannelMapValue.h>
15
16/* ROOT headers. */
17#include <TFile.h>
18#include <TH1.h>
19#include <TTree.h>
20
21using namespace Belle2;
22
23REG_MODULE(KLMChannelStatusCollector);
24
27 m_ElementNumbers(&(KLMElementNumbers::Instance())),
29 m_HitMap("KLMChannelMapHits", DataStore::c_Persistent)
30{
31 setDescription("Module for KLM channel status calibration (data collection).");
33}
34
38
40{
41 m_KLMDigits.isRequired();
42 m_HitMap.registerInDataStore();
43 m_HitMap.create();
44 m_HitMap->setDataAllChannels(0);
45 TTree* calibrationData = new TTree("calibration_data", "");
46 registerObject<TTree>("calibration_data", calibrationData);
47}
48
50{
51 for (KLMDigit& digit : m_KLMDigits) {
52 /* Multi-strip digits do not correspond to individual channels. */
53 if (digit.isMultiStrip())
54 continue;
55 KLMChannelNumber channel =
56 m_ElementNumbers->channelNumber(
57 digit.getSubdetector(), digit.getSection(), digit.getSector(),
58 digit.getLayer(), digit.getPlane(), digit.getStrip());
59 m_HitMap->setChannelData(channel, m_HitMap->getChannelData(channel) + 1);
60 }
61}
62
64{
65 KLMChannelNumber channel;
66 unsigned int hits;
67 TTree* calibrationData = getObjectPtr<TTree>("calibration_data");
68 calibrationData->Branch("channel", &channel, "channel/s");
69 calibrationData->Branch("hits", &hits, "hits/i");
70 KLMChannelIndex klmChannels;
71 for (KLMChannelIndex& klmChannel : klmChannels) {
72 channel = klmChannel.getKLMChannelNumber();
73 hits = m_HitMap->getChannelData(channel);
74 calibrationData->Fill();
75 }
76 /* Clear data for case of multiple runs. */
77 m_HitMap->setDataAllChannels(0);
78}
79
81 const char* dqmFile)
82{
83 std::string histogramName;
84 TFile* f = new TFile(dqmFile);
86 for (KLMChannelIndex& klmSector : klmSectors) {
87 int nHistograms;
88 if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
89 nHistograms = 2;
90 else
91 nHistograms = 3;
92 for (int j = 0; j < nHistograms; ++j) {
93 histogramName = "KLM/strip_hits_subdetector_" +
94 std::to_string(klmSector.getSubdetector()) +
95 "_section_" + std::to_string(klmSector.getSection()) +
96 "_sector_" + std::to_string(klmSector.getSector()) +
97 "_" + std::to_string(j);
98 TH1* histogram = (TH1*)f->Get(histogramName.c_str());
99 if (histogram == nullptr) {
100 B2ERROR("KLM DQM histogram " << histogramName << " is not found.");
101 continue;
102 }
103 int n = histogram->GetXaxis()->GetNbins();
104 for (int i = 1; i <= n; ++i) {
105 KLMChannelNumber channelIndex = std::round(histogram->GetBinCenter(i));
106 KLMChannelNumber channelNumber = m_ChannelArrayIndex->getNumber(channelIndex);
107 m_HitMap->setChannelData(channelNumber, histogram->GetBinContent(i));
108 }
109 }
110 }
111}
void registerObject(const std::string &name, T *obj)
Register object with a name, takes ownership, do not access the pointer beyond prepare()
T * getObjectPtr(const std::string &name)
Calls the CalibObjManager to get the requested stored collector data.
CalibrationCollectorModule()
Constructor. Sets the default prefix for calibration dataobjects.
In the store you can park objects that have to be accessed by various modules.
Definition DataStore.h:51
KLM channel index.
StoreObjPtr< KLMChannelMapValue< unsigned int > > m_HitMap
Hit map.
const KLMElementNumbers * m_ElementNumbers
Element numbers.
void collectFromDQM(const char *dqmFile)
Collection of data from DQM modules.
void collect() override
This method is called for each event.
void closeRun() override
This method is called at the end of run.
const KLMChannelArrayIndex * m_ChannelArrayIndex
KLM channel array index.
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition KLMDigit.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
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition Module.h:649
uint16_t KLMChannelNumber
Channel number.
Abstract base class for different kinds of events.