]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCclustererMI.cxx
Improved fitting of gating pulse (Marian)
[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 /* $Id$ */
17
18 //-------------------------------------------------------
19 //          Implementation of the TPC clusterer
20 //
21 //   Origin: Marian Ivanov 
22 //-------------------------------------------------------
23
24 #include "AliTPCReconstructor.h"
25 #include "AliTPCclustererMI.h"
26 #include "AliTPCclusterMI.h"
27 #include <TObjArray.h>
28 #include <TFile.h>
29 #include "TGraph.h"
30 #include "TF1.h"
31 #include "TRandom.h"
32 #include "AliMathBase.h"
33
34 #include "AliTPCClustersArray.h"
35 #include "AliTPCClustersRow.h"
36 #include "AliDigits.h"
37 #include "AliSimDigits.h"
38 #include "AliTPCParam.h"
39 #include "AliTPCRecoParam.h"
40 #include "AliRawReader.h"
41 #include "AliTPCRawStream.h"
42 #include "AliRunLoader.h"
43 #include "AliLoader.h"
44 #include "Riostream.h"
45 #include <TTree.h>
46 #include "AliTPCcalibDB.h"
47 #include "AliTPCCalPad.h"
48 #include "AliTPCCalROC.h"
49 #include "TTreeStream.h"
50 #include "AliLog.h"
51
52
53 ClassImp(AliTPCclustererMI)
54
55
56
57   AliTPCclustererMI::AliTPCclustererMI(const AliTPCParam* par, const AliTPCRecoParam * recoParam)
58 {
59   fIsOldRCUFormat = kFALSE;
60   fInput =0;
61   fOutput=0;
62   fParam = par;
63   if (recoParam) {
64     fRecoParam = recoParam;
65   }else{
66     //set default parameters if not specified
67     fRecoParam = AliTPCReconstructor::GetRecoParam();
68     if (!fRecoParam)  fRecoParam = AliTPCRecoParam::GetLowFluxParam();
69   }
70   fDebugStreamer = new TTreeSRedirector("TPCsignal.root");
71   fAmplitudeHisto = 0;
72 }
73
74 AliTPCclustererMI::~AliTPCclustererMI(){
75   DumpHistos();
76   if (fAmplitudeHisto) delete fAmplitudeHisto;
77   if (fDebugStreamer) delete fDebugStreamer;
78 }
79
80 void AliTPCclustererMI::SetInput(TTree * tree)
81 {
82   //
83   // set input tree with digits
84   //
85   fInput = tree;  
86   if  (!fInput->GetBranch("Segment")){
87     cerr<<"AliTPC::Digits2Clusters(): no porper input tree !\n";
88     fInput=0;
89     return;
90   }
91 }
92
93 void AliTPCclustererMI::SetOutput(TTree * tree) 
94 {
95   //
96   //
97   fOutput= tree;  
98   AliTPCClustersRow clrow;
99   AliTPCClustersRow *pclrow=&clrow;  
100   clrow.SetClass("AliTPCclusterMI");
101   clrow.SetArray(1); // to make Clones array
102   fOutput->Branch("Segment","AliTPCClustersRow",&pclrow,32000,200);    
103 }
104
105
106 Float_t  AliTPCclustererMI::GetSigmaY2(Int_t iz){
107   // sigma y2 = in digits  - we don't know the angle
108   Float_t z = iz*fParam->GetZWidth()+fParam->GetNTBinsL1()*fParam->GetZWidth();
109   Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/
110     (fPadWidth*fPadWidth);
111   Float_t sres = 0.25;
112   Float_t res = sd2+sres;
113   return res;
114 }
115
116
117 Float_t  AliTPCclustererMI::GetSigmaZ2(Int_t iz){
118   //sigma z2 = in digits - angle estimated supposing vertex constraint
119   Float_t z = iz*fZWidth+fParam->GetNTBinsL1()*fParam->GetZWidth();
120   Float_t sd2 = (z*fParam->GetDiffL()*fParam->GetDiffL())/(fZWidth*fZWidth);
121   Float_t angular = fPadLength*(fParam->GetZLength()-z)/(fRx*fZWidth);
122   angular*=angular;
123   angular/=12.;
124   Float_t sres = fParam->GetZSigma()/fZWidth;
125   sres *=sres;
126   Float_t res = angular +sd2+sres;
127   return res;
128 }
129
130 void AliTPCclustererMI::MakeCluster(Int_t k,Int_t max,Float_t *bins, UInt_t /*m*/,
131 AliTPCclusterMI &c) 
132 {
133   Int_t i0=k/max;  //central pad
134   Int_t j0=k%max;  //central time bin
135
136   // set pointers to data
137   //Int_t dummy[5] ={0,0,0,0,0};
138   Float_t * matrix[5]; //5x5 matrix with digits  - indexing i = 0 ..4  j = -2..2
139   Float_t * resmatrix[5];
140   for (Int_t di=-2;di<=2;di++){
141     matrix[di+2]  =  &bins[k+di*max];
142     resmatrix[di+2]  =  &fResBins[k+di*max];
143   }
144   //build matrix with virtual charge
145   Float_t sigmay2= GetSigmaY2(j0);
146   Float_t sigmaz2= GetSigmaZ2(j0);
147
148   Float_t vmatrix[5][5];
149   vmatrix[2][2] = matrix[2][0];
150   c.SetType(0);
151   c.SetMax((UShort_t)(vmatrix[2][2])); // write maximal amplitude
152   for (Int_t di =-1;di <=1;di++)
153     for (Int_t dj =-1;dj <=1;dj++){
154       Float_t amp = matrix[di+2][dj];
155       if ( (amp<2) && (fLoop<2)){
156         // if under threshold  - calculate virtual charge
157         Float_t ratio = TMath::Exp(-1.2*TMath::Abs(di)/sigmay2)*TMath::Exp(-1.2*TMath::Abs(dj)/sigmaz2);
158         amp = ((matrix[2][0]-2)*(matrix[2][0]-2)/(matrix[-di+2][-dj]+2))*ratio;
159         if (amp>2) amp = 2;
160         vmatrix[2+di][2+dj]=amp;
161         vmatrix[2+2*di][2+2*dj]=0;
162         if ( (di*dj)!=0){       
163           //DIAGONAL ELEMENTS
164           vmatrix[2+2*di][2+dj] =0;
165           vmatrix[2+di][2+2*dj] =0;
166         }
167         continue;
168       }
169       if (amp<4){
170         //if small  amplitude - below  2 x threshold  - don't consider other one        
171         vmatrix[2+di][2+dj]=amp;
172         vmatrix[2+2*di][2+2*dj]=0;  // don't take to the account next bin
173         if ( (di*dj)!=0){       
174           //DIAGONAL ELEMENTS
175           vmatrix[2+2*di][2+dj] =0;
176           vmatrix[2+di][2+2*dj] =0;
177         }
178         continue;
179       }
180       //if bigger then take everything
181       vmatrix[2+di][2+dj]=amp;
182       vmatrix[2+2*di][2+2*dj]= matrix[2*di+2][2*dj] ;      
183       if ( (di*dj)!=0){       
184           //DIAGONAL ELEMENTS
185           vmatrix[2+2*di][2+dj] = matrix[2*di+2][dj];
186           vmatrix[2+di][2+2*dj] = matrix[2+di][dj*2];
187         }      
188     }
189
190
191   
192   Float_t sumw=0;
193   Float_t sumiw=0;
194   Float_t sumi2w=0;
195   Float_t sumjw=0;
196   Float_t sumj2w=0;
197   //
198   for (Int_t i=-2;i<=2;i++)
199     for (Int_t j=-2;j<=2;j++){
200       Float_t amp = vmatrix[i+2][j+2];
201
202       sumw    += amp;
203       sumiw   += i*amp;
204       sumi2w  += i*i*amp;
205       sumjw   += j*amp;
206       sumj2w  += j*j*amp;
207     }    
208   //   
209   Float_t meani = sumiw/sumw;
210   Float_t mi2   = sumi2w/sumw-meani*meani;
211   Float_t meanj = sumjw/sumw;
212   Float_t mj2   = sumj2w/sumw-meanj*meanj;
213   //
214   Float_t ry = mi2/sigmay2;
215   Float_t rz = mj2/sigmaz2;
216   
217   //
218   if ( ( (ry<0.6) || (rz<0.6) ) && fLoop==2) return;
219   if ( (ry <1.2) && (rz<1.2) || (!fRecoParam->GetDoUnfold())) {
220     //
221     //if cluster looks like expected or Unfolding not switched on
222     //standard COG is used
223     //+1.2 deviation from expected sigma accepted
224     //    c.fMax = FitMax(vmatrix,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
225
226     meani +=i0;
227     meanj +=j0;
228     //set cluster parameters
229     c.SetQ(sumw);
230     c.SetY(meani*fPadWidth); 
231     c.SetZ(meanj*fZWidth); 
232     c.SetSigmaY2(mi2);
233     c.SetSigmaZ2(mj2);
234     AddCluster(c);
235     //remove cluster data from data
236     for (Int_t di=-2;di<=2;di++)
237       for (Int_t dj=-2;dj<=2;dj++){
238         resmatrix[di+2][dj] -= vmatrix[di+2][dj+2];
239         if (resmatrix[di+2][dj]<0) resmatrix[di+2][dj]=0;
240       }
241     resmatrix[2][0] =0;
242     return;     
243   }
244   //
245   //unfolding when neccessary  
246   //
247   
248   Float_t * matrix2[7]; //7x7 matrix with digits  - indexing i = 0 ..6  j = -3..3
249   Float_t dummy[7]={0,0,0,0,0,0};
250   for (Int_t di=-3;di<=3;di++){
251     matrix2[di+3] =  &bins[k+di*max];
252     if ((k+di*max)<3)  matrix2[di+3] = &dummy[3];
253     if ((k+di*max)>fMaxBin-3)  matrix2[di+3] = &dummy[3];
254   }
255   Float_t vmatrix2[5][5];
256   Float_t sumu;
257   Float_t overlap;
258   UnfoldCluster(matrix2,vmatrix2,meani,meanj,sumu,overlap);
259   //
260   //  c.fMax = FitMax(vmatrix2,meani,meanj,TMath::Sqrt(sigmay2),TMath::Sqrt(sigmaz2));
261   meani +=i0;
262   meanj +=j0;
263   //set cluster parameters
264   c.SetQ(sumu);
265   c.SetY(meani*fPadWidth); 
266   c.SetZ(meanj*fZWidth); 
267   c.SetSigmaY2(mi2);
268   c.SetSigmaZ2(mj2);
269   c.SetType(Char_t(overlap)+1);
270   AddCluster(c);
271
272   //unfolding 2
273   meani-=i0;
274   meanj-=j0;
275   if (gDebug>4)
276     printf("%f\t%f\n", vmatrix2[2][2], vmatrix[2][2]);
277 }
278
279
280
281 void AliTPCclustererMI::UnfoldCluster(Float_t * matrix2[7], Float_t recmatrix[5][5], Float_t & meani, Float_t & meanj, 
282                                       Float_t & sumu, Float_t & overlap )
283 {
284   //
285   //unfold cluster from input matrix
286   //data corresponding to cluster writen in recmatrix
287   //output meani and meanj
288
289   //take separatelly y and z
290
291   Float_t sum3i[7] = {0,0,0,0,0,0,0};
292   Float_t sum3j[7] = {0,0,0,0,0,0,0};
293
294   for (Int_t k =0;k<7;k++)
295     for (Int_t l = -1; l<=1;l++){
296       sum3i[k]+=matrix2[k][l];
297       sum3j[k]+=matrix2[l+3][k-3];
298     }
299   Float_t mratio[3][3]={{1,1,1},{1,1,1},{1,1,1}};
300   //
301   //unfold  y 
302   Float_t sum3wi    = 0;  //charge minus overlap
303   Float_t sum3wio   = 0;  //full charge
304   Float_t sum3iw    = 0;  //sum for mean value
305   for (Int_t dk=-1;dk<=1;dk++){
306     sum3wio+=sum3i[dk+3];
307     if (dk==0){
308       sum3wi+=sum3i[dk+3];     
309     }
310     else{
311       Float_t ratio =1;
312       if (  ( ((sum3i[dk+3]+3)/(sum3i[3]-3))+1 < (sum3i[2*dk+3]-3)/(sum3i[dk+3]+3))||
313            sum3i[dk+3]<=sum3i[2*dk+3] && sum3i[dk+3]>2 ){
314         Float_t xm2 = sum3i[-dk+3];
315         Float_t xm1 = sum3i[+3];
316         Float_t x1  = sum3i[2*dk+3];
317         Float_t x2  = sum3i[3*dk+3];    
318         Float_t w11   = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);      
319         Float_t w12   = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
320         ratio = w11/(w11+w12);   
321         for (Int_t dl=-1;dl<=1;dl++)
322           mratio[dk+1][dl+1] *= ratio;
323       }
324       Float_t amp = sum3i[dk+3]*ratio;
325       sum3wi+=amp;
326       sum3iw+= dk*amp;      
327     }
328   }
329   meani = sum3iw/sum3wi;
330   Float_t overlapi = (sum3wio-sum3wi)/sum3wio;
331
332
333
334   //unfold  z 
335   Float_t sum3wj    = 0;  //charge minus overlap
336   Float_t sum3wjo   = 0;  //full charge
337   Float_t sum3jw    = 0;  //sum for mean value
338   for (Int_t dk=-1;dk<=1;dk++){
339     sum3wjo+=sum3j[dk+3];
340     if (dk==0){
341       sum3wj+=sum3j[dk+3];     
342     }
343     else{
344       Float_t ratio =1;
345       if ( ( ((sum3j[dk+3]+3)/(sum3j[3]-3))+1 < (sum3j[2*dk+3]-3)/(sum3j[dk+3]+3)) ||
346            (sum3j[dk+3]<=sum3j[2*dk+3] && sum3j[dk+3]>2)){
347         Float_t xm2 = sum3j[-dk+3];
348         Float_t xm1 = sum3j[+3];
349         Float_t x1  = sum3j[2*dk+3];
350         Float_t x2  = sum3j[3*dk+3];    
351         Float_t w11   = TMath::Max((Float_t)(4.*xm1-xm2),(Float_t)0.000001);      
352         Float_t w12   = TMath::Max((Float_t)(4 *x1 -x2),(Float_t)0.);
353         ratio = w11/(w11+w12);   
354         for (Int_t dl=-1;dl<=1;dl++)
355           mratio[dl+1][dk+1] *= ratio;
356       }
357       Float_t amp = sum3j[dk+3]*ratio;
358       sum3wj+=amp;
359       sum3jw+= dk*amp;      
360     }
361   }
362   meanj = sum3jw/sum3wj;
363   Float_t overlapj = (sum3wjo-sum3wj)/sum3wjo;  
364   overlap = Int_t(100*TMath::Max(overlapi,overlapj)+3);  
365   sumu = (sum3wj+sum3wi)/2.;
366   
367   if (overlap ==3) {
368     //if not overlap detected remove everything
369     for (Int_t di =-2; di<=2;di++)
370       for (Int_t dj =-2; dj<=2;dj++){
371         recmatrix[di+2][dj+2] = matrix2[3+di][dj];
372       }
373   }
374   else{
375     for (Int_t di =-1; di<=1;di++)
376       for (Int_t dj =-1; dj<=1;dj++){
377         Float_t ratio =1;
378         if (mratio[di+1][dj+1]==1){
379           recmatrix[di+2][dj+2]     = matrix2[3+di][dj];
380           if (TMath::Abs(di)+TMath::Abs(dj)>1){
381             recmatrix[2*di+2][dj+2] = matrix2[3+2*di][dj];
382             recmatrix[di+2][2*dj+2] = matrix2[3+di][2*dj];
383           }       
384           recmatrix[2*di+2][2*dj+2] = matrix2[3+2*di][2*dj];
385         }
386         else
387           {
388             //if we have overlap in direction
389             recmatrix[di+2][dj+2] = mratio[di+1][dj+1]* matrix2[3+di][dj];    
390             if (TMath::Abs(di)+TMath::Abs(dj)>1){
391               ratio =  TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+0*di][1*dj]+1)),(Float_t)1.);
392               recmatrix[2*di+2][dj+2] = ratio*recmatrix[di+2][dj+2];
393               //
394               ratio =  TMath::Min((Float_t)(recmatrix[di+2][dj+2]/(matrix2[3+1*di][0*dj]+1)),(Float_t)1.);
395               recmatrix[di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
396             }
397             else{
398               ratio =  recmatrix[di+2][dj+2]/matrix2[3][0];
399               recmatrix[2*di+2][2*dj+2] = ratio*recmatrix[di+2][dj+2];
400             }
401           }
402       }
403   }
404   if (gDebug>4) 
405     printf("%f\n", recmatrix[2][2]);
406   
407 }
408
409 Float_t AliTPCclustererMI::FitMax(Float_t vmatrix[5][5], Float_t y, Float_t z, Float_t sigmay, Float_t sigmaz)
410 {
411   //
412   // estimate max
413   Float_t sumteor= 0;
414   Float_t sumamp = 0;
415
416   for (Int_t di = -1;di<=1;di++)
417     for (Int_t dj = -1;dj<=1;dj++){
418       if (vmatrix[2+di][2+dj]>2){
419         Float_t teor = TMath::Gaus(di,y,sigmay*1.2)*TMath::Gaus(dj,z,sigmaz*1.2);
420         sumteor += teor*vmatrix[2+di][2+dj];
421         sumamp  += vmatrix[2+di][2+dj]*vmatrix[2+di][2+dj];
422       }
423     }    
424   Float_t max = sumamp/sumteor;
425   return max;
426 }
427
428 void AliTPCclustererMI::AddCluster(AliTPCclusterMI &c){
429   //
430   // transform cluster to the global coordinata
431   // add the cluster to the array
432   //
433   Float_t meani = c.GetY()/fPadWidth;
434   Float_t meanj = c.GetZ()/fZWidth;
435
436   Int_t ki = TMath::Nint(meani-3);
437   if (ki<0) ki=0;
438   if (ki>=fMaxPad) ki = fMaxPad-1;
439   Int_t kj = TMath::Nint(meanj-3);
440   if (kj<0) kj=0;
441   if (kj>=fMaxTime-3) kj=fMaxTime-4;
442   // ki and kj shifted to "real" coordinata
443   if (fRowDig) {
444     c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,0)-2,0);
445     c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,1)-2,1);
446     c.SetLabel(fRowDig->GetTrackIDFast(kj,ki,2)-2,2);
447   }
448   
449   
450   Float_t s2 = c.GetSigmaY2();
451   Float_t w=fParam->GetPadPitchWidth(fSector);
452   
453   c.SetSigmaY2(s2*w*w);
454   s2 = c.GetSigmaZ2(); 
455   w=fZWidth;
456   c.SetSigmaZ2(s2*w*w);
457   c.SetY((meani - 2.5 - 0.5*fMaxPad)*fParam->GetPadPitchWidth(fSector));
458   c.SetZ(fZWidth*(meanj-3)); 
459   c.SetZ(c.GetZ() - 3.*fParam->GetZSigma() + fParam->GetNTBinsL1()*fParam->GetZWidth()); // PASA delay + L1 delay
460   c.SetZ(fSign*(fParam->GetZLength() - c.GetZ()));
461   c.SetX(fRx);
462   c.SetDetector(fSector);
463   c.SetRow(fRow);
464
465   if (ki<=1 || ki>=fMaxPad-1 || kj==1 || kj==fMaxTime-2) {
466     //c.SetSigmaY2(c.GetSigmaY2()*25.);
467     //c.SetSigmaZ2(c.GetSigmaZ2()*4.);
468     c.SetType(-(c.GetType()+3));  //edge clusters
469   }
470   if (fLoop==2) c.SetType(100);
471
472   TClonesArray * arr = fRowCl->GetArray();
473   // AliTPCclusterMI * cl = 
474   new ((*arr)[fNcluster]) AliTPCclusterMI(c);
475
476   fNcluster++;
477 }
478
479
480 //_____________________________________________________________________________
481 void AliTPCclustererMI::Digits2Clusters()
482 {
483   //-----------------------------------------------------------------
484   // This is a simple cluster finder.
485   //-----------------------------------------------------------------
486
487   if (!fInput) { 
488     Error("Digits2Clusters", "input tree not initialised");
489     return;
490   }
491  
492   if (!fOutput) {
493     Error("Digits2Clusters", "output tree not initialised");
494     return;
495   }
496
497   AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
498
499   AliSimDigits digarr, *dummy=&digarr;
500   fRowDig = dummy;
501   fInput->GetBranch("Segment")->SetAddress(&dummy);
502   Stat_t nentries = fInput->GetEntries();
503   
504   fMaxTime=fParam->GetMaxTBin()+6; // add 3 virtual time bins before and 3   after
505     
506   Int_t nclusters  = 0;
507
508   for (Int_t n=0; n<nentries; n++) {
509     fInput->GetEvent(n);
510     if (!fParam->AdjustSectorRow(digarr.GetID(),fSector,fRow)) {
511       cerr<<"AliTPC warning: invalid segment ID ! "<<digarr.GetID()<<endl;
512       continue;
513     }
514     Int_t row = fRow;
515     AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector);  // pad gains per given sector
516     
517     //
518     AliTPCClustersRow *clrow= new AliTPCClustersRow();
519     fRowCl = clrow;
520     clrow->SetClass("AliTPCclusterMI");
521     clrow->SetArray(1);
522
523     clrow->SetID(digarr.GetID());
524     fOutput->GetBranch("Segment")->SetAddress(&clrow);
525     fRx=fParam->GetPadRowRadii(fSector,row);
526     
527     
528     const Int_t kNIS=fParam->GetNInnerSector(), kNOS=fParam->GetNOuterSector();
529     fZWidth = fParam->GetZWidth();
530     if (fSector < kNIS) {
531       fMaxPad = fParam->GetNPadsLow(row);
532       fSign = (fSector < kNIS/2) ? 1 : -1;
533       fPadLength = fParam->GetPadPitchLength(fSector,row);
534       fPadWidth = fParam->GetPadPitchWidth();
535     } else {
536       fMaxPad = fParam->GetNPadsUp(row);
537       fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
538       fPadLength = fParam->GetPadPitchLength(fSector,row);
539       fPadWidth  = fParam->GetPadPitchWidth();
540     }
541     
542     
543     fMaxBin=fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
544     fBins    =new Float_t[fMaxBin];
545     fResBins =new Float_t[fMaxBin];  //fBins with residuals after 1 finder loop 
546     memset(fBins,0,sizeof(Float_t)*fMaxBin);
547     memset(fResBins,0,sizeof(Float_t)*fMaxBin);
548     
549     if (digarr.First()) //MI change
550       do {
551         Float_t dig=digarr.CurrentDigit();
552         if (dig<=fParam->GetZeroSup()) continue;
553         Int_t j=digarr.CurrentRow()+3, i=digarr.CurrentColumn()+3;
554         Float_t gain = gainROC->GetValue(row,digarr.CurrentColumn());
555         fBins[i*fMaxTime+j]=dig/gain;
556       } while (digarr.Next());
557     digarr.ExpandTrackBuffer();
558
559     FindClusters();
560
561     fOutput->Fill();
562     delete clrow;    
563     nclusters+=fNcluster;    
564     delete[] fBins;      
565     delete[] fResBins;  
566   }  
567
568   Info("Digits2Clusters", "Number of found clusters : %d", nclusters);
569 }
570
571 void AliTPCclustererMI::Digits2Clusters(AliRawReader* rawReader)
572 {
573 //-----------------------------------------------------------------
574 // This is a cluster finder for the TPC raw data.
575 // The method assumes NO ordering of the altro channels.
576 // The pedestal subtraction can be switched on and off
577 // using an option of the TPC reconstructor
578 //-----------------------------------------------------------------
579
580   if (!fOutput) {
581     Error("Digits2Clusters", "output tree not initialised");
582     return;
583   }
584
585   fRowDig = NULL;
586
587   AliTPCCalPad * gainTPC = AliTPCcalibDB::Instance()->GetPadGainFactor();
588
589   Int_t nclusters  = 0;
590   
591   fMaxTime = fParam->GetMaxTBin() + 6; // add 3 virtual time bins before and 3 after
592   const Int_t kNIS = fParam->GetNInnerSector();
593   const Int_t kNOS = fParam->GetNOuterSector();
594   const Int_t kNS = kNIS + kNOS;
595   fZWidth = fParam->GetZWidth();
596   Int_t zeroSup = fParam->GetZeroSup();
597
598   Float_t** allBins = NULL;
599   Float_t** allBinsRes = NULL;
600
601   // Loop over sectors
602   for(fSector = 0; fSector < kNS; fSector++) {
603
604     AliTPCCalROC * gainROC = gainTPC->GetCalROC(fSector);  // pad gains per given sector
605  
606     Int_t nRows = 0;
607     Int_t nDDLs = 0, indexDDL = 0;
608     if (fSector < kNIS) {
609       nRows = fParam->GetNRowLow();
610       fSign = (fSector < kNIS/2) ? 1 : -1;
611       nDDLs = 2;
612       indexDDL = fSector * 2;
613     }
614     else {
615       nRows = fParam->GetNRowUp();
616       fSign = ((fSector-kNIS) < kNOS/2) ? 1 : -1;
617       nDDLs = 4;
618       indexDDL = (fSector-kNIS) * 4 + kNIS * 2;
619     }
620
621     allBins = new Float_t*[nRows];
622     allBinsRes = new Float_t*[nRows];
623
624     for (Int_t iRow = 0; iRow < nRows; iRow++) {
625       Int_t maxPad;
626       if (fSector < kNIS)
627         maxPad = fParam->GetNPadsLow(iRow);
628       else
629         maxPad = fParam->GetNPadsUp(iRow);
630       
631       Int_t maxBin = fMaxTime*(maxPad+6);  // add 3 virtual pads  before and 3 after
632       allBins[iRow] = new Float_t[maxBin];
633       allBinsRes[iRow] = new Float_t[maxBin];
634       memset(allBins[iRow],0,sizeof(Float_t)*maxBin);
635       memset(allBinsRes[iRow],0,sizeof(Float_t)*maxBin);
636     }
637     
638     // Loas the raw data for corresponding DDLs
639     rawReader->Reset();
640     AliTPCRawStream input(rawReader);
641     input.SetOldRCUFormat(fIsOldRCUFormat);
642     rawReader->Select("TPC",indexDDL,indexDDL+nDDLs-1);
643     
644     // Begin loop over altro data
645     while (input.Next()) {
646
647       if (input.GetSector() != fSector)
648         AliFatal(Form("Sector index mismatch ! Expected (%d), but got (%d) !",fSector,input.GetSector()));
649
650       
651       Int_t iRow = input.GetRow();
652       if (iRow < 0 || iRow >= nRows)
653         AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
654                       iRow, 0, nRows -1));
655
656       Int_t iPad = input.GetPad() + 3;
657
658       Int_t maxPad;
659       if (fSector < kNIS)
660         maxPad = fParam->GetNPadsLow(iRow);
661       else
662         maxPad = fParam->GetNPadsUp(iRow);
663
664       if (input.GetPad() < 0 || input.GetPad() >= maxPad)
665         AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
666                       input.GetPad(), 0, maxPad -1));
667
668       Int_t iTimeBin = input.GetTime() + 3;
669       if ( input.GetTime() < 0 || input.GetTime() >= fParam->GetMaxTBin())
670         AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
671                       input.GetTime(), 0, fParam->GetMaxTBin() -1));
672
673       Int_t maxBin = fMaxTime*(maxPad+6);  // add 3 virtual pads  before and 3 after
674
675       if (((iPad*fMaxTime+iTimeBin) >= maxBin) ||
676           ((iPad*fMaxTime+iTimeBin) < 0))
677         AliFatal(Form("Index outside the allowed range"
678                       " Sector=%d Row=%d Pad=%d Timebin=%d"
679                       " (Max.index=%d)",fSector,iRow,iPad,iTimeBin,maxBin));
680
681       Float_t signal = input.GetSignal();
682       //      if (!fPedSubtraction && signal <= zeroSup) continue;
683       if (!fRecoParam->GetCalcPedestal() && signal <= zeroSup) continue;
684
685       Float_t gain = gainROC->GetValue(iRow,input.GetPad());
686       allBins[iRow][iPad*fMaxTime+iTimeBin] = signal/gain;
687       allBins[iRow][iPad*fMaxTime+0] = 1;  // pad with signal
688     } // End of the loop over altro data
689
690     // Now loop over rows and perform pedestal subtraction
691     //    if (fPedSubtraction) {
692     if (fRecoParam->GetCalcPedestal()) {
693       for (Int_t iRow = 0; iRow < nRows; iRow++) {
694         Int_t maxPad;
695         if (fSector < kNIS)
696           maxPad = fParam->GetNPadsLow(iRow);
697         else
698           maxPad = fParam->GetNPadsUp(iRow);
699
700         for (Int_t iPad = 0; iPad < maxPad + 6; iPad++) {
701           if (allBins[iRow][iPad*fMaxTime+0] !=1) continue;  // no data
702           Float_t *p = &allBins[iRow][iPad*fMaxTime+3];
703           //Float_t pedestal = TMath::Median(fMaxTime, p);      
704           Int_t id[3] = {fSector, iRow, iPad-3};
705           Float_t pedestal = ProcesSignal(p, fMaxTime, id);
706           for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
707             allBins[iRow][iPad*fMaxTime+iTimeBin] -= pedestal;
708             if (iTimeBin < AliTPCReconstructor::GetRecoParam()->GetFirstBin())  
709               allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
710             if (iTimeBin > AliTPCReconstructor::GetRecoParam()->GetLastBin())  
711               allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
712             if (allBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
713               allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
714           }
715         }
716       }
717     }
718
719     // Now loop over rows and find clusters
720     for (fRow = 0; fRow < nRows; fRow++) {
721       fRowCl = new AliTPCClustersRow;
722       fRowCl->SetClass("AliTPCclusterMI");
723       fRowCl->SetArray(1);
724       fRowCl->SetID(fParam->GetIndex(fSector, fRow));
725       fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
726
727       fRx = fParam->GetPadRowRadii(fSector, fRow);
728       fPadLength = fParam->GetPadPitchLength(fSector, fRow);
729       fPadWidth  = fParam->GetPadPitchWidth();
730       if (fSector < kNIS)
731         fMaxPad = fParam->GetNPadsLow(fRow);
732       else
733         fMaxPad = fParam->GetNPadsUp(fRow);
734       fMaxBin = fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
735
736       fBins = allBins[fRow];
737       fResBins = allBinsRes[fRow];
738
739       FindClusters();
740
741       fOutput->Fill();
742       delete fRowCl;    
743       nclusters += fNcluster;    
744     } // End of loop to find clusters
745
746  
747     for (Int_t iRow = 0; iRow < nRows; iRow++) {
748       delete [] allBins[iRow];
749       delete [] allBinsRes[iRow];
750     }
751
752     delete [] allBins;
753     delete [] allBinsRes;
754
755   } // End of loop over sectors
756
757   Info("Digits2Clusters", "Number of found clusters : %d\n", nclusters);
758
759 }
760
761 void AliTPCclustererMI::FindClusters()
762 {
763   //add virtual charge at the edge   
764   for (Int_t i=0; i<fMaxTime; i++){
765     Float_t amp1 = fBins[i+3*fMaxTime]; 
766     Float_t amp0 =0;
767     if (amp1>0){
768       Float_t amp2 = fBins[i+4*fMaxTime];
769       if (amp2==0) amp2=0.5;
770       Float_t sigma2 = GetSigmaY2(i);           
771       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);   
772       if (gDebug>4) printf("\n%f\n",amp0);
773     }  
774     fBins[i+2*fMaxTime] = amp0;    
775     amp0 = 0;
776     amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
777     if (amp1>0){
778       Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
779       if (amp2==0) amp2=0.5;
780       Float_t sigma2 = GetSigmaY2(i);           
781       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);           
782       if (gDebug>4) printf("\n%f\n",amp0);
783     }        
784     fBins[(fMaxPad+3)*fMaxTime+i] = amp0;           
785   }
786
787 //  memcpy(fResBins,fBins, fMaxBin*2);
788   memcpy(fResBins,fBins, fMaxBin);
789   //
790   fNcluster=0;
791   //first loop - for "gold cluster" 
792   fLoop=1;
793   Float_t *b=&fBins[-1]+2*fMaxTime;
794   Int_t crtime = Int_t((fParam->GetZLength()-fRecoParam->GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
795
796   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
797     b++;
798     if (*b<8) continue;   //threshold form maxima
799     if (i%fMaxTime<crtime) {
800       Int_t delta = -(i%fMaxTime)+crtime;
801       b+=delta;
802       i+=delta;
803       continue; 
804     }
805      
806     if (!IsMaximum(*b,fMaxTime,b)) continue;
807     AliTPCclusterMI c;
808     Int_t dummy=0;
809     MakeCluster(i, fMaxTime, fBins, dummy,c);
810     //}
811   }
812   //memcpy(fBins,fResBins, fMaxBin*2);
813   //second  loop - for rest cluster 
814   /*        
815   fLoop=2;
816   b=&fResBins[-1]+2*fMaxTime;
817   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
818     b++;
819     if (*b<25) continue;   // bigger threshold for maxima
820     if (!IsMaximum(*b,fMaxTime,b)) continue;
821     AliTPCclusterMI c;
822     Int_t dummy;
823     MakeCluster(i, fMaxTime, fResBins, dummy,c);
824     //}
825   }
826   */
827 }
828
829
830 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3]){
831   //
832   // process signal on given pad - + streaming of additional information in special mode
833   //
834   // id[0] - sector
835   // id[1] - row
836   // id[2] - pad  
837   Int_t offset =100;
838   Float_t kMin =fRecoParam->GetDumpAmplitudeMin();   // minimal signal to be dumped
839   Double_t median = TMath::Median(nchannels-offset, &(signal[offset]));
840   if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
841   //
842
843   Double_t mean   = TMath::Mean(nchannels-offset, &(signal[offset]));
844   Double_t rms    = TMath::RMS(nchannels-offset, &(signal[offset]));
845   Double_t *dsignal = new Double_t[nchannels];
846   Double_t *dtime   = new Double_t[nchannels];
847   Float_t max    =  0;
848   Float_t maxPos =  0;
849   for (Int_t i=0; i<nchannels; i++){
850     dtime[i] = i;
851     dsignal[i] = signal[i];
852     if (i<offset) continue;
853     if (signal[i]>max && i <fMaxTime-100) {  // temporary remove spike signals at the end
854       max = signal[i];
855       maxPos = i;
856     }    
857   }
858   
859   Double_t mean06=0, mean09=0;
860   Double_t rms06=0, rms09=0; 
861   AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean06, rms06, int(0.6*(nchannels-offset)));
862   AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean09, rms09, int(0.9*(nchannels-offset)));
863   //
864   //
865   UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
866   if (uid[0]< AliTPCROC::Instance()->GetNSectors() 
867       && uid[1]<  AliTPCROC::Instance()->GetNRows(uid[0])  && 
868       uid[2] < AliTPCROC::Instance()->GetNPads(uid[0], uid[1])){
869     if (!fAmplitudeHisto){
870       fAmplitudeHisto = new TObjArray(72);
871     }
872     TObjArray  * sectorArray = (TObjArray*)fAmplitudeHisto->UncheckedAt(uid[0]);
873     if (!sectorArray){
874       Int_t npads = AliTPCROC::Instance()->GetNChannels(uid[0]);
875       sectorArray = new TObjArray(npads);
876       fAmplitudeHisto->AddAt(sectorArray, uid[0]);
877     }
878     Int_t position =  uid[2]+ AliTPCROC::Instance()->GetRowIndexes(uid[0])[uid[1]];
879     TH1F * histo = (TH1F*)sectorArray->UncheckedAt(position);
880     if (!histo){
881       char hname[100];
882       sprintf(hname,"Amp_%d_%d_%d",uid[0],uid[1],uid[2]);
883       TFile * backup = gFile;
884       fDebugStreamer->GetFile()->cd();
885       histo = new TH1F(hname, hname, 100, 5,100);
886       //histo->SetDirectory(0);     // histogram not connected to directory -(File)
887       sectorArray->AddAt(histo, position);
888       if (backup) backup->cd();
889     }
890     for (Int_t i=0; i<nchannels; i++){
891       if (signal[i]>0) histo->Fill(signal[i]);
892     }
893   }
894   //
895   TGraph * graph;
896   Bool_t random = (gRandom->Rndm()<0.0001);
897   if (max-median>kMin || rms06>2.*fParam->GetZeroSup() || random){
898     graph =new TGraph(nchannels, dtime, dsignal);
899     if (rms06>2.*fParam->GetZeroSup() || random)
900       (*fDebugStreamer)<<"SignalN"<<    //noise pads - or random sample of pads
901         "Sector="<<uid[0]<<
902         "Row="<<uid[1]<<
903         "Pad="<<uid[2]<<
904         "Graph.="<<graph<<
905         "Max="<<max<<
906         "MaxPos="<<maxPos<<
907         //
908         "Median="<<median<<
909         "Mean="<<mean<<
910         "RMS="<<rms<<      
911         "Mean06="<<mean06<<
912         "RMS06="<<rms06<<
913         "Mean09="<<mean09<<
914         "RMS09="<<rms09<<
915         "\n";
916     if (max-median>kMin) 
917       (*fDebugStreamer)<<"SignalB"<<     // pads with signal
918         "Sector="<<uid[0]<<
919         "Row="<<uid[1]<<
920         "Pad="<<uid[2]<<
921         "Graph.="<<graph<<
922         "Max="<<max<<
923         "MaxPos="<<maxPos<<     
924         //
925         "Median="<<median<<
926         "Mean="<<mean<<
927         "RMS="<<rms<<      
928         "Mean06="<<mean06<<
929         "RMS06="<<rms06<<
930         "Mean09="<<mean09<<
931         "RMS09="<<rms09<<
932         "\n";
933     delete graph;
934   }
935   
936   (*fDebugStreamer)<<"Signal"<<
937     "Sector="<<uid[0]<<
938     "Row="<<uid[1]<<
939     "Pad="<<uid[2]<<
940     "Max="<<max<<
941     "MaxPos="<<maxPos<<
942     //
943     "Median="<<median<<
944     "Mean="<<mean<<
945     "RMS="<<rms<<      
946     "Mean06="<<mean06<<
947     "RMS06="<<rms06<<
948     "Mean09="<<mean09<<
949     "RMS09="<<rms09<<
950     "\n";
951   //
952   //
953   //  Central Electrode signal analysis  
954   //
955   Double_t ceQmax  =0, ceQsum=0, ceTime=0;
956   Double_t cemean  = mean06, cerms=rms06 ;
957   Int_t    cemaxpos= 0;
958   Double_t ceThreshold=5.*cerms;
959   Double_t ceSumThreshold=8.*cerms;
960   const Int_t    cemin=5;  // range for the analysis of the ce signal +- channels from the peak
961   const Int_t    cemax=5;
962   for (Int_t i=nchannels-2; i>1; i--){
963     if ( (dsignal[i]-mean06)>ceThreshold && dsignal[i]>=dsignal[i+1] && dsignal[i]>=dsignal[i-1] ){
964       cemaxpos=i;
965       break;
966     }
967   }
968   if (cemaxpos!=0){
969       for (Int_t i=cemaxpos-cemin; i<cemaxpos+cemax; i++){
970           if (i>0 && i<nchannels&&dsignal[i]- cemean>0){
971               Double_t val=dsignal[i]- cemean;
972               ceTime+=val*dtime[i];
973               ceQsum+=val;
974               if (val>ceQmax) ceQmax=val;
975           }
976       }
977       if (ceQmax&&ceQsum>ceSumThreshold) {
978           ceTime/=ceQsum;
979           (*fDebugStreamer)<<"Signalce"<<
980               "Sector="<<uid[0]<<
981               "Row="<<uid[1]<<
982               "Pad="<<uid[2]<<
983               "Max="<<ceQmax<<
984               "Qsum="<<ceQsum<<
985               "Time="<<ceTime<<
986               "RMS06="<<rms06<<
987               //
988               "\n";
989       }
990   }
991   // end of ce signal analysis
992   //
993
994   //
995   //  Gating grid signal analysis  
996   //
997   Double_t ggQmax  =0, ggQsum=0, ggTime=0;
998   Double_t ggmean  = mean06, ggrms=rms06 ;
999   Int_t    ggmaxpos= 0;
1000   Double_t ggThreshold=5.*ggrms;
1001   Double_t ggSumThreshold=8.*ggrms;
1002
1003   for (Int_t i=1; i<nchannels-1; i++){
1004     if ( (dsignal[i]-mean06)>ggThreshold && dsignal[i]>=dsignal[i+1] && dsignal[i]>=dsignal[i-1] &&
1005          (dsignal[i]+dsignal[i+1]+dsignal[i-1]-3*mean06)>ggSumThreshold){
1006       ggmaxpos=i;
1007       if (dsignal[i-1]>dsignal[i+1]) ggmaxpos=i-1;
1008       break;
1009     }
1010   }
1011   if (ggmaxpos!=0){
1012       for (Int_t i=ggmaxpos-1; i<ggmaxpos+3; i++){       
1013           if (i>0 && i<nchannels && dsignal[i]-ggmean>0){
1014               Double_t val=dsignal[i]- ggmean;
1015               ggTime+=val*dtime[i];
1016               ggQsum+=val;
1017               if (val>ggQmax) ggQmax=val;
1018           }
1019       }
1020       if (ggQmax&&ggQsum>ggSumThreshold) {
1021           ggTime/=ggQsum;
1022           (*fDebugStreamer)<<"Signalgg"<<
1023               "Sector="<<uid[0]<<
1024               "Row="<<uid[1]<<
1025               "Pad="<<uid[2]<<
1026               "Max="<<ggQmax<<
1027               "Qsum="<<ggQsum<<
1028               "Time="<<ggTime<<
1029               "RMS06="<<rms06<<
1030               //
1031               "\n";
1032       }
1033   }
1034   // end of gg signal analysis
1035       
1036
1037   delete [] dsignal;
1038   delete [] dtime;
1039   if (rms06>fRecoParam->GetMaxNoise()) return 1024+median; // sign noisy channel in debug mode
1040   return median;
1041 }
1042
1043
1044
1045 void AliTPCclustererMI::DumpHistos(){
1046   if (!fAmplitudeHisto) return;
1047   for (UInt_t isector=0; isector<AliTPCROC::Instance()->GetNSectors(); isector++){
1048     TObjArray * array = (TObjArray*)fAmplitudeHisto->UncheckedAt(isector);
1049     if (!array) continue;
1050     for (UInt_t ipad = 0; ipad <(UInt_t)array->GetEntriesFast(); ipad++){
1051       TH1F * histo = (TH1F*) array->UncheckedAt(ipad);
1052       if (!histo) continue;
1053       if (histo->GetEntries()<100) continue;
1054       histo->Fit("gaus","q");
1055       Float_t mean =  histo->GetMean();
1056       Float_t rms  =  histo->GetRMS();
1057       Float_t gmean = histo->GetFunction("gaus")->GetParameter(1);
1058       Float_t gsigma = histo->GetFunction("gaus")->GetParameter(2);
1059       Float_t max = histo->GetFunction("gaus")->GetParameter(0);
1060
1061       // get pad number
1062       UInt_t row=0, pad =0;
1063       const UInt_t *indexes = AliTPCROC::Instance()->GetRowIndexes(isector);
1064       for (UInt_t irow=0; irow< AliTPCROC::Instance()->GetNRows(isector); irow++){
1065         if (indexes[irow]<ipad){
1066           row = irow;
1067           pad = ipad-indexes[irow];
1068         }
1069       }      
1070       Int_t rpad = pad - (AliTPCROC::Instance()->GetNPads(isector,row))/2;
1071       //
1072       (*fDebugStreamer)<<"Fit"<<
1073         "Sector="<<isector<<
1074         "Row="<<row<<
1075         "Pad="<<pad<<
1076         "RPad="<<rpad<<
1077         "Max="<<max<<
1078         "Mean="<<mean<<
1079         "RMS="<<rms<<      
1080         "GMean="<<gmean<<
1081         "GSigma="<<gsigma<<
1082         "\n";
1083       if (array->UncheckedAt(ipad)) fDebugStreamer->StoreObject(array->UncheckedAt(ipad));
1084     }
1085   }
1086 }