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