Project
Loading...
Searching...
No Matches
Digitizer.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 <TRandom.h>
13#include <cmath>
14#include <numeric>
17#include "FV0Base/Geometry.h"
18#include "FV0Base/Constants.h"
19#include "TF1Convolution.h"
20
22
23using namespace o2::math_utils;
24using namespace o2::fv0;
25
27{
28 mEventId = -1;
29 mSrcId = -1;
30 for (auto& analogSignal : mPmtChargeVsTime) {
31 std::fill_n(std::begin(analogSignal), analogSignal.size(), 0);
32 }
33 mLastBCCache.clear();
34 mCfdStartIndex.fill(0);
35}
36
37//_______________________________________________________________________
39{
40 LOG(info) << "init";
41 mNBins = FV0DigParam::Instance().waveformNbins; //Will be computed using detector set-up from CDB
42 mBinSize = FV0DigParam::Instance().waveformBinWidth; //Will be set-up from CDB
43 mNTimeBinsPerBC = std::lround(o2::constants::lhc::LHCBunchSpacingNS / mBinSize); // 1920 bins/BC
44
45 for (Int_t detID = 0; detID < Constants::nFv0Channels; detID++) {
46 mPmtChargeVsTime[detID].resize(mNBins);
47 mLastBCCache.mPmtChargeVsTime[detID].resize(mNBins);
48 }
49
51 TF1Convolution convolutionRingA1ToA4("expo", "landau", 5.e-09, 90.e-09, false);
52 TF1 convolutionRingA1ToA4Fn("convolutionFn", convolutionRingA1ToA4, 5.e-09, 90.e-09, convolutionRingA1ToA4.GetNpar());
53 convolutionRingA1ToA4Fn.SetParameters(FV0DigParam::Instance().constRingA1ToA4, FV0DigParam::Instance().slopeRingA1ToA4,
54 FV0DigParam::Instance().mpvRingA1ToA4, FV0DigParam::Instance().sigmaRingA1ToA4);
55
57 TF1Convolution convolutionRing5("expo", "landau", 5.e-09, 90.e-09, false);
58 TF1 convolutionRing5Fn("convolutionFn", convolutionRing5, 5.e-09, 90.e-09, convolutionRing5.GetNpar());
59 convolutionRing5Fn.SetParameters(FV0DigParam::Instance().constRing5, FV0DigParam::Instance().slopeRing5,
60 FV0DigParam::Instance().mpvRing5, FV0DigParam::Instance().sigmaRing5);
62 mPmtResponseGlobalRingA1ToA4.resize(mNBins);
63 const float binSizeInNs = mBinSize * 1.e-09; // to convert ns into sec
64 double x = (binSizeInNs) / 2.0;
65 for (auto& y : mPmtResponseGlobalRingA1ToA4) {
66 y = FV0DigParam::Instance().getNormRingA1ToA4() // normalisation to have MIP adc at 16
67 * convolutionRingA1ToA4Fn.Eval(x + FV0DigParam::Instance().offsetRingA1ToA4); // offset to adjust mean position of waveform
68 x += binSizeInNs;
69 }
71 mPmtResponseGlobalRing5.resize(mNBins);
72 x = (binSizeInNs) / 2.0;
73 for (auto& y : mPmtResponseGlobalRing5) {
74 y = FV0DigParam::Instance().getNormRing5() // normalisation to have MIP adc at 16
75 * convolutionRing5Fn.Eval(x + FV0DigParam::Instance().offsetRing5); // offset to adjust mean position of waveform
76 x += binSizeInNs;
77 }
78 mLastBCCache.clear();
79 mCfdStartIndex.fill(0);
80 LOG(info) << "init -> finished";
81}
82
83void Digitizer::process(const std::vector<o2::fv0::Hit>& hits,
84 std::vector<o2::fv0::Digit>& digitsBC,
85 std::vector<o2::fv0::ChannelData>& digitsCh,
86 std::vector<o2::fv0::DetTrigInput>& digitsTrig,
88{
89 LOG(debug) << "Begin with " << hits.size() << " hits";
90 flush(digitsBC, digitsCh, digitsTrig, labels); // flush cached signal which cannot be affect by new event
91
92 std::vector<int> hitIdx(hits.size());
93 std::iota(std::begin(hitIdx), std::end(hitIdx), 0);
94 std::sort(std::begin(hitIdx), std::end(hitIdx),
95 [&hits](int a, int b) { return hits[a].GetTrackID() < hits[b].GetTrackID(); });
96
97 // use ordered hits
98 for (auto ids : hitIdx) {
99 const auto& hit = hits[ids];
100 Int_t detId = hit.GetDetectorID();
101
102 if (mDeadChannelMap && !mDeadChannelMap->isChannelAlive(detId)) {
103 continue;
104 }
105
106 Double_t hitEdep = hit.GetHitValue() * 1e3; // convert to MeV
107 Float_t const hitTime = hit.GetTime() * 1e9; // convert to ns
108 // TODO: check how big is inaccuracy if more than 1 'below-threshold' particles hit the same detector cell
110 continue;
111 }
112 float distanceFromXc = 0;
113 if (Geometry::instance()->isRing5(detId)) {
114 distanceFromXc = getDistFromCellCenter(detId, hit.GetX(), hit.GetY());
115 }
116
117 int iChannelPerCell = 0;
118 while (iChannelPerCell < 2) { // loop over 2 channels, into which signal from each cell in ring 5 is split
119 if (Geometry::instance()->isRing5(detId)) {
120 // The first channel number is located counter-clockwise from the cell center
121 // and remains identical to the detector number, the second one is clockwise and incremented by 8
122 if (iChannelPerCell == 1) {
123 detId += 8;
124 }
125 // Split signal magnitude to fractions depending on the distance of the hit from the cell center
126 hitEdep = (hit.GetHitValue() * 1e3) * getSignalFraction(distanceFromXc, iChannelPerCell == 0);
127 // LOG(info) << " detId: " << detId << "-" << iChannelPerCell << " hitEdep: " << hitEdep << " distanceFromXc: " << distanceFromXc;
128 ++iChannelPerCell;
129 } else {
130 iChannelPerCell = 2; // not a ring 5 cell -> don't repeat the loop
131 }
132 Double_t const nPhotons = hitEdep * DP::N_PHOTONS_PER_MEV;
133 float const nPhE = SimulateLightYield(detId, nPhotons);
134 float const mipFraction = float(nPhE / FV0DigParam::Instance().avgNumberPhElectronPerMip);
135 Long64_t timeHit = hitTime;
136 timeHit += mIntRecord.getTimeNS();
137 o2::InteractionTimeRecord const irHit(timeHit);
138 std::array<o2::InteractionRecord, NBC2Cache> cachedIR;
139 int nCachedIR = 0;
140 for (int i = BCCacheMin; i < BCCacheMax + 1; i++) {
141 double const tNS = timeHit + o2::constants::lhc::LHCBunchSpacingNS * i;
142 cachedIR[nCachedIR].setFromNS(tNS);
143 if (tNS < 0 && cachedIR[nCachedIR] > irHit) {
144 continue; // don't go to negative BC/orbit (it will wrap)
145 }
146 // ensure existence of cached container
147 setBCCache(cachedIR[nCachedIR++]);
148 } // BCCache loop
149
150 createPulse(mipFraction, hit.GetTrackID(), hitTime, hit.GetPos().R(), cachedIR, nCachedIR, detId);
151
152 } //while loop
153 } //hitloop
154}
155
156void Digitizer::createPulse(float mipFraction, int parID, const double hitTime, const float hitR,
157 std::array<o2::InteractionRecord, NBC2Cache> const& cachedIR, int nCachedIR, const int detId)
158{
159
160 std::array<bool, NBC2Cache> added;
161 added.fill(false);
162
163 for (int ir = 0; ir < NBC2Cache; ir++) {
164 auto bcCache = getBCCache(cachedIR[ir]);
165 for (int ich = 0; ich < Constants::nFv0Channels; ich++) {
166 (*bcCache).mPmtChargeVsTime[ich].resize(mNTimeBinsPerBC);
167 }
168 }
169
170 // Subtract time-of-flight from hit time
171 const float timeOfFlight = hitR / o2::constants::physics::LightSpeedCm2NS;
172 Int_t const NBinShift = std::lround((hitTime - timeOfFlight + FV0DigParam::Instance().hitTimeOffset) / FV0DigParam::Instance().waveformBinWidth);
173
174 if (NBinShift >= 0 && NBinShift < FV0DigParam::Instance().waveformNbins) {
175 mPmtResponseTemp.resize(FV0DigParam::Instance().waveformNbins, 0.);
176 if (isRing5(detId)) {
177 std::memcpy(&mPmtResponseTemp[NBinShift], &mPmtResponseGlobalRing5[0],
178 sizeof(double) * (FV0DigParam::Instance().waveformNbins - NBinShift));
179 } else {
180 std::memcpy(&mPmtResponseTemp[NBinShift], &mPmtResponseGlobalRingA1ToA4[0],
181 sizeof(double) * (FV0DigParam::Instance().waveformNbins - NBinShift));
182 }
183 } else {
184 if (isRing5(detId)) {
185 mPmtResponseTemp = mPmtResponseGlobalRing5;
186 mPmtResponseTemp.erase(mPmtResponseTemp.begin(), mPmtResponseTemp.begin() + abs(NBinShift));
187 } else {
188 mPmtResponseTemp = mPmtResponseGlobalRingA1ToA4;
189 mPmtResponseTemp.erase(mPmtResponseTemp.begin(), mPmtResponseTemp.begin() + abs(NBinShift));
190 }
191
192 mPmtResponseTemp.resize(FV0DigParam::Instance().waveformNbins);
193 }
194
195 for (int ir = 0; ir < int(mPmtResponseTemp.size() / mNTimeBinsPerBC); ir++) {
196 auto bcCache = getBCCache(cachedIR[ir]);
197
198 for (int iBin = 0; iBin < mNTimeBinsPerBC; iBin++) {
199 (*bcCache).mPmtChargeVsTime[detId][iBin] += (mPmtResponseTemp[ir * mNTimeBinsPerBC + iBin] * mipFraction);
200 }
201 added[ir] = true;
202 }
204 for (int ir = 0; ir < nCachedIR; ir++) {
205 if (added[ir]) {
206 auto bcCache = getBCCache(cachedIR[ir]);
207 (*bcCache).labels.emplace_back(parID, mEventId, mSrcId, detId);
208 }
209 }
210}
211
212void Digitizer::flush(std::vector<o2::fv0::Digit>& digitsBC,
213 std::vector<o2::fv0::ChannelData>& digitsCh,
214 std::vector<o2::fv0::DetTrigInput>& digitsTrig,
216{
217 ++mEventId;
218 while (!mCache.empty()) {
219 auto const& bc = mCache.front();
220 if (mIntRecord.differenceInBC(bc) > NBC2Cache) { // Build events that are separated by NBC2Cache BCs from current BC
221 storeBC(bc, digitsBC, digitsCh, digitsTrig, labels);
222 mCache.pop_front();
223 } else {
224 return;
225 }
226 }
227}
228
229void Digitizer::storeBC(const BCCache& bc,
230 std::vector<o2::fv0::Digit>& digitsBC,
231 std::vector<o2::fv0::ChannelData>& digitsCh,
232 std::vector<o2::fv0::DetTrigInput>& digitsTrig,
234
235{
236 size_t const nBC = digitsBC.size(); // save before digitsBC is being modified
237 size_t const first = digitsCh.size(); // save before digitsCh is being modified
238 int8_t nTotFiredCells = 0;
239 int8_t nTrgFiredCells = 0; // number of fired cells, that follow additional trigger conditions (time gate)
240 int totalChargeAllRing = 0;
241 int32_t avgTime = 0;
242 double nSignalInner = 0;
243 double nSignalOuter = 0;
244
245 if (mLastBCCache.differenceInBC(bc) != 1) { // if the last buffered BC is not the one before the current BC
246 mLastBCCache.clear(); // clear the bufffer (mPmtChargeVsTime set to 0s)
247 mCfdStartIndex.fill(0); // reset all start indices to 0, i.e., to the beginning of the BC
248 }
249
250 for (int iPmt = 0; iPmt < Constants::nFv0Channels; iPmt++) {
251 // run the CFD: this updates the start index for the next BC in case the CFD dead time ends in the next BC
252 double cfdWithOffset = SimulateTimeCfd(mCfdStartIndex[iPmt], mLastBCCache.mPmtChargeVsTime[iPmt], bc.mPmtChargeVsTime[iPmt]);
253 double cfdZero = cfdWithOffset - FV0DigParam::Instance().avgCfdTimeForMip;
254
255 // Conditions to sum charge are: all participating channels must have time within +/- 2.5 ns, AND
256 // at least one channel must follow more strict conditions (see below)
257 if (cfdZero < -FV0DigParam::Instance().cfdCheckWindow || cfdZero > FV0DigParam::Instance().cfdCheckWindow) {
258 continue;
259 }
260
261 int iTotalCharge = std::lround(IntegrateCharge(bc.mPmtChargeVsTime[iPmt]) * DP::INV_CHARGE_PER_ADC); // convert Coulomb to adc;
262
263 uint8_t channelBits = FV0DigParam::Instance().defaultChainQtc;
264 if (std::rand() % 2) {
266 }
267 if (iTotalCharge > (FV0DigParam::Instance().maxCountInAdc) && FV0DigParam::Instance().useMaxChInAdc) {
268 iTotalCharge = FV0DigParam::Instance().maxCountInAdc; // max adc channel for one PMT
270 }
271
272 if (iTotalCharge < FV0DigParam::Instance().getCFDTrshInAdc()) {
273 continue;
274 }
275
276 int iCfdZero = std::lround(cfdZero * DP::INV_TIME_PER_TDCCHANNEL);
277 digitsCh.emplace_back(iPmt, iCfdZero, iTotalCharge, channelBits);
278 ++nTotFiredCells;
279
280 int triggerGate = FV0DigParam::Instance().mTime_trg_gate;
281 if (std::abs(iCfdZero) < triggerGate) {
282 ++nTrgFiredCells;
283 //---trigger---
284 totalChargeAllRing += iTotalCharge;
285 avgTime += iCfdZero;
286 if (iPmt < 24) {
287 nSignalInner++;
288 } else {
289 nSignalOuter++;
290 }
291 }
292 }
293 // save BC information for the CFD detector
294 mLastBCCache = bc;
295 if (nTotFiredCells < 1) {
296 return;
297 }
298 if (nTrgFiredCells > 0) {
299 avgTime /= nTrgFiredCells;
300 } else {
302 }
304 bool isA, isAIn, isAOut, isCen, isSCen;
305 isA = nTrgFiredCells > 0;
306 isAIn = nSignalInner > 0; // ring 1,2 and 3
307 isAOut = nSignalOuter > 0; // ring 4 and 5
308 isCen = totalChargeAllRing > FV0DigParam::Instance().adcChargeCenThr;
309 isSCen = totalChargeAllRing > FV0DigParam::Instance().adcChargeSCenThr;
310
311 Triggers triggers;
312 const int unusedCharge = o2::fit::Triggers::DEFAULT_AMP;
313 const int unusedTime = o2::fit::Triggers::DEFAULT_TIME;
314 const int unusedZero = o2::fit::Triggers::DEFAULT_ZERO;
315 const bool unusedBitsInSim = false; // bits related to laser and data validity
316 const bool bitDataIsValid = true;
317 triggers.setTriggers(isA, isAIn, isAOut, isCen, isSCen, nTrgFiredCells, (int8_t)unusedZero,
318 (int32_t)(0.125 * totalChargeAllRing), (int32_t)unusedCharge, (int16_t)avgTime, (int16_t)unusedTime, unusedBitsInSim, unusedBitsInSim, bitDataIsValid);
319 digitsBC.emplace_back(first, nTotFiredCells, bc, triggers, mEventId - 1);
320 digitsTrig.emplace_back(bc, isA, isAIn, isAOut, isCen, isSCen);
321 for (auto const& lbl : bc.labels) {
322 labels.addElement(nBC, lbl);
323 }
324}
325
326// -------------------------------------------------------------------------------
327// --- Internal helper methods related to conversion of energy-deposition into ---
328// --- photons -> photoelectrons -> electrical signal ---
329// -------------------------------------------------------------------------------
330Int_t Digitizer::SimulateLightYield(Int_t pmt, Int_t nPhot) const
331{
332 const Float_t epsilon = 0.0001f;
334 if ((fabs(1.0f - p) < epsilon) || nPhot == 0) {
335 return nPhot;
336 }
337 const Int_t n = Int_t(nPhot < 100
338 ? gRandom->Binomial(nPhot, p)
339 : gRandom->Gaus((p * nPhot) + 0.5, TMath::Sqrt(p * (1. - p) * nPhot)));
340 return n;
341}
342//---------------------------------------------------------------------------
343Float_t Digitizer::IntegrateCharge(const ChannelDigitF& pulse) const
344{
345 int const chargeIntMin = FV0DigParam::Instance().isIntegrateFull ? 0 : (FV0DigParam::Instance().avgCfdTimeForMip - 6.0) / mBinSize; //Charge integration offset (cfd mean time - 6 ns)
346 int const chargeIntMax = FV0DigParam::Instance().isIntegrateFull ? mNTimeBinsPerBC : (FV0DigParam::Instance().avgCfdTimeForMip + 14.0) / mBinSize; //Charge integration offset (cfd mean time + 14 ns)
347 if (chargeIntMin < 0 || chargeIntMin > mNTimeBinsPerBC || chargeIntMax > mNTimeBinsPerBC) {
348 LOG(fatal) << "invalid indicess: chargeInMin=" << chargeIntMin << " chargeIntMax=" << chargeIntMax;
349 }
350 Float_t totalCharge = 0.0f;
351 for (int iTimeBin = chargeIntMin; iTimeBin < chargeIntMax; iTimeBin++) {
352 totalCharge += pulse[iTimeBin];
353 }
354 return totalCharge;
355}
356//---------------------------------------------------------------------------
357Float_t Digitizer::SimulateTimeCfd(int& startIndex, const ChannelDigitF& pulseLast, const ChannelDigitF& pulse) const
358{
359 Float_t timeCfd = -1024.0f;
360
361 if (pulse.empty()) {
362 startIndex = 0;
363 return timeCfd;
364 }
365
366 Float_t const cfdThrInCoulomb = FV0DigParam::Instance().mCFD_trsh * 1e-3 / 50 * mBinSize * 1e-9; // convert mV into Coulomb assuming 50 Ohm
367
368 Int_t const binShift = TMath::Nint(FV0DigParam::Instance().timeShiftCfd / mBinSize);
369 Float_t sigPrev = 5 * pulseLast[mNTimeBinsPerBC - binShift - 1] - pulseLast[mNTimeBinsPerBC - 1]; // CFD output from the last bin of the last BC
370 for (Int_t iTimeBin = 0; iTimeBin < mNTimeBinsPerBC; ++iTimeBin) {
371 Float_t const sigCurrent = 5.0f * (iTimeBin >= binShift ? pulse[iTimeBin - binShift] : pulseLast[mNTimeBinsPerBC - binShift + iTimeBin]) - pulse[iTimeBin];
372 if (iTimeBin >= startIndex && std::abs(pulse[iTimeBin]) > cfdThrInCoulomb) { // enable
373 if (sigPrev < 0.0f && sigCurrent >= 0.0f) { // test for zero-crossing
374 timeCfd = Float_t(iTimeBin) * mBinSize;
375 startIndex = iTimeBin + std::lround(FV0DigParam::Instance().mCfdDeadTime / mBinSize); // update startIndex (CFD dead time)
376 if (startIndex < mNTimeBinsPerBC) {
377 startIndex = 0; // dead-time ends in same BC: no impact on the following BC
378 } else {
379 startIndex -= mNTimeBinsPerBC;
380 }
381 if (startIndex > mNTimeBinsPerBC) {
382 LOG(fatal) << "CFD dead-time was set to > 25 ns";
383 }
384 break; // only detects the 1st zero-crossing in the BC
385 }
386 }
387 sigPrev = sigCurrent;
388 }
389 return timeCfd;
390}
391
392float Digitizer::getDistFromCellCenter(UInt_t cellId, double hitx, double hity)
393{
395
396 // Parametrize the line (ax+by+c=0) that crosses the detector center and the cell's middle point
397 Point3Dsimple* pCell = &geo->getCellCenter(cellId);
398 float x0, y0, z0;
399 geo->getGlobalPosition(x0, y0, z0);
400 double a = -(y0 - pCell->y) / (x0 - pCell->x);
401 double b = 1;
402 double c = -(y0 - a * x0);
403 //Return the distance from hit to this line
404 return (a * hitx + b * hity + c) / TMath::Sqrt(a * a + b * b);
405}
406
407float Digitizer::getSignalFraction(float distanceFromXc, bool isFirstChannel)
408{
409 float const fraction = sigmoidPmtRing5(distanceFromXc);
410 if (distanceFromXc > 0) {
411 return isFirstChannel ? fraction : (1. - fraction);
412 } else {
413 return isFirstChannel ? (1. - fraction) : fraction;
414 }
415}
416
417//_____________________________________________________________________________
418o2::fv0::Digitizer::BCCache& Digitizer::setBCCache(const o2::InteractionRecord& ir)
419{
420 if (mCache.empty() || mCache.back() < ir) {
421 mCache.emplace_back();
422 auto& cb = mCache.back();
423 cb = ir;
424 return cb;
425 }
426 if (mCache.front() > ir) {
427 mCache.emplace_front();
428 auto& cb = mCache.front();
429 cb = ir;
430 return cb;
431 }
432 for (auto cb = mCache.begin(); cb != mCache.end(); cb++) {
433 if ((*cb) == ir) {
434 return *cb;
435 }
436 if (ir < (*cb)) {
437 auto cbnew = mCache.emplace(cb); // insert new element before cb
438 (*cbnew) = ir;
439 return (*cbnew);
440 }
441 }
442 return mCache.front();
443}
444//_____________________________________________________________________________
445o2::fv0::Digitizer::BCCache* Digitizer::getBCCache(const o2::InteractionRecord& ir)
446{
447 // get pointer on existing cache
448 for (auto cb = mCache.begin(); cb != mCache.end(); cb++) {
449 if ((*cb) == ir) {
450 return &(*cb);
451 }
452 }
453 return nullptr;
454}
455
456bool Digitizer::isRing5(int detID)
457{
458 if (detID > 31) {
459 return true;
460 } else {
461 return false;
462 }
463}
464
std::vector< std::string > labels
#define O2ParamImpl(classname)
General constants in FV0.
std::ostringstream debug
Base definition of FV0 geometry.
ClassImp(o2::fv0::Digitizer)
uint64_t bc
Definition RawEventData.h:5
Configurable digitization parameters.
int32_t i
A container to hold and manage MC truth information/labels.
static const int16_t DEFAULT_AMP
Definition Triggers.h:51
static const int16_t DEFAULT_TIME
Definition Triggers.h:50
void setTriggers(uint8_t trgsig, uint8_t chanA, uint8_t chanC, int32_t aamplA, int32_t aamplC, int16_t atimeA, int16_t atimeC)
Definition Triggers.h:103
static const int16_t DEFAULT_ZERO
Definition Triggers.h:52
void flush(std::vector< o2::fv0::Digit > &digitsBC, std::vector< o2::fv0::ChannelData > &digitsCh, std::vector< o2::fv0::DetTrigInput > &digitsTrig, o2::dataformats::MCTruthContainer< o2::fv0::MCLabel > &labels)
void process(const std::vector< o2::fv0::Hit > &hits, std::vector< o2::fv0::Digit > &digitsBC, std::vector< o2::fv0::ChannelData > &digitsCh, std::vector< o2::fv0::DetTrigInput > &digitsTrig, o2::dataformats::MCTruthContainer< o2::fv0::MCLabel > &labels)
Definition Digitizer.cxx:83
FV0 Geometry.
Definition Geometry.h:51
static Geometry * instance(EGeoType initType=eUninitialized)
Point3Dsimple & getCellCenter(UInt_t cellId)
Definition Geometry.cxx:125
bool isRing5(UInt_t cellId)
Definition Geometry.cxx:135
void getGlobalPosition(float &x, float &y, float &z)
Utility functions to be accessed externally.
Definition Geometry.cxx:118
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
GLuint * ids
Definition glcorearb.h:647
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLuint GLfloat x0
Definition glcorearb.h:5034
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLuint GLfloat GLfloat y0
Definition glcorearb.h:5034
constexpr double LHCBunchSpacingNS
constexpr float LightSpeedCm2NS
float sigmoidPmtRing5(float x)
Definition Digitizer.h:152
float float & c
Definition Utils.h:111
int64_t differenceInBC(const InteractionRecord &other) const
double getTimeNS() const
get time in ns from orbit=0/bc=0
const bool isChannelAlive(const uint8_t &chId) const
void setFlag(uint8_t flag)
Definition ChannelData.h:55
static constexpr int nFv0Channels
Definition Constants.h:32
static constexpr float N_PHOTONS_PER_MEV
static constexpr float INV_CHARGE_PER_ADC
static constexpr float INV_TIME_PER_TDCCHANNEL
std::array< ChannelDigitF, Constants::nFv0Channels > mPmtChargeVsTime
Definition Digitizer.h:75
float getNormRingA1ToA4() const
Definition FV0DigParam.h:37
float getNormRing5() const
Definition FV0DigParam.h:45
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"
o2::InteractionRecord ir(0, 0)