1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
7 * Permission to use, copy, modify and distribute this software and its *
8 * documentation strictly for non-commercial purposes is hereby granted *
9 * without fee, provided that the above copyright notice appears in all *
10 * copies and that both the copyright notice and this permission notice *
11 * appear in the supporting documentation. The authors make no claims *
12 * about the suitability of this software for any purpose. It is *
13 * provided "as is" without express or implied warranty. *
14 **************************************************************************/
18 //-------------------------------------------------------
19 // Implementation of the TPC clusterer
21 // Origin: Marian Ivanov
22 //-------------------------------------------------------
24 #include "AliTPCReconstructor.h"
25 #include "AliTPCclustererMI.h"
26 #include "AliTPCclusterMI.h"
27 #include <TObjArray.h>
32 #include "AliMathBase.h"
34 #include "AliTPCClustersArray.h"
35 #include "AliTPCClustersRow.h"
36 #include "AliDigits.h"
37 #include "AliSimDigits.h"
38 #include "AliTPCParam.h"
39 #include "AliRawReader.h"
40 #include "AliTPCRawStream.h"
41 #include "AliRunLoader.h"
42 #include "AliLoader.h"
43 #include "Riostream.h"
45 #include "AliTPCcalibDB.h"
46 #include "AliTPCCalPad.h"
47 #include "AliTPCCalROC.h"
48 #include "TTreeStream.h"
52 ClassImp(AliTPCclustererMI)
56 AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par)
58 fPedSubtraction = kFALSE;
59 fIsOldRCUFormat = kFALSE;
63 fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
67 AliTPCclustererMI::~AliTPCclustererMI(){
69 if (fAmplitudeHisto) delete fAmplitudeHisto;
70 if (fDebugStreamer) delete fDebugStreamer;
73 void AliTPCclustererMI::SetInput(TTree * tree)
76 // set input tree with digits
79 if (!fInput->GetBranch("Segment")){
80 cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
86 void AliTPCclustererMI::SetOutput(TTree * tree)
91 AliTPCClustersRow clrow;
92 AliTPCClustersRow *pclrow=&clrow;
93 clrow.SetClass("AliTPCclusterMI");
94 clrow.SetArray(1); // to make Clones array
95 fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);
99 Float_t AliTPCclustererMI::GetSigmaY2(Int_t iz){
100 // sigma y2 = in digits - we don't know the angle
101 Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
102 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
103 (fPadWidth*fPadWidth);
105 Float_t res = sd2+sres;
110 Float_t AliTPCclustererMI::GetSigmaZ2(Int_t iz){
111 //sigma z2 = in digits - angle estimated supposing vertex constraint
112 Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
113 Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
114 Float_t angular = fPadLength*(fParam->GetZLength()-z)/(fRx*fZWidth);
117 Float_t sres = fParam->GetZSigma()/fZWidth;
119 Float_t res = angular +sd2+sres;
123 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
126 Int_t i0=k/max; //central pad
127 Int_t j0=k%max; //central time bin
129 // set pointers to data
130 //Int_t dummy[5] ={0,0,0,0,0};
131 Float_t * matrix[5]; //5x5 matrix with digits - indexing i = 0 ..4 j = -2..2
132 Float_t * resmatrix[5];
133 for (Int_t di=-2;di<=2;di++){
134 matrix[di+2] = &bins[k+di*max];
135 resmatrix[di+2] = &fResBins[k+di*max];
137 //build matrix with virtual charge
138 Float_t sigmay2= GetSigmaY2(j0);
139 Float_t sigmaz2= GetSigmaZ2(j0);
141 Float_t vmatrix[5][5];
142 vmatrix[2][2] = matrix[2][0];
144 c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
145 for (Int_t di =-1;di <=1;di++)
146 for (Int_t dj =-1;dj <=1;dj++){
147 Float_t amp = matrix[di+2][dj];
148 if ( (amp<2) && (fLoop<2)){
149 // if under threshold - calculate virtual charge
150 Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
151 amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
153 vmatrix[2+di][2+dj]=amp;
154 vmatrix[2+2*di][2+2*dj]=0;
157 vmatrix[2+2*di][2+dj] =0;
158 vmatrix[2+di][2+2*dj] =0;
163 //if small amplitude - below 2 x threshold - don't consider other one
164 vmatrix[2+di][2+dj]=amp;
165 vmatrix[2+2*di][2+2*dj]=0; // don't take to the account next bin
168 vmatrix[2+2*di][2+dj] =0;
169 vmatrix[2+di][2+2*dj] =0;
173 //if bigger then take everything
174 vmatrix[2+di][2+dj]=amp;
175 vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;
178 vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
179 vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
191 for (Int_t i=-2;i<=2;i++)
192 for (Int_t j=-2;j<=2;j++){
193 Float_t amp = vmatrix[i+2][j+2];
202 Float_t meani = sumiw/sumw;
203 Float_t mi2 = sumi2w/sumw-meani*meani;
204 Float_t meanj = sumjw/sumw;
205 Float_t mj2 = sumj2w/sumw-meanj*meanj;
207 Float_t ry = mi2/sigmay2;
208 Float_t rz = mj2/sigmaz2;
211 if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
212 if ( (ry <1.2) && (rz<1.2) ) {
213 //if cluster looks like expected
214 //+1.2 deviation from expected sigma accepted
215 // c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
219 //set cluster parameters
221 c.SetY(meani*fPadWidth);
222 c.SetZ(meanj*fZWidth);
226 //remove cluster data from data
227 for (Int_t di=-2;di<=2;di++)
228 for (Int_t dj=-2;dj<=2;dj++){
229 resmatrix[di+2][dj] -= vmatrix[di+2][dj+2];
230 if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
236 //unfolding when neccessary
239 Float_t * matrix2[7]; //7x7 matrix with digits - indexing i = 0 ..6 j = -3..3
240 Float_t dummy[7]={0,0,0,0,0,0};
241 for (Int_t di=-3;di<=3;di++){
242 matrix2[di+3] = &bins[k+di*max];
243 if ((k+di*max)<3) matrix2[di+3] = &dummy[3];
244 if ((k+di*max)>fMaxBin-3) matrix2[di+3] = &dummy[3];
246 Float_t vmatrix2[5][5];
249 UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
251 // c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
254 //set cluster parameters
256 c.SetY(meani*fPadWidth);
257 c.SetZ(meanj*fZWidth);
260 c.SetType(Char_t(overlap)+1);
267 printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
272 void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj,
273 Float_t & sumu, Float_t & overlap )
276 //unfold cluster from input matrix
277 //data corresponding to cluster writen in recmatrix
278 //output meani and meanj
280 //take separatelly y and z
282 Float_t sum3i[7] = {0,0,0,0,0,0,0};
283 Float_t sum3j[7] = {0,0,0,0,0,0,0};
285 for (Int_t k =0;k<7;k++)
286 for (Int_t l = -1; l<=1;l++){
287 sum3i[k]+=matrix2[k][l];
288 sum3j[k]+=matrix2[l+3][k-3];
290 Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
293 Float_t sum3wi = 0; //charge minus overlap
294 Float_t sum3wio = 0; //full charge
295 Float_t sum3iw = 0; //sum for mean value
296 for (Int_t dk=-1;dk<=1;dk++){
297 sum3wio+=sum3i[dk+3];
303 if ( ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
304 sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
305 Float_t xm2 = sum3i[-dk+3];
306 Float_t xm1 = sum3i[+3];
307 Float_t x1 = sum3i[2*dk+3];
308 Float_t x2 = sum3i[3*dk+3];
309 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
310 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
311 ratio = w11/(w11+w12);
312 for (Int_t dl=-1;dl<=1;dl++)
313 mratio[dk+1][dl+1] *= ratio;
315 Float_t amp = sum3i[dk+3]*ratio;
320 meani = sum3iw/sum3wi;
321 Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
326 Float_t sum3wj = 0; //charge minus overlap
327 Float_t sum3wjo = 0; //full charge
328 Float_t sum3jw = 0; //sum for mean value
329 for (Int_t dk=-1;dk<=1;dk++){
330 sum3wjo+=sum3j[dk+3];
336 if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
337 (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
338 Float_t xm2 = sum3j[-dk+3];
339 Float_t xm1 = sum3j[+3];
340 Float_t x1 = sum3j[2*dk+3];
341 Float_t x2 = sum3j[3*dk+3];
342 Float_t w11 = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);
343 Float_t w12 = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
344 ratio = w11/(w11+w12);
345 for (Int_t dl=-1;dl<=1;dl++)
346 mratio[dl+1][dk+1] *= ratio;
348 Float_t amp = sum3j[dk+3]*ratio;
353 meanj = sum3jw/sum3wj;
354 Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;
355 overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);
356 sumu = (sum3wj+sum3wi)/2.;
359 //if not overlap detected remove everything
360 for (Int_t di =-2; di<=2;di++)
361 for (Int_t dj =-2; dj<=2;dj++){
362 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
366 for (Int_t di =-1; di<=1;di++)
367 for (Int_t dj =-1; dj<=1;dj++){
369 if (mratio[di+1][dj+1]==1){
370 recmatrix[di+2][dj+2] = matrix2[3+di][dj];
371 if (TMath::Abs(di)+TMath::Abs(dj)>1){
372 recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
373 recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
375 recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
379 //if we have overlap in direction
380 recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];
381 if (TMath::Abs(di)+TMath::Abs(dj)>1){
382 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
383 recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
385 ratio = TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
386 recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
389 ratio = recmatrix[di+2][dj+2]/matrix2[3][0];
390 recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
396 printf("%f\n", recmatrix[2][2]);
400 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
407 for (Int_t di = -1;di<=1;di++)
408 for (Int_t dj = -1;dj<=1;dj++){
409 if (vmatrix[2+di][2+dj]>2){
410 Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
411 sumteor += teor*vmatrix[2+di][2+dj];
412 sumamp += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
415 Float_t max = sumamp/sumteor;
419 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c){
421 // transform cluster to the global coordinata
422 // add the cluster to the array
424 Float_t meani = c.GetY()/fPadWidth;
425 Float_t meanj = c.GetZ()/fZWidth;
427 Int_t ki = TMath::Nint(meani-3);
429 if (ki>=fMaxPad) ki = fMaxPad-1;
430 Int_t kj = TMath::Nint(meanj-3);
432 if (kj>=fMaxTime-3) kj=fMaxTime-4;
433 // ki and kj shifted to "real" coordinata
435 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
436 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
437 c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
441 Float_t s2 = c.GetSigmaY2();
442 Float_t w=fParam->GetPadPitchWidth(fSector);
444 c.SetSigmaY2(s2*w*w);
447 c.SetSigmaZ2(s2*w*w);
448 c.SetY((meani - 2.5 - 0.5*fMaxPad)*fParam->GetPadPitchWidth(fSector));
449 c.SetZ(fZWidth*(meanj-3));
450 c.SetZ(c.GetZ() - 3.*fParam->GetZSigma() + fParam->GetNTBinsL1()*fParam->GetZWidth()); // PASA delay + L1 delay
451 c.SetZ(fSign*(fParam->GetZLength() - c.GetZ()));
453 c.SetDetector(fSector);
456 if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
457 //c.SetSigmaY2(c.GetSigmaY2()*25.);
458 //c.SetSigmaZ2(c.GetSigmaZ2()*4.);
459 c.SetType(-(c.GetType()+3)); //edge clusters
461 if (fLoop==2) c.SetType(100);
463 TClonesArray * arr = fRowCl->GetArray();
464 // AliTPCclusterMI * cl =
465 new ((*arr)[fNcluster]) AliTPCclusterMI(c);
471 //_____________________________________________________________________________
472 void AliTPCclustererMI::Digits2Clusters()
474 //-----------------------------------------------------------------
475 // This is a simple cluster finder.
476 //-----------------------------------------------------------------
479 Error("Digits2Clusters", "input tree not initialised");
484 Error("Digits2Clusters", "output tree not initialised");
488 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
490 AliSimDigits digarr, *dummy=&digarr;
492 fInput->GetBranch("Segment")->SetAddress(&dummy);
493 Stat_t nentries = fInput->GetEntries();
495 fMaxTime=fParam->GetMaxTBin()+6; // add 3 virtual time bins before and 3 after
499 for (Int_t n=0; n<nentries; n++) {
501 if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
502 cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
506 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
509 AliTPCClustersRow *clrow= new AliTPCClustersRow();
511 clrow->SetClass("AliTPCclusterMI");
514 clrow->SetID(digarr.GetID());
515 fOutput->GetBranch("Segment")->SetAddress(&clrow);
516 fRx=fParam->GetPadRowRadii(fSector,row);
519 const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
520 fZWidth = fParam->GetZWidth();
521 if (fSector < kNIS) {
522 fMaxPad = fParam->GetNPadsLow(row);
523 fSign = (fSector < kNIS/2) ? 1 : -1;
524 fPadLength = fParam->GetPadPitchLength(fSector,row);
525 fPadWidth = fParam->GetPadPitchWidth();
527 fMaxPad = fParam->GetNPadsUp(row);
528 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
529 fPadLength = fParam->GetPadPitchLength(fSector,row);
530 fPadWidth = fParam->GetPadPitchWidth();
534 fMaxBin=fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
535 fBins =new Float_t[fMaxBin];
536 fResBins =new Float_t[fMaxBin]; //fBins with residuals after 1 finder loop
537 memset(fBins,0,sizeof(Float_t)*fMaxBin);
538 memset(fResBins,0,sizeof(Float_t)*fMaxBin);
540 if (digarr.First()) //MI change
542 Float_t dig=digarr.CurrentDigit();
543 if (dig<=fParam->GetZeroSup()) continue;
544 Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
545 Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
546 fBins[i*fMaxTime+j]=dig/gain;
547 } while (digarr.Next());
548 digarr.ExpandTrackBuffer();
554 nclusters+=fNcluster;
559 Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
562 void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
564 //-----------------------------------------------------------------
565 // This is a cluster finder for the TPC raw data.
566 // The method assumes NO ordering of the altro channels.
567 // The pedestal subtraction can be switched on and off
568 // using an option of the TPC reconstructor
569 //-----------------------------------------------------------------
572 Error("Digits2Clusters", "output tree not initialised");
578 AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
582 fMaxTime = fParam->GetMaxTBin() + 6; // add 3 virtual time bins before and 3 after
583 const Int_t kNIS = fParam->GetNInnerSector();
584 const Int_t kNOS = fParam->GetNOuterSector();
585 const Int_t kNS = kNIS + kNOS;
586 fZWidth = fParam->GetZWidth();
587 Int_t zeroSup = fParam->GetZeroSup();
589 Float_t** allBins = NULL;
590 Float_t** allBinsRes = NULL;
593 for(fSector = 0; fSector < kNS; fSector++) {
595 AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector); // pad gains per given sector
598 Int_t nDDLs = 0, indexDDL = 0;
599 if (fSector < kNIS) {
600 nRows = fParam->GetNRowLow();
601 fSign = (fSector < kNIS/2) ? 1 : -1;
603 indexDDL = fSector * 2;
606 nRows = fParam->GetNRowUp();
607 fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
609 indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
612 allBins = new Float_t*[nRows];
613 allBinsRes = new Float_t*[nRows];
615 for (Int_t iRow = 0; iRow < nRows; iRow++) {
618 maxPad = fParam->GetNPadsLow(iRow);
620 maxPad = fParam->GetNPadsUp(iRow);
622 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
623 allBins[iRow] = new Float_t[maxBin];
624 allBinsRes[iRow] = new Float_t[maxBin];
625 memset(allBins[iRow],0,sizeof(Float_t)*maxBin);
626 memset(allBinsRes[iRow],0,sizeof(Float_t)*maxBin);
629 // Loas the raw data for corresponding DDLs
631 AliTPCRawStream input(rawReader);
632 input.SetOldRCUFormat(fIsOldRCUFormat);
633 rawReader->Select(0,indexDDL,indexDDL+nDDLs-1);
635 // Begin loop over altro data
636 while (input.Next()) {
638 if (input.GetSector() != fSector)
639 AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
641 if (input.GetTime() < 40) continue;
643 Int_t iRow = input.GetRow();
644 if (iRow < 0 || iRow >= nRows)
645 AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
648 Int_t iPad = input.GetPad() + 3;
652 maxPad = fParam->GetNPadsLow(iRow);
654 maxPad = fParam->GetNPadsUp(iRow);
656 if (input.GetPad() < 0 || input.GetPad() >= maxPad)
657 AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
658 input.GetPad(), 0, maxPad -1));
660 Int_t iTimeBin = input.GetTime() + 3;
661 if ( input.GetTime() < 0 || input.GetTime() >= fParam->GetMaxTBin())
662 AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
663 input.GetTime(), 0, fParam->GetMaxTBin() -1));
665 Int_t maxBin = fMaxTime*(maxPad+6); // add 3 virtual pads before and 3 after
667 if (((iPad*fMaxTime+iTimeBin) >= maxBin) ||
668 ((iPad*fMaxTime+iTimeBin) < 0))
669 AliFatal(Form("Index outside the allowed range"
670 " Sector=%d Row=%d Pad=%d Timebin=%d"
671 " (Max.index=%d)",fSector,iRow,iPad,iTimeBin,maxBin));
673 Float_t signal = input.GetSignal();
674 if (!fPedSubtraction && signal <= zeroSup) continue;
676 Float_t gain = gainROC->GetValue(iRow,input.GetPad());
677 allBins[iRow][iPad*fMaxTime+iTimeBin] = signal/gain;
678 allBins[iRow][iPad*fMaxTime+0] = 1; // pad with signal
679 } // End of the loop over altro data
681 // Now loop over rows and perform pedestal subtraction
682 if (fPedSubtraction) {
683 for (Int_t iRow = 0; iRow < nRows; iRow++) {
686 maxPad = fParam->GetNPadsLow(iRow);
688 maxPad = fParam->GetNPadsUp(iRow);
690 for (Int_t iPad = 0; iPad < maxPad + 6; iPad++) {
691 if (allBins[iRow][iPad*fMaxTime+0] !=1) continue; // no data
692 Float_t *p = &allBins[iRow][iPad*fMaxTime+3];
693 //Float_t pedestal = TMath::Median(fMaxTime, p);
694 Int_t id[3] = {fSector, iRow, iPad-3};
695 Float_t pedestal = ProcesSignal(p, fMaxTime, id);
696 for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
697 allBins[iRow][iPad*fMaxTime+iTimeBin] -= pedestal;
698 if (allBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
699 allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
705 // Now loop over rows and find clusters
706 for (fRow = 0; fRow < nRows; fRow++) {
707 fRowCl = new AliTPCClustersRow;
708 fRowCl->SetClass("AliTPCclusterMI");
710 fRowCl->SetID(fParam->GetIndex(fSector, fRow));
711 fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
713 fRx = fParam->GetPadRowRadii(fSector, fRow);
714 fPadLength = fParam->GetPadPitchLength(fSector, fRow);
715 fPadWidth = fParam->GetPadPitchWidth();
717 fMaxPad = fParam->GetNPadsLow(fRow);
719 fMaxPad = fParam->GetNPadsUp(fRow);
720 fMaxBin = fMaxTime*(fMaxPad+6); // add 3 virtual pads before and 3 after
722 fBins = allBins[fRow];
723 fResBins = allBinsRes[fRow];
729 nclusters += fNcluster;
730 } // End of loop to find clusters
733 for (Int_t iRow = 0; iRow < nRows; iRow++) {
734 delete [] allBins[iRow];
735 delete [] allBinsRes[iRow];
739 delete [] allBinsRes;
741 } // End of loop over sectors
743 Info("Digits2Clusters", "Number of found clusters : %d\n", nclusters);
747 void AliTPCclustererMI::FindClusters()
749 //add virtual charge at the edge
750 for (Int_t i=0; i<fMaxTime; i++){
751 Float_t amp1 = fBins[i+3*fMaxTime];
754 Float_t amp2 = fBins[i+4*fMaxTime];
755 if (amp2==0) amp2=0.5;
756 Float_t sigma2 = GetSigmaY2(i);
757 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
758 if (gDebug>4) printf("\n%f\n",amp0);
760 fBins[i+2*fMaxTime] = amp0;
762 amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
764 Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
765 if (amp2==0) amp2=0.5;
766 Float_t sigma2 = GetSigmaY2(i);
767 amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);
768 if (gDebug>4) printf("\n%f\n",amp0);
770 fBins[(fMaxPad+3)*fMaxTime+i] = amp0;
773 // memcpy(fResBins,fBins, fMaxBin*2);
774 memcpy(fResBins,fBins, fMaxBin);
777 //first loop - for "gold cluster"
779 Float_t *b=&fBins[-1]+2*fMaxTime;
780 Int_t crtime = Int_t((fParam->GetZLength()-AliTPCReconstructor::GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
782 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
784 if (*b<8) continue; //threshold form maxima
785 if (i%fMaxTime<crtime) {
786 Int_t delta = -(i%fMaxTime)+crtime;
792 if (!IsMaximum(*b,fMaxTime,b)) continue;
795 MakeCluster(i, fMaxTime, fBins, dummy,c);
798 //memcpy(fBins,fResBins, fMaxBin*2);
799 //second loop - for rest cluster
802 b=&fResBins[-1]+2*fMaxTime;
803 for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
805 if (*b<25) continue; // bigger threshold for maxima
806 if (!IsMaximum(*b,fMaxTime,b)) continue;
809 MakeCluster(i, fMaxTime, fResBins, dummy,c);
816 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3]){
818 // process signal on given pad - + streaming of additional information in special mode
825 Float_t kMaxNoise = 3;
826 Double_t median = TMath::Median(nchannels-offset, &(signal[offset]));
827 if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
830 Double_t mean = TMath::Mean(nchannels-offset, &(signal[offset]));
831 Double_t rms = TMath::RMS(nchannels-offset, &(signal[offset]));
832 Double_t *dsignal = new Double_t[nchannels];
833 Double_t *dtime = new Double_t[nchannels];
836 for (Int_t i=0; i<nchannels; i++){
838 dsignal[i] = signal[i];
839 if (signal[i]>max && i <fMaxTime-100) { // temporary remove spike signals at the end
845 Double_t mean06=0,mean07=0, mean08=0, mean09=0;
846 Double_t rms06=0, rms07=0, rms08, rms09=0;
847 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean06, rms06, int(0.6*(nchannels-offset)));
848 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean07, rms07, int(0.7*(nchannels-offset)));
849 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean08, rms08, int(0.8*(nchannels-offset)));
850 AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean09, rms09, int(0.9*(nchannels-offset)));
853 UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
854 if (uid[0]< AliTPCROC::Instance()->GetNSectors()
855 && uid[1]< AliTPCROC::Instance()->GetNRows(uid[0]) &&
856 uid[2] < AliTPCROC::Instance()->GetNPads(uid[0], uid[1])){
857 if (!fAmplitudeHisto){
858 fAmplitudeHisto = new TObjArray(72);
860 TObjArray * sectorArray = (TObjArray*)fAmplitudeHisto->UncheckedAt(uid[0]);
862 Int_t npads = AliTPCROC::Instance()->GetNChannels(uid[0]);
863 sectorArray = new TObjArray(npads);
864 fAmplitudeHisto->AddAt(sectorArray, uid[0]);
866 Int_t position = uid[2]+ AliTPCROC::Instance()->GetRowIndexes(uid[0])[uid[1]];
867 TH1F * histo = (TH1F*)sectorArray->UncheckedAt(position);
870 sprintf(hname,"Amp_%d_%d_%d",uid[0],uid[1],uid[2]);
871 TFile * backup = gFile;
872 fDebugStreamer->GetFile()->cd();
873 histo = new TH1F(hname, hname, 100, 1,100);
874 //histo->SetDirectory(0); // histogram not connected to directory -(File)
875 sectorArray->AddAt(histo, position);
876 if (backup) backup->cd();
878 for (Int_t i=0; i<nchannels; i++){
879 if (signal[i]>0) histo->Fill(signal[i]);
884 Bool_t random = (gRandom->Rndm()<0.0001);
885 if (max-median>kMin || rms06>2.*fParam->GetZeroSup() || random){
886 graph =new TGraph(nchannels, dtime, dsignal);
887 if (rms06>2.*fParam->GetZeroSup() || random)
888 (*fDebugStreamer)<<"SignalN"<< //noise pads - or random sample of pads
909 (*fDebugStreamer)<<"SignalB"<< // pads with signal
932 (*fDebugStreamer)<<"Signal"<<
953 if (rms06>kMaxNoise) return 1024+median; // sign noisy channel in debug mode
959 void AliTPCclustererMI::DumpHistos(){
960 if (!fAmplitudeHisto) return;
961 for (UInt_t isector=0; isector<AliTPCROC::Instance()->GetNSectors(); isector++){
962 TObjArray * array = (TObjArray*)fAmplitudeHisto->UncheckedAt(isector);
963 if (!array) continue;
964 for (UInt_t ipad = 0; ipad <(UInt_t)array->GetEntriesFast(); ipad++){
965 TH1F * histo = (TH1F*) array->UncheckedAt(ipad);
966 if (!histo) continue;
967 if (histo->GetEntries()<100) continue;
968 histo->Fit("gaus","q");
969 Float_t mean = histo->GetMean();
970 Float_t rms = histo->GetRMS();
971 Float_t gmean = histo->GetFunction("gaus")->GetParameter(1);
972 Float_t gsigma = histo->GetFunction("gaus")->GetParameter(2);
973 Float_t max = histo->GetFunction("gaus")->GetParameter(0);
976 UInt_t row=0, pad =0;
977 const UInt_t *indexes = AliTPCROC::Instance()->GetRowIndexes(isector);
978 for (UInt_t irow=0; irow< AliTPCROC::Instance()->GetNRows(isector); irow++){
979 if (indexes[irow]<ipad){
981 pad = ipad-indexes[irow];
984 Int_t rpad = pad - (AliTPCROC::Instance()->GetNPads(isector,row))/2;
986 (*fDebugStreamer)<<"Fit"<<
997 if (array->UncheckedAt(ipad)) fDebugStreamer->StoreObject(array->UncheckedAt(ipad));