CbmRoot
CbmMvdSensorClusterfinderTask.cxx
Go to the documentation of this file.
1 // -------------------------------------------------------------------------
2 // ----- CbmMvdSensorClusterfinderTask source file -----
3 // ----- Created 03.12.2014 by P. Sitzmann -----
4 // -------------------------------------------------------------------------
5 
7 #include "FairLogger.h"
8 #include "TClonesArray.h"
9 #include "TObjArray.h"
10 
11 using std::cout;
12 using std::endl;
13 using std::pair;
14 using std::vector;
15 
16 // ----- Default constructor -------------------------------------------
19 // -------------------------------------------------------------------------
20 
21 // -------------------------------------------------------------------------
23 // -------------------------------------------------------------------------
24 
25 // -------------------------------------------------------------------------
27  Int_t iVerbose)
29  , fAdcDynamic(200)
30  , fAdcOffset(0)
31  , fAdcBits(1)
32  , fAdcSteps(-1)
33  , fAddress(0)
34  , fAdcStepSize(-1.)
35  , fDigis(NULL)
36  , fPixelChargeHistos(NULL)
37  , fTotalChargeInNpixelsArray(NULL)
38  , fResolutionHistoX(NULL)
39  , fResolutionHistoY(NULL)
40  , fResolutionHistoCleanX(NULL)
41  , fResolutionHistoCleanY(NULL)
42  , fResolutionHistoMergedX(NULL)
43  , fResolutionHistoMergedY(NULL)
44  , fBadHitHisto(NULL)
45  , fGausArray(NULL)
46  , fGausArrayIt(-1)
47  , fGausArrayLimit(5000)
48  , fDigiMap()
49  , fDigiMapIt()
50  , h(NULL)
51  , h3(NULL)
52  , h1(NULL)
53  , h2(NULL)
54  , Qseed(NULL)
55  , fFullClusterHisto(NULL)
56  , c1(NULL)
57  , fNEvent(0)
58  , fMode(iMode)
59  , fCounter(0)
60  , fVerbose(iVerbose)
61  , fSigmaNoise(15.)
62  , fSeedThreshold(1.)
63  , fNeighThreshold(1.)
64  , fUseMCInfo(kFALSE)
65  , inputSet(kFALSE)
66  , ftempPixelMap()
67  , fLayerRadius(0.)
68  , fLayerRadiusInner(0.)
69  , fLayerPosZ(0.)
70  , fHitPosX(0.)
71  , fHitPosY(0.)
72  , fHitPosZ(0.)
73  , fHitPosErrX(0.0005)
74  , fHitPosErrY(0.0005)
75  , fHitPosErrZ(0.0)
76  , fBranchName("MvdHit")
77  , fAddNoise(kFALSE) {}
78 // -------------------------------------------------------------------------
79 
80 // ----- Virtual private method Init ----------------------------------
82 
83 
84  fSensor = mysensor;
85  //cout << "-Start- " << GetName() << ": Initialisation of sensor " << fSensor->GetName() << endl;
86  fInputBuffer = new TClonesArray("CbmMvdDigi", 10000);
87  fOutputBuffer = new TClonesArray("CbmMvdCluster", 10000);
88 
89 
90  //Add charge collection histograms
91  fPixelChargeHistos = new TObjArray();
92 
93  fTotalChargeInNpixelsArray = new TObjArray();
94 
95  fAdcSteps = (Int_t) TMath::Power(2, fAdcBits);
97 
99 
100  initialized = kTRUE;
101 
102  if (fShowDebugHistos) {
103  fResolutionHistoX = new TH1F("SinglePointResolution_X",
104  "SinglePointResolution_X",
105  10000,
106  -0.0100,
107  0.0100);
108  fResolutionHistoY = new TH1F("SinglePointResolution_Y",
109  "SinglePointResolution_Y",
110  10000,
111  -0.0100,
112  0.0100);
113  fResolutionHistoCleanX = new TH1F("SinglePointResolution_X_Clean",
114  "SinglePointResolution_X_Clean",
115  10000,
116  -0.0100,
117  0.0100);
118  fResolutionHistoCleanY = new TH1F("SinglePointResolution_Y_Clean",
119  "SinglePointResolution_Y_Clean",
120  10000,
121  -0.0100,
122  0.0100);
123  fResolutionHistoMergedX = new TH1F("SinglePointResolution_X_Merged",
124  "SinglePointResolution_X_Merged",
125  10000,
126  -0.0100,
127  0.0100);
128  fResolutionHistoMergedY = new TH1F("SinglePointResolution_Y_Merged",
129  "SinglePointResolution_Y_Merged",
130  10000,
131  -0.0100,
132  0.0100);
133  fBadHitHisto = new TH2F(
134  "BadHits", "Hits above 0.003cm", 1000, -2.5, 2.5, 1000, -2.5, 2.5);
136  new TH1F("ChargeOfAllPixels", "ChargeOfAllPixels", 12000, 0, 12000);
137  //}
138 
139  TH1F* histo;
140  TH1F* histoTotalCharge;
141  char* histoName = new char[20];
142  char* histoTotalChargeName = new char[50];
143 
144  //Add charge collection histograms
145  fPixelChargeHistos = new TObjArray();
146  for (Int_t i = 0; i < 49; i++) {
147  sprintf(histoName, "ChargePixel%i", i + 1);
148  histo = new TH1F(histoName, histoName, 200, 0, 200);
149  fPixelChargeHistos->AddLast(histo);
150  };
151 
152  fTotalChargeInNpixelsArray = new TObjArray();
153  for (Int_t i = 0; i < 49; i++) {
154  sprintf(histoTotalChargeName, "totalChargeInNPixels%i", i + 1);
155  histoTotalCharge =
156  new TH1F(histoTotalChargeName, histoTotalChargeName, 12000, 0, 12000);
157  fTotalChargeInNpixelsArray->AddLast(histoTotalCharge);
158  };
159 
160  //Number 49
161  histo = new TH1F("ChargePixelSeed", "ChargePixelSeed", 200, 0, 14000);
162  fPixelChargeHistos->AddLast(histo);
163  //Number 50
164  histo = new TH1F("ChargePixel9of49", "ChargePixel 9 Of 49", 200, 0, 14000);
165  fPixelChargeHistos->AddLast(histo);
166  //Number 51
167  histo =
168  new TH1F("ChargePixel25of49", "ChargePixel 25 Of 49", 200, 0, 14000);
169  fPixelChargeHistos->AddLast(histo);
170  //Number 52
171  histo =
172  new TH1F("ChargePixel49of49", "ChargePixel 49 Of 49", 200, 0, 14000);
173  fPixelChargeHistos->AddLast(histo);
174  //Number 53
175  histo = new TH1F(
176  "ChargePixel9of49Sorted", "ChargePixel 9 Of 49 Sorted", 200, 0, 14000);
177  fPixelChargeHistos->AddLast(histo);
178  //Number 54
179  histo = new TH1F(
180  "ChargePixel25of49Sorted", "ChargePixel 25 Of 49 Sorted", 200, 0, 14000);
181  fPixelChargeHistos->AddLast(histo);
182  //Number 55
183  histo = new TH1F(
184  "ChargePixel49of49Sorted", "ChargePixel 49 Of 49 Sorted", 49, 0.5, 49.5);
185  fPixelChargeHistos->AddLast(histo);
186  }
187  //cout << "-Finished- " << GetName() << ": Initialisation of sensor " << fSensor->GetName() << endl;
188 }
189 // -------------------------------------------------------------------------
190 
191 // ----- Virtual public method Reinit ----------------------------------
193  cout << "-I- "
194  << "CbmMvdSensorClusterfinderTask::ReInt---------------" << endl;
195  return kTRUE;
196 }
197 // -------------------------------------------------------------------------
198 
199 // ----- Virtual public method ExecChain --------------
201 // -------------------------------------------------------------------------
202 
203 // ----- Public method Exec --------------
205  if (fInputBuffer->GetEntriesFast() > 0) {
206  fOutputBuffer->Delete();
207  inputSet = kFALSE;
208  vector<Int_t>* clusterArray = new vector<Int_t>;
209 
210  CbmMvdDigi* digi;
211 
212  Int_t iDigi = 0;
213  digi = (CbmMvdDigi*) fInputBuffer->At(iDigi);
214 
215 
216  if (!digi) {
217  cout << "-E- : CbmMvdSensorFindHitTask - Fatal: No Digits found in this "
218  "event."
219  << endl;
220  }
221 
222  Int_t nDigis = fInputBuffer->GetEntriesFast();
223 
224 
225  nDigis = fInputBuffer->GetEntriesFast();
226  TArrayS* pixelUsed = new TArrayS(nDigis);
227 
228  for (iDigi = 0; iDigi < nDigis; iDigi++) {
229  pixelUsed->AddAt(0, iDigi);
230  }
231 
232  fDigiMap.clear();
233  Int_t refId;
234  for (Int_t k = 0; k < nDigis; k++) {
235 
236  digi = (CbmMvdDigi*) fInputBuffer->At(k);
237  refId = digi->GetRefId();
238  if (refId < 0) {
239  LOG(fatal) << "RefID of this digi is -1 this should not happend ";
240  }
241  //apply fNeighThreshold
242 
243  if (GetAdcCharge(digi->GetCharge()) < fNeighThreshold) continue;
244 
245  pair<Int_t, Int_t> a(digi->GetPixelX(), digi->GetPixelY());
246  //cout << endl << "registerde pixel x:" << digi->GetPixelX() << " y:" << digi->GetPixelY() << endl;
247  fDigiMap[a] = k;
248  };
249 
250 
251  if (gDebug > 0) {
252  cout << "\n-I- " << GetName() << ": VolumeId " << fSensor->GetVolumeId()
253  << endl;
254  }
255 
256  for (iDigi = 0; iDigi < nDigis; iDigi++) {
257 
258  if (gDebug > 0 && iDigi % 10000 == 0) {
259  cout << "-I- " << GetName() << " Digi:" << iDigi << endl;
260  };
261 
262  digi = (CbmMvdDigi*) fInputBuffer->At(iDigi);
263  //cout << endl << "working with pixel x:" << digi->GetPixelX() << " y:" << digi->GetPixelY() << endl;
264 
265 
266  /*
267  ---------------------------------------------------------
268  check if digi is above threshold (define seed pixel)
269  then check for neighbours.
270  Once the cluster is created (seed and neighbours)
271  calculate the position of the hit
272  using center of gravity (CoG) method.
273  ---------------------------------------------------------
274  */
275 
276  if (gDebug > 0) {
277  cout << "-I- "
278  << "CbmMvdSensorFindHitTask: Checking for seed pixels..." << endl;
279  }
280 
281  if ((GetAdcCharge(digi->GetCharge()) >= fSeedThreshold)
282  && (pixelUsed->At(iDigi) == kFALSE)) {
283  clusterArray->clear();
284  clusterArray->push_back(iDigi);
285 
286  pixelUsed->AddAt(1, iDigi);
287 
288  pair<Int_t, Int_t> a(digi->GetPixelX(), digi->GetPixelY());
289  fDigiMapIt = fDigiMap.find(a);
290  fDigiMap.erase(fDigiMapIt);
291 
292  for (ULong64_t iCluster = 0; iCluster < clusterArray->size();
293  iCluster++) {
294 
295  if (gDebug > 0) {
296  cout << "-I- "
297  << " CbmMvdSensorClusterfinderTask: Calling method "
298  "CheckForNeighbours()..."
299  << endl;
300  }
301 
302  CheckForNeighbours(clusterArray, iCluster, pixelUsed);
303  //cout << endl << "checked for neighbours, create cluster" << endl;
304  }
305 
306  Int_t i = 0;
307  Int_t pixelCharge;
308 
309 
310  pair<Int_t, Int_t> pixelCoords;
311  Int_t clusterSize = clusterArray->size();
312  Int_t nClusters = fOutputBuffer->GetEntriesFast();
313  //cout << endl << "new cluster: " << nClusters << endl;
314  CbmMvdCluster* clusterNew =
315  new ((*fOutputBuffer)[nClusters]) CbmMvdCluster();
316  clusterNew->SetAddress(fAddress);
317 
318  for (i = 0; i < clusterSize; i++) {
319  CbmMvdDigi* digiInCluster =
320  (CbmMvdDigi*) fInputBuffer->At(clusterArray->at(i));
321  clusterNew->AddDigi(digiInCluster->GetRefId());
322  pixelCoords = std::make_pair(digiInCluster->GetPixelX(),
323  digiInCluster->GetPixelY());
324  pixelCharge = digiInCluster->GetCharge();
325  ftempPixelMap[pixelCoords] = pixelCharge;
326  }
327  clusterNew->SetPixelMap(ftempPixelMap);
328  ftempPixelMap.clear();
329 
330  if (fShowDebugHistos) UpdateDebugHistos(clusterNew);
331 
332  } // if AdcCharge>threshold
333  else { //cout << endl << "pixel is with " << digi->GetCharge() << " under Threshold or used" << endl;
334  }
335  } // loop on digis
336 
337 
338  delete pixelUsed;
339  clusterArray->clear();
340  delete clusterArray;
341  fInputBuffer->Delete();
342 
343  fDigiMap.clear();
344  } else { //cout << endl << "No input found." << endl;
345  }
346 }
347 // -------------------------------------------------------------------------
348 
349 // -------------------------------------------------------------------------
351  vector<Int_t>* clusterArray,
352  Int_t clusterDigi,
353  TArrayS* pixelUsed) {
354  CbmMvdDigi* seed =
355  (CbmMvdDigi*) fInputBuffer->At(clusterArray->at(clusterDigi));
356  //cout << endl << "pixel nr. " << clusterDigi << " is seed" << endl ;
357 
358 
359  // Remove Seed Pixel from list of non-used pixels
360  Int_t channelX = seed->GetPixelX();
361  Int_t channelY = seed->GetPixelY();
362  pair<Int_t, Int_t> a(channelX, channelY);
363 
364  // Find first neighbour
365 
366  a = std::make_pair(channelX + 1, channelY);
367  fDigiMapIt = fDigiMap.find(a);
368 
369  if (!(fDigiMapIt == fDigiMap.end())) {
370  Int_t i = fDigiMap[a];
371  //cout << endl << "pixel nr. " << i << " is used" << endl ;
372  // Only digis depassing fNeighThreshold are in the map, no cut required
373  clusterArray->push_back(i);
374 
375  pixelUsed->AddAt(1, i); // block pixel for the seed pixel scanner
376  fDigiMap.erase(fDigiMapIt); // block pixel for the neighbour pixel scanner
377  }
378 
379  a = std::make_pair(channelX - 1, channelY);
380  fDigiMapIt = fDigiMap.find(a);
381 
382  if (!(fDigiMapIt == fDigiMap.end())) {
383  Int_t i = fDigiMap[a];
384  //cout << endl << "pixel nr. " << i << " is used" << endl ;
385  // Only digits depassing fNeighThreshold are in the map, no cut required
386  clusterArray->push_back(i);
387  pixelUsed->AddAt(1, i); // block pixel for the seed pixel scanner
388  fDigiMap.erase(fDigiMapIt); // block pixel for the neighbour pixel scanner
389  }
390 
391  a = std::make_pair(channelX, channelY - 1);
392  fDigiMapIt = fDigiMap.find(a);
393  if (!(fDigiMapIt == fDigiMap.end())) {
394  Int_t i = fDigiMap[a];
395  // Only digits depassing fNeighThreshold are in the map, no cut required
396  //cout << endl << "pixel nr. " << i << " is used" << endl ;
397  clusterArray->push_back(i);
398  pixelUsed->AddAt(1, i); // block pixel for the seed pixel scanner
399  fDigiMap.erase(fDigiMapIt); // block pixel for the neighbour pixel scanner
400  }
401 
402  a = std::make_pair(channelX, channelY + 1);
403  fDigiMapIt = fDigiMap.find(a);
404 
405  if (!(fDigiMapIt == fDigiMap.end())) {
406  Int_t i = fDigiMap[a];
407  //cout << endl << "pixel nr. " << i << " is used" << endl ;
408  // Only digis depassing fNeighThreshold are in the map, no cut required
409  clusterArray->push_back(i);
410  pixelUsed->AddAt(1, i); // block pixel for the seed pixel scanner
411  fDigiMap.erase(fDigiMapIt); // block pixel for the neighbour pixel scanner
412  }
413 }
414 // -------------------------------------------------------------------------
415 
416 // -------------------------------------------------------------------------
418 
419  Int_t adcCharge;
420 
421  if (charge < fAdcOffset) { return 0; };
422 
423  adcCharge = int((charge - fAdcOffset) / fAdcStepSize);
424  if (adcCharge > fAdcSteps - 1) { adcCharge = fAdcSteps - 1; }
425 
426  return adcCharge;
427 }
428 // -------------------------------------------------------------------------
429 
431  /************************************************************
432  Algorithm for cluster shapes
433 
434  ************************************************************/
435  Float_t chargeArray3D[fChargeArraySize][fChargeArraySize];
436  Float_t chargeArray[fChargeArraySize * fChargeArraySize];
437  Short_t seedPixelOffset = fChargeArraySize / 2; // 3 for 7, 2 for 5
438  Int_t seedIndexX = 0, seedIndexY = 0;
439  Float_t seedCharge = 0.;
440  Float_t clusterCharge = cluster->GetClusterCharge();
441  std::map<std::pair<Int_t, Int_t>, Int_t> clusterMap = cluster->GetPixelMap();
442 
443  for (Int_t k = 0; k < fChargeArraySize; k++) {
444  for (Int_t j = 0; j < fChargeArraySize; j++) {
445  chargeArray3D[k][j] = gRandom->Gaus(0, fSigmaNoise);
446  }
447  }
448  for (std::map<std::pair<Int_t, Int_t>, Int_t>::iterator iter =
449  clusterMap.begin();
450  iter != clusterMap.end();
451  iter++) {
452  if (iter->second > seedCharge) {
453  seedCharge = iter->second;
454  seedIndexX = iter->first.first;
455  seedIndexY = iter->first.second;
456  }
457  }
458  //cout << endl << "seed pixel with " << seedCharge << " charge" << endl;
459  for (std::map<std::pair<Int_t, Int_t>, Int_t>::iterator iter =
460  clusterMap.begin();
461  iter != clusterMap.end();
462  iter++) {
463 
464  Int_t relativeX = iter->first.first + seedPixelOffset - seedIndexX;
465  Int_t relativeY = iter->first.second + seedPixelOffset - seedIndexY;
466 
467  if (fVerbose > 1)
468  cout << relativeX << " " << relativeY << " " << iter->first.first << " "
469  << seedIndexX << endl;
470 
471 
472  if (relativeX >= 0 && relativeX < fChargeArraySize && relativeY >= 0
473  && relativeY < fChargeArraySize) {
474  chargeArray3D[relativeX][relativeY] = iter->second;
475  }
476 
477  if ((relativeX - seedPixelOffset == 0)
478  && (relativeY - seedPixelOffset == 0)) { //seed digiArray
479  }
480  }
481 
482  if (fVerbose > 1) {
483  for (Int_t i = 0; i < fChargeArraySize; i++) {
484  for (Int_t j = 0; j < fChargeArraySize; j++) {
485  cout << chargeArray3D[i][j] << " ";
486  }
487  cout << endl;
488  }
489  }
490  fFullClusterHisto->Fill(clusterCharge);
491 
492  for (Int_t k = 0; k < fChargeArraySize; k++) {
493  for (Int_t j = 0; j < fChargeArraySize; j++) {
494  chargeArray[fChargeArraySize * k + j] = chargeArray3D[k][j];
495  }
496  }
497 
498  Int_t qSeed = chargeArray3D[seedPixelOffset][seedPixelOffset];
499  Int_t q9 = 0;
500 
501  for (Int_t k = seedPixelOffset - 1; k < seedPixelOffset + 1; k++) {
502  for (Int_t j = seedPixelOffset - 1; j < seedPixelOffset + 1; j++) {
503  q9 = q9 + chargeArray3D[k][j];
504  }
505  };
506 
507 
508  if (fChargeArraySize <= 7) {
509  for (Int_t i = 0; i < (fChargeArraySize * fChargeArraySize); i++) {
510  ((TH1F*) fPixelChargeHistos->At(i))->Fill(chargeArray[i]);
511  //cout << counter++<<" Charge: " << chargeArray[i]<< endl;
512  };
513  };
514 
515  //cout << "End of Cluster: "<<fChargeArraySize*fChargeArraySize << endl;
516 
517  Int_t q25 = 0;
518  Int_t q49 = 0;
519 
520 
521  for (Int_t k = seedPixelOffset - 2; k < seedPixelOffset + 2; k++) {
522  for (Int_t j = seedPixelOffset - 2; j < seedPixelOffset + 2; j++) {
523  q25 = q25 + chargeArray3D[k][j];
524  }
525  };
526 
527  if (fChargeArraySize >= 7) {
528  for (Int_t k = seedPixelOffset - 3; k < seedPixelOffset + 3; k++) {
529  for (Int_t j = seedPixelOffset - 3; j < seedPixelOffset + 3; j++) {
530  q49 = q49 + chargeArray3D[k][j];
531  }
532  }
533  }
534 
535  ((TH1F*) fPixelChargeHistos->At(49))->Fill(qSeed);
536  ((TH1F*) fPixelChargeHistos->At(50))->Fill(q9);
537  ((TH1F*) fPixelChargeHistos->At(51))->Fill(q25);
538  ((TH1F*) fPixelChargeHistos->At(52))->Fill(q49);
539 
540 
541  //Prepare selection of crowns for charge bow histograms
542 
543 
544  Int_t orderArray[fChargeArraySize * fChargeArraySize];
545 
546  TMath::Sort(
547  fChargeArraySize * fChargeArraySize, chargeArray, orderArray, kTRUE);
548 
549  Float_t qSort = 0;
550  for (Int_t i = 0; i < 9; i++) {
551  qSort += chargeArray[orderArray[i]];
552  };
553  ((TH1F*) fPixelChargeHistos->At(53))->Fill(qSort);
554 
555  for (Int_t i = 9; i < 25; i++) {
556  qSort += chargeArray[orderArray[i]];
557  };
558  ((TH1F*) fPixelChargeHistos->At(54))->Fill(qSort);
559 
560  TH1F* histoTotalCharge;
561  qSort = 0;
562  for (Int_t i = 0; i < fChargeArraySize * fChargeArraySize; i++) {
563  qSort += chargeArray[orderArray[i]];
564  ((TH1F*) fPixelChargeHistos->At(55))->Fill(i + 1, qSort);
565  histoTotalCharge = (TH1F*) fTotalChargeInNpixelsArray->At(i);
566  histoTotalCharge->Fill(qSort);
567  }
568 }
569 
570 //--------------------------------------------------------------------------
572 
573  if (fShowDebugHistos) {
574  cout << "\n============================================================"
575  << endl;
576  cout << "-I- " << GetName()
577  << "::Finish: Total events skipped: " << fCounter << endl;
578  cout << "============================================================"
579  << endl;
580  cout << "-I- Parameters used" << endl;
581  cout << "Gaussian noise [electrons] : " << fSigmaNoise << endl;
582  cout << "Noise simulated [Bool] : " << fAddNoise << endl;
583  cout << "Threshold seed [ADC] : " << fSeedThreshold << endl;
584  cout << "Threshold neighbours [ADC] : " << fNeighThreshold << endl;
585  cout << "ADC - Bits : " << fAdcBits << endl;
586  cout << "ADC - Dynamic [electrons] : " << fAdcDynamic << endl;
587  cout << "ADC - Offset [electrons] : " << fAdcOffset << endl;
588  cout << "============================================================"
589  << endl;
590 
591 
592  TH1F* histo;
593  TH2F* clusterShapeHistogram;
594  TCanvas* canvas2 = new TCanvas("HitFinderCharge", "HitFinderCharge");
595  canvas2->Divide(2, 2);
596  canvas2->cd(1);
597  if (fChargeArraySize <= 7) {
598  clusterShapeHistogram = new TH2F("MvdClusterShape",
599  "MvdClusterShape",
601  0,
604  0,
606  for (Int_t i = 0; i < fChargeArraySize * fChargeArraySize; i++) {
607  histo = (TH1F*) fPixelChargeHistos->At(i);
608  Float_t charge = histo->GetMean();
609  //cout <<i << " Charge " << charge << " xCluster: " << i%fChargeArraySize << " yCluster: " << i/fChargeArraySize << endl;
610  //histo->Fit("landau");
611  //TF1* fitFunction= histo->GetFunction("landau");
612  //Double_t MPV=fitFunction->GetParameter(1);
613  clusterShapeHistogram->Fill(
614  i % fChargeArraySize, i / fChargeArraySize, charge);
615  }
616  }
617  clusterShapeHistogram->Draw("Lego2");
618  canvas2->cd(2);
619  histo = (TH1F*) fPixelChargeHistos->At(50);
620  histo->Draw();
621  canvas2->cd(3);
622  histo = (TH1F*) fPixelChargeHistos->At(51);
623  histo->Draw();
624  canvas2->cd(4);
625  //fFullClusterHisto->Draw();
626  }
627 }
628 //--------------------------------------------------------------------------
CbmMvdSensorClusterfinderTask::ExecChain
void ExecChain()
Definition: CbmMvdSensorClusterfinderTask.cxx:200
CbmMvdSensorClusterfinderTask::fDigiMapIt
std::map< std::pair< Int_t, Int_t >, Int_t >::iterator fDigiMapIt
Definition: CbmMvdSensorClusterfinderTask.h:122
CbmMvdDigi::GetPixelY
Int_t GetPixelY()
Definition: CbmMvdDigi.cxx:143
CbmMvdSensorClusterfinderTask::fResolutionHistoMergedX
TH1F * fResolutionHistoMergedX
Definition: CbmMvdSensorClusterfinderTask.h:114
CbmMvdCluster::GetClusterCharge
Float_t GetClusterCharge()
Definition: CbmMvdCluster.h:55
CbmMvdDigi::GetCharge
Double_t GetCharge() const
Definition: CbmMvdDigi.h:49
CbmMvdSensorTask::fSensor
CbmMvdSensor * fSensor
Definition: CbmMvdSensorTask.h:59
CbmMvdSensorClusterfinderTask::fAddNoise
Bool_t fAddNoise
Definition: CbmMvdSensorClusterfinderTask.h:164
CbmMvdSensorClusterfinderTask::inputSet
Bool_t inputSet
Definition: CbmMvdSensorClusterfinderTask.h:143
CbmMvdSensorClusterfinderTask::fAdcSteps
Int_t fAdcSteps
Definition: CbmMvdSensorClusterfinderTask.h:101
CbmMvdSensorClusterfinderTask::fNeighThreshold
Double_t fNeighThreshold
Definition: CbmMvdSensorClusterfinderTask.h:141
CbmMvdSensorClusterfinderTask::InitTask
void InitTask(CbmMvdSensor *mySensor)
Definition: CbmMvdSensorClusterfinderTask.cxx:81
CbmMvdSensorClusterfinderTask.h
CbmMvdSensorClusterfinderTask::fAddress
Int_t fAddress
Definition: CbmMvdSensorClusterfinderTask.h:102
CbmCluster::AddDigi
void AddDigi(Int_t index)
Add digi to cluster.
Definition: CbmCluster.h:47
CbmMvdSensorClusterfinderTask::fResolutionHistoY
TH1F * fResolutionHistoY
Definition: CbmMvdSensorClusterfinderTask.h:111
CbmMvdSensorClusterfinderTask::fAdcOffset
Int_t fAdcOffset
Definition: CbmMvdSensorClusterfinderTask.h:99
CbmMvdSensorClusterfinderTask::fTotalChargeInNpixelsArray
TObjArray * fTotalChargeInNpixelsArray
Definition: CbmMvdSensorClusterfinderTask.h:107
CbmMvdSensorClusterfinderTask::fResolutionHistoMergedY
TH1F * fResolutionHistoMergedY
Definition: CbmMvdSensorClusterfinderTask.h:115
i
int i
Definition: L1/vectors/P4_F32vec4.h:25
CbmMvdCluster::GetPixelMap
std::map< std::pair< Int_t, Int_t >, Int_t > GetPixelMap()
Definition: CbmMvdCluster.h:47
CbmMvdSensorClusterfinderTask::fVerbose
Int_t fVerbose
Definition: CbmMvdSensorClusterfinderTask.h:138
CbmMvdSensorClusterfinderTask::fAdcDynamic
Int_t fAdcDynamic
Definition: CbmMvdSensorClusterfinderTask.h:98
CbmMvdSensorPlugin::fShowDebugHistos
Bool_t fShowDebugHistos
Definition: CbmMvdSensorPlugin.h:72
CbmMvdSensorClusterfinderTask::Exec
void Exec()
Definition: CbmMvdSensorClusterfinderTask.cxx:204
CbmMvdSensorClusterfinderTask
Definition: CbmMvdSensorClusterfinderTask.h:38
CbmMvdCluster
Definition: CbmMvdCluster.h:27
CbmMvdSensorClusterfinderTask::ReInit
Bool_t ReInit()
Definition: CbmMvdSensorClusterfinderTask.cxx:192
CbmMvdSensor::GetStationNr
Int_t GetStationNr() const
Definition: CbmMvdSensor.h:61
CbmMvdSensorClusterfinderTask::CbmMvdSensorClusterfinderTask
CbmMvdSensorClusterfinderTask()
Definition: CbmMvdSensorClusterfinderTask.cxx:17
CbmMvdSensorClusterfinderTask::Finish
void Finish()
Definition: CbmMvdSensorClusterfinderTask.cxx:571
CbmMvdSensorTask::fOutputBuffer
TClonesArray * fOutputBuffer
Definition: CbmMvdSensorTask.h:58
CbmMvdSensorClusterfinderTask::~CbmMvdSensorClusterfinderTask
virtual ~CbmMvdSensorClusterfinderTask()
Definition: CbmMvdSensorClusterfinderTask.cxx:22
CbmMvdSensorClusterfinderTask::fResolutionHistoCleanX
TH1F * fResolutionHistoCleanX
Definition: CbmMvdSensorClusterfinderTask.h:112
CbmMvdSensorTask
Definition: CbmMvdSensorTask.h:26
CbmMvdSensorClusterfinderTask::fDigiMap
std::map< std::pair< Int_t, Int_t >, Int_t > fDigiMap
Definition: CbmMvdSensorClusterfinderTask.h:121
CbmMvdCluster::SetPixelMap
void SetPixelMap(std::map< std::pair< Int_t, Int_t >, Int_t > PixelMap)
Definition: CbmMvdCluster.cxx:30
CbmMvdSensorPlugin::initialized
Bool_t initialized
Definition: CbmMvdSensorPlugin.h:71
CbmMvdSensor
Definition: CbmMvdSensor.h:40
CbmMvdSensorClusterfinderTask::fAdcStepSize
Float_t fAdcStepSize
Definition: CbmMvdSensorClusterfinderTask.h:103
h
Data class with information on a STS local track.
CbmMvdSensorClusterfinderTask::ftempPixelMap
std::map< std::pair< Int_t, Int_t >, Int_t > ftempPixelMap
Definition: CbmMvdSensorClusterfinderTask.h:146
CbmMvdSensorClusterfinderTask::fFullClusterHisto
TH1F * fFullClusterHisto
Definition: CbmMvdSensorClusterfinderTask.h:130
CbmMvdSensorClusterfinderTask::fSeedThreshold
Double_t fSeedThreshold
Definition: CbmMvdSensorClusterfinderTask.h:140
CbmMvdSensorClusterfinderTask::fSigmaNoise
Double_t fSigmaNoise
Definition: CbmMvdSensorClusterfinderTask.h:139
CbmMvdSensorPlugin::GetName
virtual const char * GetName() const
Definition: CbmMvdSensorPlugin.h:62
CbmMvdDigi::GetRefId
Int_t GetRefId() const
Definition: CbmMvdDigi.h:63
CbmMvdSensorClusterfinderTask::CheckForNeighbours
void CheckForNeighbours(std::vector< Int_t > *clusterArray, Int_t clusterDigi, TArrayS *pixelUsed)
Definition: CbmMvdSensorClusterfinderTask.cxx:350
ClassImp
ClassImp(CbmConverterManager) InitStatus CbmConverterManager
Definition: CbmConverterManager.cxx:12
CbmMvdSensorClusterfinderTask::fAdcBits
Int_t fAdcBits
Definition: CbmMvdSensorClusterfinderTask.h:100
CbmMvdSensorClusterfinderTask::fChargeArraySize
static const Short_t fChargeArraySize
Definition: CbmMvdSensorClusterfinderTask.h:161
fCounter
constexpr const Int_t fCounter(0)
CbmCluster::SetAddress
void SetAddress(Int_t address)
Definition: CbmCluster.h:94
CbmMvdSensorTask::fInputBuffer
TClonesArray * fInputBuffer
Definition: CbmMvdSensorTask.h:53
CbmMvdSensorClusterfinderTask::GetAdcCharge
Int_t GetAdcCharge(Float_t charge)
Definition: CbmMvdSensorClusterfinderTask.cxx:417
CbmMvdSensorClusterfinderTask::UpdateDebugHistos
void UpdateDebugHistos(CbmMvdCluster *cluster)
Definition: CbmMvdSensorClusterfinderTask.cxx:430
CbmMvdSensorClusterfinderTask::fBadHitHisto
TH2F * fBadHitHisto
Definition: CbmMvdSensorClusterfinderTask.h:116
CbmMvdDigi
Definition: CbmMvdDigi.h:21
CbmMvdDigi::GetPixelX
Int_t GetPixelX()
Definition: CbmMvdDigi.cxx:141
CbmMvdSensorClusterfinderTask::fResolutionHistoCleanY
TH1F * fResolutionHistoCleanY
Definition: CbmMvdSensorClusterfinderTask.h:113
CbmMvdSensor::GetSensorNr
Int_t GetSensorNr() const
Definition: CbmMvdSensor.h:64
CbmMvdSensorClusterfinderTask::fPixelChargeHistos
TObjArray * fPixelChargeHistos
Definition: CbmMvdSensorClusterfinderTask.h:106
CbmMvdSensorClusterfinderTask::fResolutionHistoX
TH1F * fResolutionHistoX
Definition: CbmMvdSensorClusterfinderTask.h:110
CbmMvdSensorClusterfinderTask::fCounter
Int_t fCounter
Definition: CbmMvdSensorClusterfinderTask.h:137
PairAnalysisStyler::Fill
static Int_t Fill[]
Definition: PairAnalysisStyleDefs.h:82
CbmMvdSensor::GetVolumeId
Int_t GetVolumeId() const
Definition: CbmMvdSensor.h:62