27#include <fairlogger/Logger.h>
39 LOG(info) <<
"Initializing digitizer";
41 mChips.resize(mNumberOfChips);
42 for (
int i = mNumberOfChips;
i--;) {
43 mChips[
i].setChipIndex(
i);
45 mChips[
i].setNoiseMap(mNoiseMap);
49 mChips[
i].setDeadChanMap(mDeadChanMap);
54 mChipSimResp = mParams.getResponse();
55 mChipSimRespVD = mChipSimResp;
56 mChipSimRespMLOT = mChipSimResp;
60 LOG(info) <<
" Depth max VD: " << mChipSimRespVD->getDepthMax();
61 LOG(info) <<
" Depth min VD: " << mChipSimRespVD->getDepthMin();
63 LOG(info) <<
" Depth max MLOT: " << mChipSimRespMLOT->getDepthMax();
64 LOG(info) <<
" Depth min MLOT: " << mChipSimRespMLOT->getDepthMin();
66 float thicknessVD = 0.0095;
69 LOG(info) <<
"Using response name: " << mRespName;
70 mSimRespOrientation =
false;
72 if (mRespName ==
"APTS") {
75 mSimRespVDShift = mChipSimRespVD->getDepthMax();
78 mSimRespOrientation =
true;
79 }
else if (mRespName ==
"ALICE3") {
82 mSimRespVDShift = mChipSimRespVD->getDepthMax();
86 LOG(fatal) <<
"Unknown response name: " << mRespName;
89 mSimRespMLOTShift = mChipSimRespMLOT->getDepthMax() - thicknessMLOT / 2.f;
94 LOGP(info,
"TRK Digitizer is initialised.");
96 LOGP(info,
"VD shift = {} ; ML/OT shift = {} = {} - {}", mSimRespVDShift, mSimRespMLOTShift, mChipSimRespMLOT->getDepthMax(), thicknessMLOT / 2.f);
97 LOGP(info,
"VD pixel scale on x = {} ; z = {}", mSimRespVDScaleX, mSimRespVDScaleZ);
98 LOGP(info,
"ML/OT pixel scale on x = {} ; z = {}", mSimRespMLOTScaleX, mSimRespMLOTScaleZ);
99 LOGP(info,
"Response orientation: {}", mSimRespOrientation ?
"flipped" :
"normal");
106 if (mGeometry->getSubDetID(chipID) == 0) {
107 return mChipSimRespVD;
110 else if (mGeometry->getSubDetID(chipID) == 1) {
111 return mChipSimRespMLOT;
121 LOG(info) <<
" Digitizing " << mGeometry->
getName() <<
" (ID: " << mGeometry->
getDetID()
122 <<
") hits of event " << evID <<
" from source " << srcID
123 <<
" 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";
147 return mGeometry->getLayerTRK((*hits)[idx].GetDetectorID()) ==
layer;
149 processHit((*hits)[
i], mROFrameMax, evID, srcID,
layer);
156 LOG(info) <<
"Setting event time to " << irt.
getTimeNS() <<
" ns after orbit 0 bc 0";
163 if (mCollisionTimeWrtROF < 0 && nbc > 0) {
169 mIsBeforeFirstRO =
true;
172 mIsBeforeFirstRO =
false;
180 if (mNewROFrame < mROFrameMin) {
181 LOG(error) <<
"New ROFrame " << mNewROFrame <<
" (" << irt <<
") precedes currently cashed " << mROFrameMin;
182 throw std::runtime_error(
"deduced ROFrame precedes already processed one");
185 if (mROFrameMax < mNewROFrame) {
186 mROFrameMax = mNewROFrame - 1;
194 if (frameLast > mROFrameMax) {
195 frameLast = mROFrameMax;
198 getExtraDigBuffer(mROFrameMax);
199 LOG(info) <<
"Filling " << mGeometry->
getName() <<
" digits output for RO frames " << mROFrameMin <<
":"
205 for (; mROFrameMin <= frameLast; mROFrameMin++) {
209 auto& extra = *(mExtraBuff.front().get());
210 for (
auto& chip : mChips) {
211 if (chip.isDisabled() || (
layer >= 0 && mGeometry->getLayerTRK(chip.getChipIndex()) !=
layer)) {
214 chip.addNoise(mROFrameMin, mROFrameMin, &mParams, mGeometry->getSubDetID(chip.getChipIndex()), mGeometry->
getLayer(chip.getChipIndex()));
215 auto&
buffer = chip.getPreDigits();
219 auto itBeg =
buffer.begin();
221 ULong64_t maxKey = chip.getOrderingKey(mROFrameMin + 1, 0, 0) - 1;
222 for (; iter !=
buffer.end(); ++iter) {
223 if (iter->first > maxKey) {
226 auto& preDig = iter->second;
228 int digID = mDigits->size();
229 mDigits->emplace_back(chip.getChipIndex(), preDig.row, preDig.col, preDig.charge);
230 LOG(
debug) <<
"Adding digit ID: " << digID <<
" with chipID: " << chip.getChipIndex() <<
", row: " << preDig.row <<
", col: " << preDig.col <<
", charge: " << preDig.charge;
231 mMCLabels->
addElement(digID, preDig.labelRef.label);
232 auto& nextRef = preDig.labelRef;
233 while (nextRef.next >= 0) {
234 nextRef = extra[nextRef.next];
239 buffer.erase(itBeg, iter);
245 mROFRecords->push_back(rcROF);
249 mExtraBuff.emplace_back(mExtraBuff.front().release());
250 mExtraBuff.pop_front();
255void Digitizer::processHit(
const o2::trk::Hit& hit, uint32_t& maxFr,
int evID,
int srcID,
int rofLayer)
258 int subDetID = mGeometry->getSubDetID(chipID);
261 int disk = mGeometry->getDisk(chipID);
264 LOG(
debug) <<
"Skipping disk " << disk;
268 LOG(
debug) <<
"Processing hit for chip " << chipID;
269 auto& chip = mChips[chipID];
270 if (chip.isDisabled()) {
271 LOG(
debug) <<
"Skipping disabled chip " << chipID;
274 float timeInROF = hit.
GetTime() * sec2ns;
275 LOG(
debug) <<
"Hit time: " << timeInROF <<
" ns";
276 if (timeInROF > 20e3) {
277 const int maxWarn = 10;
278 static int warnNo = 0;
279 if (warnNo < maxWarn) {
280 LOG(warning) <<
"Ignoring hit with time_in_event = " << timeInROF <<
" ns"
281 << ((++warnNo < maxWarn) ?
"" :
" (suppressing further warnings)");
285 timeInROF += mCollisionTimeWrtROF;
286 if (mIsBeforeFirstRO && timeInROF < 0) {
288 LOG(
debug) <<
"Ignoring hit with timeInROF = " << timeInROF;
301 int nFrames = roFrameRelMax + 1 - roFrameRel;
302 uint32_t roFrameMax = mNewROFrame + roFrameRelMax;
303 if (roFrameMax > maxFr) {
322 LOG(
debug) <<
"Called curved to flat: " << xyzLocS.x() <<
" -> " << xyFlatS.x() <<
", " << xyzLocS.y() <<
" -> " << xyFlatS.y();
324 xyzLocS.SetXYZ(xyFlatS.x(), xyFlatS.y(), xyzLocS.Z());
325 xyzLocE.SetXYZ(xyFlatE.x(), xyFlatE.y(), xyzLocE.Z());
343 LOG(
debug) <<
"Step into the sensitive volume: " <<
step <<
". Number of steps: " << nSteps;
344 int rowS = -1, colS = -1, rowE = -1, colE = -1, nSkip = 0;
349 if (++nSkip >= nSteps) {
350 LOG(
debug) <<
"Did not enter to sensitive matrix, " << nSkip <<
" >= " << nSteps;
358 if (++nSkip >= nSteps) {
359 LOG(
debug) <<
"Did not enter to sensitive matrix, " << nSkip <<
" >= " << nSteps;
370 std::swap(rowS, rowE);
373 std::swap(colS, colE);
391 int rowSpan = rowE - rowS + 1, colSpan = colE - colS + 1;
393 float respMatrix[rowSpan][colSpan];
394 std::fill(&respMatrix[0][0], &respMatrix[0][0] + rowSpan * colSpan, 0.f);
397 nElectrons *= nStepsInv;
402 int rowPrev = -1, colPrev = -1,
row,
col;
403 float cRowPix = 0.f, cColPix = 0.f;
414 xyzLocS.SetY(xyzLocS.Y() + ((subDetID == 0) ? mSimRespVDShift : mSimRespMLOTShift));
418 for (
int iStep = nSteps; iStep--;) {
421 if (
row != rowPrev ||
col != colPrev) {
428 bool flipCol =
false, flipRow =
false;
430 float rowMax{}, colMax{};
433 rowMax = 0.5f * Segmentation::PitchRowVD * mSimRespVDScaleX;
434 colMax = 0.5f * Segmentation::PitchColVD * mSimRespVDScaleZ;
435 rspmat = resp->
getResponse(mSimRespVDScaleX * (xyzLocS.X() - cRowPix), mSimRespVDScaleZ * (xyzLocS.Z() - cColPix), xyzLocS.Y(), flipRow, flipCol, rowMax, colMax);
437 rowMax = 0.5f * Segmentation::PitchRowMLOT * mSimRespMLOTScaleX;
438 colMax = 0.5f * Segmentation::PitchColMLOT * mSimRespMLOTScaleZ;
439 rspmat = resp->
getResponse(mSimRespMLOTScaleX * (xyzLocS.X() - cRowPix), mSimRespMLOTScaleZ * (xyzLocS.Z() - cColPix), xyzLocS.Y(), flipRow, flipCol, rowMax, colMax);
444 if (rspmat ==
nullptr) {
445 LOG(
debug) <<
"Error in rspmat for step " << iStep <<
" / " << nSteps;
454 if (rowDest < 0 || rowDest >= rowSpan) {
459 if (colDest < 0 || colDest >= colSpan) {
462 respMatrix[rowDest][colDest] += rspmat->getValue(irow, icol, mSimRespOrientation ? !flipRow : flipRow, flipCol);
469 auto roFrameAbs = mNewROFrame + roFrameRel;
470 LOG(
debug) <<
"\nSpanning through rows and columns; rowspan = " << rowSpan <<
" colspan = " << colSpan <<
" = " << colE <<
" - " << colS <<
" +1 ";
471 for (
int irow = rowSpan; irow--;) {
472 uint16_t rowIS = irow + rowS;
473 for (
int icol = colSpan; icol--;) {
474 float nEleResp = respMatrix[irow][icol];
475 if (nEleResp <= 1.e-36) {
478 LOG(
debug) <<
"nEleResp: value " << nEleResp <<
" for pixel " << irow <<
" " << icol;
479 int nEle = gRandom->Poisson(nElectrons * nEleResp);
480 LOG(
debug) <<
"Charge detected in the pixel: " << nEle <<
" for pixel " << irow <<
" " << icol;
483 LOG(
debug) <<
"Ignoring pixel with nEle = " << nEle <<
" < min charge to account "
488 uint16_t colIS = icol + colS;
489 if (mNoiseMap && mNoiseMap->
isNoisy(chipID, rowIS, colIS)) {
492 if (mDeadChanMap && mDeadChanMap->
isNoisy(chipID, rowIS, colIS)) {
495 registerDigits(chip, roFrameAbs, timeInROF, nFrames, rowIS, colIS, nEle, lbl, rofLayer);
507 LOG(
debug) <<
"Registering digits for chip " << chip.
getChipIndex() <<
" at ROFrame " << roFrame
508 <<
" row " <<
row <<
" col " <<
col <<
" nEle " << nEle <<
" label " << lbl;
510 for (
int i = 0;
i < nROF;
i++) {
511 uint32_t roFr = roFrame +
i;
519 if (roFr > mEventROFrameMax) {
520 mEventROFrameMax = roFr;
522 if (roFr < mEventROFrameMin) {
523 mEventROFrameMin = roFr;
529 LOG(
debug) <<
"Added digit with key: " <<
key <<
" ROF: " << roFr <<
" row: " <<
row <<
" col: " <<
col <<
" charge: " << nEleROF;
531 LOG(
debug) <<
"Added to pre-digit with key: " <<
key <<
" ROF: " << roFr <<
" row: " <<
row <<
" col: " <<
col <<
" charge: " << nEleROF;
536 ExtraDig* extra = getExtraDigBuffer(roFr);
540 if ((*extra)[nxt].
label == lbl) {
544 nxt = (*extra)[nxt].next;
551 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
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)
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
GLint GLint GLint GLint GLint GLint GLint GLbitfield GLenum filter
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"