Project
Loading...
Searching...
No Matches
ITSMFTDigitizerSpec.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
12#include "ITSMFTDigitizerSpec.h"
17#include "Framework/Lifetime.h"
18#include "Framework/Task.h"
36#include "MFTBase/GeometryTGeo.h"
37#include <TChain.h>
38#include <TStopwatch.h>
39#include <string>
40#include <format>
41
42using namespace o2::framework;
44
45namespace o2::itsmft
46{
47
48using namespace o2::base;
49template <int N>
51{
52 public:
56
58
60 {
61 mDisableQED = ic.options().get<bool>("disable-qed");
62 }
63
65 {
66 if (mFinished) {
67 return;
68 }
69
71 const o2::InteractionRecord firstIR(0, mFirstOrbitTF);
73
74 TStopwatch timer;
75 timer.Start();
76 LOG(info) << " CALLING ITS DIGITIZATION ";
77
78 // read collision context from input
79 auto context = pc.inputs().get<o2::steer::DigitizationContext*>("collisioncontext");
80 context->initSimChains(ID, mSimChains);
81 const bool withQED = context->isQEDProvided() && !mDisableQED;
82 auto& timesview = context->getEventRecords(withQED);
83 LOG(info) << "GOT " << timesview.size() << " COLLISSION TIMES";
84 LOG(info) << "SIMCHAINS: " << mSimChains.size();
85
86 // if there is nothing to do ... return
87 if (timesview.size() == 0) {
88 return;
89 }
90
91 uint64_t nDigits{0};
92 constexpr uint32_t nLayers = (DPLAlpideParam<N>::supportsStaggering()) ? NLayers : 1;
93 for (uint32_t iLayer = 0; iLayer < nLayers; ++iLayer) {
94 const int layer = (DPLAlpideParam<N>::supportsStaggering()) ? iLayer : -1;
99
100 // digits are directly put into DPL owned resource
101 auto& digitsAccum = pc.outputs().make<std::vector<itsmft::Digit>>(Output{Origin, "DIGITS", iLayer});
102
103 // rofs are accumulated first and the copied
104 const int nROFsPerOrbit = o2::constants::lhc::LHCMaxBunches / DPLAlpideParam<N>::Instance().getROFLengthInBC(iLayer);
105 const int nROFsTF = nROFsPerOrbit * raw::HBFUtils::Instance().getNOrbitsPerTF();
106 mROFRecordsAccum[iLayer].reserve(nROFsTF);
107
108 auto accumulate = [this, &digitsAccum, &iLayer]() {
109 // accumulate result of single event processing on a specific layer, called after processing every event supplied
110 // AND after the final flushing via digitizer::fillOutputContainer
111 if (!mDigits[iLayer].size()) {
112 return; // no digits were flushed, nothing to accumulate
113 }
114 auto ndigAcc = digitsAccum.size();
115 std::copy(mDigits[iLayer].begin(), mDigits[iLayer].end(), std::back_inserter(digitsAccum));
116
117 // fix ROFrecords references on ROF entries
118 auto nROFRecsOld = mROFRecordsAccum[iLayer].size();
119
120 for (int i = 0; i < mROFRecords[iLayer].size(); i++) {
121 auto& rof = mROFRecords[iLayer][i];
122 rof.setFirstEntry(ndigAcc + rof.getFirstEntry());
123 rof.print();
124
125 if (mFixMC2ROF[iLayer] < mMC2ROFRecordsAccum[iLayer].size()) { // fix ROFRecord entry in MC2ROF records
126 for (int m2rid = mFixMC2ROF[iLayer]; m2rid < mMC2ROFRecordsAccum[iLayer].size(); m2rid++) {
127 // need to register the ROFRecors entry for MC event starting from this entry
128 auto& mc2rof = mMC2ROFRecordsAccum[iLayer][m2rid];
129 if (rof.getROFrame() == mc2rof.minROF) {
130 mFixMC2ROF[iLayer]++;
131 mc2rof.rofRecordID = nROFRecsOld + i;
132 mc2rof.print();
133 }
134 }
135 }
136 }
137
138 std::copy(mROFRecords[iLayer].begin(), mROFRecords[iLayer].end(), std::back_inserter(mROFRecordsAccum[iLayer]));
139 if (mWithMCTruth) {
140 mLabelsAccum[iLayer].mergeAtBack(mLabels[iLayer]);
141 }
142 LOG(info) << "Added " << mDigits[iLayer].size() << " digits:" << iLayer;
143 // clean containers from already accumulated stuff
144 mLabels[iLayer].clear();
145 mDigits[iLayer].clear();
146 mROFRecords[iLayer].clear();
147 }; // and accumulate lambda
148
149 const auto& eventParts = context->getEventParts(withQED);
150 const int64_t bcShift = mDigitizer.getParams().getROFrameBiasInBC(layer); // this accounts the misalignment and the opt. imposed rof delay
151 // loop over all composite collisions given from context (aka loop over all the interaction records)
152 for (int collID = 0; collID < timesview.size(); ++collID) {
153 auto irt = timesview[collID];
154 if (irt.toLong() < bcShift) { // due to the ROF misalignment (+opt. delay) the collision would go to negative ROF ID, discard
155 continue;
156 }
157 irt -= bcShift; // account for the ROF start shift
158
160 mDigitizer.resetEventROFrames(); // to estimate min/max ROF for this collID
161 // for each collision, loop over the constituents event and source IDs
162 // (background signal merging is basically taking place here)
163 for (const auto& part : eventParts[collID]) {
164
165 // get the hits for this event and this source
166 mHits.clear();
167 context->retrieveHits(mSimChains, o2::detectors::SimTraits::DETECTORBRANCHNAMES[ID][0].c_str(), part.sourceID, part.entryID, &mHits);
168
169 if (mHits.size() > 0) {
170 LOG(debug) << "For collision " << collID << " eventID " << part.entryID << " found " << mHits.size() << " hits ";
171 mDigitizer.process(&mHits, part.entryID, part.sourceID, layer); // call actual digitization procedure
172 }
173 }
175 accumulate();
176 }
178 accumulate();
179 nDigits += digitsAccum.size();
180
181 // here we have all digits and labels and we can send them to consumer (aka snapshot it onto output)
182 // ensure that the rof output is continuous
183 if (nROFsTF != mROFRecordsAccum[iLayer].size()) {
184 // it can happen that in the digitization rofs without contributing hits are skipped
185 // however downstream consumers of the clusters cannot know apriori the time structure
186 // the cluster rofs do not account for the bias so it will start always at BC=0
187 // also have to account for spillage into next TF
188 const size_t nROFsLayer = std::max((size_t)nROFsTF, mROFRecordsAccum[iLayer].size());
189 std::vector<o2::itsmft::ROFRecord> expDigitRofVec(nROFsLayer);
190 for (int iROF{0}; iROF < nROFsLayer; ++iROF) {
191 auto& rof = expDigitRofVec[iROF];
192 int orb = iROF * DPLAlpideParam<N>::Instance().getROFLengthInBC(iLayer) / o2::constants::lhc::LHCMaxBunches + mFirstOrbitTF;
193 int bc = iROF * DPLAlpideParam<N>::Instance().getROFLengthInBC(iLayer) % o2::constants::lhc::LHCMaxBunches;
195 rof.setBCData(ir);
196 rof.setROFrame(iROF);
197 rof.setNEntries(0);
198 rof.setFirstEntry(-1);
199 }
200 uint32_t prevEntry{0};
201 for (const auto& rof : mROFRecordsAccum[iLayer]) {
202 const auto& ir = rof.getBCData();
203 const auto irToFirst = ir - firstIR;
204 const int irROF = irToFirst.toLong() / DPLAlpideParam<N>::Instance().getROFLengthInBC(iLayer);
205 auto& expROF = expDigitRofVec[irROF];
206 expROF.setFirstEntry(rof.getFirstEntry());
207 expROF.setNEntries(rof.getNEntries());
208 if (expROF.getBCData() != rof.getBCData()) {
209 LOGP(fatal, "detected mismatch between expected {} and received {}", expROF.asString(), rof.asString());
210 }
211 }
212 int prevFirst{0};
213 for (auto& rof : expDigitRofVec) {
214 if (rof.getFirstEntry() < 0) {
215 rof.setFirstEntry(prevFirst);
216 }
217 prevFirst = rof.getFirstEntry();
218 }
219 // if more rofs where accumulated than ROFs possible in the TF, cut them away
220 // by construction expDigitRofVec is at least nROFsTF long
221 expDigitRofVec.resize(nROFsTF);
222 pc.outputs().snapshot(Output{Origin, "DIGITSROF", iLayer}, expDigitRofVec);
223 } else {
224 pc.outputs().snapshot(Output{Origin, "DIGITSROF", iLayer}, mROFRecordsAccum[iLayer]);
225 }
226 if (mWithMCTruth) {
227 pc.outputs().snapshot(Output{Origin, "DIGITSMC2ROF", iLayer}, mMC2ROFRecordsAccum[iLayer]);
228 auto& sharedlabels = pc.outputs().make<o2::dataformats::ConstMCTruthContainer<o2::MCCompLabel>>(Output{Origin, "DIGITSMCTR", iLayer});
229 mLabelsAccum[iLayer].flatten_to(sharedlabels);
230 // free space of existing label containers
231 mLabels[iLayer].clear_andfreememory();
232 mLabelsAccum[iLayer].clear_andfreememory();
233 }
234 }
235
236 LOG(info) << ID.getName() << ": Sending ROMode= " << mROMode << " to GRPUpdater";
237 pc.outputs().snapshot(Output{Origin, "ROMode", 0}, mROMode);
238
239 timer.Stop();
240 LOG(info) << "Digitization took " << timer.CpuTime() << "s";
241 LOG(info) << "Produced " << nDigits << " digits";
242
243 // we should be only called once; tell DPL that this process is ready to exit
244 pc.services().get<ControlService>().readyToQuit(QuitRequest::Me);
245
246 mFinished = true;
247 }
248
249 void finaliseCCDB(ConcreteDataMatcher& matcher, void* obj)
250 {
251 if (matcher == ConcreteDataMatcher(Origin, "NOISEMAP", 0)) {
252 LOG(info) << ID.getName() << " noise map updated";
254 return;
255 }
256 if (matcher == ConcreteDataMatcher(Origin, "DEADMAP", 0)) {
257 LOG(info) << ID.getName() << " static dead map updated";
260 return;
261 }
262 if (matcher == ConcreteDataMatcher(Origin, "TimeDeadMap", 0)) {
264 if (!timedeadmap->isDefault()) {
265 timedeadmap->decodeMap(mFirstOrbitTF, *mDeadMap, true);
267 LOGP(fatal, "Attempt to add time-dependent map to already modified static map");
268 }
269 mTimeDeadMapUpdated = true;
271 LOG(info) << ID.getName() << " time-dependent dead map updated";
272 } else {
273 LOG(info) << ID.getName() << " time-dependent dead map is default/empty";
274 }
275
276 return;
277 }
278 if (matcher == ConcreteDataMatcher(Origin, "ALPIDEPARAM", 0)) {
279 LOG(info) << ID.getName() << " Alpide param updated";
280 const auto& par = o2::itsmft::DPLAlpideParam<N>::Instance();
281 par.printKeyValues();
282 return;
283 }
284 if (matcher == ConcreteDataMatcher(Origin, "ALPIDERESPVbb0", 0)) {
285 LOG(info) << ID.getName() << " loaded AlpideResponseData for Vbb=0V";
287 }
288 if (matcher == ConcreteDataMatcher(Origin, "ALPIDERESPVbbM3", 0)) {
289 LOG(info) << ID.getName() << " loaded AlpideResponseData for Vbb=-3V";
291 }
292 }
293
294 protected:
295 ITSMFTDPLDigitizerTask(bool mctruth = true) : BaseDPLDigitizer(InitServices::FIELD | InitServices::GEOM), mWithMCTruth(mctruth) {}
296
298 {
299 std::string detstr(o2::detectors::DetID::getName(ID));
300 pc.inputs().get<o2::itsmft::NoiseMap*>(detstr + "_noise");
301 pc.inputs().get<o2::itsmft::NoiseMap*>(detstr + "_dead");
302 // TODO: the code should run even if this object does not exist. Or: create default object
303 pc.inputs().get<o2::itsmft::TimeDeadMap*>(detstr + "_time_dead");
304 pc.inputs().get<o2::itsmft::DPLAlpideParam<N>*>(detstr + "_alppar");
305 pc.inputs().get<o2::itsmft::AlpideSimResponse*>(detstr + "_alpiderespvbb0");
306 pc.inputs().get<o2::itsmft::AlpideSimResponse*>(detstr + "_alpiderespvbbm3");
307
310 auto& digipar = mDigitizer.getParams();
311 digipar.setContinuous(dopt.continuous);
312 digipar.setROFrameBiasInBC(aopt.roFrameBiasInBC);
313 if (dopt.continuous) {
314 auto frameNS = aopt.roFrameLengthInBC * o2::constants::lhc::LHCBunchSpacingNS;
315 digipar.setROFrameLengthInBC(aopt.roFrameLengthInBC);
316 digipar.setROFrameLength(frameNS); // RO frame in ns
317 digipar.setStrobeDelay(aopt.strobeDelay); // Strobe delay wrt beginning of the RO frame, in ns
318 digipar.setStrobeLength(aopt.strobeLengthCont > 0 ? aopt.strobeLengthCont : frameNS - aopt.strobeDelay); // Strobe length in ns
319 } else {
320 digipar.setROFrameLength(aopt.roFrameLengthTrig); // RO frame in ns
321 digipar.setStrobeDelay(aopt.strobeDelay); // Strobe delay wrt beginning of the RO frame, in ns
322 digipar.setStrobeLength(aopt.strobeLengthTrig); // Strobe length in ns
323 }
324 // parameters of signal time response: flat-top duration, max rise time and q @ which rise time is 0
325 digipar.getSignalShape().setParameters(dopt.strobeFlatTop, dopt.strobeMaxRiseTime, dopt.strobeQRiseTime0);
326 digipar.setChargeThreshold(dopt.chargeThreshold); // charge threshold in electrons
327 digipar.setNoisePerPixel(dopt.noisePerPixel); // noise level
328 digipar.setTimeOffset(dopt.timeOffset);
329 digipar.setNSimSteps(dopt.nSimSteps);
330 digipar.setIBVbb(dopt.IBVbb);
331 digipar.setOBVbb(dopt.OBVbb);
332 digipar.setVbb(dopt.Vbb);
333 // staggering parameters
335 const bool withStag = aopt.withStaggering();
336 for (int iLayer{0}; iLayer < o2::itsmft::DPLAlpideParam<N>::getNLayers(); ++iLayer) {
337 const int nLayer = (withStag) ? iLayer : -1;
338 auto frameNS = aopt.getROFLengthInBC(nLayer) * o2::constants::lhc::LHCBunchSpacingNS;
339 digipar.addROFrameLayerLengthInBC(aopt.getROFLengthInBC(nLayer));
340 // NOTE: the rof delay looks from the digitizer like an additional bias
341 digipar.addROFrameLayerBiasInBC(aopt.getROFBiasInBC(nLayer) + aopt.getROFDelayInBC(nLayer));
342 digipar.addStrobeDelay(aopt.strobeDelay);
343 digipar.addStrobeLength(aopt.strobeLengthCont > 0 ? aopt.strobeLengthCont : frameNS - aopt.strobeDelay);
344 digipar.setROFrameLength(aopt.getROFLengthInBC(nLayer) * o2::constants::lhc::LHCBunchSpacingNS, iLayer);
345 }
346 }
347
349 LOG(info) << detstr << " simulated in "
350 << ((mROMode == o2::parameters::GRPObject::CONTINUOUS) ? "CONTINUOUS" : "TRIGGERED")
351 << " RO mode";
352
353 // configure digitizer
354 o2::itsmft::GeometryTGeo* geom = nullptr;
355 if constexpr (N == o2::detectors::DetID::ITS) {
357 } else {
359 }
360 geom->fillMatrixCache(o2::math_utils::bit2Mask(o2::math_utils::TransformType::L2G)); // make sure L2G matrices are loaded
363 }
364
365 bool mWithMCTruth = true;
366 bool mFinished = false;
367 bool mDisableQED = false;
368 unsigned long mFirstOrbitTF = 0x0;
370 std::array<std::vector<o2::itsmft::Digit>, NLayers> mDigits;
371 std::array<std::vector<o2::itsmft::ROFRecord>, NLayers> mROFRecords;
372 std::array<std::vector<o2::itsmft::ROFRecord>, NLayers> mROFRecordsAccum;
373 std::vector<o2::itsmft::Hit> mHits;
374 std::vector<o2::itsmft::Hit>* mHitsP = &mHits;
375 std::array<o2::dataformats::MCTruthContainer<o2::MCCompLabel>, NLayers> mLabels;
376 std::array<o2::dataformats::MCTruthContainer<o2::MCCompLabel>, NLayers> mLabelsAccum;
377 std::array<std::vector<o2::itsmft::MC2ROFRecord>, NLayers> mMC2ROFRecordsAccum;
378 std::vector<TChain*> mSimChains;
380
381 std::array<int, NLayers> mFixMC2ROF{}; // 1st entry in mc2rofRecordsAccum to be fixed for ROFRecordID
384};
385
386//_______________________________________________
387class ITSDPLDigitizerTask : public ITSMFTDPLDigitizerTask<o2::detectors::DetID::ITS>
388{
389 public:
390 ITSDPLDigitizerTask(bool mctruth = true) : ITSMFTDPLDigitizerTask<o2::detectors::DetID::ITS>(mctruth) {}
391};
392
393//_______________________________________________
394class MFTDPLDigitizerTask : public ITSMFTDPLDigitizerTask<o2::detectors::DetID::MFT>
395{
396 public:
397 MFTDPLDigitizerTask(bool mctruth = true) : ITSMFTDPLDigitizerTask<o2::detectors::DetID::MFT>(mctruth) {}
398};
399
400namespace
401{
402template <int N>
403std::vector<OutputSpec> makeOutChannels(o2::header::DataOrigin detOrig, bool mctruth)
404{
405 std::vector<OutputSpec> outputs;
406 constexpr uint32_t nLayers = (DPLAlpideParam<N>::supportsStaggering()) ? DPLAlpideParam<N>::getNLayers() : 1;
407 for (uint32_t iLayer = 0; iLayer < nLayers; ++iLayer) {
408 outputs.emplace_back(detOrig, "DIGITS", iLayer, Lifetime::Timeframe);
409 outputs.emplace_back(detOrig, "DIGITSROF", iLayer, Lifetime::Timeframe);
410 if (mctruth) {
411 outputs.emplace_back(detOrig, "DIGITSMC2ROF", iLayer, Lifetime::Timeframe);
412 outputs.emplace_back(detOrig, "DIGITSMCTR", iLayer, Lifetime::Timeframe);
413 }
414 }
415 outputs.emplace_back(detOrig, "ROMode", 0, Lifetime::Timeframe);
416 return outputs;
417}
418} // namespace
419
420DataProcessorSpec getITSDigitizerSpec(int channel, bool mctruth)
421{
423 auto detOrig = ITSDPLDigitizerTask::Origin;
424 std::vector<InputSpec> inputs;
425 inputs.emplace_back("collisioncontext", "SIM", "COLLISIONCONTEXT", static_cast<SubSpecificationType>(channel), Lifetime::Timeframe);
426 inputs.emplace_back("ITS_noise", "ITS", "NOISEMAP", 0, Lifetime::Condition, ccdbParamSpec("ITS/Calib/NoiseMap"));
427 inputs.emplace_back("ITS_dead", "ITS", "DEADMAP", 0, Lifetime::Condition, ccdbParamSpec("ITS/Calib/DeadMap"));
428 inputs.emplace_back("ITS_time_dead", "ITS", "TimeDeadMap", 0, Lifetime::Condition, ccdbParamSpec("ITS/Calib/TimeDeadMap"));
429 inputs.emplace_back("ITS_alppar", "ITS", "ALPIDEPARAM", 0, Lifetime::Condition, ccdbParamSpec("ITS/Config/AlpideParam"));
430 inputs.emplace_back("ITS_alpiderespvbb0", "ITS", "ALPIDERESPVbb0", 0, Lifetime::Condition, ccdbParamSpec("ITSMFT/Calib/ALPIDEResponseVbb0"));
431 inputs.emplace_back("ITS_alpiderespvbbm3", "ITS", "ALPIDERESPVbbM3", 0, Lifetime::Condition, ccdbParamSpec("ITSMFT/Calib/ALPIDEResponseVbbM3"));
432 return DataProcessorSpec{.name = detStr + "Digitizer",
433 .inputs = inputs,
434 .outputs = makeOutChannels<o2::detectors::DetID::ITS>(detOrig, mctruth),
435 .algorithm = AlgorithmSpec{adaptFromTask<ITSDPLDigitizerTask>(mctruth)},
436 .options = Options{
437 {"disable-qed", o2::framework::VariantType::Bool, false, {"disable QED handling"}}}};
438}
439
440DataProcessorSpec getMFTDigitizerSpec(int channel, bool mctruth)
441{
443 auto detOrig = MFTDPLDigitizerTask::Origin;
444 std::vector<InputSpec> inputs;
445 inputs.emplace_back("collisioncontext", "SIM", "COLLISIONCONTEXT", static_cast<SubSpecificationType>(channel), Lifetime::Timeframe);
446 inputs.emplace_back("MFT_noise", "MFT", "NOISEMAP", 0, Lifetime::Condition, ccdbParamSpec("MFT/Calib/NoiseMap"));
447 inputs.emplace_back("MFT_dead", "MFT", "DEADMAP", 0, Lifetime::Condition, ccdbParamSpec("MFT/Calib/DeadMap"));
448 inputs.emplace_back("MFT_time_dead", "MFT", "TimeDeadMap", 0, Lifetime::Condition, ccdbParamSpec("MFT/Calib/TimeDeadMap"));
449 inputs.emplace_back("MFT_alppar", "MFT", "ALPIDEPARAM", 0, Lifetime::Condition, ccdbParamSpec("MFT/Config/AlpideParam"));
450 inputs.emplace_back("MFT_alpiderespvbb0", "MFT", "ALPIDERESPVbb0", 0, Lifetime::Condition, ccdbParamSpec("ITSMFT/Calib/ALPIDEResponseVbb0"));
451 inputs.emplace_back("MFT_alpiderespvbbm3", "MFT", "ALPIDERESPVbbM3", 0, Lifetime::Condition, ccdbParamSpec("ITSMFT/Calib/ALPIDEResponseVbbM3"));
452 return DataProcessorSpec{.name = detStr + "Digitizer",
453 .inputs = inputs,
454 .outputs = makeOutChannels<o2::detectors::DetID::MFT>(detOrig, mctruth),
455 .algorithm = AlgorithmSpec{adaptFromTask<MFTDPLDigitizerTask>(mctruth)},
456 .options = Options{{"disable-qed", o2::framework::VariantType::Bool, false, {"disable QED handling"}}}};
457}
458
459} // namespace o2::itsmft
460 // end namespace o2
Definition of the base digitizer task class.
A const (ready only) version of MCTruthContainer.
Definition of the ITSMFT digit.
Definition of the Names Generator class.
o2::framework::DataAllocator::SubSpecificationType SubSpecificationType
std::ostringstream debug
uint64_t bc
Definition RawEventData.h:5
int32_t i
Header of the General Run Parameters object.
Definition of the GeometryTGeo class.
Definition of the ITSMFT ROFrame (trigger) record.
Definition of the ITS digitizer.
Definition of the ITSMFT NoiseMap.
Definition of the ITSMFT time-dependend dead map.
virtual void init(o2::framework::InitContext &) final
A read-only version of MCTruthContainer allowing for storage optimisation.
Static class with identifiers, bitmasks and names for ALICE detectors.
Definition DetID.h:58
static constexpr const char * getName(ID id)
names of defined detectors
Definition DetID.h:146
static constexpr ID ITS
Definition DetID.h:63
static constexpr ID MFT
Definition DetID.h:71
static const std::array< std::vector< std::string >, DetID::nDetectors > DETECTORBRANCHNAMES
Definition SimTraits.h:39
void snapshot(const Output &spec, T const &object)
o2::header::DataHeader::SubSpecificationType SubSpecificationType
decltype(auto) make(const Output &spec, Args... args)
ConfigParamRegistry const & options()
Definition InitContext.h:33
decltype(auto) get(R binding, int part=0) const
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.
static GeometryTGeo * Instance()
int getROFrameBiasInBC(int layer=-1) const
Definition DigiParams.h:73
void setContinuous(bool v)
Definition DigiParams.h:53
void setDeadChannelsMap(const o2::itsmft::NoiseMap *mp)
Definition Digitizer.h:64
void fillOutputContainer(uint32_t maxFrame=0xffffffff, int layer=-1)
void setGeometry(const o2::itsmft::GeometryTGeo *gm)
Definition Digitizer.h:88
void setMCLabels(o2::dataformats::MCTruthContainer< o2::MCCompLabel > *mclb)
Definition Digitizer.h:59
void setAlpideResponse(const o2::itsmft::AlpideSimResponse *resp, int i)
Definition Digitizer.h:67
void setNoiseMap(const o2::itsmft::NoiseMap *mp)
Definition Digitizer.h:63
void setROFRecords(std::vector< o2::itsmft::ROFRecord > *rec)
Definition Digitizer.h:60
void setDigits(std::vector< o2::itsmft::Digit > *dig)
Definition Digitizer.h:58
uint32_t getEventROFrameMax() const
Definition Digitizer.h:91
o2::itsmft::DigiParams & getParams()
Definition Digitizer.h:61
void setEventTime(const o2::InteractionTimeRecord &irt, int layer=-1)
void process(const std::vector< Hit > *hits, int evID, int srcID, int layer=-1)
Steer conversion of hits to digits.
uint32_t getEventROFrameMin() const
Definition Digitizer.h:90
std::vector< o2::itsmft::Hit > * mHitsP
void updateTimeDependentParams(ProcessingContext &pc)
o2::parameters::GRPObject::ROMode mROMode
static constexpr o2::detectors::DetID ID
std::array< o2::dataformats::MCTruthContainer< o2::MCCompLabel >, NLayers > mLabels
std::array< std::vector< o2::itsmft::ROFRecord >, NLayers > mROFRecordsAccum
std::array< std::vector< o2::itsmft::MC2ROFRecord >, NLayers > mMC2ROFRecordsAccum
void finaliseCCDB(ConcreteDataMatcher &matcher, void *obj)
void initDigitizerTask(framework::InitContext &ic) override
std::array< std::vector< o2::itsmft::ROFRecord >, NLayers > mROFRecords
std::vector< o2::itsmft::Hit > mHits
std::array< o2::dataformats::MCTruthContainer< o2::MCCompLabel >, NLayers > mLabelsAccum
static constexpr o2::header::DataOrigin Origin
std::array< std::vector< o2::itsmft::Digit >, NLayers > mDigits
void run(framework::ProcessingContext &pc)
NoiseMap class for the ITS and MFT.
Definition NoiseMap.h:39
void decodeMap(NoiseMap &noisemap) const
static GeometryTGeo * Instance()
bool initSimChains(o2::detectors::DetID detid, std::vector< TChain * > &simchains) const
GLsizeiptr size
Definition glcorearb.h:659
GLuint GLuint end
Definition glcorearb.h:469
GLenum GLuint GLint GLint layer
Definition glcorearb.h:1310
constexpr o2::header::DataOrigin gDataOriginMFT
Definition DataHeader.h:572
constexpr o2::header::DataOrigin gDataOriginITS
Definition DataHeader.h:570
constexpr int LHCMaxBunches
constexpr double LHCBunchSpacingNS
Defining PrimaryVertex explicitly as messageable.
std::vector< ConfigParamSpec > ccdbParamSpec(std::string const &path, int runDependent, std::vector< CCDBMetadata > metadata={}, int qrate=0)
std::vector< ConfigParamSpec > Options
header::DataHeader::SubSpecificationType SubSpecificationType
DataProcessorSpec getMFTDigitizerSpec(int channel, bool mctruth)
DataProcessorSpec getITSDigitizerSpec(int channel, bool mctruth)
a couple of static helper functions to create timestamp values for CCDB queries or override obsolete ...
class listing possible services
static constexpr int getNLayers()
static constexpr bool supportsStaggering() noexcept
static constexpr int L2G
Definition Cartesian.h:54
int getNOrbitsPerTF() const
get IR corresponding to start of the HBF
Definition HBFUtils.h:49
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
o2::InteractionRecord ir(0, 0)