Belle II Software  release-05-02-19
KLMUnpackerModule.cc
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2016 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Petr Katrenko *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 /* Own header. */
12 #include <klm/modules/KLMUnpacker/KLMUnpackerModule.h>
13 
14 /* KLM headers. */
15 #include <klm/dataobjects/bklm/BKLMElementNumbers.h>
16 #include <klm/dataobjects/KLMScintillatorFirmwareFitResult.h>
17 #include <klm/rawdata/RawData.h>
18 
19 /* Belle 2 headers. */
20 #include <framework/logging/Logger.h>
21 
22 /* C++ headers. */
23 #include <cstdint>
24 
25 using namespace std;
26 using namespace Belle2;
27 using namespace Belle2::KLM;
28 
29 REG_MODULE(KLMUnpacker)
30 
32  m_ElementNumbers(&(KLMElementNumbers::Instance())),
33  m_triggerCTimeOfPreviousEvent(0),
34  m_eklmElementNumbers(&(EKLMElementNumbers::Instance()))
35 {
36  setDescription("KLM unpacker (creates KLMDigits from RawKLM).");
37  setPropertyFlags(c_ParallelProcessingCertified);
38  addParam("outputKLMDigitsName", m_outputKLMDigitsName,
39  "Name of KLMDigit store array.", string(""));
40  addParam("WriteDigitRaws", m_WriteDigitRaws,
41  "Record raw data in dataobject format (e.g. for debugging).", false);
42  addParam("WriteWrongHits", m_WriteWrongHits,
43  "Record wrong hits (e.g. for debugging).", false);
44  addParam("DebugElectronicsMap", m_DebugElectronicsMap,
45  "Debug electronics map (record DAQ channel instead of strip).",
46  false);
47  addParam("DAQChannelBKLMScintillators", m_DAQChannelBKLMScintillators,
48  "Record DAQ channel for BKLM scintillators.", false);
49  addParam("DAQChannelModule", m_DAQChannelModule,
50  "Record DAQ channel for specific module.", -1);
51  addParam("IgnoreWrongHits", m_IgnoreWrongHits,
52  "Ignore wrong hits (i.e. no B2ERROR).", false);
53  addParam("IgnoreStrip0", m_IgnoreStrip0,
54  "Ignore hits with strip = 0 (normally expected for certain firmware "
55  "versions).", true);
56  addParam("keepEvenPackages", m_keepEvenPackages,
57  "Keep packages that have even length normally indicating that "
58  "data was corrupted ", false);
59  addParam("SciThreshold", m_scintThreshold,
60  "Scintillator strip hits with charge lower this value will be "
61  "marked as bad.", double(140.0));
62  addParam("loadThresholdFromDB", m_loadThresholdFromDB,
63  "Load threshold from database (true) or not (false)", true);
64 }
65 
66 KLMUnpackerModule::~KLMUnpackerModule()
67 {
68 }
69 
70 void KLMUnpackerModule::initialize()
71 {
72  m_RawKLMs.isRequired();
73  /* Digits. */
74  m_Digits.registerInDataStore(m_outputKLMDigitsName);
75  m_klmDigitsOutOfRange.registerInDataStore("KLMDigitsOutOfRange");
76  /* Event information. */
77  m_DigitEventInfos.registerInDataStore();
78  m_Digits.registerRelationTo(m_DigitEventInfos);
79  m_klmDigitsOutOfRange.registerRelationTo(m_DigitEventInfos);
80  /* Raw data in dataobject format. */
81  if (m_WriteDigitRaws) {
82  m_klmDigitRaws.registerInDataStore();
83  m_Digits.registerRelationTo(m_klmDigitRaws);
84  m_klmDigitsOutOfRange.registerRelationTo(m_klmDigitRaws);
85  }
86 }
87 
88 void KLMUnpackerModule::beginRun()
89 {
90  if (!m_ElectronicsMap.isValid())
91  B2FATAL("KLM electronics map is not available.");
92  if (!m_TimeConversion.isValid())
93  B2FATAL("EKLM time conversion parameters are not available.");
94  if (!m_eklmChannels.isValid())
95  B2FATAL("EKLM channel data are not available.");
96  if (m_loadThresholdFromDB) {
97  if (!m_bklmADCParams.isValid())
98  B2FATAL("BKLM ADC threshold paramenters are not available.");
99  m_scintThreshold = m_bklmADCParams->getADCThreshold();
100  }
101  m_triggerCTimeOfPreviousEvent = 0;
102 }
103 
104 void KLMUnpackerModule::createDigit(
105  const KLM::RawData* raw, const KLMDigitRaw* klmDigitRaw,
106  KLMDigitEventInfo* klmDigitEventInfo, int subdetector, int section,
107  int sector, int layer, int plane, int strip, int lastStrip)
108 {
109  KLMDigit* klmDigit = m_Digits.appendNew();
110  klmDigit->addRelationTo(klmDigitEventInfo);
111  if (m_WriteDigitRaws)
112  klmDigit->addRelationTo(klmDigitRaw);
113  bool isRPC = (subdetector == KLMElementNumbers::c_BKLM) &&
114  (layer >= BKLMElementNumbers::c_FirstRPCLayer);
115  if (isRPC) {
116  /*
117  * For RPC hits, digitize both the coarse (ctime) and fine (tdc) times
118  * relative to the revo9 trigger time rather than the event header's
119  * TriggerCTime. For the fine-time (tdc) measurement (11 bits), shift
120  * the revo9Trig time by 10 ticks to align the new prompt-time peak
121  * with the TriggerCtime-relative peak.
122  */
123  klmDigitEventInfo->increaseRPCHits();
124  float triggerTime = klmDigitEventInfo->getRevo9TriggerWord();
125  std::pair<int, double> rpcTimes =
126  m_TimeConversion->getRPCTimes(raw->getCTime(), raw->getTDC(),
127  triggerTime);
128  klmDigit->setTime(rpcTimes.second);
129  } else {
130  /*
131  * For scintillator hits, store the ctime relative to the event header's
132  * trigger ctime.
133  */
134  klmDigitEventInfo->increaseSciHits();
135  double time = m_TimeConversion->getScintillatorTime(
136  raw->getCTime(), klmDigitEventInfo->getTriggerCTime());
137  klmDigit->setTime(time);
138  if (subdetector == KLMElementNumbers::c_BKLM) {
139  if (raw->getCharge() < m_scintThreshold)
140  klmDigit->setFitStatus(KLM::c_ScintillatorFirmwareSuccessfulFit);
141  else
142  klmDigit->setFitStatus(KLM::c_ScintillatorFirmwareNoSignal);
143  } else {
144  int stripGlobal = m_eklmElementNumbers->stripNumber(
145  section, layer, sector, plane, strip);
146  const EKLMChannelData* channelData =
147  m_eklmChannels->getChannelData(stripGlobal);
148  if (channelData == nullptr)
149  B2FATAL("Incomplete EKLM channel data.");
150  if (raw->getCharge() < channelData->getThreshold())
151  klmDigit->setFitStatus(KLM::c_ScintillatorFirmwareSuccessfulFit);
152  else
153  klmDigit->setFitStatus(KLM::c_ScintillatorFirmwareNoSignal);
154  }
155  }
156  klmDigit->setSubdetector(subdetector);
157  klmDigit->setSection(section);
158  klmDigit->setLayer(layer);
159  klmDigit->setSector(sector);
160  klmDigit->setPlane(plane);
161  klmDigit->setStrip(strip);
162  if (lastStrip > 0)
163  klmDigit->setLastStrip(lastStrip);
164  klmDigit->setCharge(raw->getCharge());
165  klmDigit->setCTime(raw->getCTime());
166  klmDigit->setTDC(raw->getTDC());
167 }
168 
169 void KLMUnpackerModule::unpackKLMDigit(
170  const int* rawData, int copper, int hslb, int daqSubdetector,
171  KLMDigitEventInfo* klmDigitEventInfo)
172 {
173  KLMDigitRaw* klmDigitRaw;
174  KLM::RawData raw(copper, hslb + 1, rawData,
175  &m_klmDigitRaws, &klmDigitRaw, m_WriteDigitRaws);
176  const uint16_t* detectorChannel;
177  int subdetector, section, sector, layer, plane, strip;
178  /* Get channel groups. */
179  std::vector<ChannelGroup> channelGroups;
180  raw.getChannelGroups(channelGroups);
181  /* Get detector channels. */
182  KLMElectronicsChannel electronicsChannel(
183  copper, hslb + 1, raw.getLane(), raw.getAxis(), raw.getChannel());
184  bool channelFound = false;
185  for (ChannelGroup& channelGroup : channelGroups) {
186  if (channelGroup.lastChannel == 0) {
187  /* Single-strip hit. */
188  detectorChannel =
189  m_ElectronicsMap->getDetectorChannel(&electronicsChannel);
190  if (detectorChannel != nullptr) {
191  /* The channel is found, get element numbers. */
192  channelFound = true;
193  m_ElementNumbers->channelNumberToElementNumbers(
194  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
195  &strip);
196  channelGroup.firstStrip = strip;
197  channelGroup.lastStrip = 0;
198  } else {
199  /* The channel is not found, print error message. */
200  if (!(m_IgnoreWrongHits || (raw.getChannel() == 0 && m_IgnoreStrip0))) {
201  if (daqSubdetector == KLMElementNumbers::c_BKLM) {
202  B2DEBUG(20, "Channel does not exist in the KLM electronics map."
203  << LogVar("Copper", electronicsChannel.getCopper())
204  << LogVar("Slot", electronicsChannel.getSlot())
205  << LogVar("Lane", electronicsChannel.getLane())
206  << LogVar("Axis", electronicsChannel.getAxis())
207  << LogVar("Channel", electronicsChannel.getChannel()));
208  } else {
209  B2ERROR("Channel does not exist in the KLM electronics map."
210  << LogVar("Copper", electronicsChannel.getCopper())
211  << LogVar("Slot", electronicsChannel.getSlot())
212  << LogVar("Lane", electronicsChannel.getLane())
213  << LogVar("Axis", electronicsChannel.getAxis())
214  << LogVar("Channel", electronicsChannel.getChannel()));
215  }
216  }
217  }
218  } else {
219  /* Do not process multiple-strip hit in the electronics map debug mode. */
220  if (m_DebugElectronicsMap)
221  return;
222  /*
223  * Multiple-strip hit. It is necessary to find matching detector channels
224  * for all DAQ channels, because all channels in the group may not
225  * be necessary connected to strips in case of BKLM.
226  */
227  bool firstMatchedChannel = true;
228  for (int channel = channelGroup.firstChannel;
229  channel <= channelGroup.lastChannel; ++channel) {
230  electronicsChannel.setChannel(channel);
231  detectorChannel =
232  m_ElectronicsMap->getDetectorChannel(&electronicsChannel);
233  /* The channel is found, get element numbers. */
234  if (detectorChannel != nullptr) {
235  channelFound = true;
236  m_ElementNumbers->channelNumberToElementNumbers(
237  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
238  &strip);
239  if (firstMatchedChannel) {
240  firstMatchedChannel = false;
241  channelGroup.firstStrip = strip;
242  channelGroup.lastStrip = 0;
243  } else {
244  channelGroup.lastStrip = strip;
245  }
246  }
247  }
248  /* No matches found for this group at all. */
249  if (firstMatchedChannel) {
250  B2DEBUG(20, "No matching channels exist in the KLM electronics map."
251  << LogVar("Copper", electronicsChannel.getCopper())
252  << LogVar("Slot", electronicsChannel.getSlot())
253  << LogVar("Lane", electronicsChannel.getLane())
254  << LogVar("Axis", electronicsChannel.getAxis())
255  << LogVar("First channel", channelGroup.firstChannel)
256  << LogVar("Last channel", channelGroup.lastChannel)
257  << LogVar("Trigger bits", raw.getTriggerBits()));
258  }
259  }
260  }
261  /* No detector channel is found. */
262  if (!channelFound) {
263  if (!(m_WriteWrongHits || m_DebugElectronicsMap))
264  return;
265  /*
266  * Try to find channel from the same plane.
267  * BKLM phi-plane channels may start from 3 or 5.
268  */
269  electronicsChannel.setChannel(5);
270  detectorChannel = m_ElectronicsMap->getDetectorChannel(&electronicsChannel);
271  if (detectorChannel == nullptr)
272  return;
273  /* The channel is found, store out-of-range digit. */
274  m_ElementNumbers->channelNumberToElementNumbers(
275  *detectorChannel, &subdetector, &section, &sector, &layer, &plane,
276  &strip);
277  if (m_WriteWrongHits) {
278  klmDigitEventInfo->increaseOutOfRangeHits();
279  KLMDigit* klmDigitOutOfRange =
280  m_klmDigitsOutOfRange.appendNew();
281  klmDigitOutOfRange->addRelationTo(klmDigitEventInfo);
282  if (m_WriteDigitRaws)
283  klmDigitOutOfRange->addRelationTo(klmDigitRaw);
284  klmDigitOutOfRange->setSubdetector(KLMElementNumbers::c_BKLM);
285  klmDigitOutOfRange->setSection(section);
286  klmDigitOutOfRange->setLayer(layer);
287  klmDigitOutOfRange->setSector(sector);
288  klmDigitOutOfRange->setPlane(plane);
289  klmDigitOutOfRange->setStrip(strip);
290  klmDigitOutOfRange->setCharge(raw.getCharge());
291  klmDigitOutOfRange->setCTime(raw.getCTime());
292  klmDigitOutOfRange->setTDC(raw.getTDC());
293  return;
294  }
295  }
296  /* Debug mode: write raw channel number to strip number. */
297  if (m_DebugElectronicsMap) {
298  if (m_DAQChannelBKLMScintillators) {
299  if ((subdetector == KLMElementNumbers::c_BKLM) &&
300  (layer < BKLMElementNumbers::c_FirstRPCLayer))
301  strip = raw.getChannel();
302  }
303  if (m_DAQChannelModule >= 0) {
304  uint16_t klmModule =
305  m_ElementNumbers->moduleNumberByChannel(*detectorChannel);
306  if (klmModule == m_DAQChannelModule)
307  strip = raw.getChannel();
308  }
309  }
310  /* Create KLM digits. */
311  for (const ChannelGroup& channelGroup : channelGroups) {
312  if (channelGroup.firstStrip != 0) {
313  createDigit(&raw, klmDigitRaw, klmDigitEventInfo, subdetector, section,
314  sector, layer, plane, channelGroup.firstStrip,
315  channelGroup.lastStrip);
316  }
317  }
318 }
319 
320 void KLMUnpackerModule::convertPCIe40ToCOPPER(int channel, unsigned int* copper, int* hslb) const
321 {
322  if (channel >= 0 && channel < 16) {
323  int id = channel / 4;
324  *copper = BKLM_ID + id + 1;
325  *hslb = channel - id * 4;
326  } else if (channel >= 16 && channel < 32) {
327  int id = (channel - 16) / 4;
328  *copper = EKLM_ID + id + 1;
329  *hslb = (channel - 16) - id * 4;
330  } else
331  B2FATAL("The PCIe40 channel is invalid."
332  << LogVar("Channel", channel));
333 }
334 
335 void KLMUnpackerModule::event()
336 {
337  /*
338  * Length of one hit in 4-byte words. This is needed to find the hits in the
339  * detector buffer.
340  */
341  const int hitLength = 2;
342  for (int i = 0; i < m_RawKLMs.getEntries(); i++) {
343  if (m_RawKLMs[i]->GetNumEvents() != 1) {
344  B2ERROR("RawKLM a wrong number of entries (should be 1)."
345  << LogVar("RawKLM index", i)
346  << LogVar("Number of entries", m_RawKLMs[i]->GetNumEvents()));
347  continue;
348  }
349  /*
350  * GetNumEntries is defined in RawDataBlock.h and gives the numberOfNodes*numberOfEvents.
351  * Number of nodes is the number of COPPER boards.
352  */
353  for (int j = 0; j < m_RawKLMs[i]->GetNumEntries(); j++) {
354  unsigned int copper = m_RawKLMs[i]->GetNodeID(j);
355  int hslb, subdetector;
356  m_RawKLMs[i]->GetBuffer(j);
357  for (int channelReadoutBoard = 0; channelReadoutBoard < m_RawKLMs[i]->GetMaxNumOfCh(j); channelReadoutBoard++) {
358  if (m_RawKLMs[i]->GetMaxNumOfCh(j) == 4) { // COPPER data
359  hslb = channelReadoutBoard;
360  if ((copper >= EKLM_ID) && (copper <= EKLM_ID + 4))
361  subdetector = KLMElementNumbers::c_EKLM;
362  else if ((copper >= BKLM_ID) && (copper <= BKLM_ID + 4))
363  subdetector = KLMElementNumbers::c_BKLM;
364  else
365  continue;
366  } else if (m_RawKLMs[i]->GetMaxNumOfCh(j) == 48) { // PCIe40 data
367  if (channelReadoutBoard >= 0 && channelReadoutBoard < 16)
368  subdetector = KLMElementNumbers::c_BKLM;
369  else if (channelReadoutBoard >= 16 && channelReadoutBoard < 32)
370  subdetector = KLMElementNumbers::c_EKLM;
371  else
372  continue;
373  convertPCIe40ToCOPPER(channelReadoutBoard, &copper, &hslb);
374  } else {
375  B2FATAL("The maximum number of channels per readout board is invalid."
376  << LogVar("Number of channels", m_RawKLMs[i]->GetMaxNumOfCh(j)));
377  }
378  KLMDigitEventInfo* klmDigitEventInfo =
379  m_DigitEventInfos.appendNew(m_RawKLMs[i], j);
380  klmDigitEventInfo->setPreviousEventTriggerCTime(
381  m_triggerCTimeOfPreviousEvent);
382  m_triggerCTimeOfPreviousEvent = klmDigitEventInfo->getTriggerCTime();
383  int numDetNwords = m_RawKLMs[i]->GetDetectorNwords(j, channelReadoutBoard);
384  int* hslbBuffer = m_RawKLMs[i]->GetDetectorBuffer(j, channelReadoutBoard);
385  int numHits = numDetNwords / hitLength;
386  if (numDetNwords % hitLength != 1 && numDetNwords != 0) {
387  B2ERROR("Incorrect number of data words."
388  << LogVar("Number of data words", numDetNwords));
389  if (!m_keepEvenPackages)
390  continue;
391  }
392  if (numDetNwords > 0) {
393  /*
394  * In the last word there are the revo9 trigger word
395  * and the the user word (both from DCs).
396  */
397  unsigned int revo9TriggerWord =
398  (hslbBuffer[numDetNwords - 1] >> 16) & 0xFFFF;
399  klmDigitEventInfo->setRevo9TriggerWord(revo9TriggerWord);
400  unsigned int userWord = hslbBuffer[numDetNwords - 1] & 0xFFFF;
401  klmDigitEventInfo->setUserWord(userWord);
402  } else {
403  klmDigitEventInfo->setRevo9TriggerWord(0);
404  klmDigitEventInfo->setUserWord(0);
405  }
406  for (int iHit = 0; iHit < numHits; iHit++) {
407  unpackKLMDigit(&hslbBuffer[iHit * hitLength], copper, hslb,
408  subdetector, klmDigitEventInfo);
409  }
410  }
411  }
412  }
413 }
414 
415 void KLMUnpackerModule::endRun()
416 {
417 }
418 
419 void KLMUnpackerModule::terminate()
420 {
421 }
Belle2::KLMDigitEventInfo::setPreviousEventTriggerCTime
void setPreviousEventTriggerCTime(unsigned int triggerCTimeOfPreviousEvent)
Set trigger CTime of previous event.
Definition: KLMDigitEventInfo.h:122
Belle2::KLMDigitEventInfo::increaseOutOfRangeHits
void increaseOutOfRangeHits()
Increase by 1 the number of outOfRange-flagged hits in the event.
Definition: KLMDigitEventInfo.h:223
Belle2::EKLMChannelData
EKLM channel data.
Definition: EKLMChannelData.h:33
Belle2::KLMElectronicsChannel::getChannel
int getChannel() const
Get channel.
Definition: KLMElectronicsChannel.h:145
Belle2::KLMDigitEventInfo::getRevo9TriggerWord
unsigned int getRevo9TriggerWord() const
Get revo9 trigger word (from DCs).
Definition: KLMDigitEventInfo.h:240
Belle2::EKLMElementNumbers
EKLM element numbers.
Definition: EKLMElementNumbers.h:34
REG_MODULE
#define REG_MODULE(moduleName)
Register the given module (without 'Module' suffix) with the framework.
Definition: Module.h:652
Belle2::KLMDigit::setSubdetector
void setSubdetector(int subdetector)
Set subdetector number.
Definition: KLMDigit.h:98
Belle2::KLMDigit::setCharge
void setCharge(uint16_t charge)
Set charge.
Definition: KLMDigit.h:248
Belle2::KLMElectronicsChannel::getSlot
int getSlot() const
Get slot.
Definition: KLMElectronicsChannel.h:94
Belle2::KLMDigit::setStrip
void setStrip(int strip)
Set strip number.
Definition: KLMDigit.h:188
Belle2::RelationsInterface::addRelationTo
void addRelationTo(const RelationsInterface< BASE > *object, float weight=1.0, const std::string &namedRelation="") const
Add a relation from this object to another object (with caching).
Definition: RelationsObject.h:144
Belle2::KLMDigitEventInfo::increaseSciHits
void increaseSciHits()
Increase by 1 the number of scintillator hits in the event.
Definition: KLMDigitEventInfo.h:198
Belle2::KLMElectronicsChannel::setChannel
void setChannel(int channel)
Set channel.
Definition: KLMElectronicsChannel.h:154
Belle2::KLMDigitEventInfo::setUserWord
void setUserWord(unsigned int userWord)
Set user word (from DCs).
Definition: KLMDigitEventInfo.h:267
Belle2::KLMDigitEventInfo::getTriggerCTime
unsigned int getTriggerCTime() const
Get trigger CTIME.
Definition: KLMDigitEventInfo.h:76
Belle2::EKLMChannelData::getThreshold
int getThreshold() const
Get threshold.
Definition: EKLMChannelData.h:82
Belle2::KLMDigit
KLM digit (class representing a digitized hit in RPCs or scintillators).
Definition: KLMDigit.h:40
Belle2::KLMDigit::setPlane
void setPlane(int plane)
Set plane number.
Definition: KLMDigit.h:170
Belle2::KLMDigit::setTime
void setTime(float time)
Set hit time.
Definition: KLMDigit.h:302
Belle2::KLMElectronicsChannel::getAxis
int getAxis() const
Get axis.
Definition: KLMElectronicsChannel.h:128
Belle2::KLMDigitEventInfo
Class to store debugging informations from the unpacker (event based).
Definition: KLMDigitEventInfo.h:34
Belle2::Module
Base class for Modules.
Definition: Module.h:74
Belle2::KLMDigit::setTDC
void setTDC(uint16_t tdc)
Set TDC.
Definition: KLMDigit.h:284
Belle2::KLMDigitEventInfo::setRevo9TriggerWord
void setRevo9TriggerWord(unsigned int revo9TriggerWord)
Set Revo9 trigger word (from DCs).
Definition: KLMDigitEventInfo.h:249
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KLMDigit::setSector
void setSector(int sector)
Set sector number.
Definition: KLMDigit.h:134
LogVar
Class to store variables with their name which were sent to the logging service.
Definition: LogVariableStream.h:24
Belle2::KLMDigitRaw
Class to store the raw words from the unpacker, digit-by-digit.
Definition: KLMDigitRaw.h:39
Belle2::KLMDigit::setLastStrip
void setLastStrip(int lastStrip)
Set last strip number (for multi-strip digits).
Definition: KLMDigit.h:206
Belle2::KLMDigit::setFitStatus
void setFitStatus(int s)
Set fit status.
Definition: KLMDigit.h:383
Belle2::KLMElectronicsChannel::getCopper
int getCopper() const
Get copper.
Definition: KLMElectronicsChannel.h:77
Belle2::KLMElectronicsChannel
BKLM electronics channel.
Definition: KLMElectronicsChannel.h:33
Belle2::KLM::ChannelGroup
Channel group.
Definition: RawData.h:38
Belle2::KLMElementNumbers
KLM element numbers.
Definition: KLMElementNumbers.h:37
Belle2::KLMElectronicsChannel::getLane
int getLane() const
Get lane.
Definition: KLMElectronicsChannel.h:111
Belle2::KLMDigit::setLayer
void setLayer(int layer)
Set layer number.
Definition: KLMDigit.h:152
Belle2::KLMDigitEventInfo::increaseRPCHits
void increaseRPCHits()
Increase by 1 the number of RPC hits in the event.
Definition: KLMDigitEventInfo.h:173
Belle2::KLM::RawData
KLM raw data.
Definition: RawData.h:57
Belle2::KLMDigit::setSection
void setSection(int section)
Set section number.
Definition: KLMDigit.h:116
Belle2::KLMDigit::setCTime
void setCTime(uint16_t ctime)
Set CTIME.
Definition: KLMDigit.h:266
Belle2::KLMUnpackerModule
KLM unpacker.
Definition: KLMUnpackerModule.h:53