Belle II Software  release-05-01-25
EKLMChannelDataImporter.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 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/eklm/calibration/EKLMChannelDataImporter.h>
13 
14 /* KLM headers. */
15 #include <klm/dataobjects/eklm/EKLMElementNumbers.h>
16 #include <klm/dataobjects/KLMElementNumbers.h>
17 #include <klm/dbobjects/eklm/EKLMChannels.h>
18 #include <klm/dbobjects/KLMElectronicsMap.h>
19 #include <klm/eklm/geometry/GeometryData.h>
20 
21 /* Belle 2 headers. */
22 #include <framework/database/IntervalOfValidity.h>
23 #include <framework/database/DBImportObjPtr.h>
24 #include <framework/database/DBObjPtr.h>
25 #include <framework/logging/Logger.h>
26 
27 /* ROOT headers. */
28 #include <TFile.h>
29 #include <TTree.h>
30 
31 using namespace Belle2;
32 
34 {
35 }
36 
38 {
39 }
40 
41 void EKLMChannelDataImporter::setIOV(int experimentLow, int runLow,
42  int experimentHigh, int runHigh)
43 {
44  m_ExperimentLow = experimentLow;
45  m_RunLow = runLow;
46  m_ExperimentHigh = experimentHigh;
47  m_RunHigh = runHigh;
48 }
49 
51 {
52  m_Channels.construct();
54  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
55  int iSection, iLayer, iSector, iPlane, iStrip, strip;
56  for (iSection = 1; iSection <= geoDat->getNSections(); iSection++) {
57  for (iLayer = 1; iLayer <= geoDat->getNDetectorLayers(iSection);
58  iLayer++) {
59  for (iSector = 1; iSector <= geoDat->getNSectors(); iSector++) {
60  for (iPlane = 1; iPlane <= geoDat->getNPlanes(); iPlane++) {
61  for (iStrip = 1; iStrip <= geoDat->getNStrips(); iStrip++) {
62  strip = elementNumbers->stripNumber(
63  iSection, iLayer, iSector, iPlane, iStrip);
64  m_Channels->setChannelData(strip, channelData);
65  }
66  }
67  }
68  }
69  }
70 }
71 
73  int section, int layer, int sector, int plane, int strip,
74  EKLMChannelData* channelData)
75 {
76  int stripGlobal;
77  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
78  stripGlobal = elementNumbers->stripNumber(section, layer, sector, plane,
79  strip);
80  m_Channels->setChannelData(stripGlobal, channelData);
81 }
82 
83 void EKLMChannelDataImporter::loadActiveChannels(const char* activeChannelsData)
84 {
85  int i, n;
86  int copper, dataConcentrator, lane, daughterCard, channel, active;
87  /* cppcheck-suppress variableScope */
88  int subdetector, section, layer, sector, plane, strip, stripGlobal;
89  /* cppcheck-suppress variableScope */
90  const uint16_t* detectorChannel;
91  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
92  const KLMElementNumbers* klmElementNumbers = &(KLMElementNumbers::Instance());
93  DBObjPtr<KLMElectronicsMap> electronicsMap;
94  KLMElectronicsChannel electronicsChannel;
95  TFile* file;
96  TTree* tree;
97  file = new TFile(activeChannelsData, "");
98  tree = (TTree*)file->Get("tree");
99  n = tree->GetEntries();
100  tree->SetBranchAddress("copper", &copper);
101  tree->SetBranchAddress("data_concentrator", &dataConcentrator);
102  tree->SetBranchAddress("lane", &lane);
103  tree->SetBranchAddress("daughter_card", &daughterCard);
104  tree->SetBranchAddress("channel", &channel);
105  tree->SetBranchAddress("active", &active);
106  for (i = 0; i < n; i++) {
107  tree->GetEntry(i);
108  electronicsChannel.setCopper(copper);
109  electronicsChannel.setSlot(dataConcentrator);
110  electronicsChannel.setLane(lane);
111  electronicsChannel.setAxis(1);
112  electronicsChannel.setChannel(1);
113  detectorChannel = electronicsMap->getDetectorChannel(&electronicsChannel);
114  if (detectorChannel == nullptr) {
115  B2FATAL("Wrong DAQ channel in calibration data: copper = " << copper <<
116  ", data_concentrator = " << dataConcentrator << ", lane = " <<
117  lane);
118  }
119  klmElementNumbers->channelNumberToElementNumbers(
120  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
121  &strip);
122  stripGlobal = elementNumbers->stripNumber(section, layer, sector, plane,
123  strip);
124  EKLMChannelData* channelData = const_cast<EKLMChannelData*>(
125  m_Channels->getChannelData(stripGlobal));
126  if (channelData == nullptr)
127  B2FATAL("Channel data are not loaded. Use loadChannelData().");
128  }
129  delete tree;
130  delete file;
131 }
132 
133 void EKLMChannelDataImporter::loadHighVoltage(const char* highVoltageData)
134 {
135  int i, n;
136  int copper, dataConcentrator, lane, daughterCard, channel;
137  float voltage;
138  /* cppcheck-suppress variableScope */
139  int subdetector, section, layer, sector, plane, strip, stripGlobal;
140  /* cppcheck-suppress variableScope */
141  const uint16_t* detectorChannel;
142  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
143  const KLMElementNumbers* klmElementNumbers = &(KLMElementNumbers::Instance());
144  DBObjPtr<KLMElectronicsMap> electronicsMap;
145  KLMElectronicsChannel electronicsChannel;
146  TFile* file;
147  TTree* tree;
148  file = new TFile(highVoltageData, "");
149  tree = (TTree*)file->Get("tree");
150  n = tree->GetEntries();
151  tree->SetBranchAddress("copper", &copper);
152  tree->SetBranchAddress("data_concentrator", &dataConcentrator);
153  tree->SetBranchAddress("lane", &lane);
154  tree->SetBranchAddress("daughter_card", &daughterCard);
155  tree->SetBranchAddress("channel", &channel);
156  tree->SetBranchAddress("voltage", &voltage);
157  for (i = 0; i < n; i++) {
158  tree->GetEntry(i);
159  electronicsChannel.setCopper(copper);
160  electronicsChannel.setSlot(dataConcentrator);
161  electronicsChannel.setLane(lane);
162  electronicsChannel.setAxis(1);
163  electronicsChannel.setChannel(1);
164  detectorChannel = electronicsMap->getDetectorChannel(&electronicsChannel);
165  if (detectorChannel == nullptr) {
166  B2FATAL("Wrong DAQ channel in calibration data: copper = " << copper <<
167  ", data_concentrator = " << dataConcentrator << ", lane = " <<
168  lane);
169  }
170  klmElementNumbers->channelNumberToElementNumbers(
171  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
172  &strip);
173  stripGlobal = elementNumbers->stripNumber(section, layer, sector, plane,
174  strip);
175  EKLMChannelData* channelData = const_cast<EKLMChannelData*>(
176  m_Channels->getChannelData(stripGlobal));
177  if (channelData == nullptr)
178  B2FATAL("Channel data are not loaded. Use loadChannelData().");
179  channelData->setVoltage(voltage);
180  }
181  delete tree;
182  delete file;
183 }
184 
185 void EKLMChannelDataImporter::loadLookbackWindow(const char* lookbackWindowData)
186 {
187  int i, n;
188  int copper, dataConcentrator, lane, daughterCard, channel;
189  int lookbackTime, lookbackWindowWidth;
190  /* cppcheck-suppress variableScope */
191  int subdetector, section, layer, sector, plane, strip, stripGlobal;
192  /* cppcheck-suppress variableScope */
193  const uint16_t* detectorChannel;
194  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
195  const KLMElementNumbers* klmElementNumbers = &(KLMElementNumbers::Instance());
196  DBObjPtr<KLMElectronicsMap> electronicsMap;
197  KLMElectronicsChannel electronicsChannel;
198  TFile* file;
199  TTree* tree;
200  file = new TFile(lookbackWindowData, "");
201  tree = (TTree*)file->Get("tree");
202  n = tree->GetEntries();
203  tree->SetBranchAddress("copper", &copper);
204  tree->SetBranchAddress("data_concentrator", &dataConcentrator);
205  tree->SetBranchAddress("lane", &lane);
206  tree->SetBranchAddress("daughter_card", &daughterCard);
207  tree->SetBranchAddress("channel", &channel);
208  tree->SetBranchAddress("lookback_time", &lookbackTime);
209  tree->SetBranchAddress("lookback_window_width", &lookbackWindowWidth);
210  for (i = 0; i < n; i++) {
211  tree->GetEntry(i);
212  electronicsChannel.setCopper(copper);
213  electronicsChannel.setSlot(dataConcentrator);
214  electronicsChannel.setLane(lane);
215  electronicsChannel.setAxis(1);
216  electronicsChannel.setChannel(1);
217  detectorChannel = electronicsMap->getDetectorChannel(&electronicsChannel);
218  if (detectorChannel == nullptr) {
219  B2FATAL("Wrong DAQ channel in calibration data: copper = " << copper <<
220  ", data_concentrator = " << dataConcentrator << ", lane = " <<
221  lane);
222  }
223  klmElementNumbers->channelNumberToElementNumbers(
224  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
225  &strip);
226  stripGlobal = elementNumbers->stripNumber(section, layer, sector, plane,
227  strip);
228  EKLMChannelData* channelData = const_cast<EKLMChannelData*>(
229  m_Channels->getChannelData(stripGlobal));
230  if (channelData == nullptr)
231  B2FATAL("Channel data are not loaded. Use loadChannelData().");
232  channelData->setLookbackTime(lookbackTime);
233  channelData->setLookbackWindowWidth(lookbackWindowWidth);
234  }
235  delete tree;
236  delete file;
237 }
238 
239 void EKLMChannelDataImporter::loadThresholds(const char* thresholdsData)
240 {
241  int i, n;
242  int copper, dataConcentrator, lane, daughterCard, channel;
243  int active, pedestalMin, threshold, adjustmentVoltage;
244  /* cppcheck-suppress variableScope */
245  int subdetector, section, layer, sector, plane, strip, stripGlobal;
246  /* cppcheck-suppress variableScope */
247  const uint16_t* detectorChannel;
248  const EKLMElementNumbers* elementNumbers = &(EKLMElementNumbers::Instance());
249  const KLMElementNumbers* klmElementNumbers = &(KLMElementNumbers::Instance());
250  DBObjPtr<KLMElectronicsMap> electronicsMap;
251  KLMElectronicsChannel electronicsChannel;
252  TFile* file;
253  TTree* tree;
254  file = new TFile(thresholdsData, "");
255  tree = (TTree*)file->Get("tree");
256  n = tree->GetEntries();
257  tree->SetBranchAddress("copper", &copper);
258  tree->SetBranchAddress("data_concentrator", &dataConcentrator);
259  tree->SetBranchAddress("lane", &lane);
260  tree->SetBranchAddress("daughter_card", &daughterCard);
261  tree->SetBranchAddress("channel", &channel);
262  tree->SetBranchAddress("active", &active);
263  tree->SetBranchAddress("pedestal_min", &pedestalMin);
264  tree->SetBranchAddress("threshold", &threshold);
265  tree->SetBranchAddress("adjustment_voltage", &adjustmentVoltage);
266  for (i = 0; i < n; i++) {
267  tree->GetEntry(i);
268  electronicsChannel.setCopper(copper);
269  electronicsChannel.setSlot(dataConcentrator);
270  electronicsChannel.setLane(lane);
271  electronicsChannel.setAxis(1);
272  electronicsChannel.setChannel(1);
273  detectorChannel = electronicsMap->getDetectorChannel(&electronicsChannel);
274  if (detectorChannel == nullptr) {
275  B2FATAL("Wrong DAQ channel in calibration data: copper = " << copper <<
276  ", data_concentrator = " << dataConcentrator << ", lane = " <<
277  lane);
278  }
279  klmElementNumbers->channelNumberToElementNumbers(
280  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
281  &strip);
282  stripGlobal = elementNumbers->stripNumber(section, layer, sector, plane,
283  strip);
284  EKLMChannelData* channelData = const_cast<EKLMChannelData*>(
285  m_Channels->getChannelData(stripGlobal));
286  if (channelData == nullptr)
287  B2FATAL("Channel data are not loaded. Use loadChannelData().");
288  channelData->setPedestal(pedestalMin);
289  channelData->setThreshold(threshold);
290  channelData->setAdjustmentVoltage(adjustmentVoltage);
291  }
292  delete tree;
293  delete file;
294 }
295 
297 {
300  m_Channels.import(iov);
301 }
302 
Belle2::IntervalOfValidity
A class that describes the interval of experiments/runs for which an object in the database is valid.
Definition: IntervalOfValidity.h:35
Belle2::EKLMChannelDataImporter::setChannelData
void setChannelData(int section, int layer, int sector, int plane, int strip, EKLMChannelData *channelData)
Set channel data.
Definition: EKLMChannelDataImporter.cc:72
Belle2::EKLMChannelData
EKLM channel data.
Definition: EKLMChannelData.h:33
Belle2::EKLMChannelDataImporter::loadChannelData
void loadChannelData(EKLMChannelData *channelData)
Load specific channel data to all channels.
Definition: EKLMChannelDataImporter.cc:50
Belle2::EKLMChannelData::setPedestal
void setPedestal(float pedestal)
Set pedestal.
Definition: EKLMChannelData.h:58
Belle2::EKLMElementNumbers
EKLM element numbers.
Definition: EKLMElementNumbers.h:34
Belle2::KLMElectronicsChannel::setCopper
void setCopper(int copper)
Set copper.
Definition: KLMElectronicsChannel.h:86
Belle2::EKLMGeometry::getNPlanes
int getNPlanes() const
Get number of planes.
Definition: EKLMGeometry.h:1717
Belle2::EKLMChannelDataImporter::loadThresholds
void loadThresholds(const char *thresholdsData)
Load thresholds and adjustment voltages from ROOT file.
Definition: EKLMChannelDataImporter.cc:239
Belle2::EKLMChannelData::setAdjustmentVoltage
void setAdjustmentVoltage(int adjustmentVoltage)
Set adjustment voltage.
Definition: EKLMChannelData.h:122
Belle2::KLMElementNumbers::Instance
static const KLMElementNumbers & Instance()
Instantiation.
Definition: KLMElementNumbers.cc:31
Belle2::EKLMChannelDataImporter::loadActiveChannels
void loadActiveChannels(const char *activeChannelsData)
Load active channels from ROOT file.
Definition: EKLMChannelDataImporter.cc:83
Belle2::EKLMChannelDataImporter::m_RunHigh
int m_RunHigh
High run.
Definition: EKLMChannelDataImporter.h:118
Belle2::KLMElectronicsChannel::setChannel
void setChannel(int channel)
Set channel.
Definition: KLMElectronicsChannel.h:154
Belle2::KLMElectronicsChannel::setAxis
void setAxis(int axis)
Set axis.
Definition: KLMElectronicsChannel.h:137
Belle2::DBObjPtr
Class for accessing objects in the database.
Definition: DBObjPtr.h:31
Belle2::EKLMGeometry::getNStrips
int getNStrips() const
Get number of strips.
Definition: EKLMGeometry.h:1741
Belle2::EKLMChannelData::setThreshold
void setThreshold(int threshold)
Set threshold.
Definition: EKLMChannelData.h:90
Belle2::EKLMChannelDataImporter::importChannelData
void importChannelData()
Import channel data.
Definition: EKLMChannelDataImporter.cc:296
Belle2::EKLMChannelDataImporter::~EKLMChannelDataImporter
~EKLMChannelDataImporter()
Destructor.
Definition: EKLMChannelDataImporter.cc:37
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::EKLMChannelDataImporter::setIOV
void setIOV(int experimentLow, int runLow, int experimentHigh, int runHigh)
Set interval of validity.
Definition: EKLMChannelDataImporter.cc:41
Belle2::EKLMElementNumbers::Instance
static const EKLMElementNumbers & Instance()
Instantiation.
Definition: EKLMElementNumbers.cc:20
Belle2::EKLMChannelData::setVoltage
void setVoltage(float voltage)
Set voltage.
Definition: EKLMChannelData.h:106
Belle2::EKLMGeometry::getNDetectorLayers
int getNDetectorLayers(int section) const
Get number of detector layers.
Definition: EKLMGeometry.cc:295
Belle2::EKLM::GeometryData
EKLM geometry data.
Definition: GeometryData.h:40
Belle2::EKLMElementNumbers::stripNumber
int stripNumber(int section, int layer, int sector, int plane, int strip) const
Get strip number.
Definition: EKLMElementNumbers.cc:212
Belle2::EKLMChannelDataImporter::EKLMChannelDataImporter
EKLMChannelDataImporter()
Constructor.
Definition: EKLMChannelDataImporter.cc:33
Belle2::EKLMChannelData::setLookbackTime
void setLookbackTime(int lookbackTime)
Set lookback time (unit is 32 TDC counts).
Definition: EKLMChannelData.h:138
Belle2::EKLMChannelDataImporter::loadLookbackWindow
void loadLookbackWindow(const char *lookbackWindowData)
Load lookback window from ROOT file.
Definition: EKLMChannelDataImporter.cc:185
Belle2::EKLMChannelDataImporter::m_Channels
DBImportObjPtr< EKLMChannels > m_Channels
Channel data.
Definition: EKLMChannelDataImporter.h:106
Belle2::EKLMChannelDataImporter::m_ExperimentHigh
int m_ExperimentHigh
High experiment.
Definition: EKLMChannelDataImporter.h:115
Belle2::KLMElectronicsChannel
BKLM electronics channel.
Definition: KLMElectronicsChannel.h:33
Belle2::EKLMChannelDataImporter::m_RunLow
int m_RunLow
Low run.
Definition: EKLMChannelDataImporter.h:112
Belle2::KLMElectronicsChannel::setLane
void setLane(int lane)
Set lane.
Definition: KLMElectronicsChannel.h:120
Belle2::KLMElementNumbers
KLM element numbers.
Definition: KLMElementNumbers.h:37
Belle2::EKLMChannelData::setLookbackWindowWidth
void setLookbackWindowWidth(int lookbackWindowWidth)
Set lookback window width (unit is 32 TDC counts).
Definition: EKLMChannelData.h:154
Belle2::EKLMGeometry::getNSections
int getNSections() const
Get number of sections.
Definition: EKLMGeometry.h:1687
Belle2::EKLMChannelDataImporter::loadHighVoltage
void loadHighVoltage(const char *highVoltageData)
Load high voltage from ROOT file.
Definition: EKLMChannelDataImporter.cc:133
Belle2::EKLMGeometry::getNSectors
int getNSectors() const
Get number of sectors.
Definition: EKLMGeometry.h:1709
Belle2::KLMElectronicsChannel::setSlot
void setSlot(int slot)
Set slot.
Definition: KLMElectronicsChannel.h:103
Belle2::EKLM::GeometryData::Instance
static const GeometryData & Instance(enum DataSource dataSource=c_Database, const GearDir *gearDir=nullptr)
Instantiation.
Definition: GeometryData.cc:35
Belle2::KLMElementNumbers::channelNumberToElementNumbers
void channelNumberToElementNumbers(uint16_t channel, int *subdetector, int *section, int *sector, int *layer, int *plane, int *strip) const
Get element numbers by channel number.
Definition: KLMElementNumbers.cc:106
Belle2::EKLMChannelDataImporter::m_ExperimentLow
int m_ExperimentLow
Low experiment.
Definition: EKLMChannelDataImporter.h:109