]> git.uio.no Git - u/mrichter/AliRoot.git/blob - RICH/AliRICHRecon.cxx
Removing extra semicolon
[u/mrichter/AliRoot.git] / RICH / AliRICHRecon.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 //                                                                      //
18 // AliRICHRecon                                                         //
19 //                                                                      //
20 // RICH class to perfom pattern recognition based on Hough transfrom    //
21 // for single chamber                                                   //
22 //////////////////////////////////////////////////////////////////////////
23
24 #include "AliRICHRecon.h"  //class header
25 #include "AliRICHCluster.h" //CkovAngle()
26 #include <TRotation.h>     //TracePhoton()
27 #include <TH1D.h>          //HoughResponse()
28 #include <TClonesArray.h>  //CkovAngle()
29
30 #include <TTree.h>         //Display()
31 #include <TFile.h>         //Display()
32 #include <AliESD.h>        //Display()
33 #include <TPolyMarker.h>   //Display()
34 #include <TLatex.h>        //Display()
35 #include <TCanvas.h>       //Display()
36
37
38 const Double_t AliRICHRecon::fgkRadThick=1.5;
39 const Double_t AliRICHRecon::fgkWinThick=0.5;
40 const Double_t AliRICHRecon::fgkGapThick=8.0;
41 const Double_t AliRICHRecon::fgkRadIdx  =1.292;
42 const Double_t AliRICHRecon::fgkWinIdx  =1.5787;
43 const Double_t AliRICHRecon::fgkGapIdx  =1.0005;
44
45
46 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 AliRICHRecon::AliRICHRecon():TTask("RichRec","RichPat"),
48   fPhotCnt(-1),
49   fCkovSigma2(0),
50   fIsWEIGHT(kFALSE),
51   fDTheta(0.001),
52   fWindowWidth(0.045),
53   fTrkDir(TVector3(0,0,1)),fTrkPos(TVector2(30,40))  
54 {
55 // main ctor
56   for (Int_t i=0; i<3000; i++) {
57     fPhotFlag[i] =  0;
58     fPhotCkov[i] = -1;
59     fPhotPhi [i] = -1;
60     fPhotWei [i] =  0;
61   }
62 }
63 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
64 Double_t AliRICHRecon::CkovAngle(TClonesArray *pCluLst,Int_t &iNaccepted)
65 {
66 // Pattern recognition method based on Hough transform
67 // Arguments: pCluLst  - list of clusters for this chamber   
68 //   Returns:          - track ckov angle, [rad], 
69   
70   if(pCluLst->GetEntries()>200) fIsWEIGHT = kTRUE; // offset to take into account bkg in reconstruction
71   else                          fIsWEIGHT = kFALSE;
72
73   // Photon Flag:  Flag = 0 initial set; Flag = 1 good candidate (charge compatible with photon); Flag = 2 photon used for the ring;
74
75   fPhotCnt=0;                                                      
76   for (Int_t iClu=0; iClu<pCluLst->GetEntriesFast();iClu++){//clusters loop
77     AliRICHCluster *pClu=(AliRICHCluster*)pCluLst->UncheckedAt(iClu);                       //get pointer to current cluster    
78     if(pClu->Q()>100) continue;                                                             //avoid MIP clusters from bkg
79     
80     fPhotCkov[fPhotCnt]=FindPhotCkov(pClu->X(),pClu->Y());                                  //find ckov angle for this  photon candidate
81     fPhotCnt++;         //increment counter of photon candidates
82   }//clusters loop
83
84   iNaccepted=FlagPhot(HoughResponse()); //flag photons according to individual theta ckov with respect to most probable track theta ckov
85   if(iNaccepted<1) return -11; 
86   else             return FindRingCkov(pCluLst->GetEntries());  //find best Theta ckov for ring i.e. track
87 }//ThetaCerenkov()
88 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
89 Double_t AliRICHRecon::FindPhotCkov(Double_t cluX,Double_t cluY)
90 {
91 // Finds Cerenkov angle  for this photon candidate
92 // Arguments: cluX,cluY - position of cadidate's cluster  
93 //   Returns: Cerenkov angle 
94
95   TVector2 pos(cluX,cluY); Double_t cluR=(pos-fTrkPos).Mod();  Double_t phi=FindPhotPhi(cluX,cluY);      
96   Printf("new dist %f phi %f",cluR,phi);
97   Double_t ckov1=0,ckov2=0.75;
98   const Double_t kTol=0.05; 
99   Int_t iIterCnt = 0;
100   while(1){
101     if(iIterCnt>=50) return -1;
102     Double_t ckov=0.5*(ckov1+ckov2);
103     Double_t dist=cluR-TracePhot(ckov,phi,pos); iIterCnt++;   //get distance between trial point and cluster position
104     Printf("New: phi %f ckov %f dist %f",phi,ckov,dist);
105     if     (dist> kTol) ckov1=ckov;                           //cluster @ larger ckov 
106     else if(dist<-kTol) ckov2=ckov;                           //cluster @ smaller ckov
107     else                return ckov;                          //precision achived         
108   }
109 }//FindPhotTheta()
110 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
111 Double_t AliRICHRecon::FindPhotPhi(Double_t cluX,Double_t cluY)
112 {
113 // Finds phi angle og photon candidate by considering the cluster's position  of this candudate w.r.t track position
114   
115   Double_t emiss=0; 
116   return fPhotPhi[fPhotCnt]=TMath::ATan2(cluY-fTrkPos.Y()-emiss*TMath::Tan(fTrkDir.Theta())*TMath::Sin(fTrkDir.Phi()),
117                                          cluX-fTrkPos.X()-emiss*TMath::Tan(fTrkDir.Theta())*TMath::Cos(fTrkDir.Phi()));
118 }
119 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
120 Double_t AliRICHRecon::FindRingArea(Double_t ckovAng)const
121 {
122 // Find area inside the cerenkov ring which lays inside PCs
123 // Arguments: ckovThe - cernkov    
124 //   Returns: area of the ring in cm^2 for given theta ckov
125    
126   
127   TVector2 pos1,pos2;
128   
129   const Int_t kN=100;
130   Double_t area=0;
131   for(Int_t i=0;i<kN;i++){
132     TracePhot(ckovAng,Double_t(TMath::TwoPi()*i    /kN),pos1);//trace this photon 
133     TracePhot(ckovAng,Double_t(TMath::TwoPi()*(i+1)/kN),pos2);//trace this photon 
134     area+=(pos1-fTrkPos)*(pos2-fTrkPos);
135       
136   }
137   return area;
138 }//FindRingArea()
139 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
140 Double_t AliRICHRecon::FindRingCkov(Int_t)
141 {
142 // Loops on all Ckov candidates and estimates the best Theta Ckov for a ring formed by those candidates. Also estimates an error for that Theat Ckov
143 // collecting errors for all single Ckov candidates thetas. (Assuming they are independent)  
144 // Arguments: iNclus- total number of clusters in chamber for background estimation
145 //    Return: best estimation of track Theta ckov
146
147   Double_t wei = 0.;
148   Double_t weightThetaCerenkov = 0.;
149
150   Double_t ckovMin=9999.,ckovMax=0.;
151   Double_t sigma2 = 0;   //to collect error squared for this ring
152   
153   for(Int_t i=0;i<fPhotCnt;i++){//candidates loop
154     if(fPhotFlag[i] == 2){
155       if(fPhotCkov[i]<ckovMin) ckovMin=fPhotCkov[i];  //find max and min Theta ckov from all candidates within probable window
156       if(fPhotCkov[i]>ckovMax) ckovMax=fPhotCkov[i]; 
157       weightThetaCerenkov += fPhotCkov[i]*fPhotWei[i];   wei += fPhotWei[i];                 //collect weight as sum of all candidate weghts   
158       
159      //Double_t phiref=(GetPhiPoint()-GetTrackPhi());
160        if(fPhotCkov[i]<=0) continue;//?????????????????Flag photos = 2 may imply CkovEta = 0?????????????? 
161                                      
162       sigma2 += 1./Sigma2(fPhotCkov[i],fPhotPhi[i]);
163     }
164   }//candidates loop
165   
166   if(sigma2>0) fCkovSigma2=1./sigma2;
167   else         fCkovSigma2=1e10;  
168   
169
170   if(wei != 0.) weightThetaCerenkov /= wei; else weightThetaCerenkov = 0.;  
171   return weightThetaCerenkov;
172 }//FindCkovRing()
173 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
174 Int_t AliRICHRecon::FlagPhot(Double_t ckov)
175 {
176 // Flag photon candidates if their individual ckov angle is inside the window around ckov angle returned by  HoughResponse()
177 // Arguments: ckov- value of most probable ckov angle for track as returned by HoughResponse()
178 //   Returns: number of photon candidates happened to be inside the window
179
180   
181   Int_t steps = (Int_t)((ckov )/ fDTheta); //how many times we need to have fDTheta to fill the distance between 0  and thetaCkovHough
182
183   Double_t tmin = (Double_t)(steps - 1)*fDTheta;
184   Double_t tmax = (Double_t)(steps)*fDTheta;
185   Double_t tavg = 0.5*(tmin+tmax);
186
187   tmin = tavg - 0.5*fWindowWidth;  tmax = tavg + 0.5*fWindowWidth;
188
189   Int_t iInsideCnt = 0; //count photons which Theta ckov inside the window
190   for(Int_t i=0;i<fPhotCnt;i++){//photon candidates loop
191     if(fPhotCkov[i] >= tmin && fPhotCkov[i] <= tmax)    { 
192       fPhotFlag[i]=2;     
193       iInsideCnt++;
194     }
195   }
196   return iInsideCnt;
197 }//FlagPhot()
198 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
199 Double_t AliRICHRecon::TracePhot(Double_t ckovThe,Double_t ckovPhi,TVector2 &pos)const
200 {
201 // Trace a single Ckov photon from emission point somewhere in radiator up to photocathode taking into account ref indexes of materials it travereses
202 // Arguments: ckovThe,ckovPhi- photon ckov angles, [rad]  (warning: not photon theta and phi)     
203 //   Returns: distance between photon point on PC and track projection  
204   TRotation mtheta;   mtheta.RotateY(fTrkDir.Theta());
205   TRotation mphi;       mphi.RotateZ(fTrkDir.Phi());  
206   TRotation mrot=mphi*mtheta;
207   
208   TVector3  posCkov(fTrkPos.X(),fTrkPos.Y(),-0.5*fgkRadThick-fgkWinThick-fgkGapThick);   //RAD: photon position is track position @ middle of RAD 
209   TVector3  dirCkov;   dirCkov.SetMagThetaPhi(1,ckovThe,ckovPhi);                        //initially photon is directed according to requested ckov angle
210                                                dirCkov=mrot*dirCkov;                     //now we know photon direction in LORS
211                        dirCkov.SetPhi(ckovPhi);   
212   if(dirCkov.Theta() > TMath::ASin(1./fgkRadIdx)) return -999;//total refraction on WIN-GAP boundary
213   
214   Propagate(dirCkov,posCkov,-fgkWinThick-fgkGapThick); //go to RAD-WIN boundary  remeber that z=0 is PC plane
215   Refract  (dirCkov,         fgkRadIdx,fgkWinIdx    ); //RAD-WIN refraction
216   Propagate(dirCkov,posCkov,-fgkGapThick           );  //go to WIN-GAP boundary
217   Refract  (dirCkov,         fgkWinIdx,fgkGapIdx    ); //WIN-GAP refraction
218   Propagate(dirCkov,posCkov,0                     );   //go to PC
219   
220   pos.Set(posCkov.X(),posCkov.Y());
221   return (pos-fTrkPos).Mod();
222 }//TracePhoton()
223 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
224 void AliRICHRecon::Propagate(const TVector3 &dir,TVector3 &pos,Double_t z)const
225 {
226 // Finds an intersection point between a line and XY plane shifted along Z.
227 // Arguments:  dir,pos   - vector along the line and any point of the line
228 //             z         - z coordinate of plain 
229 //   Returns:  none
230 //   On exit:  pos is the position if this intesection if any
231   static TVector3 nrm(0,0,1); 
232          TVector3 pnt(0,0,z);
233   
234   TVector3 diff=pnt-pos;
235   Double_t sint=(nrm*diff)/(nrm*dir);
236   pos+=sint*dir;
237 }//Propagate()
238 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
239 void AliRICHRecon::Refract(TVector3 &dir,Double_t n1,Double_t n2)const
240 {
241 // Refract direction vector according to Snell law
242 // Arguments: 
243 //            n1 - ref idx of first substance
244 //            n2 - ref idx of second substance
245 //   Returns: none
246 //   On exit: dir is new direction
247   Double_t sinref=(n1/n2)*TMath::Sin(dir.Theta());
248   if(sinref>1.)    dir.SetXYZ(-999,-999,-999);
249   else             dir.SetTheta(TMath::ASin(sinref));
250 }//Refract()
251 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
252 Double_t AliRICHRecon::HoughResponse()
253 {
254 //
255 //
256 //       
257   Double_t kThetaMax=0.75;
258   Int_t nChannels = (Int_t)(kThetaMax/fDTheta+0.5);
259   TH1D *phots   = new TH1D("Rphot"  ,"phots"         ,nChannels,0,kThetaMax);
260   TH1D *photsw  = new TH1D("RphotWeighted" ,"photsw" ,nChannels,0,kThetaMax);
261   TH1D *resultw = new TH1D("resultw","resultw"       ,nChannels,0,kThetaMax);
262   Int_t nBin = (Int_t)(kThetaMax/fDTheta);
263   Int_t nCorrBand = (Int_t)(fWindowWidth/(2*fDTheta));
264   
265   for (Int_t i=0; i< fPhotCnt; i++){//photon cadidates loop
266     Double_t angle = fPhotCkov[i];  if(angle<0||angle>kThetaMax) continue;
267     phots->Fill(angle);
268     Int_t bin = (Int_t)(0.5+angle/(fDTheta));
269     Double_t weight=1.;
270     if(fIsWEIGHT){
271       Double_t lowerlimit = ((Double_t)bin)*fDTheta - 0.5*fDTheta;  Double_t upperlimit = ((Double_t)bin)*fDTheta + 0.5*fDTheta;   
272       Double_t diffArea = FindRingArea(upperlimit)-FindRingArea(lowerlimit);
273       if(diffArea>0) weight = 1./diffArea;
274     }
275     photsw->Fill(angle,weight);
276     fPhotWei[i]=weight;
277   }//photon candidates loop 
278    
279   for (Int_t i=1; i<=nBin;i++){
280     Int_t bin1= i-nCorrBand;
281     Int_t bin2= i+nCorrBand;
282     if(bin1<1) bin1=1;
283     if(bin2>nBin)bin2=nBin;
284     Double_t sumPhots=phots->Integral(bin1,bin2);
285     if(sumPhots<3) continue;                            // if less then 3 photons don't trust to this ring
286     Double_t sumPhotsw=photsw->Integral(bin1,bin2);
287     resultw->Fill((Double_t)((i+0.5)*fDTheta),sumPhotsw);
288   } 
289 // evaluate the "BEST" theta ckov as the maximum value of histogramm
290   Double_t *pVec = resultw->GetArray();
291   Int_t locMax = TMath::LocMax(nBin,pVec);
292   phots->Delete();photsw->Delete();resultw->Delete(); // Reset and delete objects
293   
294   return (Double_t)(locMax*fDTheta+0.5*fDTheta); //final most probable track theta ckov   
295 }//HoughResponse()
296 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
297 Double_t AliRICHRecon::Sigma2(Double_t ckovTh, Double_t ckovPh)const
298 {
299 // Analithical calculation of total error (as a sum of localization, geometrical and chromatic errors) on Cerenkov angle for a given Cerenkov photon 
300 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
301 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
302 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
303 //            MIP beta
304 //   Returns: absolute error on Cerenkov angle, [radians]    
305   
306   TVector3 v(-999,-999,-999);
307   Double_t trkBeta = 1./(TMath::Cos(ckovTh)*fgkRadIdx);
308
309   v.SetX(SigLoc (ckovTh,ckovPh,trkBeta));
310   v.SetY(SigGeom(ckovTh,ckovPh,trkBeta));
311   v.SetZ(SigCrom(ckovTh,ckovPh,trkBeta));
312
313   return v.Mag2();
314 }
315 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
316 Double_t AliRICHRecon::SigLoc(Double_t thetaC, Double_t phiC,Double_t betaM)const
317 {
318 // Analithical calculation of localization error (due to finite segmentation of PC) on Cerenkov angle for a given Cerenkov photon 
319 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
320 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
321 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
322 //            MIP beta
323 //   Returns: absolute error on Cerenkov angle, [radians]    
324   Double_t phiDelta = phiC - fTrkDir.Phi();
325
326   Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
327   Double_t k = 1.-fgkRadIdx*fgkRadIdx+alpha*alpha/(betaM*betaM);
328   if (k<0) return 1e10;
329
330   Double_t mu =TMath::Sin(fTrkDir.Theta())*TMath::Sin(fTrkDir.Phi())+TMath::Tan(thetaC)*(TMath::Cos(fTrkDir.Theta())*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Phi())+TMath::Sin(phiDelta)*TMath::Cos(fTrkDir.Phi()));
331   Double_t e  =TMath::Sin(fTrkDir.Theta())*TMath::Cos(fTrkDir.Phi())+TMath::Tan(thetaC)*(TMath::Cos(fTrkDir.Theta())*TMath::Cos(phiDelta)*TMath::Cos(fTrkDir.Phi())-TMath::Sin(phiDelta)*TMath::Sin(fTrkDir.Phi()));
332
333   Double_t kk = betaM*TMath::Sqrt(k)/(8*alpha);
334   Double_t dtdxc = kk*(k*(TMath::Cos(phiDelta)*TMath::Cos(fTrkDir.Phi())-TMath::Cos(fTrkDir.Theta())*TMath::Sin(phiDelta)*TMath::Sin(fTrkDir.Phi()))-(alpha*mu/(betaM*betaM))*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiDelta));
335   Double_t dtdyc = kk*(k*(TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Phi())+TMath::Cos(fTrkDir.Theta())*TMath::Sin(phiDelta)*TMath::Cos(fTrkDir.Phi()))+(alpha* e/(betaM*betaM))*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiDelta));
336
337   return  TMath::Sqrt(0.2*0.2*dtdxc*dtdxc + 0.25*0.25*dtdyc*dtdyc);
338 }
339 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
340 Double_t AliRICHRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)const
341 {
342 // Analithical calculation of chromatic error (due to lack of knowledge of Cerenkov photon energy) on Cerenkov angle for a given Cerenkov photon 
343 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
344 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
345 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
346 //            MIP beta
347 //   Returns: absolute error on Cerenkov angle, [radians]    
348   Double_t phiDelta = phiC - fTrkDir.Phi();
349   Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
350
351   Double_t dtdn = TMath::Cos(fTrkDir.Theta())*fgkRadIdx*betaM*betaM/(alpha*TMath::Tan(thetaC));
352             
353   Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
354
355   return f*dtdn;
356 }//SigCrom()
357 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
358 Double_t AliRICHRecon::SigGeom(Double_t thetaC, Double_t phiC,Double_t betaM)const
359 {
360 // Analithical calculation of geometric error (due to lack of knowledge of creation point in radiator) on Cerenkov angle for a given Cerenkov photon 
361 // created by a given MIP. Formulae according to CERN-EP-2000-058 
362 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
363 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
364 //            MIP beta
365 //   Returns: absolute error on Cerenkov angle, [radians]    
366
367   Double_t phiDelta = phiC - fTrkDir.Phi();
368   Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
369
370   Double_t k = 1.-fgkRadIdx*fgkRadIdx+alpha*alpha/(betaM*betaM);
371   if (k<0) return 1e10;
372
373   Double_t eTr = 0.5*1.5*betaM*TMath::Sqrt(k)/(8*alpha);
374   Double_t lambda = 1.-TMath::Sin(fTrkDir.Theta())*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiC)*TMath::Sin(phiC);
375
376   Double_t c = 1./(1.+ eTr*k/(alpha*alpha*TMath::Cos(thetaC)*TMath::Cos(thetaC)));
377   Double_t i = betaM*TMath::Tan(thetaC)*lambda*TMath::Power(k,1.5);
378   Double_t ii = 1.+eTr*betaM*i;
379
380   Double_t err = c * (i/(alpha*alpha*8) +  ii*(1.-lambda) / ( alpha*alpha*8*betaM*(1.+eTr)) );
381   Double_t trErr = 1.5/(TMath::Sqrt(12.)*TMath::Cos(fTrkDir.Theta()));
382
383   return trErr*err;
384 }//SigGeom()
385 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
386 void AliRICHRecon::Display() 
387 {
388 // Display digits, reconstructed tracks intersections and RICH rings if available 
389 // Arguments: none
390 //   Returns: none    
391   TFile *pEsdFl=TFile::Open("AliESDs.root");     if(!pEsdFl || !pEsdFl->IsOpen()) return;//open AliESDs.root                                                                    
392   TTree *pEsdTr=(TTree*) pEsdFl->Get("esdTree"); if(!pEsdTr)                      return;//get ESD tree
393                                                                  
394   AliESD *pEsd=new AliESD;  pEsdTr->SetBranchAddress("ESD", &pEsd);
395   
396   TPolyMarker  *pDigMap[7]; //digits map
397   TPolyMarker  *pTrkMap[7]; Int_t aTrkCnt[7]; //TRKxPC intersection map
398   
399   for(Int_t i=0;i<7;i++){
400                   pDigMap[i]=new TPolyMarker(); pDigMap[i]->SetMarkerStyle(25); pDigMap[i]->SetMarkerSize(0.5); pDigMap[i]->SetMarkerColor(kGreen); 
401     aTrkCnt[i]=0; pTrkMap[i]=new TPolyMarker(); pTrkMap[i]->SetMarkerStyle(4);  pTrkMap[i]->SetMarkerSize(0.5); pTrkMap[i]->SetMarkerColor(kRed); 
402   }
403
404   AliRICHRecon rec;
405   
406   TLatex t;
407   TCanvas *pC = new TCanvas("RICHDisplay","RICH Display",0,0,1226,900);  pC->Divide(3,3);
408   
409   for(Int_t iEvt=0;iEvt<pEsdTr->GetEntries();iEvt++) {                //events loop
410     pC->cd(3);  t.DrawText(0.2,0.4,Form("Event %i",iEvt));        //print current event number
411     pEsdTr->GetEntry(iEvt);                                       //get ESD for this event   
412     for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){//ESD tracks loop
413       AliESDtrack *pTrk = pEsd->GetTrack(iTrk);             //
414       Float_t th,ph,x,y; pTrk->GetRICHtrk(x,y,th,ph); if(x<0) continue;
415       Int_t ch=pTrk->GetRICHcluIdx()/1000000; Printf("ch=%i",ch);
416       pTrkMap[ch]->SetPoint(aTrkCnt[ch]++,x,y);
417     }//ESD tracks loop
418     
419 //     al->GetEvent(iEvt);   rl->TreeD()->GetEntry(0); //get digits list
420     for(Int_t iCh=0;iCh<7;iCh++) {//chambers loop
421 //       for(Int_t iDig=0;iDig < r->DigLst(iCh)->GetEntries();iDig++) {      //digits loop
422 //         AliRICHDigit *pDig = (AliRICHDigit*)r->DigLst(iCh)->At(iDig);     
423 //         pDigMap[iCh]->SetPoint(iDig,pDig->LorsX(),pDig->LorsY());
424 //       }                                                             //digits loop
425 // 
426 //       
427       if(iCh==6) pC->cd(1); if(iCh==5) pC->cd(2);
428       if(iCh==4) pC->cd(4); if(iCh==3) pC->cd(5); if(iCh==2) pC->cd(6);
429                             if(iCh==1) pC->cd(8); if(iCh==0) pC->cd(9);
430                           
431       AliRICHDigit::DrawPc();  pTrkMap[iCh]->Draw(); pDigMap[iCh]->Draw();
432     }//chambers loop
433 //    pC->Update();
434 //    pC->Modified();
435 //    if(iEvt<iEvtTo) {gPad->WaitPrimitive();pC->Clear();}
436     
437     
438     
439   }//events loop
440   delete pEsd;  pEsdFl->Close();//close AliESDs.root
441 //  rl->UnloadDigits();
442 }//Display()