]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCclustererMI.cxx
Clusterization of raw data (T.Kuhr)
[u/mrichter/AliRoot.git] / TPC / AliTPCclustererMI.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.4  2003/03/05 11:16:15  kowal2
19 Logs added
20
21 */
22
23
24
25
26 //-------------------------------------------------------
27 //          Implementation of the TPC clusterer
28 //
29 //   Origin: Marian Ivanov 
30 //-------------------------------------------------------
31
32 #include "AliTPCclustererMI.h"
33 #include "AliTPCclusterMI.h"
34 #include <TObjArray.h>
35 #include <TFile.h>
36 #include "AliTPCClustersArray.h"
37 #include "AliTPCClustersRow.h"
38 #include "AliDigits.h"
39 #include "AliSimDigits.h"
40 #include "AliTPCParam.h"
41 #include "AliTPCBuffer160.h"
42 #include "Riostream.h"
43 #include <TTree.h>
44
45 ClassImp(AliTPCclustererMI)
46
47
48
49 AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par)
50 {
51   fInput =0;
52   fOutput=0;
53   fParam = par;
54 }
55 void AliTPCclustererMI::SetInput(TTree * tree)
56 {
57   //
58   // set input tree with digits
59   //
60   fInput = tree;  
61   if  (!fInput->GetBranch("Segment")){
62     cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
63     fInput=0;
64     return;
65   }
66 }
67
68 void AliTPCclustererMI::SetOutput(TTree * tree) 
69 {
70   //
71   //
72   fOutput= tree;  
73   AliTPCClustersRow clrow;
74   AliTPCClustersRow *pclrow=&clrow;  
75   clrow.SetClass("AliTPCclusterMI");
76   clrow.SetArray(1); // to make Clones array
77   fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);    
78 }
79
80
81 Float_t  AliTPCclustererMI::GetSigmaY2(Int_t iz){
82   // sigma y2 = in digits  - we don't know the angle
83   Float_t z = iz*fParam->GetZWidth();
84   Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
85     (fPadWidth*fPadWidth);
86   Float_t sres = 0.25;
87   Float_t res = sd2+sres;
88   return res;
89 }
90
91
92 Float_t  AliTPCclustererMI::GetSigmaZ2(Int_t iz){
93   //sigma z2 = in digits - angle estimated supposing vertex constraint
94   Float_t z = iz*fZWidth;
95   Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
96   Float_t angular = fPadLength*(fParam->GetZLength()-z)/(fRx*fZWidth);
97   angular*=angular;
98   angular/=12.;
99   Float_t sres = fParam->GetZSigma()/fZWidth;
100   sres *=sres;
101   Float_t res = angular +sd2+sres;
102   return res;
103 }
104
105 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Int_t *bins, UInt_t m,
106 AliTPCclusterMI &c) 
107 {
108   Int_t i0=k/max;  //central pad
109   Int_t j0=k%max;  //central time bin
110
111   // set pointers to data
112   //Int_t dummy[5] ={0,0,0,0,0};
113   Int_t * matrix[5]; //5x5 matrix with digits  - indexing i = 0 ..4  j = -2..2
114   Int_t * resmatrix[5];
115   for (Int_t di=-2;di<=2;di++){
116     matrix[di+2]  =  &bins[k+di*max];
117     resmatrix[di+2]  =  &fResBins[k+di*max];
118   }
119   //build matrix with virtual charge
120   Float_t sigmay2= GetSigmaY2(j0);
121   Float_t sigmaz2= GetSigmaZ2(j0);
122
123   Float_t vmatrix[5][5];
124   vmatrix[2][2] = matrix[2][0];
125   c.SetType(0);
126   c.SetMax(Short_t(vmatrix[2][2])); // write maximal amplitude
127   for (Int_t di =-1;di <=1;di++)
128     for (Int_t dj =-1;dj <=1;dj++){
129       Float_t amp = matrix[di+2][dj];
130       if ( (amp<2) && (fLoop<2)){
131         // if under threshold  - calculate virtual charge
132         Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
133         amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
134         if (amp>2) amp = 2;
135         vmatrix[2+di][2+dj]=amp;
136         vmatrix[2+2*di][2+2*dj]=0;
137         if ( (di*dj)!=0){       
138           //DIAGONAL ELEMENTS
139           vmatrix[2+2*di][2+dj] =0;
140           vmatrix[2+di][2+2*dj] =0;
141         }
142         continue;
143       }
144       if (amp<4){
145         //if small  amplitude - below  2 x threshold  - don't consider other one        
146         vmatrix[2+di][2+dj]=amp;
147         vmatrix[2+2*di][2+2*dj]=0;  // don't take to the account next bin
148         if ( (di*dj)!=0){       
149           //DIAGONAL ELEMENTS
150           vmatrix[2+2*di][2+dj] =0;
151           vmatrix[2+di][2+2*dj] =0;
152         }
153         continue;
154       }
155       //if bigger then take everything
156       vmatrix[2+di][2+dj]=amp;
157       vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;      
158       if ( (di*dj)!=0){       
159           //DIAGONAL ELEMENTS
160           vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
161           vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
162         }      
163     }
164
165
166   
167   Float_t sumw=0;
168   Float_t sumiw=0;
169   Float_t sumi2w=0;
170   Float_t sumjw=0;
171   Float_t sumj2w=0;
172   //
173   for (Int_t i=-2;i<=2;i++)
174     for (Int_t j=-2;j<=2;j++){
175       Float_t amp = vmatrix[i+2][j+2];
176
177       sumw    += amp;
178       sumiw   += i*amp;
179       sumi2w  += i*i*amp;
180       sumjw   += j*amp;
181       sumj2w  += j*j*amp;
182     }    
183   //   
184   Float_t meani = sumiw/sumw;
185   Float_t mi2   = sumi2w/sumw-meani*meani;
186   Float_t meanj = sumjw/sumw;
187   Float_t mj2   = sumj2w/sumw-meanj*meanj;
188   //
189   Float_t ry = mi2/sigmay2;
190   Float_t rz = mj2/sigmaz2;
191   
192   //
193   if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
194   if ( (ry <1.2) && (rz<1.2) ) {
195     //if cluster looks like expected 
196     //+1.2 deviation from expected sigma accepted
197     //    c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
198
199     meani +=i0;
200     meanj +=j0;
201     //set cluster parameters
202     c.SetQ(sumw);
203     c.SetY(meani*fPadWidth); 
204     c.SetZ(meanj*fZWidth); 
205     c.SetSigmaY2(mi2);
206     c.SetSigmaZ2(mj2);
207     AddCluster(c);
208     //remove cluster data from data
209     for (Int_t di=-2;di<=2;di++)
210       for (Int_t dj=-2;dj<=2;dj++){
211         resmatrix[di+2][dj] -= Int_t(vmatrix[di+2][dj+2]);
212         if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
213       }
214     resmatrix[2][0] =0;
215     return;     
216   }
217   //
218   //unfolding when neccessary  
219   //
220   
221   Int_t * matrix2[7]; //7x7 matrix with digits  - indexing i = 0 ..6  j = -3..3
222   Int_t dummy[7]={0,0,0,0,0,0};
223   for (Int_t di=-3;di<=3;di++){
224     matrix2[di+3] =  &bins[k+di*max];
225     if ((k+di*max)<3)  matrix2[di+3] = &dummy[3];
226     if ((k+di*max)>fMaxBin-3)  matrix2[di+3] = &dummy[3];
227   }
228   Float_t vmatrix2[5][5];
229   Float_t sumu;
230   Float_t overlap;
231   UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
232   //
233   //  c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
234   meani +=i0;
235   meanj +=j0;
236   //set cluster parameters
237   c.SetQ(sumu);
238   c.SetY(meani*fPadWidth); 
239   c.SetZ(meanj*fZWidth); 
240   c.SetSigmaY2(mi2);
241   c.SetSigmaZ2(mj2);
242   c.SetType(Char_t(overlap)+1);
243   AddCluster(c);
244
245   //unfolding 2
246   meani-=i0;
247   meanj-=j0;
248   if (gDebug>4)
249     printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
250 }
251
252
253
254 void AliTPCclustererMI::UnfoldCluster(Int_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj, 
255                                       Float_t & sumu, Float_t & overlap )
256 {
257   //
258   //unfold cluster from input matrix
259   //data corresponding to cluster writen in recmatrix
260   //output meani and meanj
261
262   //take separatelly y and z
263
264   Float_t sum3i[7] = {0,0,0,0,0,0,0};
265   Float_t sum3j[7] = {0,0,0,0,0,0,0};
266
267   for (Int_t k =0;k<7;k++)
268     for (Int_t l = -1; l<=1;l++){
269       sum3i[k]+=matrix2[k][l];
270       sum3j[k]+=matrix2[l+3][k-3];
271     }
272   Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
273   //
274   //unfold  y 
275   Float_t sum3wi    = 0;  //charge minus overlap
276   Float_t sum3wio   = 0;  //full charge
277   Float_t sum3iw    = 0;  //sum for mean value
278   for (Int_t dk=-1;dk<=1;dk++){
279     sum3wio+=sum3i[dk+3];
280     if (dk==0){
281       sum3wi+=sum3i[dk+3];     
282     }
283     else{
284       Float_t ratio =1;
285       if (  ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
286            sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
287         Float_t xm2 = sum3i[-dk+3];
288         Float_t xm1 = sum3i[+3];
289         Float_t x1  = sum3i[2*dk+3];
290         Float_t x2  = sum3i[3*dk+3];    
291         Float_t w11   = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);      
292         Float_t w12   = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
293         ratio = w11/(w11+w12);   
294         for (Int_t dl=-1;dl<=1;dl++)
295           mratio[dk+1][dl+1] *= ratio;
296       }
297       Float_t amp = sum3i[dk+3]*ratio;
298       sum3wi+=amp;
299       sum3iw+= dk*amp;      
300     }
301   }
302   meani = sum3iw/sum3wi;
303   Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
304
305
306
307   //unfold  z 
308   Float_t sum3wj    = 0;  //charge minus overlap
309   Float_t sum3wjo   = 0;  //full charge
310   Float_t sum3jw    = 0;  //sum for mean value
311   for (Int_t dk=-1;dk<=1;dk++){
312     sum3wjo+=sum3j[dk+3];
313     if (dk==0){
314       sum3wj+=sum3j[dk+3];     
315     }
316     else{
317       Float_t ratio =1;
318       if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
319            (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
320         Float_t xm2 = sum3j[-dk+3];
321         Float_t xm1 = sum3j[+3];
322         Float_t x1  = sum3j[2*dk+3];
323         Float_t x2  = sum3j[3*dk+3];    
324         Float_t w11   = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);      
325         Float_t w12   = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
326         ratio = w11/(w11+w12);   
327         for (Int_t dl=-1;dl<=1;dl++)
328           mratio[dl+1][dk+1] *= ratio;
329       }
330       Float_t amp = sum3j[dk+3]*ratio;
331       sum3wj+=amp;
332       sum3jw+= dk*amp;      
333     }
334   }
335   meanj = sum3jw/sum3wj;
336   Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;  
337   overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);  
338   sumu = (sum3wj+sum3wi)/2.;
339   
340   if (overlap ==3) {
341     //if not overlap detected remove everything
342     for (Int_t di =-2; di<=2;di++)
343       for (Int_t dj =-2; dj<=2;dj++){
344         recmatrix[di+2][dj+2] = matrix2[3+di][dj];
345       }
346   }
347   else{
348     for (Int_t di =-1; di<=1;di++)
349       for (Int_t dj =-1; dj<=1;dj++){
350         Float_t ratio =1;
351         if (mratio[di+1][dj+1]==1){
352           recmatrix[di+2][dj+2]     = matrix2[3+di][dj];
353           if (TMath::Abs(di)+TMath::Abs(dj)>1){
354             recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
355             recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
356           }       
357           recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
358         }
359         else
360           {
361             //if we have overlap in direction
362             recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];    
363             if (TMath::Abs(di)+TMath::Abs(dj)>1){
364               ratio =  TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
365               recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
366               //
367               ratio =  TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
368               recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
369             }
370             else{
371               ratio =  recmatrix[di+2][dj+2]/matrix2[3][0];
372               recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
373             }
374           }
375       }
376   }
377   if (gDebug>4) 
378     printf("%f\n", recmatrix[2][2]);
379   
380 }
381
382 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
383 {
384   //
385   // estimate max
386   Float_t sumteor= 0;
387   Float_t sumamp = 0;
388
389   for (Int_t di = -1;di<=1;di++)
390     for (Int_t dj = -1;dj<=1;dj++){
391       if (vmatrix[2+di][2+dj]>2){
392         Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
393         sumteor += teor*vmatrix[2+di][2+dj];
394         sumamp  += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
395       }
396     }    
397   Float_t max = sumamp/sumteor;
398   return max;
399 }
400
401 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c){
402   //
403   // transform cluster to the global coordinata
404   // add the cluster to the array
405   //
406   Float_t meani = c.GetY()/fPadWidth;
407   Float_t meanj = c.GetZ()/fZWidth;
408
409   Int_t ki = TMath::Nint(meani-3);
410   if (ki<0) ki=0;
411   if (ki>=fMaxPad) ki = fMaxPad-1;
412   Int_t kj = TMath::Nint(meanj-3);
413   if (kj<0) kj=0;
414   if (kj>=fMaxTime-3) kj=fMaxTime-4;
415   // ki and kj shifted to "real" coordinata
416   if (fRowDig) {
417     c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
418     c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
419     c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
420   }
421   
422   
423   Float_t s2 = c.GetSigmaY2();
424   Float_t w=fParam->GetPadPitchWidth(fSector);
425   
426   c.SetSigmaY2(s2*w*w);
427   s2 = c.GetSigmaZ2(); 
428   w=fZWidth;
429   c.SetSigmaZ2(s2*w*w);
430   c.SetY((meani - 2.5 - 0.5*fMaxPad)*fParam->GetPadPitchWidth(fSector));
431   c.SetZ(fZWidth*(meanj-3)); 
432   c.SetZ(c.GetZ()   - 3.*fParam->GetZSigma()); // PASA delay 
433   c.SetZ(fSign*(fParam->GetZLength() - c.GetZ()));
434   
435   if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
436     //c.SetSigmaY2(c.GetSigmaY2()*25.);
437     //c.SetSigmaZ2(c.GetSigmaZ2()*4.);
438     c.SetType(-(c.GetType()+3));  //edge clusters
439   }
440   if (fLoop==2) c.SetType(100);
441
442   TClonesArray * arr = fRowCl->GetArray();
443   // AliTPCclusterMI * cl = 
444   new ((*arr)[fNcluster]) AliTPCclusterMI(c);
445
446   fNcluster++;
447 }
448
449
450 //_____________________________________________________________________________
451 void AliTPCclustererMI::Digits2Clusters(const AliTPCParam *par, Int_t eventn)
452 {
453   //-----------------------------------------------------------------
454   // This is a simple cluster finder.
455   //-----------------------------------------------------------------
456   TDirectory *savedir=gDirectory; 
457
458  if (fInput==0){ 
459     cerr<<"AliTPC::Digits2Clusters(): input tree not initialised !\n";
460     return;
461   }
462  
463   if (fOutput==0) {
464      cerr<<"AliTPC::Digits2Clusters(): output tree not initialised !\n";
465      return;
466   }
467
468   AliSimDigits digarr, *dummy=&digarr;
469   fRowDig = dummy;
470   fInput->GetBranch("Segment")->SetAddress(&dummy);
471   Stat_t nentries = fInput->GetEntries();
472   
473   fMaxTime=par->GetMaxTBin()+6; // add 3 virtual time bins before and 3   after
474     
475   fParam = par;
476   ((AliTPCParam*)par)->Write(par->GetTitle());
477
478   Int_t nclusters  = 0;
479   
480   for (Int_t n=0; n<nentries; n++) {
481     fInput->GetEvent(n);
482     Int_t row;
483     if (!par->AdjustSectorRow(digarr.GetID(),fSector,row)) {
484       cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
485       continue;
486     }
487         
488     AliTPCClustersRow *clrow= new AliTPCClustersRow();
489     fRowCl = clrow;
490     clrow->SetClass("AliTPCclusterMI");
491     clrow->SetArray(1);
492
493     clrow->SetID(digarr.GetID());
494     fOutput->GetBranch("Segment")->SetAddress(&clrow);
495     fRx=par->GetPadRowRadii(fSector,row);
496     
497     
498     const Int_t kNIS=par->GetNInnerSector(), kNOS=par->GetNOuterSector();
499     fZWidth = fParam->GetZWidth();
500     if (fSector < kNIS) {
501       fMaxPad = par->GetNPadsLow(row);
502       fSign = (fSector < kNIS/2) ? 1 : -1;
503       fPadLength = par->GetPadPitchLength(fSector,row);
504       fPadWidth = par->GetPadPitchWidth();
505     } else {
506       fMaxPad = par->GetNPadsUp(row);
507       fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
508       fPadLength = par->GetPadPitchLength(fSector,row);
509       fPadWidth  = par->GetPadPitchWidth();
510     }
511     
512     
513     fMaxBin=fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
514     fBins    =new Int_t[fMaxBin];
515     fResBins =new Int_t[fMaxBin];  //fBins with residuals after 1 finder loop 
516     memset(fBins,0,sizeof(Int_t)*fMaxBin);
517     
518     if (digarr.First()) //MI change
519       do {
520         Short_t dig=digarr.CurrentDigit();
521         if (dig<=par->GetZeroSup()) continue;
522         Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
523         fBins[i*fMaxTime+j]=dig;
524       } while (digarr.Next());
525     digarr.ExpandTrackBuffer();
526
527     FindClusters();
528
529     fOutput->Fill();
530     delete clrow;    
531     nclusters+=fNcluster;    
532     delete[] fBins;      
533     delete[] fResBins;  
534   }  
535   cerr<<"Number of found clusters : "<<nclusters<<"                        \n";
536   
537   fOutput->Write();
538   savedir->cd();
539 }
540
541 void AliTPCclustererMI::Digits2Clusters(const char* fileName)
542 {
543 //-----------------------------------------------------------------
544 // This is a cluster finder for raw data.
545 //
546 // It assumes, that the input file contains uncompressed data 
547 // in TPC Altro format without mini headers.
548 // It also assumes, that the data is sorted in a such a way that
549 // data from different pad rows is not mixed. The order of pads
550 // in a pad row can be random.
551 //-----------------------------------------------------------------
552
553   TDirectory* savedir = gDirectory; 
554
555   AliTPCBuffer160 buffer(fileName, 0);  // buffer for reading raw data
556   fRowDig = NULL;
557   if (!fOutput) {
558     Error("Digits2Clusters", "output tree not initialised !\n");
559     return;
560   }
561   fOutput->GetCurrentFile()->cd();
562   const_cast<AliTPCParam*>(fParam)->Write();
563
564   Int_t nclusters  = 0;
565   
566   fMaxTime = fParam->GetMaxTBin() + 6; // add 3 virtual time bins before and 3 after
567   const Int_t kNIS = fParam->GetNInnerSector();
568   const Int_t kNOS = fParam->GetNOuterSector();
569   fZWidth = fParam->GetZWidth();
570   Int_t zeroSup = fParam->GetZeroSup();
571   Int_t offset = 1;
572
573   Int_t nWords, iPad, iRow;
574   Int_t prevSector = -1;
575   Int_t prevRow = -1;
576   fBins = NULL;
577   while (buffer.ReadTrailerBackward(nWords, iPad, iRow, fSector) == 0) {
578     Int_t nSkip = (4 - (nWords % 4)) % 4;
579
580     // when the sector or row number has changed ...
581     if ((fSector != prevSector) || (iRow != prevRow)) {
582
583       // ... find clusters in the previous pad row, and ...
584       if (fBins) {
585         FindClusters();
586
587         fOutput->Fill();
588         delete fRowCl;    
589         nclusters += fNcluster;    
590         delete[] fBins;      
591         delete[] fResBins;  
592       }
593
594       // ... prepare for the next pad row
595       fRowCl = new AliTPCClustersRow;
596       fRowCl->SetClass("AliTPCclusterMI");
597       fRowCl->SetArray(1);
598
599       fRowCl->SetID(fParam->GetIndex(fSector, iRow));
600       fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
601       fRx = fParam->GetPadRowRadii(fSector, iRow);
602     
603       if (fSector < kNIS) {
604         fMaxPad = fParam->GetNPadsLow(iRow);
605         fSign = (fSector < kNIS/2) ? 1 : -1;
606       } else {
607         fMaxPad = fParam->GetNPadsUp(iRow);
608         fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
609       }
610       fPadLength = fParam->GetPadPitchLength(fSector, iRow);
611       fPadWidth  = fParam->GetPadPitchWidth();
612     
613       fMaxBin = fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
614       fBins    = new Int_t[fMaxBin];
615       fResBins = new Int_t[fMaxBin];  //fBins with residuals after 1 finder loop 
616       memset(fBins, 0, sizeof(Int_t)*fMaxBin);
617
618       prevSector = fSector;
619       prevRow = iRow;
620     }
621
622     // fill fBins with digits data
623     for (Int_t iWord = 0; iWord < nSkip; iWord++) {
624       if (buffer.GetNextBackWord() != 0x2AA) {
625         Error("Digits2Clusters", "could not read dummy data\n");
626         break;
627       }
628     }
629     while (nWords > 0) {
630       Int_t bunchLength = buffer.GetNextBackWord() - 2;
631       if (bunchLength < 0) {
632         Error("Digits2Clusters", "could not read bunch length\n");
633         break;
634       }
635       Int_t timeBin = buffer.GetNextBackWord();
636       if (timeBin < 0) {
637         Error("Digits2Clusters", "could not read time bin\n");
638         break;
639       }
640       for (Int_t iTimeBin = timeBin; iTimeBin > timeBin-bunchLength; iTimeBin--) {
641         Int_t dig = buffer.GetNextBackWord();
642         if (dig < 0) {
643           Error("Digits2Clusters", "could not read sample amplitude\n");
644           break;
645         }
646         dig += offset;
647         if (dig <= zeroSup) continue;
648         Int_t i = iPad + 3;
649         Int_t j = iTimeBin + 3;
650         fBins[i*fMaxTime+j] = dig;
651       }
652       nWords -= 2+bunchLength;
653     }
654   }  
655
656   // find clusters in the last pad row
657   if (fBins) {
658     FindClusters();
659
660     fOutput->Fill();
661     delete fRowCl;    
662     nclusters += fNcluster;    
663     delete[] fBins;      
664     delete[] fResBins;  
665   }
666
667   Info("Digits2Clusters", "Number of found clusters : %d\n", nclusters);
668
669   fOutput->Write();
670   savedir->cd();
671 }
672
673 void AliTPCclustererMI::FindClusters()
674 {
675   //add virtual charge at the edge   
676   for (Int_t i=0; i<fMaxTime; i++){
677     Float_t amp1 = fBins[i+3*fMaxTime]; 
678     Float_t amp0 =0;
679     if (amp1>0){
680       Float_t amp2 = fBins[i+4*fMaxTime];
681       if (amp2==0) amp2=0.5;
682       Float_t sigma2 = GetSigmaY2(i);           
683       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);   
684       if (gDebug>4) printf("\n%f\n",amp0);
685     }  
686     fBins[i+2*fMaxTime] = Int_t(amp0);    
687     amp0 = 0;
688     amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
689     if (amp1>0){
690       Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
691       if (amp2==0) amp2=0.5;
692       Float_t sigma2 = GetSigmaY2(i);           
693       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);           
694       if (gDebug>4) printf("\n%f\n",amp0);
695     }        
696     fBins[(fMaxPad+3)*fMaxTime+i] = Int_t(amp0);           
697   }
698
699 //  memcpy(fResBins,fBins, fMaxBin*2);
700   memcpy(fResBins,fBins, fMaxBin);
701   //
702   fNcluster=0;
703   //first loop - for "gold cluster" 
704   fLoop=1;
705   Int_t *b=&fBins[-1]+2*fMaxTime;
706   Int_t crtime = Int_t((fParam->GetZLength()-1.05*fRx)/fZWidth-5);
707
708   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
709     b++;
710     if (*b<8) continue;   //threshold form maxima
711     if (i%fMaxTime<crtime) {
712       Int_t delta = -(i%fMaxTime)+crtime;
713       b+=delta;
714       i+=delta;
715       continue; 
716     }
717      
718     if (!IsMaximum(*b,fMaxTime,b)) continue;
719     AliTPCclusterMI c;
720     Int_t dummy=0;
721     MakeCluster(i, fMaxTime, fBins, dummy,c);
722     //}
723   }
724   //memcpy(fBins,fResBins, fMaxBin*2);
725   //second  loop - for rest cluster 
726   /*        
727   fLoop=2;
728   b=&fResBins[-1]+2*fMaxTime;
729   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
730     b++;
731     if (*b<25) continue;   // bigger threshold for maxima
732     if (!IsMaximum(*b,fMaxTime,b)) continue;
733     AliTPCclusterMI c;
734     Int_t dummy;
735     MakeCluster(i, fMaxTime, fResBins, dummy,c);
736     //}
737   }
738   */
739 }