]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHClusterFinder.cxx
Will be recrated later in AliRICH.h after AliRICHRecon will be tested
[u/mrichter/AliRoot.git] / RICH / AliRICHClusterFinder.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 #include "AliRICHClusterFinder.h"
18 #include "AliRun.h"
19 #include "AliRICH.h"
20 #include "AliRICHMap.h"
21 #include "AliRICHSDigit.h"
22 #include "AliRICHDigit.h"
23 #include "AliRICHRawCluster.h"
24 #include "AliRICHParam.h"
25 #include <TTree.h>
26 #include <TCanvas.h>
27 #include <TH1.h>
28 #include <TF1.h>
29 #include <TPad.h>
30 #include <TGraph.h> 
31 #include <TMinuit.h>
32
33 static AliSegmentation     *gSegmentation;
34 static AliRICHResponse*     gResponse;
35 static Int_t                gix[500];
36 static Int_t                giy[500];
37 static Float_t              gCharge[500];
38 static Int_t                gNbins;
39 static Bool_t               gFirst=kTRUE;
40 static TMinuit *gMyMinuit ;
41 void fcn(Int_t &npar, Double_t *gin, Double_t &f, Double_t *par, Int_t);
42 static Int_t                gChargeTot;
43
44 ClassImp(AliRICHClusterFinder)
45 //__________________________________________________________________________________________________
46 AliRICHClusterFinder::AliRICHClusterFinder(AliRICH *pRICH)   
47 {//main ctor
48   Info("main ctor","Start.");
49   
50   fRICH=pRICH;
51   
52   fSegmentation=Rich()->C(1)->GetSegmentationModel();
53   fResponse    =Rich()->C(1)->GetResponseModel();
54     
55   fDigits=0;    fNdigits=0;
56   fChamber=0;
57   fHitMap=0;  
58   
59   fCogCorr = 0;
60   SetNperMax();
61   SetClusterSize();
62   fNPeaks=-1;
63 }//main ctor
64 //__________________________________________________________________________________________________
65 void AliRICHClusterFinder::Decluster(AliRICHRawCluster *cluster)
66 {// Decluster algorithm
67   Info("Decluster","Start.");    
68   if(cluster->fMultiplicity==1||cluster->fMultiplicity==2){//Nothing special for 1- and 2-clusters
69     if(fNPeaks != 0) {cluster->fNcluster[0]=fNPeaks; cluster->fNcluster[1]=0;} 
70     AddRawCluster(*cluster); 
71     fNPeaks++;
72   }else if(cluster->fMultiplicity==3){// 3-cluster, check topology
73       Centered(cluster);// ok, cluster is centered and added in Centered()
74   }else{//4-and more-pad clusters
75     if(cluster->fMultiplicity<= fMaxClusterSize){
76         SplitByLocalMaxima(cluster);
77     }//if <= fClusterSize
78   }//if multiplicity 
79 }//Decluster()
80 //__________________________________________________________________________________________________
81 Bool_t AliRICHClusterFinder::Centered(AliRICHRawCluster *cluster)
82 {//Is the cluster centered?
83
84   AliRICHDigit* dig;
85   dig= (AliRICHDigit*)fDigits->UncheckedAt(cluster->fIndexMap[0]);
86   Int_t x[kMaxNeighbours], y[kMaxNeighbours], xN[kMaxNeighbours], yN[kMaxNeighbours];
87   Int_t nn=Rich()->Param()->Neighbours(dig->PadX(),dig->PadY(),x,y);
88     
89   
90   Int_t nd=0;
91   for (Int_t i=0; i<nn; i++){//neighbours loop
92     if(fHitMap->TestHit(x[i],y[i]) == kUsed){
93       xN[nd]=x[i];
94       yN[nd]=y[i];
95       nd++;
96     }
97   }//neighbours loop
98     
99   if(nd==2){// cluster is centered !
100         if (fNPeaks != 0) {
101             cluster->fNcluster[0]=fNPeaks;
102             cluster->fNcluster[1]=0;
103         }  
104         cluster->fCtype=0;
105         AddRawCluster(*cluster);
106         fNPeaks++;
107         return kTRUE;
108     } else if (nd ==1) {
109 // Highest signal on an edge, split cluster into 2+1
110 // who is the neighbour ?            
111         Int_t nind=fHitMap->GetHitIndex(xN[0], yN[0]);
112         Int_t i1= (nind==cluster->fIndexMap[1]) ? 1:2;
113         Int_t i2= (nind==cluster->fIndexMap[1]) ? 2:1;    
114 // 2-cluster
115         AliRICHRawCluster cnew;
116         if (fNPeaks == 0) {
117             cnew.fNcluster[0]=-1;
118             cnew.fNcluster[1]=fNRawClusters;
119         } else {
120             cnew.fNcluster[0]=fNPeaks;
121             cnew.fNcluster[1]=0;
122         }
123         cnew.fMultiplicity=2;
124         cnew.fIndexMap[0]=cluster->fIndexMap[0];
125         cnew.fIndexMap[1]=cluster->fIndexMap[i1];
126         FillCluster(&cnew);
127         cnew.fClusterType=cnew.PhysicsContribution();
128         AddRawCluster(cnew);
129         fNPeaks++;
130 // 1-cluster
131         cluster->fMultiplicity=1;
132         cluster->fIndexMap[0]=cluster->fIndexMap[i2];
133         cluster->fIndexMap[1]=0;
134         cluster->fIndexMap[2]=0;        
135         FillCluster(cluster);
136         if (fNPeaks != 0) {
137             cluster->fNcluster[0]=fNPeaks;
138             cluster->fNcluster[1]=0;
139         }  
140         cluster->fClusterType=cluster->PhysicsContribution();
141         AddRawCluster(*cluster);
142         fNPeaks++;
143         return kFALSE;
144     } else {
145         Warning("Centered","\n Completely screwed up %d !! \n",nd);
146         
147     }    
148   return kFALSE;
149 }//Centered()
150 //__________________________________________________________________________________________________
151 void AliRICHClusterFinder::SplitByLocalMaxima(AliRICHRawCluster *c)
152 {// Split the cluster according to the number of maxima inside
153   Info("SplitbyLocalMaxima","Start.");
154   
155       AliRICHDigit* dig[100], *digt;
156     Int_t ix[100], iy[100], q[100];
157     Float_t x[100], y[100];
158     Int_t i; // loops over digits
159     Int_t j; // loops over local maxima
160     Int_t mul=c->fMultiplicity;
161 //  dump digit information into arrays
162   for (i=0; i<mul; i++){
163     dig[i]= (AliRICHDigit*)fDigits->UncheckedAt(c->fIndexMap[i]);
164     ix[i]= dig[i]->PadX();
165     iy[i]= dig[i]->PadY();
166     q[i] = dig[i]->Signal();
167     AliRICHParam::Pad2Local(ix[i], iy[i], x[i], y[i]);
168   }
169 //  Find local maxima
170     Bool_t isLocal[100];
171     Int_t nLocal=0;
172     Int_t associatePeak[100];
173     Int_t indLocal[100];
174     Int_t nn;
175     Int_t xNei[kMaxNeighbours], yNei[kMaxNeighbours];
176     for (i=0; i<mul; i++) {
177         fSegmentation->Neighbours(ix[i], iy[i], &nn, xNei, yNei);
178         isLocal[i]=kTRUE;
179         for (j=0; j<nn; j++) {
180             if (fHitMap->TestHit(xNei[j], yNei[j])==kEmpty) continue;
181             digt=(AliRICHDigit*) fHitMap->GetHit(xNei[j], yNei[j]);
182             if (digt->Signal() > q[i]) {
183                 isLocal[i]=kFALSE;
184                 break;
185 // handle special case of neighbouring pads with equal signal
186             } else if (digt->Signal() == q[i]) {
187                 if (nLocal >0) {
188                     for (Int_t k=0; k<nLocal; k++) {
189                         if (xNei[j]==ix[indLocal[k]] && yNei[j]==iy[indLocal[k]]){
190                             isLocal[i]=kFALSE;
191                         }
192                     }
193                 }
194             } 
195         } // loop over next neighbours
196         // Maxima should not be on the edge
197         if (isLocal[i]) {
198             indLocal[nLocal]=i;
199             nLocal++;
200         } 
201     } // loop over all digits
202 // If only one local maximum found but multiplicity is high take global maximum from the list of digits.    
203     if (nLocal==1 && mul>5) {
204         Int_t nnew=0;
205         for (i=0; i<mul; i++) {
206             if (!isLocal[i]) {
207                 indLocal[nLocal]=i;
208                 isLocal[i]=kTRUE;
209                 nLocal++;
210                 nnew++;
211             }
212             if (nnew==1) break;
213         }
214     }
215     if(nLocal==2) {// If number of local maxima is 2 try to fit a double gaussian
216
217 //  Initialise global variables for fit
218         gFirst=kTRUE;
219         gSegmentation=fSegmentation;
220         gResponse    =fResponse;
221         gNbins=mul;
222         
223         for (i=0; i<mul; i++) {
224             gix[i]=ix[i];
225             giy[i]=iy[i];
226             gCharge[i]=Float_t(q[i]);
227         }
228         if (gFirst)    gMyMinuit = new TMinuit(5);
229         
230         gMyMinuit->SetFCN(fcn);
231         gMyMinuit->mninit(5,10,7);
232         Double_t arglist[20];
233         arglist[0]=1;
234 // Set starting values 
235         static Double_t vstart[5];
236         vstart[0]=x[indLocal[0]];
237         vstart[1]=y[indLocal[0]];       
238         vstart[2]=x[indLocal[1]];
239         vstart[3]=y[indLocal[1]];       
240         vstart[4]=Float_t(q[indLocal[0]])/Float_t(q[indLocal[0]]+q[indLocal[1]]);
241 // lower and upper limits
242         static Double_t lower[5], upper[5];
243         lower[0]=vstart[0]-AliRICHParam::PadSizeX()/2;
244         lower[1]=vstart[1]-AliRICHParam::PadSizeY()/2;
245         
246         upper[0]=vstart[0]+AliRICHParam::PadSizeX()/2;
247         upper[1]=vstart[1]+AliRICHParam::PadSizeY()/2;
248         
249         lower[2]=vstart[2]-AliRICHParam::PadSizeX()/2;
250         lower[3]=vstart[3]-AliRICHParam::PadSizeY()/2;
251         
252         upper[2]=vstart[2]+AliRICHParam::PadSizeX()/2;
253         upper[3]=vstart[3]+AliRICHParam::PadSizeY()/2;
254         
255         lower[4]=0.;
256         upper[4]=1.;
257 // step sizes
258         static Double_t step[5]={0.005, 0.03, 0.005, 0.03, 0.01};
259         Int_t iErr;
260         
261         gMyMinuit->mnparm(0,"x1",vstart[0],step[0],lower[0],upper[0],iErr);
262         gMyMinuit->mnparm(1,"y1",vstart[1],step[1],lower[1],upper[1],iErr);
263         gMyMinuit->mnparm(2,"x2",vstart[2],step[2],lower[2],upper[2],iErr);
264         gMyMinuit->mnparm(3,"y2",vstart[3],step[3],lower[3],upper[3],iErr);
265         gMyMinuit->mnparm(4,"a0",vstart[4],step[4],lower[4],upper[4],iErr);
266 // ready for minimisation       
267         gMyMinuit->SetPrintLevel(-1);
268         gMyMinuit->mnexcm("SET OUT", arglist, 0, iErr);
269         arglist[0]= -1;
270         arglist[1]= 0;
271         
272         gMyMinuit->mnexcm("SET NOGR", arglist, 0, iErr);
273         gMyMinuit->mnexcm("SIMPLEX", arglist, 0, iErr);
274         gMyMinuit->mnexcm("MIGRAD", arglist, 0, iErr);
275         gMyMinuit->mnexcm("EXIT" , arglist, 0, iErr);
276
277         Double_t xrec[2], yrec[2], qfrac;
278         TString chname;
279         Double_t epxz, b1, b2;
280         gMyMinuit->mnpout(0, chname, xrec[0], epxz, b1, b2, iErr);      
281         gMyMinuit->mnpout(1, chname, yrec[0], epxz, b1, b2, iErr);      
282         gMyMinuit->mnpout(2, chname, xrec[1], epxz, b1, b2, iErr);      
283         gMyMinuit->mnpout(3, chname, yrec[1], epxz, b1, b2, iErr);      
284         gMyMinuit->mnpout(4, chname, qfrac,   epxz, b1, b2, iErr);      
285         
286         cout<<"xrex[0]="<<xrec[0]<<"yrec[0]="<<yrec[0]<<"xrec[1]="<<xrec[1]<<"yrec[1]="<<yrec[1]<<"qfrac="<<qfrac<<endl;
287         for (j=0; j<2; j++) { // One cluster for each maximum
288             AliRICHRawCluster cnew;
289             if (fNPeaks == 0) {
290                 cnew.fNcluster[0]=-1;
291                 cnew.fNcluster[1]=fNRawClusters;
292             } else {
293                 cnew.fNcluster[0]=fNPeaks;
294                 cnew.fNcluster[1]=0;
295             }
296             cnew.fMultiplicity=0;
297             cnew.fX=Float_t(xrec[j]);
298             cnew.fY=Float_t(yrec[j]);
299             if (j==0) {
300                 cnew.fQ=Int_t(gChargeTot*qfrac);
301             } else {
302                 cnew.fQ=Int_t(gChargeTot*(1-qfrac));
303             }
304             for (i=0; i<mul; i++) {
305                 cnew.fIndexMap[cnew.fMultiplicity]=c->fIndexMap[i];
306                 cnew.fContMap[cnew.fMultiplicity]=AliRICHParam::AssignChargeToPad(xrec[j],yrec[j],gix[i], giy[i]);
307                 cnew.fMultiplicity++;
308             }
309             FillCluster(&cnew,0);
310             cnew.fClusterType=cnew.PhysicsContribution();
311             AddRawCluster(cnew);
312             fNPeaks++;
313         }
314     }//if 2 maximum in cluster
315     Bool_t fitted=kTRUE;
316
317     if (nLocal >2 || !fitted) {
318         // Check if enough local clusters have been found, if not add global maxima to the list 
319         Int_t nPerMax;
320         if (nLocal!=0) {
321             nPerMax=mul/nLocal;
322         } else {
323             Warning("SplitByLocalMaxima","no local maximum found");
324             nPerMax=fNperMax+1;
325         }
326         
327         if (nPerMax > fNperMax) {
328             Int_t nGlob=mul/fNperMax-nLocal+1;
329             if (nGlob > 0) {
330                 Int_t nnew=0;
331                 for (i=0; i<mul; i++) {
332                     if (!isLocal[i]) {
333                         indLocal[nLocal]=i;
334                         isLocal[i]=kTRUE;
335                         nLocal++;
336                         nnew++;
337                     }
338                     if (nnew==nGlob) break;
339                 }
340             }
341         }
342         for (i=0; i<mul; i++) { // Associate hits to peaks
343             Float_t dmin=1.E10;
344             Float_t qmax=0;
345             if (isLocal[i]) continue;
346             for (j=0; j<nLocal; j++) {
347                 Int_t il=indLocal[j];
348                 Float_t d=TMath::Sqrt((x[i]-x[il])*(x[i]-x[il])
349                                       +(y[i]-y[il])*(y[i]-y[il]));
350                 Float_t ql=q[il];
351                 if (d<dmin) {           // Select nearest peak
352                     dmin=d;
353                     qmax=ql;
354                     associatePeak[i]=j;
355                 } else if (d==dmin) {               // If more than one take highest peak
356                     if (ql>qmax) {
357                         dmin=d;
358                         qmax=ql;
359                         associatePeak[i]=j;
360                     }
361                 }
362             }
363         }       
364  // One cluster for each maximum
365         for (j=0; j<nLocal; j++) {
366             AliRICHRawCluster cnew;
367             if (fNPeaks == 0) {
368                 cnew.fNcluster[0]=-1;
369                 cnew.fNcluster[1]=fNRawClusters;
370             } else {
371                 cnew.fNcluster[0]=fNPeaks;
372                 cnew.fNcluster[1]=0;
373             }
374             cnew.fIndexMap[0]=c->fIndexMap[indLocal[j]];
375             cnew.fMultiplicity=1;
376             for (i=0; i<mul; i++) {
377                 if (isLocal[i]) continue;
378                 if (associatePeak[i]==j) {
379                     cnew.fIndexMap[cnew.fMultiplicity]=c->fIndexMap[i];
380                     cnew.fMultiplicity++;
381                 }
382             }
383             FillCluster(&cnew);
384             cnew.fClusterType=cnew.PhysicsContribution();
385             AddRawCluster(cnew);
386             fNPeaks++;
387         }
388     }
389 }//SplitByLocalMaxima(AliRICHRawCluster *c)
390 //__________________________________________________________________________________________________
391 void  AliRICHClusterFinder::FillCluster(AliRICHRawCluster* c, Int_t flag) 
392 {//  Completes cluster information starting from list of digits
393     AliRICHDigit* dig;
394     Float_t x, y;
395     Int_t  ix, iy;
396     Float_t fraction=0;
397     
398     c->fPeakSignal=0;
399     if (flag) {
400         c->fX=0;
401         c->fY=0;
402         c->fQ=0;
403     }
404  
405
406     for (Int_t i=0; i<c->fMultiplicity; i++){
407         dig= (AliRICHDigit*)fDigits->UncheckedAt(c->fIndexMap[i]);
408         ix=dig->PadX();
409         iy=dig->PadY();
410         Int_t q=dig->Signal();
411         if (dig->Physics() >= dig->Signal()) {
412           c->fPhysicsMap[i]=2;
413         } else if (dig->Physics() == 0) {
414           c->fPhysicsMap[i]=0;
415         } else  c->fPhysicsMap[i]=1;
416 // peak signal and track list
417         if (flag) {
418            if (q>c->fPeakSignal) {
419               c->fPeakSignal=q;
420             c->fTracks[0]=dig->Hit();
421             c->fTracks[1]=dig->Track(0);
422             c->fTracks[2]=dig->Track(1);
423            }
424         } else {
425            if (c->fContMap[i] > fraction) {
426               fraction=c->fContMap[i];
427               c->fPeakSignal=q;
428             c->fTracks[0]=dig->Hit();
429             c->fTracks[1]=dig->Track(0);
430             c->fTracks[2]=dig->Track(1);
431            }
432         }
433         if (flag) {
434             AliRICHParam::Pad2Local(ix,iy,x,y);
435             c->fX += q*x;
436             c->fY += q*y;
437             c->fQ += q;
438         }
439
440     } // loop over digits
441
442  if (flag) {
443      
444      c->fX/=c->fQ;
445      c->fX=fSegmentation->GetAnod(c->fX);
446      c->fY/=c->fQ; 
447 //  apply correction to the coordinate along the anode wire
448      x=c->fX;   
449      y=c->fY;
450      AliRICHParam::Local2Pad(x,y,ix,iy);
451      AliRICHParam::Pad2Local(ix,iy,x,y);
452      Int_t isec=fSegmentation->Sector(ix,iy);
453      TF1* cogCorr = fSegmentation->CorrFunc(isec-1);
454      
455      if (cogCorr) {
456          Float_t yOnPad=(c->fY-y)/fSegmentation->Dpy(isec);
457          c->fY=c->fY-cogCorr->Eval(yOnPad, 0, 0);
458      }
459  }
460 }//FillCluster() 
461 //__________________________________________________________________________________________________
462 void  AliRICHClusterFinder::AddDigit2Cluster(Int_t i, Int_t j, AliRICHRawCluster &c)
463 {//Find clusters Add i,j as element of the cluster  
464   Info("AddDigit2Cluster","Start with digit(%i,%i)",i,j);
465   
466   Int_t idx = fHitMap->GetHitIndex(i,j);
467   AliRICHDigit* dig = (AliRICHDigit*) fHitMap->GetHit(i,j);
468   Int_t q=dig->Signal();
469   if(q>TMath::Abs(c.fPeakSignal)){
470         c.fPeakSignal=q;
471         c.fTracks[0]=dig->Hit();
472         c.fTracks[1]=dig->Track(0);
473         c.fTracks[2]=dig->Track(1);
474     }
475 //  Make sure that list of digits is ordered 
476     Int_t mu=c.fMultiplicity;
477     c.fIndexMap[mu]=idx;
478
479     if (dig->Physics() >= dig->Signal()) {
480         c.fPhysicsMap[mu]=2;
481     } else if (dig->Physics() == 0) {
482         c.fPhysicsMap[mu]=0;
483     } else  c.fPhysicsMap[mu]=1;
484
485     if (mu > 0) {
486         for (Int_t ind=mu-1; ind>=0; ind--) {
487             Int_t ist=(c.fIndexMap)[ind];
488             Int_t ql=((AliRICHDigit*)fDigits->UncheckedAt(ist))->Signal();
489             if (q>ql) {
490                 c.fIndexMap[ind]=idx;
491                 c.fIndexMap[ind+1]=ist;
492             } else {
493                 break;
494             }
495         }
496     }
497     
498   c.fMultiplicity++;    
499     if (c.fMultiplicity >= 50 ) {
500         Info("AddDigit2CLuster","multiplicity >50  %d \n",c.fMultiplicity);
501         c.fMultiplicity=49;
502     }
503   Float_t x,y;// Prepare center of gravity calculation
504   AliRICHParam::Pad2Local(i,j,x,y);
505   c.fX+=q*x;    c.fY+=q*y;    c.fQ += q;
506   fHitMap->FlagHit(i,j);// Flag hit as taken  
507
508
509   Int_t xList[4], yList[4];    //  Now look recursively for all neighbours
510   for (Int_t iNei=0;iNei<Rich()->Param()->Neighbours(i,j,xList,yList);iNei++)
511     if(fHitMap->TestHit(xList[iNei],yList[iNei])==kUnused) AddDigit2Cluster(xList[iNei],yList[iNei],c);    
512 }//AddDigit2Cluster()
513 //__________________________________________________________________________________________________
514 void AliRICHClusterFinder::FindRawClusters()
515 {//finds neighbours and fill the tree with raw clusters
516   Info("FindRawClusters","Start for Chamber %i.",fChamber);
517   
518   if(!fNdigits)return;
519
520   fHitMap=new AliRICHMap(fDigits);
521
522   for(Int_t iDigN=0;iDigN<fNdigits;iDigN++){//digits loop
523     AliRICHDigit *dig=(AliRICHDigit*)fDigits->UncheckedAt(iDigN);
524     Int_t i=dig->PadX();   Int_t j=dig->PadY();
525     if(fHitMap->TestHit(i,j)==kUsed||fHitMap->TestHit(i,j)==kEmpty) continue;
526         
527     AliRICHRawCluster c;
528     c.fMultiplicity=0; c.fPeakSignal=dig->Signal();
529     c.fTracks[0]=dig->Hit();c.fTracks[1]=dig->Track(0);c.fTracks[2]=dig->Track(1);        
530     c.fNcluster[0]=-1;// tag the beginning of cluster list in a raw cluster
531     
532     AddDigit2Cluster(i,j,c);//form initial cluster
533         
534     c.fX /= c.fQ;       // center of gravity
535     //c.fX=fSegmentation->GetAnod(c.fX);
536     c.fY /= c.fQ;
537     //AddRawCluster(c);
538     
539 //    Int_t ix,iy;//  apply correction to the coordinate along the anode wire
540 //    Float_t x=c.fX, y=c.fY;   
541 //    Rich()->Param()->Local2Pad(x,y,ix,iy);
542 //    Rich()->Param()->Pad2Local(ix,iy,x,y);
543 //    Int_t isec=fSegmentation->Sector(ix,iy);
544 //    TF1* cogCorr=fSegmentation->CorrFunc(isec-1);
545 //    if(cogCorr){
546 //      Float_t yOnPad=(c.fY-y)/fSegmentation->Dpy(isec);
547 //      c.fY=c.fY-cogCorr->Eval(yOnPad,0,0);
548 //    }
549
550     c.fNcluster[1]=fNRawClusters; c.fClusterType=c.PhysicsContribution();
551     
552     Decluster(&c);
553     
554     fNPeaks=0;
555
556     c.fMultiplicity=0; for(int k=0;k<c.fMultiplicity;k++) c.fIndexMap[k]=0;//reset cluster object    
557   }//digits loop
558   delete fHitMap;
559   Info("FindRawClusters","Stop.");
560 }//FindRawClusters()
561 //__________________________________________________________________________________________________
562 void AliRICHClusterFinder::CalibrateCOG()
563 {// Calibration
564
565     Float_t x[5];
566     Float_t y[5];
567     Int_t n, i;
568     if (fSegmentation) {
569         TF1 *func;
570         fSegmentation->GiveTestPoints(n, x, y);
571         for (i=0; i<n; i++) {
572             func = 0;
573             Float_t xtest=x[i];
574             Float_t ytest=y[i];     
575             SinoidalFit(xtest, ytest, func);
576             if (func) fSegmentation->SetCorrFunc(i, new TF1(*func));
577         }
578     }
579 }//CalibrateCOG()
580 //__________________________________________________________________________________________________
581 void AliRICHClusterFinder::SinoidalFit(Float_t x, Float_t y, TF1 *func)
582 {//Sinoidal fit
583   static Int_t count=0;
584     
585     count++;
586
587     const Int_t kNs=101;
588     Float_t xg[kNs], yg[kNs], xrg[kNs], yrg[kNs];
589     Float_t xsig[kNs], ysig[kNs];
590    
591     Int_t ix,iy;
592     AliRICHParam::Local2Pad(x,y,ix,iy);   
593     AliRICHParam::Pad2Local(ix,iy,x,y);   
594     Int_t isec=fSegmentation->Sector(ix,iy);
595 // Pad Limits    
596     Float_t xmin = x-Rich()->Param()->PadSizeX()/2;
597     Float_t ymin = y-Rich()->Param()->PadSizeY()/2;
598 //              
599 //      Integration Limits
600     Float_t dxI=Rich()->Param()->SigmaIntegration()*Rich()->Param()->ChargeSpreadX();
601     Float_t dyI=Rich()->Param()->SigmaIntegration()*Rich()->Param()->ChargeSpreadY();
602
603 //
604 //  Scanning
605 //
606     Int_t i;
607     Float_t qp=0;
608
609 //  y-position
610     Float_t yscan=ymin;
611     Float_t dy=Rich()->Param()->PadSizeY()/(kNs-1);
612
613     for (i=0; i<kNs; i++) {//      Pad Loop
614         Float_t sum=0;
615         Float_t qcheck=0;
616         fSegmentation->SigGenInit(x, yscan, 0);
617         
618         for (fSegmentation->FirstPad(x, yscan,0, dxI, dyI); 
619              fSegmentation->MorePads(); 
620              fSegmentation->NextPad()) 
621         {
622             qp=fResponse->IntXY(fSegmentation);
623             qp=TMath::Abs(qp);
624             if (qp > 1.e-4) {
625                 qcheck+=qp;
626                 Int_t ixs=fSegmentation->Ix();
627                 Int_t iys=fSegmentation->Iy();
628                 Float_t xs,ys;
629                 AliRICHParam::Pad2Local(ixs,iys,xs,ys);
630                 sum+=qp*ys;
631             }
632         } // Pad loop
633         Float_t ycog=sum/qcheck;
634         yg[i]=(yscan-y)/fSegmentation->Dpy(isec);
635         yrg[i]=(ycog-y)/fSegmentation->Dpy(isec);
636         ysig[i]=ycog-yscan;
637         yscan+=dy;
638     } // scan loop
639 //  x-position
640     Float_t xscan=xmin;
641     Float_t dx=fSegmentation->Dpx(isec)/(kNs-1);
642
643     for (i=0; i<kNs; i++) {//      Pad Loop
644         Float_t sum=0;
645         Float_t qcheck=0;
646         fSegmentation->SigGenInit(xscan, y, 0);
647         
648         for (fSegmentation->FirstPad(xscan, y, 0, dxI, dyI); 
649              fSegmentation->MorePads(); 
650              fSegmentation->NextPad()) 
651         {
652             qp=fResponse->IntXY(fSegmentation);
653             qp=TMath::Abs(qp);
654             if (qp > 1.e-2) {
655                 qcheck+=qp;
656                 Int_t ixs=fSegmentation->Ix();
657                 Int_t iys=fSegmentation->Iy();
658                 Float_t xs,ys;
659                 AliRICHParam::Pad2Local(ixs,iys,xs,ys);
660                 sum+=qp*xs;
661             }
662         } // Pad loop
663         Float_t xcog=sum/qcheck;
664         xcog=fSegmentation->GetAnod(xcog);
665         
666         xg[i]=(xscan-x)/fSegmentation->Dpx(isec);
667         xrg[i]=(xcog-x)/fSegmentation->Dpx(isec);
668         xsig[i]=xcog-xscan;
669         xscan+=dx;
670     }
671 // Creates a Root function based on function sinoid above and perform the fit
672     TGraph *graphyr= new TGraph(kNs,yrg,ysig);
673     Double_t sinoid(Double_t *x, Double_t *par);
674     new TF1("sinoidf",sinoid,0.5,0.5,5);
675     graphyr->Fit("sinoidf","Q");
676     func = (TF1*)graphyr->GetListOfFunctions()->At(0);
677 }//SinoidalFit()
678 //__________________________________________________________________________________________________
679 Double_t sinoid(Double_t *x, Double_t *par)
680 {// Sinoid function
681
682     Double_t arg = -2*TMath::Pi()*x[0];
683     Double_t fitval= par[0]*TMath::Sin(arg)+
684         par[1]*TMath::Sin(2*arg)+
685         par[2]*TMath::Sin(3*arg)+
686         par[3]*TMath::Sin(4*arg)+
687         par[4]*TMath::Sin(5*arg);
688     return fitval;
689 }//sinoid()
690 //__________________________________________________________________________________________________
691 Double_t DoubleGauss(Double_t *x, Double_t *par)
692 {//Double gaussian function
693   Double_t arg1 = (x[0]-par[1])/0.18;
694   Double_t arg2 = (x[0]-par[3])/0.18;
695   return par[0]*TMath::Exp(-arg1*arg1/2)+par[2]*TMath::Exp(-arg2*arg2/2);
696 }
697 //__________________________________________________________________________________________________
698 Float_t DiscrCharge(Int_t i,Double_t *par) 
699 {
700 // par[0]    x-position of first  cluster
701 // par[1]    y-position of first  cluster
702 // par[2]    x-position of second cluster
703 // par[3]    y-position of second cluster
704 // par[4]    charge fraction of first  cluster
705 // 1-par[4]  charge fraction of second cluster
706
707     static Float_t qtot;
708     if (gFirst) {
709         qtot=0;
710         for (Int_t jbin=0; jbin<gNbins; jbin++) {
711             qtot+=gCharge[jbin];
712         }
713         gFirst=kFALSE;
714         gChargeTot=Int_t(qtot);
715         
716     }
717     Float_t q1=AliRICHParam::AssignChargeToPad(par[0],par[1],gix[i],giy[i]);
718     
719     Float_t q2=AliRICHParam::AssignChargeToPad(par[2],par[3],gix[i],giy[i]);
720 //    cout<<"qtot="<<gChargeTot<<" q1="<<q1<<" q2="<<q2<<" px="<<gix[i]<<" py="<<giy[i]<<endl;
721     
722     Float_t value = qtot*(par[4]*q1+(1.-par[4])*q2);
723     return value;
724 }//DiscrCharge(Int_t i,Double_t *par) 
725 //__________________________________________________________________________________________________
726 void fcn(Int_t &npar, Double_t */*gin*/, Double_t &f, Double_t *par, Int_t)
727 {// Minimisation function
728   npar=1;
729     Int_t i;
730     Float_t delta;
731     Float_t chisq=0;
732     Float_t qcont=0;
733     Float_t qtot=0;
734     
735     for (i=0; i<gNbins; i++) {
736         Float_t q0=gCharge[i];
737         Float_t q1=DiscrCharge(i,par);
738         delta=(q0-q1)/TMath::Sqrt(q0);
739         chisq+=delta*delta;
740         qcont+=q1;
741         qtot+=q0;
742     }
743 //    chisq=chisq+=(qtot-qcont)*(qtot-qcont)*0.5;
744     f=chisq;
745 }//
746 //__________________________________________________________________________________________________
747 void AliRICHClusterFinder::Exec()
748 {
749   Info("Exec","Start.");
750   
751   Rich()->GetLoader()->LoadDigits(); 
752   
753   for(Int_t iEventN=0;iEventN<gAlice->GetEventsPerRun();iEventN++){//events loop
754     gAlice->GetRunLoader()->GetEvent(iEventN);
755     
756     Rich()->GetLoader()->MakeTree("R");  Rich()->MakeBranch("R");
757     Rich()->ResetDigitsOld();  Rich()->ResetRawClusters();
758     
759     Rich()->GetLoader()->TreeD()->GetEntry(0);
760     for(fChamber=1;fChamber<=kNCH;fChamber++){//chambers loop
761       fDigits=Rich()->DigitsOld(fChamber); fNdigits=fDigits->GetEntries();
762       
763       FindRawClusters();
764         
765     }//chambers loop
766     
767     Rich()->GetLoader()->TreeR()->Fill();
768     Rich()->GetLoader()->WriteRecPoints("OVERWRITE");
769   }//events loop  
770   Rich()->GetLoader()->UnloadDigits(); Rich()->GetLoader()->UnloadRecPoints();  
771   Rich()->ResetDigitsOld();  Rich()->ResetRawClusters();
772   Info("Exec","Stop.");      
773 }//Exec()
774 //__________________________________________________________________________________________________