]> git.uio.no Git - u/mrichter/AliRoot.git/blob - TPC/AliTPCclustererMI.cxx
Set Tracking parameters (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       if (input.GetTime() < AliTPCReconstructor::GetRecoParam()->GetFirstBin()) continue;
651       if (input.GetTime() > AliTPCReconstructor::GetRecoParam()->GetLastBin()) continue;
652       
653       Int_t iRow = input.GetRow();
654       if (iRow < 0 || iRow >= nRows)
655         AliFatal(Form("Pad-row index (%d) outside the range (%d -> %d) !",
656                       iRow, 0, nRows -1));
657
658       Int_t iPad = input.GetPad() + 3;
659
660       Int_t maxPad;
661       if (fSector < kNIS)
662         maxPad = fParam->GetNPadsLow(iRow);
663       else
664         maxPad = fParam->GetNPadsUp(iRow);
665
666       if (input.GetPad() < 0 || input.GetPad() >= maxPad)
667         AliFatal(Form("Pad index (%d) outside the range (%d -> %d) !",
668                       input.GetPad(), 0, maxPad -1));
669
670       Int_t iTimeBin = input.GetTime() + 3;
671       if ( input.GetTime() < 0 || input.GetTime() >= fParam->GetMaxTBin())
672         AliFatal(Form("Timebin index (%d) outside the range (%d -> %d) !",
673                       input.GetTime(), 0, fParam->GetMaxTBin() -1));
674
675       Int_t maxBin = fMaxTime*(maxPad+6);  // add 3 virtual pads  before and 3 after
676
677       if (((iPad*fMaxTime+iTimeBin) >= maxBin) ||
678           ((iPad*fMaxTime+iTimeBin) < 0))
679         AliFatal(Form("Index outside the allowed range"
680                       " Sector=%d Row=%d Pad=%d Timebin=%d"
681                       " (Max.index=%d)",fSector,iRow,iPad,iTimeBin,maxBin));
682
683       Float_t signal = input.GetSignal();
684       //      if (!fPedSubtraction && signal <= zeroSup) continue;
685       if (!fRecoParam->GetCalcPedestal() && signal <= zeroSup) continue;
686
687       Float_t gain = gainROC->GetValue(iRow,input.GetPad());
688       allBins[iRow][iPad*fMaxTime+iTimeBin] = signal/gain;
689       allBins[iRow][iPad*fMaxTime+0] = 1;  // pad with signal
690     } // End of the loop over altro data
691
692     // Now loop over rows and perform pedestal subtraction
693     //    if (fPedSubtraction) {
694     if (fRecoParam->GetCalcPedestal()) {
695       for (Int_t iRow = 0; iRow < nRows; iRow++) {
696         Int_t maxPad;
697         if (fSector < kNIS)
698           maxPad = fParam->GetNPadsLow(iRow);
699         else
700           maxPad = fParam->GetNPadsUp(iRow);
701
702         for (Int_t iPad = 0; iPad < maxPad + 6; iPad++) {
703           if (allBins[iRow][iPad*fMaxTime+0] !=1) continue;  // no data
704           Float_t *p = &allBins[iRow][iPad*fMaxTime+3];
705           //Float_t pedestal = TMath::Median(fMaxTime, p);      
706           Int_t id[3] = {fSector, iRow, iPad-3};
707           Float_t pedestal = ProcesSignal(p, fMaxTime, id);
708           for (Int_t iTimeBin = 0; iTimeBin < fMaxTime; iTimeBin++) {
709             allBins[iRow][iPad*fMaxTime+iTimeBin] -= pedestal;
710             if (allBins[iRow][iPad*fMaxTime+iTimeBin] < zeroSup)
711               allBins[iRow][iPad*fMaxTime+iTimeBin] = 0;
712           }
713         }
714       }
715     }
716
717     // Now loop over rows and find clusters
718     for (fRow = 0; fRow < nRows; fRow++) {
719       fRowCl = new AliTPCClustersRow;
720       fRowCl->SetClass("AliTPCclusterMI");
721       fRowCl->SetArray(1);
722       fRowCl->SetID(fParam->GetIndex(fSector, fRow));
723       fOutput->GetBranch("Segment")->SetAddress(&fRowCl);
724
725       fRx = fParam->GetPadRowRadii(fSector, fRow);
726       fPadLength = fParam->GetPadPitchLength(fSector, fRow);
727       fPadWidth  = fParam->GetPadPitchWidth();
728       if (fSector < kNIS)
729         fMaxPad = fParam->GetNPadsLow(fRow);
730       else
731         fMaxPad = fParam->GetNPadsUp(fRow);
732       fMaxBin = fMaxTime*(fMaxPad+6);  // add 3 virtual pads  before and 3 after
733
734       fBins = allBins[fRow];
735       fResBins = allBinsRes[fRow];
736
737       FindClusters();
738
739       fOutput->Fill();
740       delete fRowCl;    
741       nclusters += fNcluster;    
742     } // End of loop to find clusters
743
744  
745     for (Int_t iRow = 0; iRow < nRows; iRow++) {
746       delete [] allBins[iRow];
747       delete [] allBinsRes[iRow];
748     }
749
750     delete [] allBins;
751     delete [] allBinsRes;
752
753   } // End of loop over sectors
754
755   Info("Digits2Clusters", "Number of found clusters : %d\n", nclusters);
756
757 }
758
759 void AliTPCclustererMI::FindClusters()
760 {
761   //add virtual charge at the edge   
762   for (Int_t i=0; i<fMaxTime; i++){
763     Float_t amp1 = fBins[i+3*fMaxTime]; 
764     Float_t amp0 =0;
765     if (amp1>0){
766       Float_t amp2 = fBins[i+4*fMaxTime];
767       if (amp2==0) amp2=0.5;
768       Float_t sigma2 = GetSigmaY2(i);           
769       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);   
770       if (gDebug>4) printf("\n%f\n",amp0);
771     }  
772     fBins[i+2*fMaxTime] = amp0;    
773     amp0 = 0;
774     amp1 = fBins[(fMaxPad+2)*fMaxTime+i];
775     if (amp1>0){
776       Float_t amp2 = fBins[i+(fMaxPad+1)*fMaxTime];
777       if (amp2==0) amp2=0.5;
778       Float_t sigma2 = GetSigmaY2(i);           
779       amp0 = (amp1*amp1/amp2)*TMath::Exp(-1./sigma2);           
780       if (gDebug>4) printf("\n%f\n",amp0);
781     }        
782     fBins[(fMaxPad+3)*fMaxTime+i] = amp0;           
783   }
784
785 //  memcpy(fResBins,fBins, fMaxBin*2);
786   memcpy(fResBins,fBins, fMaxBin);
787   //
788   fNcluster=0;
789   //first loop - for "gold cluster" 
790   fLoop=1;
791   Float_t *b=&fBins[-1]+2*fMaxTime;
792   Int_t crtime = Int_t((fParam->GetZLength()-fRecoParam->GetCtgRange()*fRx)/fZWidth-fParam->GetNTBinsL1()-5);
793
794   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
795     b++;
796     if (*b<8) continue;   //threshold form maxima
797     if (i%fMaxTime<crtime) {
798       Int_t delta = -(i%fMaxTime)+crtime;
799       b+=delta;
800       i+=delta;
801       continue; 
802     }
803      
804     if (!IsMaximum(*b,fMaxTime,b)) continue;
805     AliTPCclusterMI c;
806     Int_t dummy=0;
807     MakeCluster(i, fMaxTime, fBins, dummy,c);
808     //}
809   }
810   //memcpy(fBins,fResBins, fMaxBin*2);
811   //second  loop - for rest cluster 
812   /*        
813   fLoop=2;
814   b=&fResBins[-1]+2*fMaxTime;
815   for (Int_t i=2*fMaxTime; i<fMaxBin-2*fMaxTime; i++) {
816     b++;
817     if (*b<25) continue;   // bigger threshold for maxima
818     if (!IsMaximum(*b,fMaxTime,b)) continue;
819     AliTPCclusterMI c;
820     Int_t dummy;
821     MakeCluster(i, fMaxTime, fResBins, dummy,c);
822     //}
823   }
824   */
825 }
826
827
828 Double_t AliTPCclustererMI::ProcesSignal(Float_t *signal, Int_t nchannels, Int_t id[3]){
829   //
830   // process signal on given pad - + streaming of additional information in special mode
831   //
832   // id[0] - sector
833   // id[1] - row
834   // id[2] - pad  
835   Int_t offset =100;
836   Float_t kMin =fRecoParam->GetDumpAmplitudeMin();   // minimal signal to be dumped
837   Double_t median = TMath::Median(nchannels-offset, &(signal[offset]));
838   if (AliLog::GetDebugLevel("","AliTPCclustererMI")==0) return median;
839   //
840
841   Double_t mean   = TMath::Mean(nchannels-offset, &(signal[offset]));
842   Double_t rms    = TMath::RMS(nchannels-offset, &(signal[offset]));
843   Double_t *dsignal = new Double_t[nchannels];
844   Double_t *dtime   = new Double_t[nchannels];
845   Float_t max    =  0;
846   Float_t maxPos =  0;
847   for (Int_t i=0; i<nchannels; i++){
848     dtime[i] = i;
849     dsignal[i] = signal[i];
850     if (signal[i]>max && i <fMaxTime-100) {  // temporary remove spike signals at the end
851       max = signal[i];
852       maxPos = i;
853     }    
854   }
855   
856   Double_t mean06=0, mean09=0;
857   Double_t rms06=0, rms09=0; 
858   AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean06, rms06, int(0.6*(nchannels-offset)));
859   AliMathBase::EvaluateUni(nchannels-offset, &(dsignal[offset]), mean09, rms09, int(0.9*(nchannels-offset)));
860   
861   //
862   UInt_t uid[3] = {UInt_t(id[0]),UInt_t(id[1]),UInt_t(id[2])};
863   if (uid[0]< AliTPCROC::Instance()->GetNSectors() 
864       && uid[1]<  AliTPCROC::Instance()->GetNRows(uid[0])  && 
865       uid[2] < AliTPCROC::Instance()->GetNPads(uid[0], uid[1])){
866     if (!fAmplitudeHisto){
867       fAmplitudeHisto = new TObjArray(72);
868     }
869     TObjArray  * sectorArray = (TObjArray*)fAmplitudeHisto->UncheckedAt(uid[0]);
870     if (!sectorArray){
871       Int_t npads = AliTPCROC::Instance()->GetNChannels(uid[0]);
872       sectorArray = new TObjArray(npads);
873       fAmplitudeHisto->AddAt(sectorArray, uid[0]);
874     }
875     Int_t position =  uid[2]+ AliTPCROC::Instance()->GetRowIndexes(uid[0])[uid[1]];
876     TH1F * histo = (TH1F*)sectorArray->UncheckedAt(position);
877     if (!histo){
878       char hname[100];
879       sprintf(hname,"Amp_%d_%d_%d",uid[0],uid[1],uid[2]);
880       TFile * backup = gFile;
881       fDebugStreamer->GetFile()->cd();
882       histo = new TH1F(hname, hname, 100, 5,100);
883       //histo->SetDirectory(0);     // histogram not connected to directory -(File)
884       sectorArray->AddAt(histo, position);
885       if (backup) backup->cd();
886     }
887     for (Int_t i=0; i<nchannels; i++){
888       if (signal[i]>0) histo->Fill(signal[i]);
889     }
890   }
891   //
892   TGraph * graph;
893   Bool_t random = (gRandom->Rndm()<0.0001);
894   if (max-median>kMin || rms06>2.*fParam->GetZeroSup() || random){
895     graph =new TGraph(nchannels, dtime, dsignal);
896     if (rms06>2.*fParam->GetZeroSup() || random)
897       (*fDebugStreamer)<<"SignalN"<<    //noise pads - or random sample of pads
898         "Sector="<<uid[0]<<
899         "Row="<<uid[1]<<
900         "Pad="<<uid[2]<<
901         "Graph.="<<graph<<
902         "Max="<<max<<
903         "MaxPos="<<maxPos<<
904         //
905         "Median="<<median<<
906         "Mean="<<mean<<
907         "RMS="<<rms<<      
908         "Mean06="<<mean06<<
909         "RMS06="<<rms06<<
910         "Mean09="<<mean09<<
911         "RMS09="<<rms09<<
912         "\n";
913     if (max-median>kMin) 
914       (*fDebugStreamer)<<"SignalB"<<     // pads with signal
915         "Sector="<<uid[0]<<
916         "Row="<<uid[1]<<
917         "Pad="<<uid[2]<<
918         "Graph.="<<graph<<
919         "Max="<<max<<
920         "MaxPos="<<maxPos<<     
921         //
922         "Median="<<median<<
923         "Mean="<<mean<<
924         "RMS="<<rms<<      
925         "Mean06="<<mean06<<
926         "RMS06="<<rms06<<
927         "Mean09="<<mean09<<
928         "RMS09="<<rms09<<
929         "\n";
930     delete graph;
931   }
932   
933   (*fDebugStreamer)<<"Signal"<<
934     "Sector="<<uid[0]<<
935     "Row="<<uid[1]<<
936     "Pad="<<uid[2]<<
937     "Max="<<max<<
938     "MaxPos="<<maxPos<<
939     //
940     "Median="<<median<<
941     "Mean="<<mean<<
942     "RMS="<<rms<<      
943     "Mean06="<<mean06<<
944     "RMS06="<<rms06<<
945     "Mean09="<<mean09<<
946     "RMS09="<<rms09<<
947     "\n";
948   delete [] dsignal;
949   delete [] dtime;
950   if (rms06>fRecoParam->GetMaxNoise()) return 1024+median; // sign noisy channel in debug mode
951   return median;
952 }
953
954
955
956 void AliTPCclustererMI::DumpHistos(){
957   if (!fAmplitudeHisto) return;
958   for (UInt_t isector=0; isector<AliTPCROC::Instance()->GetNSectors(); isector++){
959     TObjArray * array = (TObjArray*)fAmplitudeHisto->UncheckedAt(isector);
960     if (!array) continue;
961     for (UInt_t ipad = 0; ipad <(UInt_t)array->GetEntriesFast(); ipad++){
962       TH1F * histo = (TH1F*) array->UncheckedAt(ipad);
963       if (!histo) continue;
964       if (histo->GetEntries()<100) continue;
965       histo->Fit("gaus","q");
966       Float_t mean =  histo->GetMean();
967       Float_t rms  =  histo->GetRMS();
968       Float_t gmean = histo->GetFunction("gaus")->GetParameter(1);
969       Float_t gsigma = histo->GetFunction("gaus")->GetParameter(2);
970       Float_t max = histo->GetFunction("gaus")->GetParameter(0);
971
972       // get pad number
973       UInt_t row=0, pad =0;
974       const UInt_t *indexes = AliTPCROC::Instance()->GetRowIndexes(isector);
975       for (UInt_t irow=0; irow< AliTPCROC::Instance()->GetNRows(isector); irow++){
976         if (indexes[irow]<ipad){
977           row = irow;
978           pad = ipad-indexes[irow];
979         }
980       }      
981       Int_t rpad = pad - (AliTPCROC::Instance()->GetNPads(isector,row))/2;
982       //
983       (*fDebugStreamer)<<"Fit"<<
984         "Sector="<<isector<<
985         "Row="<<row<<
986         "Pad="<<pad<<
987         "RPad="<<rpad<<
988         "Max="<<max<<
989         "Mean="<<mean<<
990         "RMS="<<rms<<      
991         "GMean="<<gmean<<
992         "GSigma="<<gsigma<<
993         "\n";
994       if (array->UncheckedAt(ipad)) fDebugStreamer->StoreObject(array->UncheckedAt(ipad));
995     }
996   }
997 }