Belle II Software  release-05-01-25
KLMChannelStatusCollectorModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2019 Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Kirill Chilikin *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/modules/KLMChannelStatusCollector/KLMChannelStatusCollectorModule.h>
13 
14 /* KLM headers. */
15 #include <klm/dataobjects/KLMChannelIndex.h>
16 #include <klm/dataobjects/KLMChannelMapValue.h>
17 
18 /* ROOT headers. */
19 #include <TFile.h>
20 #include <TH1.h>
21 #include <TTree.h>
22 
23 using namespace Belle2;
24 
25 REG_MODULE(KLMChannelStatusCollector)
26 
29  m_ElementNumbers(&(KLMElementNumbers::Instance())),
30  m_ChannelArrayIndex(&(KLMChannelArrayIndex::Instance())),
31  m_HitMap("KLMChannelMapHits", DataStore::c_Persistent)
32 {
33  setDescription("Module for KLM channel status calibration (data collection).");
34  setPropertyFlags(c_ParallelProcessingCertified);
35 }
36 
38 {
39 }
40 
42 {
43  m_KLMDigits.isRequired();
44  m_HitMap.registerInDataStore();
45  m_HitMap.create();
46  m_HitMap->setDataAllChannels(0);
47  TTree* calibrationData = new TTree("calibration_data", "");
48  registerObject<TTree>("calibration_data", calibrationData);
49 }
50 
52 {
53  for (KLMDigit& digit : m_KLMDigits) {
54  /* Multi-strip digits do not correspond to individual channels. */
55  if (digit.isMultiStrip())
56  continue;
57  uint16_t channel =
59  digit.getSubdetector(), digit.getSection(), digit.getSector(),
60  digit.getLayer(), digit.getPlane(), digit.getStrip());
61  m_HitMap->setChannelData(channel, m_HitMap->getChannelData(channel) + 1);
62  }
63 }
64 
66 {
67  uint16_t channel;
68  unsigned int hits;
69  TTree* calibrationData = getObjectPtr<TTree>("calibration_data");
70  calibrationData->Branch("channel", &channel, "channel/s");
71  calibrationData->Branch("hits", &hits, "hits/i");
72  KLMChannelIndex klmChannels;
73  for (KLMChannelIndex& klmChannel : klmChannels) {
74  channel = klmChannel.getKLMChannelNumber();
75  hits = m_HitMap->getChannelData(channel);
76  calibrationData->Fill();
77  }
78  /* Clear data for case of multiple runs. */
79  m_HitMap->setDataAllChannels(0);
80 }
81 
83  const char* dqmFile)
84 {
85  std::string histogramName;
86  TFile* f = new TFile(dqmFile);
88  for (KLMChannelIndex& klmSector : klmSectors) {
89  int nHistograms;
90  if (klmSector.getSubdetector() == KLMElementNumbers::c_BKLM)
91  nHistograms = 2;
92  else
93  nHistograms = 3;
94  for (int j = 0; j < nHistograms; ++j) {
95  histogramName = "KLM/strip_hits_subdetector_" +
96  std::to_string(klmSector.getSubdetector()) +
97  "_section_" + std::to_string(klmSector.getSection()) +
98  "_sector_" + std::to_string(klmSector.getSector()) +
99  "_" + std::to_string(j);
100  TH1* histogram = (TH1*)f->Get(histogramName.c_str());
101  if (histogram == nullptr) {
102  B2ERROR("KLM DQM histogram " << histogramName << " is not found.");
103  continue;
104  }
105  int n = histogram->GetXaxis()->GetNbins();
106  for (int i = 1; i <= n; ++i) {
107  uint16_t channelIndex = std::round(histogram->GetBinCenter(i));
108  uint16_t channelNumber = m_ChannelArrayIndex->getNumber(channelIndex);
109  m_HitMap->setChannelData(channelNumber, histogram->GetBinContent(i));
110  }
111  }
112  }
113 }
Belle2::CalibrationCollectorModule
Calibration collector module base class.
Definition: CalibrationCollectorModule.h:44
Belle2::KLMChannelStatusCollectorModule::collect
void collect() override
This method is called for each event.
Definition: KLMChannelStatusCollectorModule.cc:51
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::KLMChannelStatusCollectorModule::collectFromDQM
void collectFromDQM(const char *dqmFile)
Collection of data from DQM modules.
Definition: KLMChannelStatusCollectorModule.cc:82
Belle2::KLMChannelStatusCollectorModule::m_HitMap
StoreObjPtr< KLMChannelMapValue< unsigned int > > m_HitMap
Hit map.
Definition: KLMChannelStatusCollectorModule.h:89
Belle2::KLMElementNumbers::channelNumber
uint16_t channelNumber(int subdetector, int section, int sector, int layer, int plane, int strip) const
Get channel number.
Definition: KLMElementNumbers.cc:37
Belle2::KLMDigit
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:40
Belle2::KLMChannelIndex::c_IndexLevelSector
@ c_IndexLevelSector
Sector.
Definition: KLMChannelIndex.h:49
Belle2::KLMElementArrayIndex::getNumber
uint16_t getNumber(uint16_t index) const
Get element number.
Definition: KLMElementArrayIndex.cc:63
Belle2::KLMChannelStatusCollectorModule
KLM channel status calibration (data collection).
Definition: KLMChannelStatusCollectorModule.h:41
Belle2::KLMChannelStatusCollectorModule::m_ElementNumbers
const KLMElementNumbers * m_ElementNumbers
Element numbers.
Definition: KLMChannelStatusCollectorModule.h:80
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KLMChannelStatusCollectorModule::closeRun
void closeRun() override
This method is called at the end of run.
Definition: KLMChannelStatusCollectorModule.cc:65
Belle2::KLMChannelStatusCollectorModule::m_ChannelArrayIndex
const KLMChannelArrayIndex * m_ChannelArrayIndex
KLM channel array index.
Definition: KLMChannelStatusCollectorModule.h:83
Belle2::KLMElementNumbers::c_BKLM
@ c_BKLM
BKLM.
Definition: KLMElementNumbers.h:47
Belle2::KLMChannelStatusCollectorModule::prepare
void prepare() override
Initializer.
Definition: KLMChannelStatusCollectorModule.cc:41
Belle2::KLMChannelStatusCollectorModule::m_KLMDigits
StoreArray< KLMDigit > m_KLMDigits
KLM digits.
Definition: KLMChannelStatusCollectorModule.h:86
Belle2::KLMChannelIndex
KLM channel index.
Definition: KLMChannelIndex.h:33
Belle2::KLMElementNumbers
KLM element numbers.
Definition: KLMElementNumbers.h:37
Belle2::KLMChannelStatusCollectorModule::~KLMChannelStatusCollectorModule
~KLMChannelStatusCollectorModule()
Destructor.
Definition: KLMChannelStatusCollectorModule.cc:37
Belle2::KLMChannelArrayIndex
KLM channel array index.
Definition: KLMChannelArrayIndex.h:33
Belle2::DataStore
In the store you can park objects that have to be accessed by various modules.
Definition: DataStore.h:53