Project
Loading...
Searching...
No Matches
RawToCellConverterSpec.cxx
Go to the documentation of this file.
1// Copyright 2019-2020 CERN and copyright holders of ALICE O2.
2// See https://alice-o2.web.cern.ch/copyright for details of the copyright holders.
3// All rights not expressly granted are reserved.
4//
5// This software is distributed under the terms of the GNU General Public
6// License v3 (GPL Version 3), copied verbatim in the file "COPYING".
7//
8// In applying this license CERN does not waive the privileges and immunities
9// granted to it by virtue of its status as an Intergovernmental Organization
10// or submit itself to any jurisdiction.
11#include <string>
12#include <iomanip>
13#include <iostream>
14#include <bitset>
15#include <set>
16
17#include <InfoLogger/InfoLogger.hxx>
18
26#include "Framework/Logger.h"
33#include "EMCALBase/Geometry.h"
34#include "EMCALBase/Mapper.h"
36#include "EMCALCalib/FeeDCS.h"
48
49using namespace o2::emcal::reco_workflow;
50
52{
53 if (mErrorMessagesSuppressed) {
54 LOG(warning) << "Suppressed further " << mErrorMessagesSuppressed << " error messages";
55 }
56}
57
59{
60 if (ctx.services().active<AliceO2::InfoLogger::InfoLoggerContext>()) {
61 auto& ilctx = ctx.services().get<AliceO2::InfoLogger::InfoLoggerContext>();
62 ilctx.setField(AliceO2::InfoLogger::InfoLoggerContext::FieldName::Detector, "EMC");
63 }
64
65 LOG(debug) << "[EMCALRawToCellConverter - init] Initialize converter ";
66 if (!mGeometry) {
67 mGeometry = Geometry::GetInstanceFromRunNumber(223409);
68 }
69 if (!mGeometry) {
70 LOG(error) << "Failure accessing geometry";
71 }
72
73 if (!mMapper) {
74 mMapper = std::unique_ptr<MappingHandler>(new o2::emcal::MappingHandler);
75 }
76 if (!mMapper) {
77 LOG(error) << "Failed to initialize mapper";
78 }
79
80 if (!mTriggerMapping) {
81 mTriggerMapping = std::make_unique<TriggerMappingV2>(mGeometry);
82 }
83
84 auto fitmethod = ctx.options().get<std::string>("fitmethod");
85 if (fitmethod == "standard") {
86 LOG(info) << "Using standard raw fitter";
87 mRawFitter = std::unique_ptr<CaloRawFitter>(new o2::emcal::CaloRawFitterStandard);
88 } else if (fitmethod == "gamma2") {
89 LOG(info) << "Using gamma2 raw fitter";
90 mRawFitter = std::unique_ptr<CaloRawFitter>(new o2::emcal::CaloRawFitterGamma2);
91 } else {
92 LOG(fatal) << "Unknown fit method" << fitmethod;
93 }
94 LOG(info) << "Creating decoding errors: " << (mCreateRawDataErrors ? "yes" : "no");
95
96 mPrintTrailer = ctx.options().get<bool>("printtrailer");
97
98 mMaxErrorMessages = ctx.options().get<int>("maxmessage");
99 LOG(info) << "Suppressing error messages after " << mMaxErrorMessages << " messages";
100
101 mMergeLGHG = !ctx.options().get<bool>("no-mergeHGLG");
102 mDisablePedestalEvaluation = ctx.options().get<bool>("no-evalpedestal");
103 mActiveLinkCheck = !ctx.options().get<bool>("no-checkactivelinks");
104
105 LOG(info) << "Running gain merging mode: " << (mMergeLGHG ? "yes" : "no");
106 LOG(info) << "Checking for active links: " << (mActiveLinkCheck ? "yes" : "no");
107 LOG(info) << "Calculate pedestals: " << (mDisablePedestalEvaluation ? "no" : "yes");
108 LOG(info) << "Using L0LM delay: " << o2::ctp::TriggerOffsetsParam::Instance().LM_L0 << " BCs";
109
110 mRawFitter->setAmpCut(mNoiseThreshold);
111 mRawFitter->setL1Phase(0.);
112}
113
115{
116 LOG(debug) << "[EMCALRawToCellConverter - run] called";
117 mCalibHandler->checkUpdates(ctx);
118 updateCalibrationObjects();
119
120 // container with BCid and feeID
121 std::unordered_map<int64_t, std::bitset<46>> bcFreq;
122
123 double timeshift = RecoParam::Instance().getCellTimeShiftNanoSec(); // subtract offset in ns in order to center the time peak around the nominal delay
124 auto maxBunchLengthRP = RecoParam::Instance().getMaxAllowedBunchLength(); // exclude bunches where either the start time or the bunch length is above the expected maximum
125 constexpr auto originEMC = o2::header::gDataOriginEMC;
126 constexpr auto descRaw = o2::header::gDataDescriptionRawData;
127
128 // reset message counter after 10 minutes
129 auto currenttime = std::chrono::system_clock::now();
130 std::chrono::duration<double> interval = mReferenceTime - currenttime;
131 if (interval.count() > 600) {
132 // Resetting error messages after 10 minutes
133 mNumErrorMessages = 0;
134 mReferenceTime = currenttime;
135 }
136
137 mOutputCells.clear();
138 mOutputTriggerRecords.clear();
139 mOutputDecoderErrors.clear();
140 mOutputTRUs.clear();
141 mOutputTRUTriggerRecords.clear();
142 mOutputPatches.clear();
143 mOutputPatchTriggerRecords.clear();
144 mOutputTimesums.clear();
145 mOutputTimesumTriggerRecords.clear();
146
147 if (isLostTimeframe(ctx)) {
148 sendData(ctx);
149 return;
150 }
151
152 mCellHandler.reset();
153
154 // Get the first orbit of the timeframe later used to check whether the corrected
155 // BC is within the timeframe
156 const auto tfOrbitFirst = ctx.services().get<o2::framework::TimingInfo>().firstTForbit;
157 auto lml0delay = o2::ctp::TriggerOffsetsParam::Instance().LM_L0;
158
159 std::vector<framework::InputSpec> filter{{"filter", framework::ConcreteDataTypeMatcher(originEMC, descRaw)}};
160 int firstEntry = 0;
161 for (const auto& rawData : framework::InputRecordWalker(ctx.inputs(), filter)) {
162 // Skip SOX headers
163 auto rdhblock = reinterpret_cast<const o2::header::RDHAny*>(rawData.payload);
164 if (o2::raw::RDHUtils::getHeaderSize(rdhblock) == static_cast<int>(o2::framework::DataRefUtils::getPayloadSize(rawData))) {
165 continue;
166 }
167
168 o2::emcal::RawReaderMemory rawreader(framework::DataRefUtils::as<const char>(rawData));
169 rawreader.setRangeSRUDDLs(0, 39);
170
171 // loop over all the DMA pages
172 while (rawreader.hasNext()) {
173 try {
174 rawreader.next();
175 } catch (RawDecodingError& e) {
176 handlePageError(e);
178 // We must break in case of header decoding as the offset to the next payload is lost
179 // consequently the parser does not know where to continue leading to an infinity loop
180 break;
181 }
182 // We must skip the page as payload is not consistent
183 // otherwise the next functions will rethrow the exceptions as
184 // the page format does not follow the expected format
185 continue;
186 }
187 for (auto& e : rawreader.getMinorErrors()) {
188 handleMinorPageError(e);
189 // For minor errors we do not need to skip the page, just print and send the error to the QC
190 }
191
192 auto& header = rawreader.getRawHeader();
193 auto triggerBC = raw::RDHUtils::getTriggerBC(header);
194 auto triggerOrbit = raw::RDHUtils::getTriggerOrbit(header);
195 auto feeID = raw::RDHUtils::getFEEID(header);
196 auto triggerbits = raw::RDHUtils::getTriggerType(header);
197
198 int correctionShiftBCmod4 = 0;
199 o2::InteractionRecord currentIR(triggerBC, triggerOrbit);
200 // Correct physics triggers for the shift of the BC due to the LM-L0 delay
201 if (triggerbits & o2::trigger::PhT) {
202 if (currentIR.differenceInBC({0, tfOrbitFirst}) >= lml0delay) {
203 currentIR -= lml0delay; // guaranteed to stay in the TF containing the collision
204 // in case we correct for the L0LM delay we need to adjust the BC mod 4, because if the L0LM delay % 4 != 0 it will change the permutation of trigger peaks
205 // we need to add back the correction we applied % 4 to the corrected BC during the correction of the cell time in order to keep the same permutation
206 correctionShiftBCmod4 = lml0delay % 4;
207 } else {
208 // discard the data associated with this IR as it was triggered before the start of timeframe
209 continue;
210 }
211 }
212
213 bcFreq[currentIR.toLong()].set(feeID, true);
214
215 // Correct the cell time for the bc mod 4 (LHC: 40 MHz clock - ALTRO: 10 MHz clock)
216 // Convention: All times shifted with respect to BC % 4 = 0 for trigger BC
217 // Attention: Correction only works for the permutation (0 1 2 3) of the BC % 4, if the permutation is
218 // different the BC for the correction has to be shifted by n BCs to obtain permutation (0 1 2 3)
219 // We apply here the following shifts:
220 // - correction for the L0-LM delay mod 4 in order to restore the original ordering of the BCs mod 4
221 // - phase shift in order to adjust for permutations different from (0 1 2 3)
222 int bcmod4 = (currentIR.bc + correctionShiftBCmod4 + RecoParam::Instance().getPhaseBCmod4()) % 4;
223 LOG(debug) << "Original BC " << triggerBC << ", L0LM corrected " << currentIR.bc;
224 LOG(debug) << "Applying correction for LM delay: " << correctionShiftBCmod4;
225 LOG(debug) << "BC mod original: " << triggerBC % 4 << ", corrected " << bcmod4;
226 LOG(debug) << "Applying time correction: " << -1 * 25 * bcmod4;
227 auto& currentEvent = mCellHandler.getEventContainer(currentIR);
228 if (!currentEvent.getTriggerBits()) {
229 currentEvent.setTriggerBits(triggerbits);
230 }
231 CellTimeCorrection timeCorrector{timeshift, bcmod4};
232
233 if (feeID >= 40) {
234 continue; // skip STU ddl
235 }
236
237 // std::cout<<rawreader.getRawHeader()<<std::endl;
238
239 // use the altro decoder to decode the raw data, and extract the RCU trailer
240 AltroDecoder decoder(rawreader);
241 if (maxBunchLengthRP) {
242 // apply user-defined max. bunch length
243 decoder.setMaxBunchLength(maxBunchLengthRP);
244 }
245 // check the words of the payload exception in altrodecoder
246 try {
247 decoder.decode();
248 } catch (AltroDecoderError& e) {
249 handleAltroError(e, feeID);
250 continue;
251 }
252 for (const auto& minorerror : decoder.getMinorDecodingErrors()) {
253 handleMinorAltroError(minorerror, feeID);
254 }
255
256 if (mPrintTrailer) {
257 // Can become very verbose, therefore must be switched on explicitly in addition
258 // to high debug level
259 LOG(debug4) << decoder.getRCUTrailer();
260 }
261 // Apply zero suppression only in case it was enabled
262 if (decoder.getRCUTrailer().hasZeroSuppression()) {
263 LOG(debug3) << "Zero suppression enabled";
264 } else {
265 LOG(debug3) << "Zero suppression disabled";
266 }
267 if (mDisablePedestalEvaluation) {
268 // auto-disable pedestal evaluation in the raw fitter
269 // treat all channels as zero-suppressed independent of
270 // what is provided from the RCU trailer
271 mRawFitter->setIsZeroSuppressed(true);
272 } else {
273 mRawFitter->setIsZeroSuppressed(decoder.getRCUTrailer().hasZeroSuppression());
274 }
275
276 try {
277
278 const auto& map = mMapper->getMappingForDDL(feeID);
279 uint16_t iSM = feeID / 2;
280
281 // Loop over all the channels
282 int nBunchesNotOK = 0;
283 for (auto& chan : decoder.getChannels()) {
284 try {
285 auto iRow = map.getRow(chan.getHardwareAddress());
286 auto iCol = map.getColumn(chan.getHardwareAddress());
287 auto chantype = map.getChannelType(chan.getHardwareAddress());
288 LocalPosition channelPosition{iSM, feeID, iCol, iRow};
289 switch (chantype) {
292 addFEEChannelToEvent(currentEvent, chan, timeCorrector, channelPosition, chantype);
293 break;
295 // Drop LEDMON reconstruction in case of physics triggers
296 if (triggerbits & o2::trigger::Cal) {
297 addFEEChannelToEvent(currentEvent, chan, timeCorrector, channelPosition, chantype);
298 }
299 break;
301 addTRUChannelToEvent(currentEvent, chan, channelPosition);
302 break;
303 default:
304 LOG(error) << "Unknown channel type for HW address " << chan.getHardwareAddress();
305 break;
306 }
308 handleAddressError(ex, feeID, chan.getHardwareAddress());
309 continue;
310 }
311 }
312 } catch (o2::emcal::MappingHandler::DDLInvalid& ddlerror) {
313 // Unable to catch mapping
314 handleDDLError(ddlerror, feeID);
315 }
316 }
317 }
318
319 std::bitset<46> bitSetActiveLinks;
320 if (mActiveLinkCheck) {
321 // build expected active mask from DCS
322 FeeDCS* feedcs = mCalibHandler->getFEEDCS();
323 auto list0 = feedcs->getDDLlist0();
324 auto list1 = feedcs->getDDLlist1();
325 // links 21 and 39 do not exist, but they are set active in DCS
326 list0.set(21, false);
327 list1.set(7, false);
328 // must be 0x307FFFDFFFFF if all links are active
329 bitSetActiveLinks = std::bitset<46>((list1.to_ullong() << 32) + list0.to_ullong());
330 }
331
332 // Loop over BCs, sort cells with increasing tower ID and write to output containers
333 RecoContainerReader eventIterator(mCellHandler);
334 while (eventIterator.hasNext()) {
335 int ncellsEvent = 0, nLEDMONsEvent = 0;
336 int eventstart = mOutputCells.size();
337 auto& currentevent = eventIterator.nextEvent();
338 const auto interaction = currentevent.getInteractionRecord();
339 if (mActiveLinkCheck) {
340 // check for current event if all links are present
341 // discard event if not all links are present
342 auto bcfreqFound = bcFreq.find(interaction.toLong());
343 if (bcfreqFound != bcFreq.end()) {
344 const auto& activelinks = bcfreqFound->second;
345 if (activelinks != bitSetActiveLinks) {
346 static int nErrors = 0;
347 if (nErrors++ < 3) {
348 LOG(error) << "Not all EMC active links contributed in global BCid=" << interaction.toLong() << ": mask=" << (activelinks ^ bitSetActiveLinks) << (nErrors == 3 ? " (not reporting further errors to avoid spamming)" : "");
349 }
350 if (mCreateRawDataErrors) {
351 for (std::size_t ilink = 0; ilink < bitSetActiveLinks.size(); ilink++) {
352 if (!bitSetActiveLinks.test(ilink)) {
353 continue;
354 }
355 if (!activelinks.test(ilink)) {
356 mOutputDecoderErrors.emplace_back(ilink, ErrorTypeFEE::ErrorSource_t::LINK_ERROR, 0, -1, -1);
357 }
358 }
359 }
360 // discard event
361 // create empty trigger record with dedicated trigger bit marking as rejected
362 mOutputTriggerRecords.emplace_back(interaction, currentevent.getTriggerBits() | o2::emcal::triggerbits::Inc, eventstart, 0);
363 continue;
364 }
365 }
366 }
367 // Add cells
368 if (currentevent.getNumberOfCells()) {
369 LOG(debug) << "Event has " << currentevent.getNumberOfCells() << " cells";
370 currentevent.sortCells(false);
371 ncellsEvent = bookEventCells(currentevent.getCells(), false);
372 }
373 // Add LEDMONs (if present)
374 if (currentevent.getNumberOfLEDMONs()) {
375 LOG(debug) << "Event has " << currentevent.getNumberOfLEDMONs() << " LEDMONs";
376 currentevent.sortCells(true);
377 nLEDMONsEvent = bookEventCells(currentevent.getLEDMons(), true);
378 }
379 LOG(debug) << "Next event [Orbit " << interaction.orbit << ", BC (" << interaction.bc << "]: Accepted " << ncellsEvent << " cells and " << nLEDMONsEvent << " LEDMONS";
380 mOutputTriggerRecords.emplace_back(interaction, currentevent.getTriggerBits(), eventstart, ncellsEvent + nLEDMONsEvent);
381
382 // Add trigger data
383 if (mDoTriggerReconstruction) {
384 auto [trus, patches] = buildL0Patches(currentevent);
385 LOG(debug) << "Found " << patches.size() << " L0 patches from " << trus.size() << " TRUs";
386 auto trusstart = mOutputTRUs.size();
387 std::copy(trus.begin(), trus.end(), std::back_inserter(mOutputTRUs));
388 mOutputTRUTriggerRecords.emplace_back(interaction, currentevent.getTriggerBits(), trusstart, trus.size());
389 auto patchesstart = mOutputPatches.size();
390 std::copy(patches.begin(), patches.end(), std::back_inserter(mOutputPatches));
391 mOutputPatchTriggerRecords.emplace_back(interaction, currentevent.getTriggerBits(), patchesstart, patches.size());
392 // For L0 timesums use fixed time, across TRUs and triggers, determined from the patch time QC
393 // average found to be - will be made configurable
394 auto timesumsstart = mOutputTimesums.size();
395 auto timesums = buildL0Timesums(currentevent, 8);
396 std::copy(timesums.begin(), timesums.end(), std::back_inserter(mOutputTimesums));
397 mOutputTimesumTriggerRecords.emplace_back(interaction, currentevent.getTriggerBits(), timesumsstart, timesums.size());
398 }
399 }
400
401 LOG(info) << "[EMCALRawToCellConverter - run] Writing " << mOutputCells.size() << " cells from " << mOutputTriggerRecords.size() << " events ...";
402 if (mDoTriggerReconstruction) {
403 LOG(info) << "[EMCALRawToCellConverter - run] Writing " << mOutputTRUs.size() << " TRU infos and " << mOutputPatches.size() << " trigger patches and " << mOutputTimesums.size() << " timesums from " << mOutputTRUTriggerRecords.size() << " events ...";
404 }
405 sendData(ctx);
406}
407
409{
410 if (mCalibHandler->finalizeCCDB(matcher, obj)) {
411 return;
412 }
413}
414
415void RawToCellConverterSpec::updateCalibrationObjects()
416{
417 if (mCalibHandler->hasUpdateRecoParam()) {
418 LOG(info) << "RecoParams updated";
419 o2::emcal::RecoParam::Instance().printKeyValues(true, true);
420 }
421 if (mCalibHandler->hasUpdateFEEDCS()) {
422 LOG(info) << "DCS params updated";
423 }
424}
425
426bool RawToCellConverterSpec::isLostTimeframe(framework::ProcessingContext& ctx) const
427{
428 constexpr auto originEMC = header::gDataOriginEMC;
429 o2::framework::InputSpec dummy{"dummy",
432 0xDEADBEEF}};
433 static size_t contDeadBeef = 0; // number of times 0xDEADBEEF was seen continuously
434 for (const auto& ref : o2::framework::InputRecordWalker(ctx.inputs(), {dummy})) {
435 const auto dh = o2::framework::DataRefUtils::getHeader<o2::header::DataHeader*>(ref);
437 if (payloadSize == 0) {
439 if (++contDeadBeef <= maxWarn) {
440 LOGP(warning, "Found input [{}/{}/{:#x}] TF#{} 1st_orbit:{} Payload {} : assuming no payload for all links in this TF{}",
441 dh->dataOrigin.str, dh->dataDescription.str, dh->subSpecification, dh->tfCounter, dh->firstTForbit, payloadSize,
442 contDeadBeef == maxWarn ? fmt::format(". {} such inputs in row received, stopping reporting", contDeadBeef) : "");
443 }
444 return true;
445 }
446 }
447 contDeadBeef = 0; // if good data, reset the counter
448 return false;
449}
450
451void RawToCellConverterSpec::addFEEChannelToEvent(o2::emcal::EventContainer& currentEvent, const o2::emcal::Channel& currentchannel, const CellTimeCorrection& timeCorrector, const LocalPosition& position, ChannelType_t chantype)
452{
453 int CellID = -1;
454 bool isLowGain = false;
455 try {
457 // high- / low-gain cell
458 CellID = getCellAbsID(position.mSupermoduleID, position.mColumn, position.mRow);
459
460 isLowGain = chantype == o2::emcal::ChannelType_t::LOW_GAIN;
461 } else {
462 CellID = geLEDMONAbsID(position.mSupermoduleID, position.mColumn); // Module index encoded in colum for LEDMONs
463 isLowGain = position.mRow == 0; // For LEDMONs gain type is encoded in the row (0 - low gain, 1 - high gain)
464 }
465 } catch (ModuleIndexException& e) {
466 handleGeometryError(e, position.mSupermoduleID, CellID, currentchannel.getHardwareAddress(), chantype);
467 return;
468 }
469
470 // define the conatiner for the fit results, and perform the raw fitting using the stadnard raw fitter
471 CaloFitResults fitResults;
472 try {
473 fitResults = mRawFitter->evaluate(currentchannel.getBunches());
474 // Prevent negative entries - we should no longer get here as the raw fit usually will end in an error state
475 if (fitResults.getAmp() < 0) {
476 fitResults.setAmp(0.);
477 }
478 if (fitResults.getTime() < 0) {
479 fitResults.setTime(0.);
480 }
481 // apply correction for bc mod 4
482 double celltime = timeCorrector.getCorrectedTime(fitResults.getTime());
483 double amp = fitResults.getAmp() * o2::emcal::constants::EMCAL_ADCENERGY;
484 if (isLowGain) {
485 amp *= o2::emcal::constants::EMCAL_HGLGFACTOR;
486 }
487 if (chantype == o2::emcal::ChannelType_t::LEDMON) {
488 // Mark LEDMONs as HIGH_GAIN/LOW_GAIN for gain type merging - will be flagged as LEDMON later when pushing to the output container
489 currentEvent.setLEDMONCell(CellID, amp, celltime, isLowGain ? o2::emcal::ChannelType_t::LOW_GAIN : o2::emcal::ChannelType_t::HIGH_GAIN, currentchannel.getHardwareAddress(), position.mFeeID, mMergeLGHG);
490 } else {
491 currentEvent.setCell(CellID, amp, celltime, chantype, currentchannel.getHardwareAddress(), position.mFeeID, mMergeLGHG);
492 }
493 } catch (CaloRawFitter::RawFitterError_t& fiterror) {
494 handleFitError(fiterror, position.mFeeID, CellID, currentchannel.getHardwareAddress());
495 }
496}
497
498void RawToCellConverterSpec::addTRUChannelToEvent(o2::emcal::EventContainer& currentEvent, const o2::emcal::Channel& currentchannel, const LocalPosition& position)
499{
500 try {
501 auto tru = mTriggerMapping->getTRUIndexFromOnlineHardareAddree(currentchannel.getHardwareAddress(), position.mFeeID, position.mSupermoduleID);
502 if (position.mColumn >= 96 && position.mColumn <= 105) {
503 auto& trudata = currentEvent.getTRUData(tru);
504 // Trigger patch information encoded columns 95-105
505 for (auto& bunch : currentchannel.getBunches()) {
506 LOG(debug) << "Found bunch of length " << static_cast<int>(bunch.getBunchLength()) << " with start time " << static_cast<int>(bunch.getStartTime()) << " (column " << static_cast<int>(position.mColumn) << ")";
507 auto l0time = bunch.getStartTime();
508 int isample = 0;
509 for (auto& adc : bunch.getADC()) {
510 // patch word might be in any of the samples, need to check all of them
511 // in case of colum 105 the first 6 bits are the patch word, the remaining 4 bits are the header word
512 if (adc == 0) {
513 isample++;
514 continue;
515 }
516 if (position.mColumn == 105) {
517 std::bitset<6> patchBits(adc & 0x3F);
518 std::bitset<4> headerbits((adc >> 6) & 0xF);
519 for (auto localindex = 0; localindex < patchBits.size(); localindex++) {
520 if (patchBits.test(localindex)) {
521 auto globalindex = (position.mColumn - 96) * 10 + localindex;
522 LOG(debug) << "Found patch with index " << globalindex << " in sample " << isample;
523 // std::cout << "Found patch with index " << globalindex << " in sample " << isample << " (" << (bunch.getStartTime() - isample) << ")" << std::endl;
524 try {
525 trudata.setPatch(globalindex, bunch.getStartTime() - isample);
527 handlePatchError(e, position.mFeeID, tru);
528 }
529 }
530 }
531 if (headerbits.test(2)) {
532 LOG(debug) << "TRU " << tru << ": Found TRU fired (" << tru << ") in sample " << isample;
533 // std::cout << "TRU " << tru << ": Found TRU fired (" << tru << ") in sample " << isample << " (" << (bunch.getStartTime() - isample) << ")" << std::endl;
534 trudata.setFired(true);
535 trudata.setL0time(bunch.getStartTime() - isample);
536 }
537 } else {
538 std::bitset<10> patchBits(adc & 0x3FF);
539 for (auto localindex = 0; localindex < patchBits.size(); localindex++) {
540 if (patchBits.test(localindex)) {
541 auto globalindex = (position.mColumn - 96) * 10 + localindex;
542 LOG(debug) << "TRU " << tru << ": Found patch with index " << globalindex << " in sample " << isample;
543 // std::cout << "TRU " << tru << ": Found patch with index " << globalindex << " in sample " << isample << " (" << (bunch.getStartTime() - isample) << ")" << std::endl;
544 try {
545 trudata.setPatch(globalindex, bunch.getStartTime() - isample);
547 handlePatchError(e, position.mFeeID, tru);
548 }
549 }
550 }
551 }
552 isample++;
553 }
554 }
555 } else {
556 try {
557 auto absFastOR = mTriggerMapping->getAbsFastORIndexFromIndexInTRU(tru, position.mColumn);
558 for (auto& bunch : currentchannel.getBunches()) {
559 // FastOR data reversed internally (positive in time direction)
560 // -> Start time marks the first timebin, consequently it must be also reversed.
561 // std::cout << "Adding non-reversed FastOR time series for FastOR " << absFastOR << " (TRU " << tru << ", index " << static_cast<int>(position.mColumn) << ") with start time " << static_cast<int>(bunch.getStartTime()) << " (reversed " << bunch.getStartTime() + 1 - bunch.getADC().size() << "): ";
562 // for (auto adc : bunch.getADC()) {
563 // std::cout << adc << ", ";
564 //}
565 // std::cout << std::endl;
566
567 try {
568 currentEvent.setFastOR(absFastOR, bunch.getStartTime(), bunch.getADC());
570 handleFastORStartTimeErrors(e, position.mFeeID, tru);
571 }
572 }
573 } catch (FastORIndexException& e) {
574 handleFastORErrors(e, position.mFeeID, tru);
575 }
576 }
577 } catch (TRUIndexException& e) {
578 handleTRUIndexError(e, position.mFeeID, currentchannel.getHardwareAddress());
579 }
580}
581
582std::tuple<RawToCellConverterSpec::TRUContainer, RawToCellConverterSpec::PatchContainer> RawToCellConverterSpec::buildL0Patches(const EventContainer& currentevent) const
583{
584 LOG(debug) << "Reconstructing patches for Orbit " << currentevent.getInteractionRecord().orbit << ", BC " << currentevent.getInteractionRecord().bc;
585 TRUContainer eventTRUs;
586 PatchContainer eventPatches;
587 auto& fastOrs = currentevent.getTimeSeriesContainer();
588 std::set<uint16_t> foundFastOrs;
589 for (auto fastor : fastOrs) {
590 foundFastOrs.insert(fastor.first);
591 }
592 for (std::size_t itru = 0; itru < TriggerMappingV2::ALLTRUS; itru++) {
593 auto& currenttru = currentevent.readTRUData(itru);
594 if (!currenttru.hasAnyPatch()) {
595 continue;
596 }
597 auto l0time = currenttru.getL0time();
598 LOG(debug) << "Found patches in TRU " << itru << ", fired: " << (currenttru.isFired() ? "yes" : "no") << ", L0 time " << static_cast<int>(l0time);
599 uint8_t npatches = 0;
600 for (auto ipatch = 0; ipatch < o2::emcal::TriggerMappingV2::PATCHESINTRU; ipatch++) {
601 if (currenttru.hasPatch(ipatch)) {
602 auto patchtime = currenttru.getPatchTime(ipatch);
603 LOG(debug) << "Found patch " << ipatch << " in TRU " << itru << " with time " << static_cast<int>(patchtime);
604 auto fastorStart = mTriggerMapping->getAbsFastORIndexFromIndexInTRU(itru, ipatch);
605 auto fastORs = mTriggerMapping->getFastORIndexFromL0Index(itru, ipatch, 4);
606 std::array<const FastORTimeSeries*, 4> fastors;
607 std::fill(fastors.begin(), fastors.end(), nullptr);
608 int indexFastorInTRU = 0;
609 for (auto fastor : fastORs) {
610 auto [truID, fastorTRU] = mTriggerMapping->getTRUFromAbsFastORIndex(fastor);
611 // std::cout << "Patch has abs FastOR " << fastor << " -> " << fastorTRU << " (in TRU)" << std::endl;
612 auto timeseriesFound = fastOrs.find(fastor);
613 if (timeseriesFound != fastOrs.end()) {
614 LOG(debug) << "Adding FastOR (" << indexFastorInTRU << ") with index " << fastor << " to patch";
615 fastors[indexFastorInTRU] = &(timeseriesFound->second);
616 indexFastorInTRU++;
617 }
618 }
619 auto [patchADC, recpatchtime] = reconstructTriggerPatch(fastors);
620 // Correct for bit shift 12->10 bits due to ALTRO format
621 patchADC = patchADC << 2;
622 LOG(debug) << "Reconstructed patch at index " << ipatch << " with peak time " << static_cast<int>(recpatchtime) << " (time sample " << static_cast<int>(patchtime) << ") and energy " << patchADC;
623 eventPatches.push_back({static_cast<uint8_t>(itru), static_cast<uint8_t>(ipatch), patchtime, patchADC});
624 }
625 }
626 eventTRUs.push_back({static_cast<uint8_t>(itru), static_cast<uint8_t>(l0time), currenttru.isFired(), npatches});
627 }
628 return std::make_tuple(eventTRUs, eventPatches);
629}
630
631std::vector<o2::emcal::CompressedL0TimeSum> RawToCellConverterSpec::buildL0Timesums(const o2::emcal::EventContainer& currentevent, uint8_t l0time) const
632{
633 std::vector<o2::emcal::CompressedL0TimeSum> timesums;
634 for (const auto& [fastorID, timeseries] : currentevent.getTimeSeriesContainer()) {
635 timesums.push_back({fastorID, static_cast<uint16_t>(timeseries.calculateL1TimeSum(l0time) << 2)});
636 }
637 return timesums;
638}
639
640std::tuple<uint16_t, uint8_t> RawToCellConverterSpec::reconstructTriggerPatch(const gsl::span<const FastORTimeSeries*> fastors) const
641{
642 constexpr size_t INTEGRATE_SAMPLES = 4,
643 MAX_SAMPLES = 12;
644 double maxpatchenergy = 0;
645 uint8_t foundtime = 0;
646 for (size_t itime = 0; itime < MAX_SAMPLES - INTEGRATE_SAMPLES; ++itime) {
647 double currenttimesum = 0;
648 for (size_t isample = 0; isample < INTEGRATE_SAMPLES; ++isample) {
649 for (auto ifastor = 0; ifastor < fastors.size(); ++ifastor) {
650 if (fastors[ifastor]) {
651 currenttimesum += fastors[ifastor]->getADCs()[itime + isample];
652 }
653 }
654 }
655 if (currenttimesum > maxpatchenergy) {
656 maxpatchenergy = currenttimesum;
657 foundtime = itime;
658 }
659 }
660
661 return std::make_tuple(maxpatchenergy, foundtime);
662}
663
664int RawToCellConverterSpec::bookEventCells(const gsl::span<const o2::emcal::RecCellInfo>& cells, bool isLELDMON)
665{
666 double noiseThresholLGnoHG = RecoParam::Instance().getNoiseThresholdLGnoHG();
667 int ncellsSelected = 0;
668 for (const auto& cell : cells) {
669 if (cell.mIsLGnoHG) {
670 // Treat error only in case the LG is above the noise threshold
671 // no HG cell found, we can assume the cell amplitude is the LG amplitude
672 int ampLG = cell.mCellData.getAmplitude() / (o2::emcal::constants::EMCAL_ADCENERGY * o2::emcal::constants::EMCAL_HGLGFACTOR);
673 // use cut at 3 sigma where sigma for the LG digitizer is 0.4 ADC counts (EMCAL-502)
674 if (ampLG > noiseThresholLGnoHG) {
675 handleGainError(reconstructionerrors::GainError_t::LGNOHG, cell.mDDLID, cell.mHWAddressLG);
676 }
677 continue;
678 }
679 if (cell.mHGOutOfRange) {
680 handleGainError(reconstructionerrors::GainError_t::HGNOLG, cell.mDDLID, cell.mHWAddressHG);
681 continue;
682 }
683 ncellsSelected++;
684 mOutputCells.push_back(cell.mCellData);
685 if (isLELDMON) {
686 // LEDMON was handled internally in the reco container as HG/LG cell for gain type merging - reflag as LEDMON
687 mOutputCells.back().setLEDMon();
688 }
689 }
690 return ncellsSelected;
691}
692
693int RawToCellConverterSpec::getCellAbsID(int supermoduleID, int column, int row)
694{
695 auto [phishift, etashift] = mGeometry->ShiftOnlineToOfflineCellIndexes(supermoduleID, row, column);
696 int cellID = mGeometry->GetAbsCellIdFromCellIndexes(supermoduleID, phishift, etashift);
697 if (cellID > 17664 || cellID < 0) {
698 throw ModuleIndexException(cellID, column, row, etashift, phishift);
699 }
700 return cellID;
701}
702
703int RawToCellConverterSpec::geLEDMONAbsID(int supermoduleID, int moduleID)
704{
705 if (moduleID >= o2::emcal::EMCAL_LEDREFS || moduleID < 0) {
706 throw ModuleIndexException(moduleID);
707 }
708 return supermoduleID * o2::emcal::EMCAL_LEDREFS + moduleID;
709}
710
711void RawToCellConverterSpec::handleAddressError(const Mapper::AddressNotFoundException& error, int feeID, int hwaddress)
712{
713 if (mNumErrorMessages < mMaxErrorMessages) {
714 LOG(warning) << "Mapping error DDL " << feeID << ": " << error.what();
715 mNumErrorMessages++;
716 if (mNumErrorMessages == mMaxErrorMessages) {
717 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
718 }
719 } else {
720 mErrorMessagesSuppressed++;
721 }
722 if (mCreateRawDataErrors) {
724 mOutputDecoderErrors.push_back(mappingError);
725 }
726}
727
728void RawToCellConverterSpec::handleAltroError(const o2::emcal::AltroDecoderError& altroerror, int ddlID)
729{
730 if (mNumErrorMessages < mMaxErrorMessages) {
731 std::string errormessage;
733 switch (altroerror.getErrorType()) {
734 case AltroErrType::RCU_TRAILER_ERROR:
735 errormessage = " RCU Trailer Error ";
736 break;
737 case AltroErrType::RCU_VERSION_ERROR:
738 errormessage = " RCU Version Error ";
739 break;
740 case AltroErrType::RCU_TRAILER_SIZE_ERROR:
741 errormessage = " RCU Trailer Size Error ";
742 break;
743 case AltroErrType::ALTRO_BUNCH_HEADER_ERROR:
744 errormessage = " ALTRO Bunch Header Error ";
745 break;
746 case AltroErrType::ALTRO_BUNCH_LENGTH_ERROR:
747 errormessage = " ALTRO Bunch Length Error ";
748 break;
749 case AltroErrType::ALTRO_PAYLOAD_ERROR:
750 errormessage = " ALTRO Payload Error ";
751 break;
752 case AltroErrType::ALTRO_MAPPING_ERROR:
753 errormessage = " ALTRO Mapping Error ";
754 break;
755 case AltroErrType::CHANNEL_ERROR:
756 errormessage = " Channel Error ";
757 break;
758 default:
759 break;
760 };
761 LOG(warning) << " EMCAL raw task: " << errormessage << " in DDL " << ddlID;
762 mNumErrorMessages++;
763 if (mNumErrorMessages == mMaxErrorMessages) {
764 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
765 }
766 } else {
767 mErrorMessagesSuppressed++;
768 }
769 if (mCreateRawDataErrors) {
770 // fill histograms with error types
772 mOutputDecoderErrors.push_back(errornum);
773 }
774}
775
776void RawToCellConverterSpec::handleMinorAltroError(const o2::emcal::MinorAltroDecodingError& minorerror, int ddlID)
777{
778 if (mNumErrorMessages < mMaxErrorMessages) {
779 LOG(warning) << " EMCAL raw task - Minor error in DDL " << ddlID << ": " << minorerror.what();
780 mNumErrorMessages++;
781 if (mNumErrorMessages == mMaxErrorMessages) {
782 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
783 }
784 } else {
785 mErrorMessagesSuppressed++;
786 }
787 if (mCreateRawDataErrors) {
788 int fecID = -1, hwaddress = -1;
789 try {
791 fecID = mMapper->getFEEForChannelInDDL(ddlID, Channel::getFecIndexFromHwAddress(hwaddress), Channel::getBranchIndexFromHwAddress(hwaddress));
793 // Unfortunately corrupted FEC IDs will not have useful information, so we need to initalize with -1
794 } catch (MappingHandler::DDLInvalid& e) {
795 // Unfortunately corrupted FEC IDs will not have useful information, so we need to initalize with -1
796 }
798 mOutputDecoderErrors.push_back(errornum);
799 }
800}
801
802void RawToCellConverterSpec::handleDDLError(const MappingHandler::DDLInvalid& error, int feeID)
803{
804 if (mNumErrorMessages < mMaxErrorMessages) {
805 LOG(error) << "Failed obtaining mapping for DDL " << error.getDDDL();
806 mNumErrorMessages++;
807 if (mNumErrorMessages == mMaxErrorMessages) {
808 LOG(error) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
809 }
810 }
811 if (mCreateRawDataErrors) {
813 }
814}
815
816void RawToCellConverterSpec::handleFitError(const o2::emcal::CaloRawFitter::RawFitterError_t& fiterror, int ddlID, int cellID, int hwaddress)
817{
819 // Display
820 if (mNumErrorMessages < mMaxErrorMessages) {
821 LOG(warning) << "Failure in raw fitting: " << CaloRawFitter::createErrorMessage(fiterror);
822 mNumErrorMessages++;
823 if (mNumErrorMessages == mMaxErrorMessages) {
824 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
825 }
826 } else {
827 mErrorMessagesSuppressed++;
828 }
829 // Exclude BUNCH_NOT_OK also from raw error objects
830 if (mCreateRawDataErrors) {
831 mOutputDecoderErrors.emplace_back(ddlID, ErrorTypeFEE::ErrorSource_t::FIT_ERROR, CaloRawFitter::getErrorNumber(fiterror), cellID, hwaddress);
832 }
833 } else {
834 LOG(debug2) << "Failure in raw fitting: " << CaloRawFitter::createErrorMessage(fiterror);
835 }
836}
837
838void RawToCellConverterSpec::handleGainError(const o2::emcal::reconstructionerrors::GainError_t& errortype, int ddlID, int hwaddress)
839{
840 int fecID = mMapper->getFEEForChannelInDDL(ddlID, Channel::getFecIndexFromHwAddress(hwaddress), Channel::getBranchIndexFromHwAddress(hwaddress));
841 if (mNumErrorMessages < mMaxErrorMessages) {
842 switch (errortype) {
844 LOG(warning) << "FEC " << fecID << ": 0x" << std::hex << hwaddress << std::dec << " (DDL " << ddlID << ") has only high-gain out-of-range";
845 break;
847 LOG(warning) << "FEC " << fecID << ": 0x" << std::hex << hwaddress << std::dec << " (DDL " << ddlID << ") has low gain but no high-gain";
848 break;
849 default:
850 LOG(warning) << "FEC " << fecID << ": 0x" << std::hex << hwaddress << std::dec << " (DDL " << ddlID << ") falsely flagged as gain error";
851 };
852 if (mNumErrorMessages == mMaxErrorMessages) {
853 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
854 }
855 } else {
856 mErrorMessagesSuppressed++;
857 }
858 if (mCreateRawDataErrors) {
859 mOutputDecoderErrors.emplace_back(ddlID, ErrorTypeFEE::GAIN_ERROR, reconstructionerrors::getErrorCodeFromGainError(errortype), fecID, hwaddress);
860 }
861}
862
863void RawToCellConverterSpec::handleGeometryError(const ModuleIndexException& error, int feeID, int cellID, int hwaddress, ChannelType_t chantype)
864{
865 if (mNumErrorMessages < mMaxErrorMessages) {
866 std::string celltypename;
867 switch (chantype) {
869 celltypename = "high gain";
870 break;
872 celltypename = "low-gain";
873 break;
875 celltypename = "TRU";
876 break;
878 celltypename = "LEDMON";
879 break;
880 };
881 int supermoduleID = feeID / 2;
882 switch (error.getModuleType()) {
883 case ModuleIndexException::ModuleType_t::CELL_MODULE:
884 LOG(warning) << "Sending invalid or negative cell ID " << error.getIndex() << " (SM " << supermoduleID << ", row " << error.getRow() << " - shift " << error.getRowShifted() << ", col " << error.getColumn() << " - shift " << error.getColumnShifted() << ") of type " << celltypename;
885 break;
886 case ModuleIndexException::ModuleType_t::LEDMON_MODULE:
887 LOG(warning) << "Sending invalid or negative LEDMON module ID " << error.getIndex() << "( SM" << supermoduleID << ")";
888 break;
889 };
890 mNumErrorMessages++;
891 if (mNumErrorMessages == mMaxErrorMessages) {
892 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
893 }
894 } else {
895 mErrorMessagesSuppressed++;
896 }
897 if (mCreateRawDataErrors) {
898 mOutputDecoderErrors.emplace_back(feeID, ErrorTypeFEE::ErrorSource_t::GEOMETRY_ERROR, reconstructionerrors::getErrorCodeFromGeometryError(cellID < 0 ? reconstructionerrors::GeometryError_t::CELL_INDEX_NEGATIVE : reconstructionerrors::GeometryError_t::CELL_RANGE_EXCEED), cellID, hwaddress); // 0 -> Cell ID out of range
899 }
900}
901
902void RawToCellConverterSpec::handlePageError(const RawDecodingError& e)
903{
904 if (mCreateRawDataErrors) {
905 mOutputDecoderErrors.emplace_back(e.getFECID(), ErrorTypeFEE::ErrorSource_t::PAGE_ERROR, RawDecodingError::ErrorTypeToInt(e.getErrorType()), -1, -1);
906 }
907 if (mNumErrorMessages < mMaxErrorMessages) {
908 LOG(warning) << " Page decoding: " << e.what() << " in FEE ID " << e.getFECID();
909 mNumErrorMessages++;
910 if (mNumErrorMessages == mMaxErrorMessages) {
911 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
912 }
913 } else {
914 mErrorMessagesSuppressed++;
915 }
916}
917
918void RawToCellConverterSpec::handleMinorPageError(const RawReaderMemory::MinorError& e)
919{
920 if (mCreateRawDataErrors) {
921 mOutputDecoderErrors.emplace_back(e.getFEEID(), ErrorTypeFEE::ErrorSource_t::PAGE_ERROR, RawDecodingError::ErrorTypeToInt(e.getErrorType()), -1, -1);
922 }
923 if (mNumErrorMessages < mMaxErrorMessages) {
924 LOG(warning) << " Page decoding: " << RawDecodingError::getErrorCodeDescription(e.getErrorType()) << " in FEE ID " << e.getFEEID();
925 mNumErrorMessages++;
926 if (mNumErrorMessages == mMaxErrorMessages) {
927 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
928 }
929 } else {
930 mErrorMessagesSuppressed++;
931 }
932}
933
934void RawToCellConverterSpec::handleFastORErrors(const FastORIndexException& e, unsigned int linkID, unsigned int indexTRU)
935{
936 if (mCreateRawDataErrors) {
938 }
939 if (mNumErrorMessages < mMaxErrorMessages) {
940 LOG(warning) << " TRU decoding: " << e.what() << " in FEE ID " << linkID << ", TRU " << indexTRU;
941 mNumErrorMessages++;
942 if (mNumErrorMessages == mMaxErrorMessages) {
943 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
944 }
945 } else {
946 mErrorMessagesSuppressed++;
947 }
948}
949
950void RawToCellConverterSpec::handleFastORStartTimeErrors(const FastOrStartTimeInvalidException& e, unsigned int linkID, unsigned int indexTRU)
951{
952 if (mCreateRawDataErrors) {
954 }
955 if (mNumErrorMessages < mMaxErrorMessages) {
956 LOG(warning) << " TRU decoding: " << e.what() << " in FEE ID " << linkID << ", TRU " << indexTRU;
957 mNumErrorMessages++;
958 if (mNumErrorMessages == mMaxErrorMessages) {
959 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
960 }
961 } else {
962 mErrorMessagesSuppressed++;
963 }
964}
965
966void RawToCellConverterSpec::handlePatchError(const TRUDataHandler::PatchIndexException& e, unsigned int linkID, unsigned int indexTRU)
967{
968 if (mCreateRawDataErrors) {
970 }
971 if (mNumErrorMessages < mMaxErrorMessages) {
972 LOG(warning) << " TRU decoding: " << e.what() << " in FEE ID " << linkID << ", TRU " << indexTRU;
973 mNumErrorMessages++;
974 if (mNumErrorMessages == mMaxErrorMessages) {
975 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
976 }
977 } else {
978 mErrorMessagesSuppressed++;
979 }
980}
981
982void RawToCellConverterSpec::handleTRUIndexError(const TRUIndexException& e, unsigned int linkID, unsigned int hwaddress)
983{
984 if (mCreateRawDataErrors) {
986 }
987 if (mNumErrorMessages < mMaxErrorMessages) {
988 LOG(warning) << " TRU decoding: " << e.what() << " in FEE ID " << linkID << ", TRU " << e.getTRUIndex() << "(hardware address: " << hwaddress << ")";
989 mNumErrorMessages++;
990 if (mNumErrorMessages == mMaxErrorMessages) {
991 LOG(warning) << "Max. amount of error messages (" << mMaxErrorMessages << " reached, further messages will be suppressed";
992 }
993 } else {
994 mErrorMessagesSuppressed++;
995 }
996}
997
998void RawToCellConverterSpec::sendData(framework::ProcessingContext& ctx) const
999{
1000 constexpr auto originEMC = o2::header::gDataOriginEMC;
1001 ctx.outputs().snapshot(framework::Output{originEMC, "CELLS", mSubspecification}, mOutputCells);
1002 ctx.outputs().snapshot(framework::Output{originEMC, "CELLSTRGR", mSubspecification}, mOutputTriggerRecords);
1003 if (mCreateRawDataErrors) {
1004 LOG(debug) << "Sending " << mOutputDecoderErrors.size() << " decoding errors";
1005 ctx.outputs().snapshot(framework::Output{originEMC, "DECODERERR", mSubspecification}, mOutputDecoderErrors);
1006 }
1007 if (mDoTriggerReconstruction) {
1008 ctx.outputs().snapshot(framework::Output{originEMC, "TRUS", mSubspecification}, mOutputTRUs);
1009 ctx.outputs().snapshot(framework::Output{originEMC, "TRUSTRGR", mSubspecification}, mOutputTRUTriggerRecords);
1010 ctx.outputs().snapshot(framework::Output{originEMC, "PATCHES", mSubspecification}, mOutputPatches);
1011 ctx.outputs().snapshot(framework::Output{originEMC, "PATCHESTRGR", mSubspecification}, mOutputPatchTriggerRecords);
1012 ctx.outputs().snapshot(framework::Output{originEMC, "FASTORS", mSubspecification}, mOutputTimesums);
1013 ctx.outputs().snapshot(framework::Output{originEMC, "FASTORSTRGR", mSubspecification}, mOutputTimesumTriggerRecords);
1014 }
1015}
1016
1017RawToCellConverterSpec::ModuleIndexException::ModuleIndexException(int moduleIndex, int column, int row, int shiftedColumn, int shiftedRow) : mModuleType(ModuleType_t::CELL_MODULE),
1018 mIndex(moduleIndex),
1019 mColumn(column),
1020 mRow(row),
1021 mColumnShifted(shiftedColumn),
1022 mRowShifted(shiftedRow) {}
1023
1024RawToCellConverterSpec::ModuleIndexException::ModuleIndexException(int moduleIndex) : mModuleType(ModuleType_t::LEDMON_MODULE), mIndex(moduleIndex) {}
1025
1026o2::framework::DataProcessorSpec o2::emcal::reco_workflow::getRawToCellConverterSpec(bool askDISTSTF, bool disableDecodingErrors, bool disableTriggerReconstruction, int subspecification)
1027{
1028 constexpr auto originEMC = o2::header::gDataOriginEMC;
1029 std::vector<o2::framework::OutputSpec> outputs;
1030
1031 outputs.emplace_back(originEMC, "CELLS", subspecification, o2::framework::Lifetime::Timeframe);
1032 outputs.emplace_back(originEMC, "CELLSTRGR", subspecification, o2::framework::Lifetime::Timeframe);
1033 if (!disableDecodingErrors) {
1034 outputs.emplace_back(originEMC, "DECODERERR", subspecification, o2::framework::Lifetime::Timeframe);
1035 }
1036 if (!disableTriggerReconstruction) {
1037 outputs.emplace_back(originEMC, "TRUS", subspecification, o2::framework::Lifetime::Timeframe);
1038 outputs.emplace_back(originEMC, "TRUSTRGR", subspecification, o2::framework::Lifetime::Timeframe);
1039 outputs.emplace_back(originEMC, "PATCHES", subspecification, o2::framework::Lifetime::Timeframe);
1040 outputs.emplace_back(originEMC, "PATCHESTRGR", subspecification, o2::framework::Lifetime::Timeframe);
1041 outputs.emplace_back(originEMC, "FASTORS", subspecification, o2::framework::Lifetime::Timeframe);
1042 outputs.emplace_back(originEMC, "FASTORSTRGR", subspecification, o2::framework::Lifetime::Timeframe);
1043 }
1044
1045 std::vector<o2::framework::InputSpec> inputs{{"stf", o2::framework::ConcreteDataTypeMatcher{originEMC, o2::header::gDataDescriptionRawData}, o2::framework::Lifetime::Timeframe}};
1046 if (askDISTSTF) {
1047 inputs.emplace_back("stdDist", "FLP", "DISTSUBTIMEFRAME", 0, o2::framework::Lifetime::Timeframe);
1048 }
1049 // CCDB objects
1050 auto calibhandler = std::make_shared<o2::emcal::CalibLoader>();
1051 calibhandler->enableRecoParams(true);
1052 calibhandler->enableFEEDCS(true);
1053 calibhandler->defineInputSpecs(inputs);
1054
1056 "EMCALRawToCellConverterSpec",
1057 inputs,
1058 outputs,
1059 o2::framework::adaptFromTask<o2::emcal::reco_workflow::RawToCellConverterSpec>(subspecification, !disableDecodingErrors, !disableTriggerReconstruction, calibhandler),
1061 {"fitmethod", o2::framework::VariantType::String, "gamma2", {"Fit method (standard or gamma2)"}},
1062 {"maxmessage", o2::framework::VariantType::Int, 100, {"Max. amout of error messages to be displayed"}},
1063 {"printtrailer", o2::framework::VariantType::Bool, false, {"Print RCU trailer (for debugging)"}},
1064 {"no-mergeHGLG", o2::framework::VariantType::Bool, false, {"Do not merge HG and LG channels for same tower"}},
1065 {"no-checkactivelinks", o2::framework::VariantType::Bool, false, {"Do not check for active links per BC"}},
1066 {"no-evalpedestal", o2::framework::VariantType::Bool, false, {"Disable pedestal evaluation"}}}};
1067}
Definition of the 32 Central Trigger System (CTS) Trigger Types defined in https://twiki....
A helper class to iteratate over all parts of all input routes.
Definition of a container to keep Monte Carlo truth external to simulation objects.
std::ostringstream debug
Errors appearing in geometry access obtaining the tower ID.
Error handling of the ALTRO Decoder.
ErrorType_t
Error codes connected with the ALTRO decoding.
@ ALTRO_MAPPING_ERROR
Incorrect ALTRO channel mapping.
static int errorTypeToInt(ErrorType_t errortype)
convert the error type from symoblic constant into int
const ErrorType_t getErrorType() const noexcept
Access to the error type connected to the erro.
Decoder of the ALTRO data in the raw page.
void decode()
Decode the ALTRO stream.
const RCUTrailer & getRCUTrailer() const
Get reference to the RCU trailer object.
const std::vector< Channel > & getChannels() const
Get the reference to the channel container.
const std::vector< MinorAltroDecodingError > & getMinorDecodingErrors() const
Get list of minor decoding errors.
void setMaxBunchLength(int maxBunchLength)
Set the max. allowed bunch length.
Container class to hold results from fitting.
Raw data fitting: Gamma-2 function.
Raw data fitting: standard TMinuit fit.
static int getErrorNumber(RawFitterError_t fiterror)
Convert error type to numeric representation.
RawFitterError_t
Error codes for failures in raw fitter procedure.
static std::string createErrorMessage(RawFitterError_t fiterror)
Create error message for a given error type.
ALTRO channel representation.
Definition Channel.h:46
static int getFecIndexFromHwAddress(int hwaddress)
Extracting FEC index in branch from the hardware address.
Definition Channel.h:163
static int getHardwareAddressFromChannelHeader(int channelheader)
Extrcting hardware address from the channel header word.
Definition Channel.h:151
uint16_t getHardwareAddress() const
Get the full hardware address.
Definition Channel.h:94
static int getBranchIndexFromHwAddress(int hwaddress)
Extracting branch index from the hardware address.
Definition Channel.h:159
const std::vector< Bunch > & getBunches() const
Get list of bunches in the channel.
Definition Channel.h:102
Errors per FEE information.
@ FIT_ERROR
Raw fit failed.
@ LINK_ERROR
Error due to missing DDL links.
@ ALTRO_ERROR
Decoding of the ALTRO payload failed.
@ GEOMETRY_ERROR
Decoded position outside EMCAL.
@ TRU_ERROR
Errors from TRU data.
@ MINOR_ALTRO_ERROR
Non-fatal error in decoding of the ALTRO payload.
@ GAIN_ERROR
Error due to gain type.
@ PAGE_ERROR
Raw page decoding failed.
Containter of cells for a given event.
const std::unordered_map< uint16_t, FastORTimeSeries > & getTimeSeriesContainer() const
Access to container with FastOR time series.
void setTriggerBits(uint64_t triggerbits)
Set trigger bits of the interaction.
void setCell(int tower, double energy, double time, ChannelType_t celltype, int hwaddress, int ddlID, bool doMergeHGLG)
Add cell information to the event container.
void setFastOR(uint16_t fastORAbsID, uint8_t starttime, const gsl::span< const uint16_t > timesamples)
Add bunch of time series to the container.
const TRUDataHandler & readTRUData(std::size_t truIndex) const
Read-only access TRU data of a given TRU.
const o2::InteractionRecord & getInteractionRecord() const
Get interaction record of the event.
TRUDataHandler & getTRUData(std::size_t truIndex)
Read and write access TRU data of a given TRU.
void setLEDMONCell(int tower, double energy, double time, ChannelType_t celltype, int hwaddress, int ddlID, bool doMergeHGLG)
Add LEDMON information to the event container.
Error handling of faulty FastOR indices.
const char * what() const noexcept final
Get error message.
Handling of error if starttime is to large (>=14). This is most likely caused by a corrupted channel ...
const char * what() const noexcept final
Access to error message.
std::bitset< 14 > getDDLlist1() const
Definition FeeDCS.h:54
std::bitset< 32 > getDDLlist0() const
Definition FeeDCS.h:53
static Geometry * GetInstanceFromRunNumber(Int_t runNumber, const std::string_view="", const std::string_view mcname="TGeant3", const std::string_view mctitle="")
Instanciate geometry depending on the run number. Mostly used in analysis and MC anchors.
Definition Geometry.cxx:219
std::tuple< int, int > ShiftOnlineToOfflineCellIndexes(Int_t supermoduleID, Int_t iphi, Int_t ieta) const
Adapt cell indices in supermodule to online indexing.
Int_t GetAbsCellIdFromCellIndexes(Int_t nSupMod, Int_t iphi, Int_t ieta) const
Transition from super module number (nSupMod) and cell indexes (ieta,iphi) to cell absolute ID number...
Definition Geometry.cxx:791
Error handling requests for unknown hardware addresses.
Definition Mapper.h:85
const char * what() const noexcept override
Access to error message.
Definition Mapper.h:103
Error handling for invalid DDL IDs (not in range for EMCAL)
Definition Mapper.h:321
int getDDDL() const
Access to DDL ID responsible for the exception.
Definition Mapper.h:334
Handler providing the correct mapping for the given DDL.
Definition Mapper.h:313
Error handling for the ALTRO decoder for non-crashing errors.
static int errorTypeToInt(ErrorType_t errortype)
convert the error type from symoblic constant into int
ErrorType_t getErrorType() const noexcept
Get the type of the error.
std::string what() const noexcept
Create and return error message for different error types.
uint32_t getChannelHeader() const noexcept
Get the header of the channel raising the error.
bool hasZeroSuppression() const
Check whether zero suppression has been applied.
Definition RCUTrailer.h:260
Error handling of the raw reader.
static const char * getErrorCodeDescription(ErrorType_t errortype)
Get description of error type.
const char * what() const noexcept override
Providing error message of the exception.
@ HEADER_DECODING
Header cannot be decoded (format incorrect)
@ HEADER_INVALID
Header in memory not belonging to requested superpage.
static int ErrorTypeToInt(RawDecodingError::ErrorType_t errortype)
Convert error type to error code number.
int getFECID() const
Get the ID of the frontend electronics responsible for the error.
ErrorType_t getErrorType() const
Get the type identifier of the error handled with this exception.
Minor (non-crashing) raw decoding errors.
RawDecodingError::ErrorType_t getErrorType() const
Get type of the error.
int getFEEID() const
Get ID of the FEE.
Reader for raw data produced by the Readout application in in-memory format.
bool hasNext() const
check if more pages are available in the raw file
const o2::header::RDHAny & getRawHeader() const
access to the raw header of the current page
void setRangeSRUDDLs(uint16_t minDDL, uint16_t maxDDL)
Set range for DDLs from SRU (for RCU trailer merging)
void next()
Read next payload from the stream.
gsl::span< const MinorError > getMinorErrors() const
Get minor (non-crashing) raw decoding errors.
Iterator over reco containers.
EventContainer & nextEvent()
Get the next event in container.
bool hasNext() const
Check whehter there are more events in the container.
EventContainer & getEventContainer(const o2::InteractionRecord &currentIR)
Get container for trigger.
void reset()
Clear container.
Handler of errors related to invalid trigger patch IDs.
const char * what() const noexcept final
Access Error message.
Error handling of faulty TRU indices.
const char * what() const noexcept final
Get error message.
unsigned int getTRUIndex() const noexcept
Get the index of the TRU raising the exception.
static constexpr unsigned int ALLTRUS
Total number of TRUs in EMCAL.
static constexpr unsigned int PATCHESINTRU
void finaliseCCDB(framework::ConcreteDataMatcher &matcher, void *obj) final
Handle objects obtained from the CCDB.
void run(framework::ProcessingContext &ctx) final
Run conversion of raw data to cells.
void init(framework::InitContext &ctx) final
Initializing the RawToCellConverterSpec.
void snapshot(const Output &spec, T const &object)
ServiceRegistryRef services()
Definition InitContext.h:34
ConfigParamRegistry const & options()
Definition InitContext.h:33
A helper class to iteratate over all parts of all input routes.
DataAllocator & outputs()
The data allocator is used to allocate memory for the output data.
InputRecord & inputs()
The inputs associated with this processing context.
ServiceRegistryRef services()
The services registry associated with this processing context.
bool active() const
Check if service of type T is currently active.
GLint GLint GLint GLint GLint GLint GLint GLbitfield GLenum filter
Definition glcorearb.h:1308
constexpr o2::header::DataDescription gDataDescriptionRawData
Definition DataHeader.h:597
constexpr o2::header::DataOrigin gDataOriginEMC
Definition DataHeader.h:565
uint8_t itsSharedClusterMap uint8_t
framework::DataProcessorSpec getRawToCellConverterSpec(bool askDISTSTF, bool disableDecodingError, bool disableTriggerReconstruction, int subspecification)
Creating DataProcessorSpec for the EMCAL Cell Converter Spec.
constexpr int getErrorCodeFromGainError(GainError_t errortype)
Convert gain error type into numberic representation.
@ LGNOHG
LG found below HG/LG transition, HG missing.
@ CELL_RANGE_EXCEED
Requested cell value exceeding allowed cell range.
@ CELL_INDEX_NEGATIVE
Requested cell value negative.
constexpr int getErrorCodeFromGeometryError(GeometryError_t errortype)
Convert geometry error type into numberic representation.
constexpr int getErrorCodeFromTRUDecodingError(TRUDecodingError_t error)
Convert TRU decoding error type into numberic representation.
@ FASTOR_STARTTIME_INVALID
FastOr stattime is larger than 14.
@ EMCAL_LEDREFS
Number of LEDs (reference/monitors) per module for EMCAL; one per StripModule.
Definition Constants.h:27
ChannelType_t
Type of a raw data channel.
Definition Constants.h:33
@ TRU
TRU channel.
Definition Constants.h:36
@ HIGH_GAIN
High gain channel.
Definition Constants.h:35
@ LOW_GAIN
Low gain channel.
Definition Constants.h:34
@ LEDMON
LED monitor channel.
Definition Constants.h:37
std::vector< ConfigParamSpec > Options
constexpr uint32_t PhT
Definition Triggers.h:30
constexpr uint32_t Cal
Definition Triggers.h:32
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
uint32_t orbit
LHC orbit.
uint16_t bc
bunch crossing ID of interaction
int64_t differenceInBC(const InteractionRecord &other) const
static o2::header::DataHeader::PayloadSizeType getPayloadSize(const DataRef &ref)
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
std::vector< Cell > cells
std::vector< int > row
ArrayADC adc