CbmRoot
Loading...
Searching...
No Matches
CbmMvdSensorHitfinderTask.cxx
Go to the documentation of this file.
1/* Copyright (C) 2013-2025 Institut fuer Kernphysik, Goethe-Universitaet Frankfurt, Frankfurt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Philipp Sitzmann [committer], Samir Amar-Youcef, Florian Uhlig */
4
5// ---------------------------------------------------------------------------------------------
6// ----- CbmMvdSensorHitfinderTask source file -----
7// ----- Created 11/09/13 by P.Sitzmann -----
8// ----- -----
9// ---------------------------------------------------------------------------------------------
11
12#include "CbmMvdCluster.h" // for CbmMvdCluster
13#include "CbmMvdDigi.h" // for CbmMvdDigi
14#include "CbmMvdHit.h" // for CbmMvdHit
15#include "CbmMvdSensor.h" // for CbmMvdSensor
16
17#include <TClonesArray.h> // for TClonesArray
18#include <TGeoMatrix.h> // for TGeoHMatrix
19#include <TMath.h> // for Power
20#include <TMathBase.h> // for Abs
21#include <TObjArray.h> // for TObjArray
22#include <TRandom.h> // for TRandom
23#include <TRandom3.h> // for gRandom
24#include <TVector3.h> // for TVector3
25
26#include <iomanip> // for setprecision, setw
27#include <iostream> // for operator<<, basic_ostream, char_traits, endl
28#include <map> // for map, __map_iterator, operator!=, operator==
29#include <vector> // for vector
30
31using std::endl;
32using std::fixed;
33using std::ios_base;
34using std::left;
35using std::map;
36using std::pair;
37using std::right;
38using std::setprecision;
39using std::setw;
40using std::vector;
41
42
43// ----- Default constructor ------------------------------------------
46 , fAdcDynamic(200)
47 , fAdcOffset(0)
48 , fAdcBits(1)
49 , fAdcSteps(-1)
50 , fAdcStepSize(-1.)
51 , fClusters(new TClonesArray("CbmMvdCluster", 10000))
52 , fPixelChargeHistos(nullptr)
54 , fResolutionHistoX(nullptr)
55 , fResolutionHistoY(nullptr)
56 , fResolutionHistoCleanX(nullptr)
57 , fResolutionHistoCleanY(nullptr)
60 , fBadHitHisto(nullptr)
61 , fGausArray(nullptr)
62 , fGausArrayIt(-1)
63 , fGausArrayLimit(5000)
64 , fDigiMap()
65 , fDigiMapIt()
66 , h(nullptr)
67 , h3(nullptr)
68 , h1(nullptr)
69 , h2(nullptr)
70 , Qseed(nullptr)
71 , fFullClusterHisto(nullptr)
72 , c1(nullptr)
73 , fNEvent(0)
74 , fMode(0)
75 , fCounter(0)
76 , fSigmaNoise(15.)
77 , fSeedThreshold(1.)
78 , fNeighThreshold(1.)
79 , fShowDebugHistos(kFALSE)
80 , fUseMCInfo(kFALSE)
81 , inputSet(kFALSE)
82 , fLayerRadius(0.)
84 , fLayerPosZ(0.)
85 , fHitPosX(0.)
86 , fHitPosY(0.)
87 , fHitPosZ(0.)
88 , fHitPosErrX(0.0005)
89 , fHitPosErrY(0.0005)
90 , fHitPosErrZ(0.0)
91 , fBranchName("MvdHit")
93 , fAddNoise(kFALSE)
94{
95 fPluginIDNumber = 300;
96}
97// -------------------------------------------------------------------------
98
99
100// ----- Standard constructor ------------------------------------------
103 , fAdcDynamic(200)
104 , fAdcOffset(0)
105 , fAdcBits(1)
106 , fAdcSteps(-1)
107 , fAdcStepSize(-1.)
108 , fClusters(new TClonesArray("CbmMvdCluster", 10000))
109 , fPixelChargeHistos(nullptr)
111 , fResolutionHistoX(nullptr)
112 , fResolutionHistoY(nullptr)
113 , fResolutionHistoCleanX(nullptr)
114 , fResolutionHistoCleanY(nullptr)
115 , fResolutionHistoMergedX(nullptr)
116 , fResolutionHistoMergedY(nullptr)
117 , fBadHitHisto(nullptr)
118 , fGausArray(nullptr)
119 , fGausArrayIt(-1)
120 , fGausArrayLimit(5000)
121 , fDigiMap()
122 , fDigiMapIt()
123 , h(nullptr)
124 , h3(nullptr)
125 , h1(nullptr)
126 , h2(nullptr)
127 , Qseed(nullptr)
128 , fFullClusterHisto(nullptr)
129 , c1(nullptr)
130 , fNEvent(0)
131 , fMode(iMode)
132 , fCounter(0)
133 , fSigmaNoise(15.)
134 , fSeedThreshold(1.)
135 , fNeighThreshold(1.)
136 , fShowDebugHistos(kFALSE)
137 , fUseMCInfo(kFALSE)
138 , inputSet(kFALSE)
139 , fLayerRadius(0.)
141 , fLayerPosZ(0.)
142 , fHitPosX(0.)
143 , fHitPosY(0.)
144 , fHitPosZ(0.)
145 , fHitPosErrX(0.0005)
146 , fHitPosErrY(0.0005)
147 , fHitPosErrZ(0.0)
148 , fBranchName("MvdHit")
149 , fDigisInCluster(0)
150 , fAddNoise(kFALSE)
151{
152 fPluginIDNumber = 300;
153}
154// -------------------------------------------------------------------------
155
156// ----- Destructor ----------------------------------------------------
158{
159
160 if (fOutputBuffer) {
161 fOutputBuffer->Delete();
162 delete fOutputBuffer;
163 }
164
165 if (fClusters) {
166 fClusters->Delete();
167 delete fClusters;
168 }
169
170 if (fInputBuffer) {
171 fInputBuffer->Delete();
172 delete fInputBuffer;
173 }
174}
175// -------------------------------------------------------------------------
176
177// ----- Virtual private method Init ----------------------------------
179{
180
181
182 fSensor = mysensor;
183 LOG(debug) << "CbmMvdSensorHitfinderTask: Initialisation of sensor " << fSensor->GetName();
184
185 fInputBuffer = new TClonesArray("CbmMvdCluster", 100);
186 fOutputBuffer = new TClonesArray("CbmMvdHit", 100);
187
188
189 //Add charge collection histograms
190 fPixelChargeHistos = new TObjArray();
191
192 fTotalChargeInNpixelsArray = new TObjArray();
193
194 fAdcSteps = (Int_t) TMath::Power(2, fAdcBits);
196
198 for (Int_t i = 0; i < fGausArrayLimit; i++) {
199 fGausArray[i] = gRandom->Gaus(0, fSigmaNoise);
200 };
201 fGausArrayIt = 0;
202
203
204 initialized = kTRUE;
205
206 //LOG(debug) << "-Finished- " << GetName() << ": Initialisation of sensor " << fSensor->GetName();
207}
208// -------------------------------------------------------------------------
209
210// ----- Virtual public method Reinit ----------------------------------
212{
213 LOG(info) << "CbmMvdSensorHitfinderTask::ReInt---------------";
214 return kSUCCESS;
215}
216// -------------------------------------------------------------------------
217
218// ----- Virtual public method ExecChain --------------
220// -------------------------------------------------------------------------
221
222// ----- Virtual public method Exec --------------
224{
225
226 // if(!inputSet)
227 // {
228 // fInputBuffer->Clear();
229 // fInputBuffer->AbsorbObjects(fPreviousPlugin->GetOutputArray());
230 // LOG(debug) << "absorbt object from previous plugin at " << fSensor->GetName() << " got " << fInputBuffer->GetEntriesFast() << " entrys";
231 // }
232 if (fInputBuffer->GetEntriesFast() > 0) {
233
234 fOutputBuffer->Clear();
235 inputSet = kFALSE;
236 for (Int_t i = 0; i < fInputBuffer->GetEntriesFast(); i++) {
237 CbmMvdCluster* cluster = (CbmMvdCluster*) fInputBuffer->At(i);
238 TVector3 pos(0, 0, 0);
239 TVector3 dpos(0, 0, 0);
240 CreateHit(cluster, pos, dpos);
241 }
242 }
243 fInputBuffer->Delete();
244}
245//--------------------------------------------------------------------------
246
247
248//--------------------------------------------------------------------------
249void CbmMvdSensorHitfinderTask::CreateHit(CbmMvdCluster* cluster, TVector3& pos, TVector3& dpos)
250{
251 // Calculate the center of gravity of the charge of a cluster
252 double startTime = 0;
253 double endTime = 0;
254 double indicatedTime = 0;
255
256 ComputeCenterOfGravity(cluster, pos, dpos, indicatedTime, startTime, endTime);
257
258
259 Int_t indexX, indexY;
260 Double_t local[3] = {pos.X(), pos.Y(), pos.Z()};
261
262 //LOG(debug)<< "found center of gravity at: " << local[0] << " , " << local[1];
263
264 fSensor->TopToPixel(local, indexX, indexY);
265
266 //LOG(debug) << "Center is on pixel: " << indexX << " , " << indexY;
267
268
269 // Save hit into array
270
271 //LOG(debug) << "adding new hit to fHits at X: " << pos.X() << " , Y: "<< pos.Y() << " , Z: " << pos.Z() << " , from cluster nr " << cluster->GetRefId();
272 Int_t nHits = fOutputBuffer->GetEntriesFast();
273
274 // FIXME: SZh 03.09.2025:
275 // Store address directly in fSensor, when the addressing map or properly named geo-nodes are available
276 uint32_t address = CbmMvdAddress::GetAddressFromSensorNrAndPixelXY(fSensor->GetSensorNr(), indexY, indexX);
277 new ((*fOutputBuffer)[nHits])
278 CbmMvdHit(address, fSensor->GetStationNr(), pos, dpos, indexX, indexY, cluster->GetRefId(),
279 0); // Bug - Should be the sensor ID
280 CbmMvdHit* currentHit = (CbmMvdHit*) fOutputBuffer->At(nHits);
281
282 currentHit->SetTime(indicatedTime);
283 currentHit->SetTimeError((endTime - startTime) / 3.64); //Assume Sqrt(12) - Error
284 currentHit->SetValidityStartTime(startTime);
285 currentHit->SetValidityEndTime(endTime);
286 currentHit->SetRefId(cluster->GetRefId());
287}
288
289//--------------------------------------------------------------------------
290
291/*void CbmMvdSensorHitfinderTask::UpdateDebugHistos(vector<Int_t>* clusterArray, Int_t seedIndexX, Int_t seedIndexY)
292{
293;
294} */
295
296
297//--------------------------------------------------------------------------
298
300 double& indicatedTime, double& startTime, double& endTime)
301{
302
303 Double_t numeratorX = 0;
304 Double_t numeratorY = 0;
305 Double_t denominator = 0;
306 Int_t charge;
309 Double_t x, y;
310 Double_t layerPosZ = fSensor->GetZ();
311 Double_t lab[3] = {0, 0, 0};
312 std::map<pair<Int_t, Int_t>, Int_t> PixelMap = cluster->GetPixelMap();
313 Int_t clusterSize = cluster->GetNofDigis();
314
315 // UInt_t shape = 0;
316 // Int_t xIndex0 = 0;
317 // Int_t yIndex0 = 0;
318
319 Double_t sigmaIn[3] = {0., 0., 0.}, sigmaOut[3] = {0., 0., 0.}, shiftIn[3] = {0., 0., 0.}, shiftOut[3] = {0., 0., 0.};
320
321 //Double_t GetFrameStartTime(Int_t frameNumber);
322 //Double_t GetFrameEndTime(Int_t frameNumber) {return GetFrameStartTime (frameNumber+1);}
323
324
325 //Scan digis for the digi with minimal frame number.
326 //For the time being, it is anticipated that this digi provides the best time estimate as it has lowest signal delay
327
328 // CbmMvdDigi* myDigi=(CbmMvdDigi*)cluster->GetDigi(0);
329 // int32_t minFrameNumber=myDigi->GetFrameNumber();
330 //
331 // for(Int_t i=0;i<clusterSize; i++){
332 // myDigi=(CbmMvdDigi*)cluster->GetDigi(i);
333 // int32_t frame=myDigi->GetFrameNumber();
334 // if(frame<minFrameNumber){minFrameNumber=frame;}
335 // }
336
337 //Estimate the highest possible delay and jitter of the analog amplifier chain from the sensor data sheet.
338
339 CbmMvdSensorDataSheet* sensorDataSheet = fSensor->GetDataSheet();
340 int32_t minFrameNumber = cluster->GetEarliestFrameNumber();
341
342
343 //Compute time walk and jitter assuming that the maximum possible signal is reasonably well represented with 10k signal electrons
344 Double_t timeWalk =
345 sensorDataSheet->GetDelay(sensorDataSheet->GetAnalogThreshold()) - sensorDataSheet->GetDelay(10000);
346 timeWalk =
347 timeWalk
348 + 2 * sensorDataSheet->GetDelaySigma(sensorDataSheet->GetAnalogThreshold()); //include pixel-to-pixel spread.
349 Double_t jitter = sensorDataSheet->GetJitter(sensorDataSheet->GetAnalogThreshold()); //include jitter.
350
351 // fixme: pixel-to-pixel spread not accounted for
352
353 startTime = fSensor->GetFrameStartTime(minFrameNumber) - timeWalk - 2 * jitter;
354 endTime = fSensor->GetFrameEndTime(minFrameNumber) + 2 * jitter;
355
356 indicatedTime = startTime + (endTime - startTime) / 2;
357
358
359 for (map<pair<Int_t, Int_t>, Int_t>::iterator it = PixelMap.begin(); it != PixelMap.end(); ++it) {
360 pair<Int_t, Int_t> pixel = it->first;
361
362 charge = GetAdcCharge(it->second);
363 xIndex = pixel.first;
364 yIndex = pixel.second;
365
366 // Determine Cluster Shape
367 //if (PixelMap.size() <= 4) {
368 // if (it == PixelMap.begin()) {
369 // xIndex0 = xIndex;
370 // yIndex0 = yIndex;
371 // }
372 // shape += TMath::Power(2, (4 * (yIndex - yIndex0 + 3)) + (xIndex - xIndex0));
373 //}
374
375 LOG(debug) << "CbmMvdSensorHitfinderTask:: iCluster= " << cluster->GetRefId() << " , clusterSize= " << clusterSize;
376 LOG(debug) << "CbmMvdSensorHitfinderTask::xIndex " << xIndex << " , yIndex " << yIndex << " , charge = " << charge;
377
378 fSensor->PixelToTop(xIndex, yIndex, lab);
379
380 x = lab[0];
381 y = lab[1];
382
383 Double_t xc = x * charge;
384 Double_t yc = y * charge;
385
386 numeratorX += xc;
387 numeratorY += yc;
388 denominator += charge;
389 }
390
391 LOG(debug) << "CbmMvdSensorHitfinderTask::=========================";
392 LOG(debug) << "CbmMvdSensorHitfinderTask::numeratorX: " << numeratorX << " , numeratorY: " << numeratorY
393 << ", denominator: " << denominator << ", indicated time: " << indicatedTime;
394
395 //Calculate x,y coordinates of the pixel in the laboratory ref frame
396 if (denominator != 0) {
399 fHitPosZ = layerPosZ;
400 }
401 else {
402 fHitPosX = 0;
403 fHitPosY = 0;
404 fHitPosZ = 0;
405 }
406
407 LOG(debug) << "CbmMvdSensorHitfinderTask::-----------------------------------";
408 LOG(debug) << "CbmMvdSensorHitfinderTask::X hit= " << fHitPosX << " Y hit= " << fHitPosY << " Z hit= " << fHitPosZ;
409 LOG(debug) << "CbmMvdSensorHitfinderTask::-----------------------------------";
410
411 // Quick fix to get some errors for the MvdHits.
412 // Proposed value 8mum
413 // See Redmine issue 3066
414 // See plots with residuals from GitLab MR 1504
415 sigmaIn[0] = 0.0008;
416 sigmaIn[1] = 0.0008;
417 sigmaIn[2] = 0.;
418
419 // // Treat Sigma/Shift of the Cluster according to the Shape
420 // // This goes back to some results of a specific MIMOSA-26 beam time. Switch off for the time being.
421 //
422 // if (shape == 12288) {
423 // sigmaIn[0] = 0.00053;
424 // sigmaIn[1] = 0.00063;
425 // sigmaIn[2] = 0.;
426 // shiftIn[0] = -0.00000;
427 // shiftIn[1] = -0.00001;
428 // shiftIn[2] = 0.;
429 // }
430 // else if (shape == 208896) {
431 // sigmaIn[0] = 0.00035;
432 // sigmaIn[1] = 0.00036;
433 // sigmaIn[2] = 0.;
434 // shiftIn[0] = -0.00000;
435 // shiftIn[1] = -0.00002;
436 // shiftIn[2] = 0.;
437 // }
438 // else if (shape == 69632) {
439 // sigmaIn[0] = 0.00028;
440 // sigmaIn[1] = 0.00028;
441 // sigmaIn[2] = 0.;
442 // shiftIn[0] = -0.00000;
443 // shiftIn[1] = -0.00002;
444 // shiftIn[2] = 0.;
445 // }
446 // else if (shape == 28672) {
447 // sigmaIn[0] = 0.00028;
448 // sigmaIn[1] = 0.00039;
449 // sigmaIn[2] = 0.;
450 // shiftIn[0] = -0.00000;
451 // shiftIn[1] = -0.00001;
452 // shiftIn[2] = 0.;
453 // }
454 // else if (shape == 143360) {
455 // sigmaIn[0] = 0.00024;
456 // sigmaIn[1] = 0.00022;
457 // sigmaIn[2] = 0.;
458 // shiftIn[0] = +0.00020;
459 // shiftIn[1] = +0.00008;
460 // shiftIn[2] = 0.;
461 // }
462 // else if (shape == 200704) {
463 // sigmaIn[0] = 0.00024;
464 // sigmaIn[1] = 0.00022;
465 // sigmaIn[2] = 0.;
466 // shiftIn[0] = -0.00020;
467 // shiftIn[1] = -0.00011;
468 // shiftIn[2] = 0.;
469 // }
470 // else if (shape == 77824) {
471 // sigmaIn[0] = 0.00024;
472 // sigmaIn[1] = 0.00022;
473 // sigmaIn[2] = 0.;
474 // shiftIn[0] = -0.00020;
475 // shiftIn[1] = +0.00008;
476 // shiftIn[2] = 0.;
477 // }
478 // else if (shape == 12800) {
479 // sigmaIn[0] = 0.00024;
480 // sigmaIn[1] = 0.00022;
481 // sigmaIn[2] = 0.;
482 // shiftIn[0] = +0.00020;
483 // shiftIn[1] = -0.00011;
484 // shiftIn[2] = 0.;
485 // }
486 // else if (shape == 4096) {
487 // sigmaIn[0] = 0.00027;
488 // sigmaIn[1] = 0.00092;
489 // sigmaIn[2] = 0.;
490 // shiftIn[0] = +0.00002;
491 // shiftIn[1] = +0.00004;
492 // shiftIn[2] = 0.;
493 // }
494 // else {
495 // sigmaIn[0] = 0.00036;
496 // sigmaIn[1] = 0.00044;
497 // sigmaIn[2] = 0.;
498 // shiftIn[0] = -0.00000;
499 // shiftIn[1] = -0.00002;
500 // shiftIn[2] = 0.;
501 // }
502 // Consider Sensor Orientation
503 TGeoHMatrix* RecoMatrix = fSensor->GetRecoMatrix();
504 TGeoHMatrix RotMatrix;
505 RotMatrix.SetRotation(RecoMatrix->GetRotationMatrix());
506
507 RotMatrix.LocalToMaster(sigmaIn, sigmaOut);
508 RotMatrix.LocalToMaster(shiftIn, shiftOut);
509
510 fHitPosErrX = TMath::Abs(sigmaOut[0]);
511 fHitPosErrY = TMath::Abs(sigmaOut[1]);
512 fHitPosErrZ = TMath::Abs(sigmaOut[2]);
513
514 fHitPosX += shiftOut[0];
515 fHitPosY += shiftOut[1];
516 fHitPosZ += shiftOut[2];
517
518 // pos = center of gravity (labframe), dpos uncertainty LOG(info)<<setw(10)<<setprecision(2)<< VolumeShape->GetDX();
519 pos.SetXYZ(fHitPosX, fHitPosY, fHitPosZ);
521}
522
523//--------------------------------------------------------------------------
524
525
526//--------------------------------------------------------------------------
528//--------------------------------------------------------------------------
529
531{
532
533 Int_t adcCharge;
534
535 if (charge < fAdcOffset) {
536 return 0;
537 };
538
539 adcCharge = int((charge - fAdcOffset) / fAdcStepSize);
540 if (adcCharge > fAdcSteps - 1) {
541 adcCharge = fAdcSteps - 1;
542 }
543
544 return adcCharge;
545}
546
547
548// ----- Private method Reset ------------------------------------------
550
551// -------------------------------------------------------------------------
552
553
Double_t denominator
TVector3 dpos
Double_t numeratorY
Double_t lab[3]
Double_t numeratorX
ClassImp(CbmMvdSensorHitfinderTask)
float Float_t
int Int_t
int32_t GetNofDigis() const
Number of digis in cluster.
Definition CbmCluster.h:69
void SetTimeError(double error)
Definition CbmHit.h:94
void SetRefId(int32_t refId)
Definition CbmHit.h:85
void SetTime(double time)
Definition CbmHit.h:88
static uint32_t GetAddressFromSensorNrAndPixelXY(int32_t sensorNr, int32_t sensorY, int32_t sensorX)
int32_t GetRefId() const
std::map< std::pair< int32_t, int32_t >, int32_t > GetPixelMap() const
int32_t GetEarliestFrameNumber() const
void SetValidityStartTime(double time)
Definition CbmMvdHit.h:82
void SetValidityEndTime(double time)
Definition CbmMvdHit.h:83
virtual Double_t GetDelaySigma(Float_t charge)
virtual Double_t GetJitter(Float_t charge)
virtual Double_t GetDelay(Float_t charge)
void InitTask(CbmMvdSensor *mySensor)
std::map< std::pair< Int_t, Int_t >, Int_t >::iterator fDigiMapIt
void ComputeCenterOfGravity(CbmMvdCluster *clusterArray, TVector3 &pos, TVector3 &dpos, double &indicatedTime, double &startTime, double &endTime)
std::map< std::pair< Int_t, Int_t >, Int_t > fDigiMap
void CreateHit(CbmMvdCluster *clusterArray, TVector3 &pos, TVector3 &dpos)
CbmMvdSensor * fSensor
TClonesArray * fOutputBuffer
TClonesArray * fInputBuffer