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
14#include "FT0Base/FT0DigParam.h"
18
20#include "FT0Base/Constants.h"
21#include <map>
22#include <array>
23#include <regex>
24#include <string>
25
26#include "TMath.h"
27#include "TRandom.h"
28#include <TH1F.h>
29#include <algorithm>
30#include <cassert>
31#include <iostream>
32#include <optional>
33#include <Vc/Vc>
34
35using namespace o2::ft0;
36
38
39namespace o2::ft0
40{
41// signal shape function
42template <typename Float>
44{
45 float p0, p1, p2, p3, p4, p5, p6, p7;
46 p0 = 1.30853;
47 p1 = 0.516807;
48 p2 = 3.36714;
49 p3 = -1.01206;
50 p4 = 1.42832;
51 p5 = 1.1589;
52 p6 = 1.22019;
53 p7 = 0.426818;
54
55 Double_t val = 0;
56
57 if (x > p3) {
58 Double_t x1 = x - p3;
59 Double_t arg1 = (log(x1) - p0) / p1;
60 val += p2 * (1.0 / (x1 * p1 * sqrt(2 * TMath::Pi()))) * exp(-0.5 * arg1 * arg1);
61 }
62
63 if (x > p7) {
64 Double_t x2 = x - p7;
65 Double_t arg2 = (log(x2) - p4) / p5;
66 val += p6 * (1.0 / (x2 * p5 * sqrt(2 * TMath::Pi()))) * exp(-0.5 * arg2 * arg2);
67 }
68
69 return val;
70};
71
72// integrated signal shape function
73inline float signalForm_integral(float x)
74{
75 float p0, p1, p2, p3, p4, p5, p6, p7;
76 p0 = 1.30853;
77 p1 = 0.516807;
78 p2 = 3.36714;
79 p3 = -1.01206;
80 p4 = 1.42832;
81 p5 = 1.1589;
82 p6 = 1.22019;
83 p7 = 0.426818;
84 Double_t val = 0;
85
86 if (x > p3) {
87 Double_t x1 = x - p3;
88 Double_t z1 = (log(x1) - p0) / (sqrt(2) * p1);
89 val += p2 * 0.5 * (1 + TMath::Erf(z1)); // norm1 * CDF1
90 }
91
92 if (x > p7) {
93 Double_t x2 = x - p7;
94 Double_t z2 = (log(x2) - p4) / (sqrt(2) * p5);
95 val += p6 * 0.5 * (1 + TMath::Erf(z2)); // norm2 * CDF2
96 }
97
98 return val;
99};
100/*
101// signal shape function
102template <typename Float>
103Float signalForm_i(Float x)
104{
105using namespace std;
106Float const a = -0.45458;
107Float const b = -0.83344945;
108return x > Float(0) ? -(exp(b * x) - exp(a * x)) / Float(7.8446501) : Float(0);
109};
110
111// integrated signal shape function
112inline float signalForm_integral(float x)
113{
114using namespace std;
115double const a = -0.45458;
116double const b = -0.83344945;
117if (x < 0) {
118 x = 0;
119}
120return -(exp(b * x) / b - exp(a * x) / a) / 7.8446501;
121};
122*/
123// SIMD version of the integrated signal shape function
124inline Vc::float_v signalForm_integralVc(Vc::float_v x)
125{
126 auto const mask = (x >= 0.0f);
127 Vc::float_v arg(0);
128 arg.assign(x, mask); // branchless if
129 Vc::float_v const a(-0.45458f);
130 Vc::float_v const b(-0.83344945f);
131 Vc::float_v result = -(Vc::exp(b * arg) / b - Vc::exp(a * arg) / a) / 7.8446501f;
132 return result;
133};
134} // namespace o2::ft0
135
136Digitizer::CFDOutput Digitizer::get_time(const std::vector<float>& times, float deadTime)
137{
138 assert(std::is_sorted(std::begin(times), std::end(times)));
139
140 // get a new batch of random values
141 for (float& n : mNoiseSamples) {
142 n = mRndGaus.getNextValue();
143 }
144
145 // this lambda function evaluates the signal shape at a give time
146 auto value_at = [&](float time) {
147 float val = 0;
148 // (1) sum over individual hits
149 Vc::float_v acc(0);
150 Vc::float_v tableVal(0);
151 const float* tp = times.data();
152 size_t m = times.size() / Vc::float_v::size();
153 for (size_t i = 0; i < m; ++i) {
154 tableVal.load(tp);
155 tp += Vc::float_v::size();
156 Vc::prefetchForOneRead(tp);
157 acc += signalFormVc(time - tableVal);
158 }
159 val += acc.sum();
160 // non-SIMD tail
161 for (size_t i = Vc::float_v::size() * m; i < times.size(); ++i, ++tp) {
162 val += signalForm(time - (*tp));
163 }
164 // (2) add noise
165 // find the right indices into the sinc table
166 int timeIndex = std::lround(time / FT0DigParam::Instance().mNoisePeriod * mSincTable.size());
167 int timeOffset = timeIndex / mSincTable.size();
168 timeIndex %= mSincTable.size();
169 if (timeOffset >= mNumNoiseSamples) { // this happens when time >= 25 ns
170 timeOffset = mNumNoiseSamples - 1;
171 LOG(debug) << "timeOffset >= mNumNoiseSamples";
172 }
173 if (timeOffset <= -mNumNoiseSamples) { // this happens when time <= -25 ns
174 timeOffset = -mNumNoiseSamples + 1;
175 LOG(debug) << "timeOffset <= -mNumNoiseSamples";
176 }
177 Vc::float_v noiseVal(0);
178 const float* np = mNoiseSamples.data();
179 tp = mSincTable[timeIndex].data() + mNumNoiseSamples - timeOffset;
180 acc = 0.0f;
181 m = mNumNoiseSamples / Vc::float_v::Size;
182 for (size_t i = 0; i < m; ++i) {
183 tableVal.load(tp);
184 tp += Vc::float_v::Size;
185 Vc::prefetchForOneRead(tp);
186 noiseVal.load(np);
187 np += Vc::float_v::Size;
188 Vc::prefetchForOneRead(np);
189 acc += noiseVal * tableVal;
190 }
191 val += acc.sum(); // horizontal sum
192 // non-SIMD tail
193 for (size_t i = Vc::float_v::Size * m; i < mNumNoiseSamples; ++i, ++tp, ++np) {
194 val += (*np) * (*tp);
195 }
196 return val;
197 };
198 auto const min_time = std::max(deadTime, *std::min_element(std::begin(times),
199 std::end(times)));
200 CFDOutput result{std::nullopt, -0.5f * FT0DigParam::Instance().mBunchWidth};
201 bool is_positive = true;
202
203 // reset the chache
204 std::fill_n(std::begin(mSignalCache), std::size(mSignalCache), -1.0f);
205 const auto& params = FT0DigParam::Instance();
206 // we need double precision for time in order to match previous behaviour
207 for (double time = min_time; time < 0.5 * params.mBunchWidth; time += DP::SIGNAL_CACHE_DT) {
208 float const val = value_at(time);
209 int const index = std::lround((time + 0.5 * params.mBunchWidth) / DP::SIGNAL_CACHE_DT);
210 if (index >= 0 && index < mSignalCache.size()) { // save the value for later use
211 mSignalCache[index] = val;
212 }
213 // look up the time-shifted signal value from the past
214 float val_prev = 0.0f;
215 int const index_prev = std::lround((time - params.mCFDShiftPos + 0.5f * params.mBunchWidth) / DP::SIGNAL_CACHE_DT);
216 val_prev = ((index_prev < 0 || index_prev >= mSignalCache.size() || mSignalCache[index_prev] < 0.0f)
217 ? value_at(time - params.mCFDShiftPos) // was not computed before
218 : mSignalCache[index_prev]); // is available in the cache
219 float const cfd_val = 5.0f * val_prev - val;
220 if (std::abs(val) > params.mCFD_trsh && !is_positive && cfd_val > 0.0f) {
221 if (!result.particle) {
222 result.particle = time;
223 }
224 result.deadTime = time + params.mCFDdeadTime;
225 time += params.mCFDdeadTime - DP::SIGNAL_CACHE_DT;
226 is_positive = true;
227 } else {
228 is_positive = cfd_val > 0.0f;
229 }
230 }
231 return result;
232}
233
234double Digitizer::measure_amplitude(const std::vector<float>& times) const
235{
236 float const from = FT0DigParam::Instance().mAmpRecordLow;
237 float const to = from + FT0DigParam::Instance().mAmpRecordUp;
238 // SIMD version has a negligible effect on the total wall time
239 Vc::float_v acc(0);
240 Vc::float_v tv(0);
241 const float* tp = times.data();
242 size_t const m = times.size() / Vc::float_v::Size;
243 for (size_t i = 0; i < m; ++i) {
244 tv.load(tp);
245 tp += Vc::float_v::Size;
246 Vc::prefetchForOneRead(tp);
247 acc += signalForm_integralVc(to - tv) - signalForm_integralVc(from - tv);
248 }
249 float result = acc.sum(); // horizontal sum
250 // non-SIMD tail
251 for (size_t i = Vc::float_v::Size * m; i < times.size(); ++i, ++tp) {
252 result += signalForm_integral(to - (*tp)) - signalForm_integral(from - (*tp));
253 }
254 return result;
255}
256
257void Digitizer::process(const std::vector<o2::ft0::HitType>* hits,
258 std::vector<o2::ft0::Digit>& digitsBC,
259 std::vector<o2::ft0::ChannelData>& digitsCh,
260 std::vector<o2::ft0::DetTrigInput>& digitsTrig,
262{
263 ;
264 // Calculating signal time, amplitude in mean_time +- time_gate --------------
265 LOG(debug) << " process firstBCinDeque " << firstBCinDeque << " mIntRecord " << mIntRecord;
266 if (firstBCinDeque != mIntRecord) {
267 flush(digitsBC, digitsCh, digitsTrig, label);
268 }
269
270 Int_t parent = -10;
271 for (auto const& hit : *hits) {
272 if (hit.GetEnergyLoss() > 0) {
273 continue;
274 }
275
276 Int_t hit_ch = hit.GetDetectorID();
277
278 // If the dead channel map is used, and the channel with ID 'hit_ch' is dead, don't process this hit.
279 if (mDeadChannelMap && !mDeadChannelMap->isChannelAlive(hit_ch)) {
280 continue;
281 }
282
283 const auto& params = FT0DigParam::Instance();
284
285 Bool_t is_A_side = (hit_ch < 4 * mGeometry.NCellsA);
286
287 // Subtract time-of-flight from hit time
288 const Float_t timeOfFlight = hit.GetPos().R() / o2::constants::physics::LightSpeedCm2NS;
289 const Float_t timeOffset = is_A_side ? params.hitTimeOffsetA : params.hitTimeOffsetC;
290 Double_t hit_time = hit.GetTime() - timeOfFlight + timeOffset + mIntRecord.getTimeOffsetWrtBC();
291
292 if (hit_time > 150) {
293 continue; // not collect very slow particles
294 }
295
296 auto relBC = o2::InteractionRecord{hit_time};
297 if (mCache.size() <= relBC.bc) {
298 mCache.resize(relBC.bc + 1);
299 }
300 mCache[relBC.bc].hits.emplace_back(BCCache::particle{hit_ch, hit_time - relBC.bc2ns()});
301 // charge particles in MCLabel
302 Int_t parentID = hit.GetTrackID();
303 if (parentID != parent) {
304 mCache[relBC.bc].labels.emplace(parentID, mEventID, mSrcID, hit_ch);
305 parent = parentID;
306 }
307 }
308}
309
310void Digitizer::storeBC(BCCache& bc,
311 std::vector<o2::ft0::Digit>& digitsBC,
312 std::vector<o2::ft0::ChannelData>& digitsCh,
313 std::vector<o2::ft0::DetTrigInput>& digitsTrig,
315{
316 if (bc.hits.empty()) {
317 return;
318 }
319 // Initialize mapping channelID -> PM hash and PM side (A/C) using FT0 LUT
320 static bool pmLutInitialized = false;
321 static std::array<uint8_t, o2::ft0::Constants::sNCHANNELS_PM> mChID2PMhash{};
322 static std::map<uint8_t, bool> mMapPMhash2isAside; // hashed PM -> is A side
323
324 if (!pmLutInitialized) {
325 std::map<std::string, uint8_t> mapFEE2hash; // module name -> hashed PM id
326 uint8_t tcmHash = 0;
327
328 const auto& lut = o2::ft0::SingleLUT::Instance().getVecMetadataFEE();
329 auto lutSorted = lut;
330 std::sort(lutSorted.begin(), lutSorted.end(),
331 [](const auto& first, const auto& second) { return first.mModuleName < second.mModuleName; });
332
333 uint8_t binPos = 0;
334 for (const auto& lutEntry : lutSorted) {
335 const auto& moduleName = lutEntry.mModuleName;
336 const auto& moduleType = lutEntry.mModuleType;
337 const auto& strChID = lutEntry.mChannelID;
338
339 auto [it, inserted] = mapFEE2hash.insert({moduleName, binPos});
340 if (inserted) {
341 if (moduleName.find("PMA") != std::string::npos) {
342 mMapPMhash2isAside.insert({binPos, true});
343 } else if (moduleName.find("PMC") != std::string::npos) {
344 mMapPMhash2isAside.insert({binPos, false});
345 }
346 ++binPos;
347 }
348
349 if (std::regex_match(strChID, std::regex("^[0-9]{1,3}$"))) {
350 int chID = std::stoi(strChID);
352 mChID2PMhash[chID] = mapFEE2hash[moduleName];
353 } else {
354 LOG(fatal) << "Incorrect LUT entry: chID " << strChID << " | " << moduleName;
355 }
356 } else if (moduleType != "TCM") {
357 LOG(fatal) << "Non-TCM module w/o numerical chID: chID " << strChID << " | " << moduleName;
358 } else { // TCM
359 tcmHash = mapFEE2hash[moduleName];
360 }
361 }
362
363 pmLutInitialized = true;
364 }
365
366 int n_hit_A = 0, n_hit_C = 0, mean_time_A = 0, mean_time_C = 0;
367 int summ_ampl_A = 0, summ_ampl_C = 0;
368 int sum_A_ampl = 0, sum_C_ampl = 0;
369 int nPMTs = mGeometry.NCellsA * 4 + mGeometry.NCellsC * 4;
370 std::vector<int> sum_ampl_ipmt(nPMTs, 0);
371 // Per-PM summed charge (like in digits2trgFT0)
372 std::map<uint8_t, int> mapPMhash2sumAmpl;
373 for (const auto& entry : mMapPMhash2isAside) {
374 mapPMhash2sumAmpl.insert({entry.first, 0});
375 }
376
377 int vertex_time;
378 const auto& params = FT0DigParam::Instance();
379 int first = digitsCh.size(), nStored = 0;
380 auto& particles = bc.hits;
381 std::sort(std::begin(particles), std::end(particles));
382 auto channel_end = particles.begin();
383 std::vector<float> channel_times;
384 for (Int_t ipmt = 0; ipmt < params.mMCPs; ++ipmt) {
385 auto channel_begin = channel_end;
386 channel_end = std::find_if(channel_begin, particles.end(),
387 [ipmt](BCCache::particle const& p) { return p.hit_ch != ipmt; });
388
389 // The hits between 'channel_begin' and 'channel_end' now contains all hits for channel 'ipmt'
390
391 if (channel_end - channel_begin < params.mAmp_trsh) {
392 continue;
393 }
394 channel_times.resize(channel_end - channel_begin);
395 std::transform(channel_begin, channel_end, channel_times.begin(), [](BCCache::particle const& p) { return p.hit_time; });
396 int chain = (std::rand() % 2) ? 1 : 0;
397 auto cfd = get_time(channel_times, mDeadTimes[ipmt].intrec.bc2ns() -
398 firstBCinDeque.bc2ns() +
399 mDeadTimes[ipmt].deadTime);
400 mDeadTimes[ipmt].intrec = firstBCinDeque;
401 mDeadTimes[ipmt].deadTime = cfd.deadTime;
402
403 if (!cfd.particle) {
404 continue;
405 }
406 // miscalibrate CFD with cahnnel offsets
407 int miscalib = 0;
408 if (mCalibOffset) {
409 miscalib = mCalibOffset->mTimeOffsets[ipmt];
410 }
411 int smeared_time = 1000. * (*cfd.particle - params.mCfdShift) * params.mChannelWidthInverse + miscalib; // + int(1000. * mIntRecord.getTimeOffsetWrtBC() * params.mChannelWidthInverse);
412 bool is_time_in_signal_gate = (smeared_time > -params.mTime_trg_gate && smeared_time < params.mTime_trg_gate);
413 float charge = measure_amplitude(channel_times) * params.mCharge2amp;
414 float amp = is_time_in_signal_gate ? params.mMV_2_Nchannels * charge : 0;
415 if (amp > 4095) {
416 amp = 4095;
417 }
418
419 LOG(debug) << mEventID << " bc " << firstBCinDeque.bc << " orbit " << firstBCinDeque.orbit << ", ipmt " << ipmt << ", smeared_time " << smeared_time << " nStored " << nStored << " offset " << miscalib;
420 if (is_time_in_signal_gate) {
423 // Sum channel charge per PM (similar logic as in digits2trgFT0)
425 mapPMhash2sumAmpl[mChID2PMhash[static_cast<uint8_t>(ipmt)]] += static_cast<int>(amp);
426 }
427 }
428 digitsCh.emplace_back(ipmt, smeared_time, int(amp), chain);
429 nStored++;
430
431 // fill triggers
432
433 Bool_t is_A_side = (ipmt < 4 * mGeometry.NCellsA);
434 if (!is_time_in_signal_gate) {
435 continue;
436 }
437
438 sum_ampl_ipmt[ipmt] += amp;
439
440 if (is_A_side) {
441 n_hit_A++;
442 summ_ampl_A += amp;
443 mean_time_A += smeared_time;
444 } else {
445 n_hit_C++;
446 summ_ampl_C += amp;
447 mean_time_C += smeared_time;
448 }
449 }
450
451 for (size_t i = 0; i < sum_ampl_ipmt.size(); i++) {
452 sum_ampl_ipmt[i] = sum_ampl_ipmt[i] >> 3;
453 if (i < 4 * mGeometry.NCellsA) {
454 sum_A_ampl += sum_ampl_ipmt[i];
455 } else {
456 sum_C_ampl += sum_ampl_ipmt[i];
457 }
458 }
459
460 // Sum over PMs (using per-PM map) for debug/monitoring
461 int sum_PM_ampl_debug = 0;
462 int sum_PM_ampl_A_debug = 0;
463 int sum_PM_ampl_C_debug = 0;
464 for (const auto& entry : mapPMhash2sumAmpl) {
465 int pmAmpl = (entry.second >> 3);
466 sum_PM_ampl_debug += pmAmpl;
467 auto itSide = mMapPMhash2isAside.find(entry.first);
468 if (itSide != mMapPMhash2isAside.end()) {
469 if (itSide->second) {
470 sum_PM_ampl_A_debug += pmAmpl;
471 } else {
472 sum_PM_ampl_C_debug += pmAmpl;
473 }
474 }
475 }
476 LOG(debug) << "Sum PM amplitude (LUT-based): total=" << sum_PM_ampl_debug
477 << " A-side=" << sum_PM_ampl_A_debug
478 << " C-side=" << sum_PM_ampl_C_debug;
479
480 Bool_t is_A, is_C, isVertex, is_Central, is_SemiCentral = 0;
481 is_A = n_hit_A > 0;
482 is_C = n_hit_C > 0;
483 is_Central = sum_PM_ampl_A_debug + sum_PM_ampl_C_debug >= 2 * params.mtrg_central_trh;
484 is_SemiCentral = sum_PM_ampl_A_debug + sum_PM_ampl_C_debug >= 2 * params.mtrg_semicentral_trh && !is_Central;
485 uint32_t amplA = is_A ? summ_ampl_A * 0.125 : -5000; // sum amplitude A side / 8 (hardware)
486 uint32_t amplC = is_C ? summ_ampl_C * 0.125 : -5000; // sum amplitude C side / 8 (hardware)
487 int timeA = is_A ? mean_time_A / n_hit_A : -5000; // average time A side
488 int timeC = is_C ? mean_time_C / n_hit_C : -5000; // average time C side
489 vertex_time = (timeC - timeA) * 0.5;
490 isVertex = is_A && is_C && (vertex_time > -params.mTime_trg_vertex_gate && vertex_time < params.mTime_trg_vertex_gate);
491 LOG(debug) << " A " << is_A << " timeA " << timeA << " mean_time_A " << mean_time_A << " n_hit_A " << n_hit_A << " C " << is_C << " timeC " << timeC << " mean_time_C " << mean_time_C << " n_hit_C " << n_hit_C << " vertex_time " << vertex_time;
492 Triggers triggers;
493 bool isLaser = false;
494 bool isOutputsAreBlocked = false;
495 bool isDataValid = true;
496 if (nStored > 0) {
497 triggers.setTriggers(is_A, is_C, isVertex, is_Central, is_SemiCentral, int8_t(n_hit_A), int8_t(n_hit_C),
498 amplA, amplC, timeA, timeC, isLaser, isOutputsAreBlocked, isDataValid);
499 digitsBC.emplace_back(first, nStored, firstBCinDeque, triggers, mEventID - 1);
500 digitsTrig.emplace_back(firstBCinDeque, is_A, is_C, isVertex, is_Central, is_SemiCentral);
501 size_t const nBC = digitsBC.size();
502 for (auto const& lbl : bc.labels) {
503 labels.addElement(nBC - 1, lbl);
504 }
505 }
506 // Debug output -------------------------------------------------------------
507
508 LOG(debug) << "Event ID: " << mEventID << ", bc " << firstBCinDeque.bc << ", N hit " << bc.hits.size();
509 LOG(debug) << "N hit A: " << int(triggers.getNChanA()) << " N hit C: " << int(triggers.getNChanC()) << " summ ampl A: " << int(triggers.getAmplA())
510 << " summ ampl C: " << int(triggers.getAmplC()) << " mean time A: " << triggers.getTimeA()
511 << " mean time C: " << triggers.getTimeC() << " nStored " << nStored;
512
513 LOG(debug) << "IS A " << triggers.getOrA() << " IsC " << triggers.getOrC() << " vertex " << triggers.getVertex() << " is Central " << triggers.getCen() << " is SemiCentral " << triggers.getSCen();
514}
515
516//------------------------------------------------------------------------
517void Digitizer::flush(std::vector<o2::ft0::Digit>& digitsBC,
518 std::vector<o2::ft0::ChannelData>& digitsCh,
519 std::vector<o2::ft0::DetTrigInput>& digitsTrig,
521{
522
523 assert(firstBCinDeque <= mIntRecord);
524
525 while (firstBCinDeque < mIntRecord && !mCache.empty()) {
526 storeBC(mCache.front(), digitsBC, digitsCh, digitsTrig, labels);
527 mCache.pop_front();
528 ++firstBCinDeque;
529 }
530 firstBCinDeque = mIntRecord;
531}
532
533void Digitizer::flush_all(std::vector<o2::ft0::Digit>& digitsBC,
534 std::vector<o2::ft0::ChannelData>& digitsCh,
535 std::vector<o2::ft0::DetTrigInput>& digitsTrig,
537{
538
539 assert(firstBCinDeque <= mIntRecord);
540 ++mEventID;
541 while (!mCache.empty()) {
542 storeBC(mCache.front(), digitsBC, digitsCh, digitsTrig, labels);
543 mCache.pop_front();
544 ++firstBCinDeque;
545 }
546}
547
549{
550 auto const sinc = [](double x) { x *= TMath::Pi(); return (std::abs(x) < 1e-12) ? 1.0 : std::sin(x) / x; };
551
552 // number of noise samples in one BC
553 const auto& params = FT0DigParam::Instance();
554 mNumNoiseSamples = std::ceil(params.mBunchWidth / params.mNoisePeriod);
555 mNoiseSamples.resize(mNumNoiseSamples);
556
557 // set up tables with sinc function values (times noiseVar)
558 for (size_t i = 0, n = mSincTable.size(); i < n; ++i) {
559 float const time = i / float(n) * params.mNoisePeriod; // [0 .. 1/params.mNoisePeriod)
560 LOG(debug) << "initParameters " << i << "/" << n << " " << time;
561 // we make a table of sinc values between -num_noise_samples and 2*num_noise_samples
562 mSincTable[i].resize(3 * mNumNoiseSamples);
563 for (int j = -mNumNoiseSamples; j < 2 * mNumNoiseSamples; ++j) {
564 mSincTable[i][mNumNoiseSamples + j] = params.mNoiseVar * sinc((time + 0.5f * params.mBunchWidth) / params.mNoisePeriod - j);
565 }
566 }
567 // set up the lookup table for the signal form
568 for (size_t i = 0; i < DP::SIGNAL_TABLE_SIZE; ++i) {
569 float const x = float(i) / float(DP::SIGNAL_TABLE_SIZE) * params.mBunchWidth;
570 mSignalTable[i] = signalForm_i(x);
571 }
572
573 // cache for signal time series used by the CFD -BC/2 .. +3BC/2
574 mSignalCache.resize(std::lround(params.mBunchWidth / DP::SIGNAL_CACHE_DT));
575}
576//_______________________________________________________________________
578{
579 LOG(info) << " @@@ Digitizer::init " << std::endl;
580 mDeadTimes.fill({InteractionRecord(0), -100.});
582}
583//_______________________________________________________________________
585{
587}
588
589//_______________________________________________________________________
591{
592 const auto& params = FT0DigParam::Instance();
593 LOG(info) << " Run Digitzation with parametrs: \n"
594 << " CFD amplitude threshold \n " << params.mCFD_trsh << " CFD signal gate in ps \n"
595 << params.mTime_trg_gate << "shift to have signal around zero after CFD trancformation \n"
596 << params.mCfdShift << "CFD distance between 0.3 of max amplitude to max \n"
597 << params.mCFDShiftPos << "MIP -> mV " << params.mMip_in_V << " Pe in MIP \n"
598 << params.mPe_in_mip << "noise level " << params.mNoiseVar << " noise frequency \n"
599 << params.mNoisePeriod << " mMCPs " << params.mMCPs;
600}
601
std::vector< std::string > labels
std::vector< unsigned long > times
#define O2ParamImpl(classname)
uint64_t exp(uint64_t base, uint8_t exp) noexcept
General constants in FT0.
std::ostringstream debug
ClassImp(Digitizer)
int64_t timeC
int16_t charge
Definition RawEventData.h:5
int64_t amplA
int64_t timeA
int64_t amplC
uint64_t bc
Definition RawEventData.h:5
int16_t time
Definition RawEventData.h:4
Configurable digitization parameters.
int32_t i
GPUChain * chain
constexpr int p2()
constexpr int p1()
constexpr to accelerate the coordinates changing
Definition of a container to keep Monte Carlo truth external to simulation objects.
Header to collect physics constants.
uint32_t j
Definition RawData.h:0
A container to hold and manage MC truth information/labels.
static SingleLUT & Instance(const Table_t *table=nullptr, long timestamp=-1)
uint8_t getNChanC() const
Definition Triggers.h:97
int32_t getAmplC() const
Definition Triggers.h:99
bool getCen() const
Definition Triggers.h:85
int32_t getAmplA() const
Definition Triggers.h:98
bool getVertex() const
Definition Triggers.h:87
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
int16_t getTimeC() const
Definition Triggers.h:101
uint8_t getNChanA() const
Definition Triggers.h:96
bool getOrA() const
Definition Triggers.h:80
bool getOrC() const
Definition Triggers.h:81
int16_t getTimeA() const
Definition Triggers.h:100
bool getSCen() const
Definition Triggers.h:83
void flush(std::vector< o2::ft0::Digit > &digitsBC, std::vector< o2::ft0::ChannelData > &digitsCh, std::vector< o2::ft0::DetTrigInput > &digitsTrig, o2::dataformats::MCTruthContainer< o2::ft0::MCLabel > &label)
void printParameters() const
CFDOutput get_time(const std::vector< float > &times, float deadTime)
VcType signalFormVc(VcType x) const
Definition Digitizer.h:121
void process(const std::vector< o2::ft0::HitType > *hits, std::vector< o2::ft0::Digit > &digitsBC, std::vector< o2::ft0::ChannelData > &digitsCh, std::vector< o2::ft0::DetTrigInput > &digitsTrig, o2::dataformats::MCTruthContainer< o2::ft0::MCLabel > &label)
void flush_all(std::vector< o2::ft0::Digit > &digitsBC, std::vector< o2::ft0::ChannelData > &digitsCh, std::vector< o2::ft0::DetTrigInput > &digitsTrig, o2::dataformats::MCTruthContainer< o2::ft0::MCLabel > &label)
double measure_amplitude(const std::vector< float > &times) const
float signalForm(float x) const
Definition Digitizer.h:106
static constexpr int NCellsA
Definition Geometry.h:52
static constexpr int NCellsC
Definition Geometry.h:53
GLdouble n
Definition glcorearb.h:1982
GLint GLenum GLint x
Definition glcorearb.h:403
const GLfloat * m
Definition glcorearb.h:4066
GLuint64EXT * result
Definition glcorearb.h:5662
GLuint entry
Definition glcorearb.h:5735
GLuint GLfloat GLfloat GLfloat x1
Definition glcorearb.h:5034
GLuint index
Definition glcorearb.h:781
GLboolean GLboolean GLboolean b
Definition glcorearb.h:1233
GLenum const GLfloat * params
Definition glcorearb.h:272
GLuint GLsizei const GLchar * label
Definition glcorearb.h:2519
GLuint GLfloat * val
Definition glcorearb.h:1582
GLboolean GLboolean GLboolean GLboolean a
Definition glcorearb.h:1233
GLint GLuint mask
Definition glcorearb.h:291
uint8_t itsSharedClusterMap uint8_t
constexpr float LightSpeedCm2NS
Float signalForm_i(Float x)
Definition Digitizer.cxx:43
Vc::float_v signalForm_integralVc(Vc::float_v x)
float signalForm_integral(float x)
Definition Digitizer.cxx:73
uint32_t orbit
LHC orbit.
uint16_t bc
bunch crossing ID of interaction
static double bc2ns(int bc, unsigned int orbit)
const bool isChannelAlive(const uint8_t &chId) const
static constexpr std::size_t sNCHANNELS_PM
Definition Constants.h:33
static constexpr double SIGNAL_CACHE_DT
std::array< int16_t, o2::ft0::Geometry::Nchannels > mTimeOffsets
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"