]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCClusterFinder.cxx
Put back the two bellows in front of the absorber.
[u/mrichter/AliRoot.git] / TPC / AliTPCClusterFinder.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
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  **************************************************************************/
15
16 /*
17 $Log$
18 Revision 1.2  2000/06/30 12:07:49  kowal2
19 Updated from the TPC-PreRelease branch
20
21 Revision 1.1.2.1  2000/06/25 08:52:51  kowal2
22 replacing AliClusterFinder
23
24 */
25
26 //-----------------------------------------------------------------------------
27 //
28 //  Implementation of class ALICLUSTERFINDER
29 // 
30 //Class for cluster finding in two dimension.
31 //In the present there are implemented two algorithm
32 //primitive recursion algorithm. (FindPeaks) 
33 //Algorithm is not working in case of overlaping clusters
34 //Maximum - minimum in direction algoritm (Find clusters)
35 //In this algoritm we suppose that each cluster has local 
36 //maximum. From this local maximum I mus see each point 
37 //of cluster.
38 //From maximum i can accept every point in radial 
39 //direction which is before border in direction
40 //Border in direction occur if we have next in
41 //direction nder threshold or response begin
42 //to increase in given radial direction
43 //-----------------------------------------------------------------------------
44
45 #include "TMinuit.h"
46 #include "AliArrayI.h"
47 #include "TClonesArray.h"
48 #include "AliTPC.h"
49 #include "TRandom.h"
50 #include "AliH2F.h"
51 #include "TMarker.h"
52 #include "AliCluster.h"
53 #include "AliTPCClusterFinder.h"
54 #include <fstream.h>
55
56 //direction constants possible direction in 8 different sectors
57 //
58
59
60 const Int_t kClStackSize =1000;
61
62
63
64
65 static AliTPCClusterFinder * gClusterFinder; //for fitting routine
66
67 void gauss(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t iflag)
68 {
69   AliArrayI * points = gClusterFinder->GetStack();
70   const Int_t nbins = gClusterFinder->GetStackIndex();
71   Int_t i;
72   //calculate chisquare
73    Double_t chisq = 0;
74    Double_t delta;
75    for (i=0;i<nbins; i++) {
76      Float_t x = points->At(i*3);
77      Float_t y = points->At(i*3+1);
78      Float_t z = points->At(i*3+2);
79      Float_t deltax2 = (x-par[1]);
80      deltax2*=deltax2;
81      deltax2*=par[3];
82      Float_t deltay2 = (y-par[2]);
83      deltay2*=deltay2;
84      deltay2*=par[4];
85      
86      delta  = z-par[0]*TMath::Exp(-deltax2-deltay2);
87      chisq += delta*delta;
88    }
89    f = chisq;   
90 }
91
92
93 ClassImp(AliTPCClusterFinder)
94   //ClassImp(AliCell)
95
96 AliTPCClusterFinder::AliTPCClusterFinder()
97 {
98   fDigits =0;
99   fDimX = 0;
100   fDimY = 0;
101   fNoiseTh = 3;
102   fMulSigma2 = 16; //4 sigma
103   fDirSigmaFac = 1.4;
104   fDirAmpFac =1.3;
105   fNType=8;
106   fThreshold = 2;
107   fStack = new AliArrayI;
108   fStack->Set(kClStackSize);  
109   fClustersArray =0;
110   SetSigmaX(1,0,0);
111   SetSigmaY(1,0,0);
112
113
114   fDetectorParam = 0;
115   ResetStatus(); 
116   fBFit = kFALSE;
117   fMinuit= new TMinuit(5);
118   fMinuit->SetFCN(gauss);
119   gClusterFinder = this;
120   
121 }
122
123
124 AliTPCClusterFinder::~AliTPCClusterFinder()
125 {
126  if (fDigits  != 0) delete fDigits;
127 }
128
129 void AliTPCClusterFinder::SetSigmaX(Float_t s0, Float_t s1x, Float_t s1y)
130 {
131   fSigmaX[0]=s0;
132   fSigmaX[1]=s1x;
133   fSigmaX[2]=s1y;
134
135 }
136 void AliTPCClusterFinder::SetSigmaY(Float_t s0, Float_t s1x, Float_t s1y)
137 {
138   fSigmaY[0]=s0;
139   fSigmaY[1]=s1x;
140   fSigmaY[2]=s1y;
141 }
142
143
144
145 Bool_t AliTPCClusterFinder::SetSigma2(Int_t i, Int_t j, Float_t & sigmax2, Float_t &sigmay2)
146 {
147   //
148   //set sigmax2 and sigma y2  accordig i and j position of cell 
149   //
150   
151   //  Float_t x[3] = {ItoX(i),JtoY(j),0};
152   Float_t x= ItoX(i);
153   Float_t y= JtoY(j);
154
155   sigmax2= fSigmaX[0]+fSigmaX[1]*x+fSigmaX[2]*y;
156   sigmay2= fSigmaY[0]+fSigmaY[1]*x+fSigmaY[2]*y;
157   return kTRUE;  
158 }
159
160 /*
161 Bool_t AliTPCClusterFinder::SetSigma2(Int_t i, Int_t j, Float_t & sigmax2, Float_t &sigmay2)
162 {
163   //
164   //set sigmax2 and sigma y2  accordig i and j position of cell 
165   //
166   if (fDetectorParam==0) {
167     sigmax2=1;
168     sigmay2=1;
169     return kFALSE;
170   }
171   Float_t x[3] = {ItoX(i),JtoY(j),0};
172   Float_t sigma[2];
173   fDetectorParam->GetClusterSize(x,fDetectorIndex,0,0,sigma);
174   sigmax2=sigma[0]*(fX2-fX1)*(fX2-fX1)/(fDimX*fDimX);
175   sigmay2=sigma[1]*(fY2-fY1)*(fY2-fY1)/(fDimY*fDimY);
176   return kTRUE;
177 }
178 */
179
180
181 void AliTPCClusterFinder::GetHisto(TH2F * his2)
182 {
183   
184   UInt_t idim =his2->GetNbinsX();
185   UInt_t jdim =his2->GetNbinsY();
186   fX1 = his2->GetXaxis()->GetXmin();
187   fX2 = his2->GetXaxis()->GetXmax();
188   fY1 = his2->GetYaxis()->GetXmin();
189   fY2 = his2->GetYaxis()->GetXmax();
190  
191   if ( (idim>0) && (jdim>0))
192     {
193       rOK = kTRUE;
194       fDimX = idim;
195       fDimY = jdim;
196       Int_t size =idim*jdim;       
197       if (fDigits !=0) delete fDigits;
198       fDigits  = (Int_t*) new Int_t[size];
199       fCells  = (AliCell*) new AliCell[size];
200
201     }  else 
202       rOK=kFALSE;
203   for (Int_t i = 0; i<(Int_t)idim;i++)    
204     for (Int_t j = 0; j<(Int_t)jdim;j++)
205       {
206         Int_t index = his2->GetBin(i+1,j+1);
207         //AliCell * cell = GetCell(i,j);
208         //if (cell!=0) cell->SetSignal(his2->GetBinContent(index));
209         SetSignal(his2->GetBinContent(index),i,j);
210       }
211    
212 }
213
214
215
216
217 void AliTPCClusterFinder::FindMaxima()
218 {
219   for (Int_t i=0; i<fDimX; i++) 
220     for (Int_t j=0;j<fDimY; j++)      
221       if (IsMaximum(i,j)) cout<<i<<"   "<<j<<"\n";                   
222 }
223
224  
225 void  AliTPCClusterFinder::Transform(AliDigitCluster * c)
226 {
227   //transform coordinata from bin coordinata to "normal coordinata"
228   //for example if we initialize finder with histogram
229   //it transform values from bin coordinata to the histogram coordinata
230   c->fX=ItoX(c->fX);
231   c->fY=JtoY(c->fY);
232   c->fMaxX=ItoX(c->fMaxX);
233   c->fMaxY=JtoY(c->fMaxY);
234
235   c->fSigmaX2=c->fSigmaX2*(fX2-fX1)*(fX2-fX1)/(fDimX*fDimX);
236   c->fSigmaY2=c->fSigmaY2*(fY2-fY1)*(fY2-fY1)/(fDimY*fDimY);  
237   c->fArea   =c->fArea*(fX2-fX1)*(fY2-fY1)/(fDimX*fDimY); 
238 }
239 void  AliTPCClusterFinder::AddToStack(Int_t i, Int_t j, Int_t signal)
240 {
241   //
242   //add digit to stack
243   //
244   if ( ((fStackIndex+2)>=kClStackSize) || (fStackIndex<0) ) return; 
245   fStack->AddAt(i,fStackIndex);
246   fStack->AddAt(j,fStackIndex+1);
247   fStack->AddAt(signal,fStackIndex+2);
248   fStackIndex+=3;  
249 }
250
251 void AliTPCClusterFinder::GetClusterStatistic(AliDigitCluster & cluster)
252 {
253   //
254   //calculate statistic of cluster 
255   //
256   Double_t sumxw,sumyw,sumx2w,sumy2w,sumxyw,sumw;
257   Int_t minx,maxx,miny,maxy;
258   sumxw=sumyw=sumx2w=sumy2w=sumxyw=sumw=0;
259   minx=fDimX;
260   maxx=-fDimX;
261   miny=fDimY;
262   maxy=-fDimY;
263   Int_t x0=fStack->At(0);
264   Int_t y0=fStack->At(1);
265   Int_t maxQx =x0;
266   Int_t maxQy =y0;  
267   Int_t maxQ=fStack->At(2);
268   
269
270   for (Int_t i = 0; i<fStackIndex;i+=3){
271     Int_t x = fStack->At(i);
272     Int_t y = fStack->At(i+1);
273     Int_t dx=x-x0;
274     Int_t dy=y-y0;
275     Int_t w = fStack->At(i+2);
276     if (w>maxQ){
277       maxQ = w;
278       maxQx = x;
279       maxQy=y;
280     }   
281     if (x<minx) minx=x;
282     if (y<miny) miny=y;
283     if (x>maxx) maxx=x;
284     if (y>maxy) maxy=y;   
285     sumxw+=dx*w;
286     sumyw+=dy*w;
287     sumx2w+=dx*dx*w;
288     sumy2w+=dy*dy*w;
289     sumxyw+=dx*dy*w;
290     sumw+=w;    
291   }
292   cluster.fQ = sumw;
293   if (sumw>0){
294     cluster.fX = sumxw/sumw;
295     cluster.fY = sumyw/sumw;
296     cluster.fQ = sumw;
297     cluster.fSigmaX2 = sumx2w/sumw-cluster.fX*cluster.fX;
298     cluster.fSigmaY2 = sumy2w/sumw-cluster.fY*cluster.fY;
299     cluster.fSigmaXY = sumxyw/sumw-cluster.fX*cluster.fY;
300     cluster.fMaxX = maxQx;
301     cluster.fMaxY = maxQy; 
302     cluster.fMax = maxQ;
303     cluster.fArea = fStackIndex/3;  
304     cluster.fNx = maxx-minx+1;
305     cluster.fNy = maxy-miny+1;
306     cluster.fX +=x0; 
307     cluster.fY +=y0;
308   }
309 }
310 void AliTPCClusterFinder::GetClusterFit(AliDigitCluster & cluster)
311 {
312   //
313   //calculate statistic of cluster 
314   //
315   Double_t arglist[10];
316   Int_t ierflg = 0;
317   
318   arglist[0] = 1;
319   fMinuit->mnexcm("SET ERR", arglist ,1,ierflg);
320   
321   //fistly find starting parameters
322   Int_t minx,maxx,miny,maxy,maxQ,maxQx,maxQy;
323   Int_t over =0;
324   Float_t sumxw,sumyw,sumw;
325   sumxw=sumyw=sumw=0;
326   minx=fDimX;
327   maxx=-fDimX;
328   miny=fDimY;
329   maxy=-fDimY;
330   maxQx=fStack->At(0);
331   maxQy=fStack->At(1);
332   maxQ=fStack->At(2);
333  
334   for (Int_t i = 0; i<fStackIndex;i+=3){
335     Int_t x = fStack->At(i);
336     Int_t y = fStack->At(i+1);
337     Int_t w = fStack->At(i+2);
338     if (w>fThreshold) {
339       over++;
340       sumw+=w;    
341       sumxw+=x*w;
342       sumyw+=y*w;
343       if (x<minx) minx=x;
344       if (y<miny) miny=y;
345       if (x>maxx) maxx=x;
346       if (y>maxy) maxy=y;
347       if (w>maxQ) {
348         maxQ=w;   
349         maxQx=x;
350         maxQy=y;
351       }    
352     }
353   }
354   Int_t nx = maxx-minx+1;
355   Int_t ny = maxy-miny+1;
356   
357   SetSigma2(maxQx,maxQy,fCurrentSigmaX2,fCurrentSigmaY2);
358   Double_t vstart[5]={maxQ,sumxw/sumw,sumyw/sumw,1/(2*fCurrentSigmaX2),1/(2*fCurrentSigmaY2)};
359   Double_t step[5]={1.,0.01,0.01,0.01,0.01};
360   fMinuit->mnparm(0, "amp", vstart[0], step[0], 0,0,ierflg);
361   fMinuit->mnparm(1, "x0", vstart[1], step[1], 0,0,ierflg);
362   fMinuit->mnparm(2, "y0", vstart[2], step[2], 0,0,ierflg);
363   fMinuit->mnparm(3, "sx2", vstart[3], step[3], 0,0,ierflg);
364   fMinuit->mnparm(4, "sy2", vstart[4], step[4], 0,0,ierflg);
365   arglist[0] = 500;
366   arglist[1] = 1.;
367
368   fMinuit->mnfree(0);  //set unfixed all parameters
369   //if we have area less then
370   if (over<=21) {  //if we dont't have more then  7  points
371     fMinuit->FixParameter(3); 
372     fMinuit->FixParameter(4);
373   }
374   else {
375     if (nx<3)  fMinuit->FixParameter(3); //fix sigma x if no data in x direction
376     if (ny<3)  fMinuit->FixParameter(4);  //fix sigma y if no data in y direction
377   }
378   fMinuit->mnexcm("MIGRAD", arglist ,2,ierflg);
379  
380   if (sumw>0){
381     Double_t  x[5];
382     Double_t  error[5];
383     fMinuit->GetParameter(0,x[0],error[0]);
384     fMinuit->GetParameter(1,x[1],error[1]);
385     fMinuit->GetParameter(2,x[2],error[2]);
386     fMinuit->GetParameter(3,x[3],error[3]);
387     fMinuit->GetParameter(4,x[4],error[4]);
388
389     cluster.fX = x[1];
390     cluster.fY = x[2];
391     cluster.fMaxX = maxQx;
392     cluster.fMaxY = maxQy;
393     
394     cluster.fQ = sumw;
395     cluster.fSigmaX2 = 1/TMath::Sqrt(2*x[3]);
396     cluster.fSigmaY2 = 1/TMath::Sqrt(2*x[4]);
397     cluster.fSigmaXY = 0;
398     cluster.fMax = x[0];
399     cluster.fArea = over;  
400     cluster.fNx = nx;
401     cluster.fNy = ny;
402   }
403 }
404
405 Bool_t   AliTPCClusterFinder::CheckIfDirBorder(Float_t x, Float_t y, 
406                                           Int_t i,Int_t j)
407 {
408   //
409   //function which control if given cell with index i, j is the 
410   //minimum in direction
411   // x and y are estimate of local maximum 
412   //direction is given by the 
413   Float_t virtualcell;
414   AliCell * cellor= GetCell(i,j);
415   Int_t     sigor = GetSignal(i,j);
416
417   //control derivation in direction
418   //if function grows up in direction then there is border
419   Float_t dx = i-x;
420   Float_t dy = j-y; 
421   Float_t dd = TMath::Sqrt(dx*dx+dy*dy);
422   Float_t ddx = TMath::Abs(dx);
423   ddx =   (ddx>0.5) ? ddx-0.5: 0;
424   ddx*=ddx;
425   Float_t ddy = TMath::Abs(dy);
426   ddy = (ddy>0.5) ? ddy-0.5: 0;
427   ddy*=ddy;
428   Float_t d2 = ddx/(2*fDirSigmaFac*fCurrentSigmaX2)+ddy/(2*fDirSigmaFac*fCurrentSigmaY2); //safety factor 
429   //I accept sigmax and sigma y bigge by factor sqrt(fDirsigmaFac)
430   Float_t amp = TMath::Exp(-d2)*fCurrentMaxAmp*fDirAmpFac; //safety factor fDirFac>1
431
432   if (sigor>amp) return kTRUE; 
433   if (dd==0) return kFALSE;
434
435   dx/=dd;
436   dy/=dd;  
437   virtualcell = GetVirtualSignal(i+dx,j+dy);
438   if (virtualcell <=fThreshold) return kFALSE;
439   if (virtualcell>sigor)
440     if (virtualcell>(sigor+fNoiseTh))
441       {cellor->SetDirBorder(fIndex); return kTRUE;}
442     else
443       {
444         virtualcell = GetVirtualSignal(i+2*dx,j+2*dy);
445         if (virtualcell>sigor)
446           { cellor->SetDirBorder(fIndex); return kTRUE;}       
447       };
448   return kFALSE;  
449 }
450
451
452
453
454
455 Bool_t  AliTPCClusterFinder::IsMaximum(Int_t i, Int_t  j)
456 {
457   //there is maximum if given digits is 1 sigma over all adjacent
458   //in 8 neighborow 
459   //or ther exist virual maximum
460   //is maximum on 24 points neighboring
461   //  Bool_t res = kFALSE;
462   Int_t over =0;
463   Int_t overth=0;
464   Int_t oversigma =0;
465   AliCell * cell = GetCell(i,j); 
466   Int_t signal = GetSignal(i,j);
467   if (cell == 0) return kFALSE;
468   for ( Int_t di=-1;di<=1;di++)
469     for ( Int_t dj=-1;dj<=1;dj++){      
470       if ( (di!=0) || (dj!=0))
471         {
472           AliCell * cell2=GetCell(i+di,j+dj);
473           Int_t signal2 = GetSignal(i+di,j+dj);
474           if (cell2 == 0) {
475             over+=1;
476             oversigma+=1;
477           }
478           else
479             {
480               if (signal2>signal) return kFALSE;
481               if (signal2>fThreshold) overth++;
482               if (signal2==signal) {
483                 if (di<0) return kFALSE; 
484                 if ( (di+dj)<0) return kFALSE;
485               }
486               //              if (signal>=signal2){
487               over+=1;
488               if (signal>fNoiseTh+signal2)           
489                 oversigma+=1;           
490               //}
491             }
492         }
493     }
494   //if I have only one neighborough over threshold 
495   if (overth<2) return kFALSE;
496   if (over<8)   return kFALSE;
497   if (oversigma==8) {
498     fCurrentMaxX = i;
499     fCurrentMaxY = j;
500     fCurrentMaxAmp =signal;
501     SetMaximum(fIndex,i,j);
502     return kTRUE;
503   }
504   //check if there exist virtual maximum
505   for (Float_t ddi=0.;(ddi<1.);ddi+=0.5)
506     for (Float_t ddj=0.;(ddj<1.);ddj+=0.5)                
507       if (IsVirtualMaximum(Float_t(i)+ddi,Float_t(j)+ddj)){
508         fCurrentMaxX = i+ddi;
509         fCurrentMaxY = j+ddj;
510         fCurrentMaxAmp =signal; 
511         SetMaximum(fIndex,i,j);
512         return kTRUE;   
513       }
514   return kFALSE;
515 }
516
517 Bool_t  AliTPCClusterFinder::IsVirtualMaximum(Float_t x, Float_t  y)
518 {
519   //there is maximum if given digits is 1 sigma over all adjacent
520   //in 8 neighborow or 
521   //is maximum on 24 points neighboring
522   Bool_t res = kFALSE;
523   Int_t over =0;
524   Int_t overth=0;
525   Int_t oversigma =0;
526   Float_t virtualcell = GetVirtualSignal(x,y); 
527   if (virtualcell < 0) return kFALSE;
528   for ( Int_t di=-1;di<=1;di++)
529     for ( Int_t dj=-1;dj<=1;dj++)
530       if ( (di!=0) || (dj!=0))
531         {
532           Float_t virtualcell2=GetVirtualSignal(x+di,y+dj);
533           if (virtualcell2 < 0) {
534             over+=1;
535             oversigma+=1;
536           }
537           else
538             {
539               if (virtualcell2>fThreshold) overth++;
540               if (virtualcell>=virtualcell2){
541                 over+=1;
542                 if (virtualcell>fNoiseTh+virtualcell2)       
543                   oversigma+=1;
544               }
545             }
546         }
547   if (overth<2) return kFALSE;
548   //if there exist only one or less neighboring above threshold
549   if (oversigma==8)  res = kTRUE;
550   else if ((over==8)&&(GetNType()==8)) res=kTRUE;
551   else if (over ==8 )
552     for ( Int_t di=-2;di<=2;di++)
553       for ( Int_t dj=-2;dj<=2;dj++)
554         if ( (di==2)||(di==-2) || (dj==2)|| (dj==-2) )
555           {
556             Float_t virtualcell2=GetVirtualSignal(x+di,y+dj);
557             if (virtualcell2 < 0) {
558               over+=1;
559               oversigma+=1;
560             }
561             else
562               {
563                 if (virtualcell>=virtualcell2) over+=1;
564               }
565           }     
566   if (over == 24) res=kTRUE;
567   return res;
568   
569 }
570
571
572 void AliTPCClusterFinder::ResetSignal()
573 {
574    //reset dignals to 0
575   Int_t size = fDimX*fDimY;
576   AliCell *dig=fCells;
577   if (rOK==kTRUE) for (Int_t i=0 ; i<size;i++) dig[i] = 0; 
578 }
579
580
581
582 void AliTPCClusterFinder::ResetStatus()
583 {
584    //reset status of signals to not used
585   Int_t size = fDimX*fDimY;
586   AliCell *dig=fCells;
587   if (rOK==kTRUE) for (Int_t i=0 ; i<size;i++) 
588       dig[i].SetStatus(0);     
589
590
591
592 AliCell  *  AliTPCClusterFinder::GetCell(Int_t i, Int_t j)
593 {
594   //return reference to the cell with index i,j 
595   if (rOK == kTRUE)
596     if ( (i>=0) && (i<fDimX) && (j>=0) && (j<fDimY) )
597       return &fCells[i+j*fDimX];
598   return 0; 
599 }
600
601 Float_t   AliTPCClusterFinder::GetVirtualSignal(Float_t ri, Float_t rj)
602 {
603   //it generate virtual cell as mean value from different cels
604   //after using it must be destructed !!!  
605   Int_t i=(Int_t)ri;
606   Int_t j=(Int_t)rj;
607   Int_t ddi = (ri>i)? 1:0;
608   Int_t ddj = (rj>j)? 1:0;
609   Float_t sum = 0;
610   Float_t sumw= 0;
611   for (Int_t di=0;di<=ddi;di++)   
612     for (Int_t dj=0;dj<=ddj;dj++)
613       {             
614         Float_t w = (ri-i-di)*(ri-i-di)+(rj-j-dj)*(rj-j-dj);    
615         if (w>0) w=1/TMath::Sqrt(w);
616         else w=9999999;
617         AliCell * cel2 =GetCell(i+di,j+dj);
618         Int_t signal2 = GetSignal(i+di,j+dj);
619         if (cel2!=0) {
620           sumw+=w;
621           sum+= signal2*w;
622         }
623       }
624   if (sumw>0)  return (sum/sumw);
625   else 
626     return -1;
627 }
628
629
630
631 void AliTPCClusterFinder::Streamer(TBuffer & R__b)
632 {
633   if (R__b.IsReading()) {
634     //      Version_t R__v = R__b.ReadVersion();
635    } else {
636       R__b.WriteVersion(AliTPCClusterFinder::IsA());    
637    } 
638 }
639
640
641
642 void AliTPCClusterFinder::SetBlockIndex(Int_t * index)
643 {
644   //
645   //calculate which indexes we must check for border
646   //
647   if (TMath::Abs(index[0])<2) index[2] = 0;
648   else {
649     index[2] = TMath::Abs(index[0])-1;
650     if (index[0]<0) index[2]*=-1;   //first x block
651   } 
652   if (TMath::Abs(index[1])<2) index[3] = 0;
653   else {
654     index[3] = TMath::Abs(index[1])-1;
655     if (index[1]<0) index[3]*=-1;   //first y block
656   } 
657   if (TMath::Abs(index[0])<TMath::Abs(index[1])){
658     index[4]=index[0];
659     index[5]=index[3];
660   }
661   else
662     if (index[0]==index[1]) {
663       index[4]=0;
664       index[5]=0;
665     }
666     else{
667       index[4]=index[2];
668       index[5]=index[1]; 
669     }
670   return;  
671 }
672
673 //***********************************************************************
674 //***********************************************************************
675
676 TClonesArray * AliTPCClusterFinder::FindPeaks1(TClonesArray *arr)
677 {
678   //find peaks and write it in form of AliTPCcluster to array
679   if (arr==0){
680     fClustersArray=new TClonesArray("AliDigitCluster",300);
681     fIndex=1;
682   }
683   else {
684     fClustersArray = arr;
685     fIndex = fClustersArray->GetEntriesFast();
686   }
687  
688   AliDigitCluster c;           
689   ResetStatus();  
690    for (Int_t i=0; i<fDimX; i++) 
691      for (Int_t j=0;j<fDimY; j++) 
692        {        
693          fStackIndex=0;          
694          fBDistType = kFALSE;
695          AliCell * cell = GetCell(i,j);
696          if (!(cell->IsChecked()))  Adjacent(i,j);
697          //if there exists more then 2 digits cluster 
698          if (fStackIndex >2 ){     
699            if (fBFit==kFALSE) GetClusterStatistic(c);
700              else GetClusterFit(c);
701            //write some important chracteristic area of cluster
702            //      
703            Transform(&c);
704            //write cluster information to array
705            TClonesArray &lclusters = *fClustersArray;
706            new (lclusters[fIndex++])  AliDigitCluster(c);
707            //             cout<<"fx="<<c.fX<<"   fy"<<c.fY<<"\n";          
708          } 
709        }
710    return fClustersArray;
711 }
712
713
714 TClonesArray * AliTPCClusterFinder::FindPeaks2(TClonesArray *arr)
715 {
716   //find peaks and write it in form of AliTPCcluster to array
717   if (arr==0){
718     fClustersArray=new TClonesArray("AliDigitCluster",300);
719     fIndex=1;
720   }
721   else {
722     fClustersArray = arr;
723     fIndex = fClustersArray->GetEntriesFast();
724   }
725
726   AliDigitCluster c;           
727   ResetStatus();  
728   
729    for (Int_t i=0; i<fDimX; i++) 
730      for (Int_t j=0;j<fDimY; j++) 
731        {
732          fStackIndex=0;  
733          if (IsMaximum(i,j) == kTRUE){
734            SetSigma2(i,j,fCurrentSigmaX2,fCurrentSigmaY2);
735            fBDistType = kTRUE;
736            Adjacent(i,j);
737            //if there exists more then 2 digits cluster 
738            if (fStackIndex >2 ){
739              if (fBFit==kFALSE) GetClusterStatistic(c);
740              else GetClusterFit(c);
741              //write some important chracteristic area of cluster
742              //    
743              Transform(&c);
744              //write cluster information to array
745              TClonesArray &lclusters = *fClustersArray;
746              new(lclusters[fIndex++]) AliDigitCluster(c);
747              //             cout<<"fx="<<c.fX<<"   fy"<<c.fY<<"\n";        
748            } 
749          }
750        }
751    return fClustersArray;
752 }
753
754
755 TClonesArray * AliTPCClusterFinder::FindPeaks3(TClonesArray *arr)
756 {
757   //find peaks and write it in form of AliTPCcluster to array
758   if (arr==0){
759     fClustersArray=new TClonesArray("AliDigitCluster",300);
760     fIndex=1;
761   }
762   else {
763     fClustersArray = arr;
764     fIndex = fClustersArray->GetEntriesFast();
765   }
766   
767   AliDigitCluster c;    
768   ResetStatus();  
769   
770   Int_t dmax=5;
771   Int_t naccepted =1;
772    for (Int_t i=0; i<fDimX; i++) 
773      for (Int_t j=0;j<fDimY; j++) 
774        {
775          fStackIndex=0;  
776          if (IsMaximum(i,j) == kTRUE){
777            SetSigma2(i,j,fCurrentSigmaX2,fCurrentSigmaY2);
778            AddToStack(i,j,GetSignal(i,j));
779            
780            //loop over different distance 
781            naccepted =1;
782            for ( Int_t dd =1;((dd<=dmax) && (naccepted>0));dd++){
783              naccepted=0; 
784              for (Int_t di = -dd;di<=dd;di++){
785                Int_t ddj = dd-abs(di);
786                Int_t sigstart = (ddj>0) ?  -1 : 0;
787                for (Int_t sig = sigstart;sig<=1;sig+=2){
788                  Int_t dj= sig*ddj; 
789                  AliCell *cell= GetCell(i+di,j+dj);
790                  Int_t signal = GetSignal(i+di,j+dj);
791                  if (cell==0) continue;
792                  Int_t index[6];
793                  index[0]=di;
794                  index[1]=dj;
795                  if (dd>2) {
796                    SetBlockIndex(index);  //adjust index to control            
797                    if ( IsBorder(fIndex,i+index[2],j+index[3]) || 
798                         IsBorder(fIndex,i+index[4],j+index[5])) {
799                      cell->SetBorder(fIndex);   
800                      continue;
801                    }
802                  }
803                  if ( signal<=fThreshold ){
804                    //if under threshold
805                    cell->SetThBorder(fIndex);
806                    if (fBFit==kTRUE)  AddToStack(i+di,j+dj,signal);
807                    continue;
808                  }
809                  naccepted++;          
810                  if (CheckIfDirBorder(fCurrentMaxX,fCurrentMaxY,i+di,j+dj) == kTRUE) {
811                    if (fBFit==kFALSE) AddToStack(i+di,j+dj,signal/2);
812                    continue; 
813                  }
814                  AddToStack(i+di,j+dj,signal);
815
816                } //loop over sig dj 
817              } //loop over di
818              
819            }//loop over dd
820          } //if there is maximum
821          //if there exists more then 2 digits cluster 
822          if (fStackIndex >2 ){
823            if (fBFit==kFALSE) GetClusterStatistic(c);
824            else GetClusterFit(c);
825            //write some important chracteristic area of cluster
826            //      
827            Transform(&c);
828            //write cluster information to array
829            TClonesArray &lclusters = *fClustersArray;
830            new(lclusters[fIndex++]) AliDigitCluster(c);
831            //             cout<<"fx="<<c.fX<<"   fy"<<c.fY<<"\n";          
832          }
833        } //lopp over all digits
834
835    return fClustersArray;
836 }
837
838
839
840
841
842
843 void AliTPCClusterFinder::Adjacent(Int_t i,Int_t j)
844 {
845   //
846   //recursive agorithm program
847   //
848   if (fBDistType==kTRUE) {
849     Float_t delta = (i-fCurrentMaxX)*(i-fCurrentMaxX)/fCurrentSigmaX2;
850     delta+=(j-fCurrentMaxY)*(j-fCurrentMaxY)/fCurrentSigmaY2;
851     if (delta > fMulSigma2) {
852        SetDirBorder(fIndex,i,j);
853       return;
854     }
855   }
856   AliCell *cell = GetCell(i,j);
857   Int_t signal = GetSignal(i,j);
858   Int_t q=signal;  
859   cell->SetChecked(fIndex);  
860   if ( (q>fThreshold) || (fBFit==kTRUE))   AddToStack(i,j,q);
861   if ( q >fThreshold )
862     {
863       
864       AliCell * newcel;      
865       newcel = GetCell(i-1,j);
866       if (newcel !=0) if (!newcel->IsChecked(fIndex) ) Adjacent(i-1,j);
867       newcel = GetCell(i,j-1);
868       if (newcel !=0) if (!newcel->IsChecked(fIndex) ) Adjacent(i,j-1);
869       newcel = GetCell(i+1,j);
870       if (newcel !=0) if (!newcel->IsChecked(fIndex) ) Adjacent(i+1,j);
871       newcel = GetCell(i,j+1);
872       if (newcel !=0) if (!newcel->IsChecked(fIndex) ) Adjacent(i,j+1);
873     }      
874   else cell->SetThBorder(fIndex);
875 }
876
877
878
879 AliH2F *  AliTPCClusterFinder::DrawHisto( const char *option=0, 
880                               Float_t x1, Float_t x2, Float_t y1, Float_t y2)
881 {
882   //
883   //draw digits in given array
884   //  
885   //make digits histo 
886   char ch[30];
887   sprintf(ch,"Cluster finder digits ");
888   if ( (fDimX<1)|| (fDimY<1)) {
889     return 0;
890   }
891   AliH2F * his  = new AliH2F(ch,ch,fDimX,fX1,fX2,fDimY,fY1,fY2);
892   //set histogram  values
893   for (Int_t i = 0; i<fDimX;i++)    
894     for (Int_t j = 0; j<fDimY;j++){
895       Float_t x = ItoX(i);
896       Float_t y= JtoY(j);
897       his->Fill(x,y,GetSignal(i,j));
898     }
899   if (x1>=0) {
900       AliH2F *h2fsub = his->GetSubrange2d(x1,x2,y1,y2);
901       delete his;
902       his=h2fsub;
903   }
904   if (his==0) return 0;
905   if (option!=0) his->Draw(option);
906   else his->Draw();
907   return his;  
908 }
909
910
911 void AliTPCClusterFinder::DrawCluster(
912                                   Int_t color, Int_t size, Int_t style)
913 {
914
915   if (fClustersArray==0) return;  
916   //draw marker for each of cluster
917   Int_t ncl=fClustersArray->GetEntriesFast();
918   for (Int_t i=0;i<ncl;i++){
919     AliCluster *cl = (AliCluster*)fClustersArray->UncheckedAt(i);
920     TMarker * marker = new TMarker;
921     marker->SetX(cl->fX);
922     marker->SetY(cl->fY);
923     marker->SetMarkerSize(size);
924     marker->SetMarkerStyle(style);
925     marker->SetMarkerColor(color);
926     marker->Draw();    
927   }
928 }
929
930
931
932 AliH2F *  AliTPCClusterFinder::DrawBorders( const char *option,  AliH2F *h, Int_t type ,
933                               Float_t x1, Float_t x2, Float_t y1, Float_t y2)
934 {
935   //
936   //draw digits in given array
937   //  
938   //make digits histo 
939   char ch[30];
940   sprintf(ch,"Cluster finder digits borders");
941   if ( (fDimX<1)|| (fDimY<1)) {
942     return 0;
943   }
944   AliH2F * his;
945   if (h!=0) his =h;
946   else his  = new AliH2F(ch,ch,fDimX,fX1,fX2,fDimY,fY1,fY2);
947   //set histogram  values
948   for (Int_t i = 0; i<fDimX;i++)    
949     for (Int_t j = 0; j<fDimY;j++){      
950       Float_t x = ItoX(i);
951       Float_t y= JtoY(j);
952       if (((type==1)||(type==0))&&IsMaximum(0,i,j)) his->Fill(x,y,16);   
953       if (((type==3)||(type==0))&&(IsDirBorder(0,i,j))) his->Fill(x,y,8);
954       if (((type==4)||(type==0))&&(IsThBorder(0,i,j))) his->Fill(x,y,4);       
955       if (((type==2)||(type==0))&&IsBorder(0,i,j)) his->Fill(x,y,1);
956
957     }
958          
959   if (x1>=0) {
960       AliH2F *h2fsub = his->GetSubrange2d(x1,x2,y1,y2);
961       delete his;
962       his=h2fsub;
963   }
964   if (his==0) return 0;
965   if (option!=0) his->Draw(option);
966   else his->Draw();
967   return his;  
968 }