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