26#include <fairlogger/Logger.h>
38 LOG(info) <<
"Initializing digitizer";
40 mChips.resize(mNumberOfChips);
41 for (
int i = mNumberOfChips;
i--;) {
42 mChips[
i].setChipIndex(
i);
44 mChips[
i].setNoiseMap(mNoiseMap);
48 mChips[
i].setDeadChanMap(mDeadChanMap);
53 mChipSimResp = mParams.getResponse();
54 mChipSimRespVD = mChipSimResp;
55 mChipSimRespMLOT = mChipSimResp;
59 LOG(info) <<
" Depth max VD: " << mChipSimRespVD->getDepthMax();
60 LOG(info) <<
" Depth min VD: " << mChipSimRespVD->getDepthMin();
62 LOG(info) <<
" Depth max MLOT: " << mChipSimRespMLOT->getDepthMax();
63 LOG(info) <<
" Depth min MLOT: " << mChipSimRespMLOT->getDepthMin();
65 float thicknessVD = 0.0095;
68 LOG(info) <<
"Using response name: " << mRespName;
69 mSimRespOrientation =
false;
71 if (mRespName ==
"APTS") {
74 mSimRespVDShift = mChipSimRespVD->getDepthMax();
77 mSimRespOrientation =
true;
78 }
else if (mRespName ==
"ALICE3") {
81 mSimRespVDShift = mChipSimRespVD->getDepthMax();
85 LOG(fatal) <<
"Unknown response name: " << mRespName;
88 mSimRespMLOTShift = mChipSimRespMLOT->getDepthMax() - thicknessMLOT / 2.f;
93 LOGP(info,
"TRK Digitizer is initialised.");
95 LOGP(info,
"VD shift = {} ; ML/OT shift = {} = {} - {}", mSimRespVDShift, mSimRespMLOTShift, mChipSimRespMLOT->getDepthMax(), thicknessMLOT / 2.f);
96 LOGP(info,
"VD pixel scale on x = {} ; z = {}", mSimRespVDScaleX, mSimRespVDScaleZ);
97 LOGP(info,
"ML/OT pixel scale on x = {} ; z = {}", mSimRespMLOTScaleX, mSimRespMLOTScaleZ);
98 LOGP(info,
"Response orientation: {}", mSimRespOrientation ?
"flipped" :
"normal");
105 if (mGeometry->getSubDetID(chipID) == 0) {
106 return mChipSimRespVD;
109 else if (mGeometry->getSubDetID(chipID) == 1) {
110 return mChipSimRespMLOT;
120 LOG(info) <<
" Digitizing " << mGeometry->
getName() <<
" (ID: " << mGeometry->
getDetID()
121 <<
") hits of event " << evID <<
" from source " << srcID
122 <<
" at time " << mEventTime.
getTimeNS() <<
" ROFrame = " << mNewROFrame
124 <<
" Min/Max ROFrames " << mROFrameMin <<
"/" << mROFrameMax;
126 std::cout <<
"Printing segmentation info: " << std::endl;
130 if (mNewROFrame > mROFrameMin) {
134 int nHits = hits->size();
135 std::vector<int> hitIdx(nHits);
136 std::iota(std::begin(hitIdx), std::end(hitIdx), 0);
138 std::sort(hitIdx.begin(), hitIdx.end(),
139 [hits](
auto lhs,
auto rhs) {
140 return (*hits)[lhs].GetDetectorID() < (*hits)[rhs].GetDetectorID();
142 LOG(info) <<
"Processing " << nHits <<
" hits";
143 for (
int i : hitIdx) {
144 processHit((*hits)[
i], mROFrameMax, evID, srcID);
158 LOG(info) <<
"Setting event time to " << irt.
getTimeNS() <<
" ns after orbit 0 bc 0";
170 if (mCollisionTimeWrtROF < 0 && nbc > 0) {
176 mIsBeforeFirstRO =
true;
179 mIsBeforeFirstRO =
false;
182 LOG(
debug) <<
" NewROFrame " << mNewROFrame <<
" = " << nbc <<
"/" << mParams.
getROFrameLengthInBC() <<
" (nbc/mParams.getROFrameLengthInBC()";
188 mIsBeforeFirstRO =
false;
191 if (mNewROFrame < mROFrameMin) {
192 LOG(error) <<
"New ROFrame " << mNewROFrame <<
" (" << irt <<
") precedes currently cashed " << mROFrameMin;
193 throw std::runtime_error(
"deduced ROFrame precedes already processed one");
196 if (mParams.
isContinuous() && mROFrameMax < mNewROFrame) {
197 mROFrameMax = mNewROFrame - 1;
205 if (frameLast > mROFrameMax) {
206 frameLast = mROFrameMax;
209 getExtraDigBuffer(mROFrameMax);
210 LOG(info) <<
"Filling " << mGeometry->
getName() <<
" digits output for RO frames " << mROFrameMin <<
":"
216 for (; mROFrameMin <= frameLast; mROFrameMin++) {
220 auto& extra = *(mExtraBuff.front().get());
221 for (
auto& chip : mChips) {
222 if (chip.isDisabled()) {
225 chip.addNoise(mROFrameMin, mROFrameMin, &mParams, mGeometry->getSubDetID(chip.getChipIndex()), mGeometry->
getLayer(chip.getChipIndex()));
226 auto&
buffer = chip.getPreDigits();
230 auto itBeg =
buffer.begin();
232 ULong64_t maxKey = chip.getOrderingKey(mROFrameMin + 1, 0, 0) - 1;
233 for (; iter !=
buffer.end(); ++iter) {
234 if (iter->first > maxKey) {
237 auto& preDig = iter->second;
239 int digID = mDigits->size();
240 mDigits->emplace_back(chip.getChipIndex(), preDig.row, preDig.col, preDig.charge);
241 LOG(
debug) <<
"Adding digit ID: " << digID <<
" with chipID: " << chip.getChipIndex() <<
", row: " << preDig.row <<
", col: " << preDig.col <<
", charge: " << preDig.charge;
242 mMCLabels->
addElement(digID, preDig.labelRef.label);
243 auto& nextRef = preDig.labelRef;
244 while (nextRef.next >= 0) {
245 nextRef = extra[nextRef.next];
250 buffer.erase(itBeg, iter);
260 mROFRecords->push_back(rcROF);
264 mExtraBuff.emplace_back(mExtraBuff.front().release());
265 mExtraBuff.pop_front();
270void Digitizer::processHit(
const o2::trk::Hit& hit, uint32_t& maxFr,
int evID,
int srcID)
273 int subDetID = mGeometry->getSubDetID(chipID);
276 int disk = mGeometry->getDisk(chipID);
279 LOG(
debug) <<
"Skipping disk " << disk;
283 LOG(
debug) <<
"Processing hit for chip " << chipID;
284 auto& chip = mChips[chipID];
285 if (chip.isDisabled()) {
286 LOG(
debug) <<
"Skipping disabled chip " << chipID;
289 float timeInROF = hit.
GetTime() * sec2ns;
290 LOG(
debug) <<
"Hit time: " << timeInROF <<
" ns";
291 if (timeInROF > 20e3) {
292 const int maxWarn = 10;
293 static int warnNo = 0;
294 if (warnNo < maxWarn) {
295 LOG(warning) <<
"Ignoring hit with time_in_event = " << timeInROF <<
" ns"
296 << ((++warnNo < maxWarn) ?
"" :
" (suppressing further warnings)");
301 timeInROF += mCollisionTimeWrtROF;
303 if (mIsBeforeFirstRO && timeInROF < 0) {
305 LOG(
debug) <<
"Ignoring hit with timeInROF = " << timeInROF;
318 int nFrames = roFrameRelMax + 1 - roFrameRel;
319 uint32_t roFrameMax = mNewROFrame + roFrameRelMax;
320 if (roFrameMax > maxFr) {
339 LOG(
debug) <<
"Called curved to flat: " << xyzLocS.x() <<
" -> " << xyFlatS.x() <<
", " << xyzLocS.y() <<
" -> " << xyFlatS.y();
341 xyzLocS.SetXYZ(xyFlatS.x(), xyFlatS.y(), xyzLocS.Z());
342 xyzLocE.SetXYZ(xyFlatE.x(), xyFlatE.y(), xyzLocE.Z());
360 LOG(
debug) <<
"Step into the sensitive volume: " <<
step <<
". Number of steps: " << nSteps;
361 int rowS = -1, colS = -1, rowE = -1, colE = -1, nSkip = 0;
366 if (++nSkip >= nSteps) {
367 LOG(
debug) <<
"Did not enter to sensitive matrix, " << nSkip <<
" >= " << nSteps;
375 if (++nSkip >= nSteps) {
376 LOG(
debug) <<
"Did not enter to sensitive matrix, " << nSkip <<
" >= " << nSteps;
387 std::swap(rowS, rowE);
390 std::swap(colS, colE);
408 int rowSpan = rowE - rowS + 1, colSpan = colE - colS + 1;
410 float respMatrix[rowSpan][colSpan];
411 std::fill(&respMatrix[0][0], &respMatrix[0][0] + rowSpan * colSpan, 0.f);
414 nElectrons *= nStepsInv;
419 int rowPrev = -1, colPrev = -1,
row,
col;
420 float cRowPix = 0.f, cColPix = 0.f;
431 xyzLocS.SetY(xyzLocS.Y() + ((subDetID == 0) ? mSimRespVDShift : mSimRespMLOTShift));
435 for (
int iStep = nSteps; iStep--;) {
438 if (
row != rowPrev ||
col != colPrev) {
445 bool flipCol =
false, flipRow =
false;
447 float rowMax{}, colMax{};
450 rowMax = 0.5f * Segmentation::PitchRowVD * mSimRespVDScaleX;
451 colMax = 0.5f * Segmentation::PitchColVD * mSimRespVDScaleZ;
452 rspmat = resp->
getResponse(mSimRespVDScaleX * (xyzLocS.X() - cRowPix), mSimRespVDScaleZ * (xyzLocS.Z() - cColPix), xyzLocS.Y(), flipRow, flipCol, rowMax, colMax);
454 rowMax = 0.5f * Segmentation::PitchRowMLOT * mSimRespMLOTScaleX;
455 colMax = 0.5f * Segmentation::PitchColMLOT * mSimRespMLOTScaleZ;
456 rspmat = resp->
getResponse(mSimRespMLOTScaleX * (xyzLocS.X() - cRowPix), mSimRespMLOTScaleZ * (xyzLocS.Z() - cColPix), xyzLocS.Y(), flipRow, flipCol, rowMax, colMax);
461 if (rspmat ==
nullptr) {
462 LOG(
debug) <<
"Error in rspmat for step " << iStep <<
" / " << nSteps;
471 if (rowDest < 0 || rowDest >= rowSpan) {
476 if (colDest < 0 || colDest >= colSpan) {
479 respMatrix[rowDest][colDest] += rspmat->getValue(irow, icol, mSimRespOrientation ? !flipRow : flipRow, flipCol);
486 auto roFrameAbs = mNewROFrame + roFrameRel;
487 LOG(
debug) <<
"\nSpanning through rows and columns; rowspan = " << rowSpan <<
" colspan = " << colSpan <<
" = " << colE <<
" - " << colS <<
" +1 ";
488 for (
int irow = rowSpan; irow--;) {
489 uint16_t rowIS = irow + rowS;
490 for (
int icol = colSpan; icol--;) {
491 float nEleResp = respMatrix[irow][icol];
492 if (nEleResp <= 1.e-36) {
495 LOG(
debug) <<
"nEleResp: value " << nEleResp <<
" for pixel " << irow <<
" " << icol;
496 int nEle = gRandom->Poisson(nElectrons * nEleResp);
497 LOG(
debug) <<
"Charge detected in the pixel: " << nEle <<
" for pixel " << irow <<
" " << icol;
500 LOG(
debug) <<
"Ignoring pixel with nEle = " << nEle <<
" < min charge to account "
505 uint16_t colIS = icol + colS;
506 if (mNoiseMap && mNoiseMap->
isNoisy(chipID, rowIS, colIS)) {
509 if (mDeadChanMap && mDeadChanMap->
isNoisy(chipID, rowIS, colIS)) {
512 registerDigits(chip, roFrameAbs, timeInROF, nFrames, rowIS, colIS, nEle, lbl);
524 LOG(
debug) <<
"Registering digits for chip " << chip.
getChipIndex() <<
" at ROFrame " << roFrame
525 <<
" row " <<
row <<
" col " <<
col <<
" nEle " << nEle <<
" label " << lbl;
527 for (
int i = 0;
i < nROF;
i++) {
528 uint32_t roFr = roFrame +
i;
536 if (roFr > mEventROFrameMax) {
537 mEventROFrameMax = roFr;
539 if (roFr < mEventROFrameMin) {
540 mEventROFrameMin = roFr;
546 LOG(
debug) <<
"Added digit with key: " <<
key <<
" ROF: " << roFr <<
" row: " <<
row <<
" col: " <<
col <<
" charge: " << nEleROF;
548 LOG(
debug) <<
"Added to pre-digit with key: " <<
key <<
" ROF: " << roFr <<
" row: " <<
row <<
" col: " <<
col <<
" charge: " << nEleROF;
553 ExtraDig* extra = getExtraDigBuffer(roFr);
557 if ((*extra)[nxt].
label == lbl) {
561 nxt = (*extra)[nxt].next;
568 extra->emplace_back(lbl);
Definition of the SegmentationChipclass.
Definition of the TRK digitizer.
math_utils::Point3D< T > GetPos() const
unsigned short GetDetectorID() const
static const DPLDigitizerParam< N > & Instance()
const char * getName() const
const o2::detectors::DetID & getDetID() const
const Mat3D & getMatrixL2G(int sensID) const
static int constexpr NPix
float getMaxDuration() const
float getCollectedCharge(float totalNEle, float tMin, float tMax) const
bool getResponse(float vRow, float vCol, float cDepth, AlpideRespSimMat &dest) const
void addDigit(ULong64_t key, UInt_t roframe, UShort_t row, UShort_t col, int charge, o2::MCCompLabel lbl)
o2::itsmft::PreDigit * findDigit(ULong64_t key)
UShort_t getChipIndex() const
int getMinChargeToAccount() const
float getStrobeDelay(int layer=-1) const
virtual void print() const
float getROFrameLengthInv(int layer=-1) const
const SignalShape & getSignalShape() const
float getStrobeLength(int layer=-1) const
float getEnergyToNElectrons() const
int getROFrameLengthInBC(int layer=-1) const
bool isContinuous() const
int getChargeThreshold() const
float getNSimStepsInv() const
float getROFrameLength(int layer=-1) const
void fillOutputContainer(uint32_t maxFrame=0xffffffff, int layer=-1)
auto getChipResponse(int chipID)
bool isContinuous() const
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.
virtual Int_t getLayer(Int_t index) const
Int_t getNumberOfChips() const
math_utils::Point3D< Float_t > GetPosStart() const
bool isFullChipMasked(int chip) const
bool isNoisy(int chip, int row, int col) const
const BCData & getBCData() const
void setFirstEntry(int idx)
int getFirstEntry() const
void setROFrame(ROFtype rof)
static bool localToDetector(float x, float z, int &iRow, int &iCol)
static bool detectorToLocal(L row, L col, T &xRow, T &zCol)
static ULong64_t getOrderingKey(UInt_t roframe, UShort_t row, UShort_t col)
Get global ordering key made of readout frame, column and row.
static constexpr float PitchColVD
static constexpr float PitchColMLOT
static constexpr float PitchRowMLOT
static void Print() noexcept
Print segmentation info.
static constexpr float PitchRowVD
static constexpr float SiliconThicknessMLOT
GLuint GLsizei const GLchar * label
GLenum GLuint GLint GLint layer
constexpr double LHCBunchSpacingNS
constexpr std::array< int, nLayers > nRows
int64_t differenceInBC(const InteractionRecord &other) const
void setFromLong(int64_t l)
double timeInBCNS
time in NANOSECONDS relative to orbit/bc
double getTimeNS() const
get time in ns from orbit=0/bc=0
int next
eventual next contribution to the same pixel
o2::MCCompLabel label
hit label
PreDigitLabelRef labelRef
label and reference to the next one
IR getFirstSampledTFIR() const
get TF and HB (abs) for this IR
LOG(info)<< "Compressed in "<< sw.CpuTime()<< " s"