]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDRecon.cxx
Magnetic field properly set in display.
[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 "AliHMPIDParam.h"   //CkovAngle()
26 #include "AliHMPIDCluster.h" //CkovAngle()
27 #include <TRotation.h>       //TracePhot()
28 #include <TH1D.h>            //HoughResponse()
29 #include <TClonesArray.h>    //CkovAngle()
30 #include <AliESDtrack.h>     //CkovAngle()
31
32 const Double_t AliHMPIDRecon::fgkRadThick=1.5;
33 const Double_t AliHMPIDRecon::fgkWinThick=0.5;
34 const Double_t AliHMPIDRecon::fgkGapThick=8.0;
35 const Double_t AliHMPIDRecon::fgkWinIdx  =1.5787;
36 const Double_t AliHMPIDRecon::fgkGapIdx  =1.0005;
37
38 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
39 AliHMPIDRecon::AliHMPIDRecon():TTask("RichRec","RichPat"),
40   fRadNmean(1.292),  
41   fPhotCnt(-1),
42   fCkovSigma2(0),
43   fIsWEIGHT(kFALSE),
44   fDTheta(0.001),
45   fWindowWidth(0.045),
46   fTrkDir(TVector3(0,0,1)),fTrkPos(TVector2(30,40))  
47 {
48 // main ctor
49   for (Int_t i=0; i<3000; i++) {
50     fPhotFlag[i] =  0;
51     fPhotCkov[i] = -1;
52     fPhotPhi [i] = -1;
53     fPhotWei [i] =  0;
54   }
55 }
56 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
57 void AliHMPIDRecon::CkovAngle(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean)
58 {
59 // Pattern recognition method based on Hough transform
60 // Arguments:   pTrk     - track for which Ckov angle is to be found
61 //              pCluLst  - list of clusters for this chamber   
62 //   Returns:            - track ckov angle, [rad], 
63     
64   AliHMPIDParam *pParam=AliHMPIDParam::Instance();
65   
66   if(pCluLst->GetEntries()>pParam->MultCut()) fIsWEIGHT = kTRUE; // offset to take into account bkg in reconstruction
67   else                                        fIsWEIGHT = kFALSE;
68
69   Float_t xRa,yRa,th,ph;       
70   pTrk->GetHMPIDtrk(xRa,yRa,th,ph);        //initialize this track: th and ph angles at middle of RAD 
71   SetTrack(xRa,yRa,th,ph);
72   
73   fRadNmean=nmean;
74
75   Float_t dMin=999,mipX=-1,mipY=-1;Int_t chId=-1,mipId=-1,mipQ=-1;                                                                           
76   fPhotCnt=0;                                                      
77   for (Int_t iClu=0; iClu<pCluLst->GetEntriesFast();iClu++){//clusters loop
78     AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu);                       //get pointer to current cluster    
79     chId=pClu->Ch();
80     if(pClu->Q()>pParam->QCut()){                                                             //charge compartible with MIP clusters      
81       Float_t dX=fPc.X()-pClu->X(),dY=fPc.Y()-pClu->Y(),d =TMath::Sqrt(dX*dX+dY*dY);          //distance between current cluster and intersection point
82       if( d < dMin) {mipId=iClu; dMin=d;mipX=pClu->X();mipY=pClu->Y();mipQ=(Int_t)pClu->Q();} //current cluster is closer, overwrite data for min cluster
83     }else{                                                                                    //charge compatible with photon cluster
84       Double_t thetaCer,phiCer;
85       if(FindPhotCkov(pClu->X(),pClu->Y(),thetaCer,phiCer)){                                  //find ckov angle for this  photon candidate
86         fPhotCkov[fPhotCnt]=thetaCer;                                                         //actual theta Cerenkov (in TRS)
87         fPhotPhi [fPhotCnt]=phiCer;                                                           //actual phi   Cerenkov (in TRS): -pi to come back to "unusual" ref system (X,Y,-Z)
88         fPhotCnt++;                                                                           //increment counter of photon candidates
89       }
90     }
91   }//clusters loop
92   Int_t iNacc=FlagPhot(HoughResponse());                                                      //flag photons according to individual theta ckov with respect to most probable
93   pTrk->SetHMPIDmip(mipX,mipY,mipQ,iNacc);                                                    //store mip info 
94
95   if(mipId==-1)              {pTrk->SetHMPIDsignal(kMipQdcCut);  return;}                     //no clusters with QDC more the threshold at all
96   if(dMin>pParam->DistCut()) {pTrk->SetHMPIDsignal(kMipDistCut); return;}                     //closest cluster with enough charge is still too far from intersection
97   pTrk->SetHMPIDcluIdx(chId,mipId);                                                           //set index of cluster
98   if(iNacc<1)    pTrk->SetHMPIDsignal(kNoPhotAccept);                                         //no photon candidates is accepted
99   else           pTrk->SetHMPIDsignal(FindRingCkov(pCluLst->GetEntries()));                   //find best Theta ckov for ring i.e. track
100   
101   pTrk->SetHMPIDchi2(fCkovSigma2);                                                            //errors squared 
102
103 }//ThetaCerenkov()
104 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
105 Bool_t AliHMPIDRecon::FindPhotCkov(Double_t cluX,Double_t cluY,Double_t &thetaCer,Double_t &phiCer)
106 {
107 // Finds Cerenkov angle  for this photon candidate
108 // Arguments: cluX,cluY - position of cadidate's cluster  
109 // Returns: Cerenkov angle 
110
111   TVector3 dirCkov;
112   
113   Double_t zRad= -0.5*fgkRadThick-0.5*fgkWinThick;                   //z position of middle of RAD
114   TVector3 rad(fTrkPos.X(),fTrkPos.Y(),zRad);                        //impact point at middle of RAD
115   TVector3  pc(cluX,cluY,0.5*fgkWinThick+fgkGapIdx);                 //mip at PC
116   Double_t cluR = TMath::Sqrt((cluX-fTrkPos.X())*(cluX-fTrkPos.X())+
117                               (cluY-fTrkPos.Y())*(cluY-fTrkPos.Y()));//ref. distance impact RAD-CLUSTER   
118   Double_t phi=(pc-rad).Phi();                                       //phi of photon
119     
120   Double_t ckov1=0;
121   Double_t ckov2=0.75+fTrkDir.Theta();                        //start to find theta cerenkov in DRS
122   const Double_t kTol=0.01;
123   Int_t iIterCnt = 0;
124   while(1){
125     if(iIterCnt>=50) return kFALSE;
126     Double_t ckov=0.5*(ckov1+ckov2);
127     dirCkov.SetMagThetaPhi(1,ckov,phi);
128     TVector2 posC=TraceForward(dirCkov);                      //trace photon with actual angles
129     Double_t dist=cluR-(posC-fTrkPos).Mod();                  //get distance between trial point and cluster position
130     if(posC.X()==-999) dist = - 999;                          //total reflection problem
131     iIterCnt++;                                               //counter step
132     if     (dist> kTol) ckov1=ckov;                           //cluster @ larger ckov
133     else if(dist<-kTol) ckov2=ckov;                           //cluster @ smaller ckov
134     else{                                                     //precision achived: ckov in DRS found
135       dirCkov.SetMagThetaPhi(1,ckov,phi);                     //
136       RecPhot(dirCkov,thetaCer,phiCer);                       //find ckov (in TRS:the effective Cherenkov angle!)
137       return kTRUE;
138     }
139   }
140 }//FindPhotTheta()
141 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
142 TVector2 AliHMPIDRecon::TraceForward(TVector3 dirCkov)const
143 {
144   //Trace forward a photon from (x,y) up to PC
145   // Arguments: dirCkov photon vector in LORS
146   //   Returns: pos of traced photon at PC
147   TVector2 pos(-999,-999);
148   Double_t thetaCer = dirCkov.Theta();
149   if(thetaCer > TMath::ASin(1./fRadNmean))  return pos;         //total refraction on WIN-GAP boundary
150   Double_t zRad= -0.5*fgkRadThick-0.5*fgkWinThick;              //z position of middle of RAD
151   TVector3  posCkov(fTrkPos.X(),fTrkPos.Y(),zRad);              //RAD: photon position is track position @ middle of RAD 
152   Propagate(dirCkov,posCkov,           -0.5*fgkWinThick);       //go to RAD-WIN boundary  
153   Refract  (dirCkov,         fRadNmean,fgkWinIdx);              //RAD-WIN refraction
154   Propagate(dirCkov,posCkov,            0.5*fgkWinThick);       //go to WIN-GAP boundary
155   Refract  (dirCkov,         fgkWinIdx,fgkGapIdx);              //WIN-GAP refraction
156   Propagate(dirCkov,posCkov,0.5*fgkWinThick+fgkGapThick);       //go to PC
157   pos.Set(posCkov.X(),posCkov.Y());
158   return pos;
159 }//TraceForward()
160 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
161 void AliHMPIDRecon::RecPhot(TVector3 dirCkov,Double_t &thetaCer,Double_t &phiCer)
162 {
163   //Theta Cerenkov reconstruction 
164   // Arguments: (x,y) of initial point in LORS, dirCkov photon vector in LORS
165   //   Returns: thetaCer theta cerenkov reconstructed
166 //  TVector3 dirTrk;
167 //  dirTrk.SetMagThetaPhi(1,fTrkDir.Theta(),fTrkDir.Phi());
168 //  Double_t thetaCer = TMath::ACos(dirCkov*dirTrk);
169   TRotation mtheta;   mtheta.RotateY(- fTrkDir.Theta());
170   TRotation mphi;       mphi.RotateZ(- fTrkDir.Phi());
171   TRotation mrot=mtheta*mphi;
172   TVector3 dirCkovTRS;
173   dirCkovTRS=mrot*dirCkov;
174   phiCer  = dirCkovTRS.Phi();                                          //actual value of the phi of the photon
175   thetaCer= dirCkovTRS.Theta();                                        //actual value of thetaCerenkov of the photon
176 }
177 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
178 Double_t AliHMPIDRecon::FindRingArea(Double_t ckovAng)const
179 {
180 // Find area inside the cerenkov ring which lays inside PCs
181 // Arguments: ckovAng - cerenkov angle    
182 //   Returns: area of the ring in cm^2 for given theta ckov
183    
184   const Int_t kN=100;
185   Double_t area=0;
186   for(Int_t i=0;i<kN;i++){
187     TVector2 pos1=TracePhot(ckovAng,Double_t(TMath::TwoPi()*i    /kN));//trace this photon 
188     TVector2 pos2=TracePhot(ckovAng,Double_t(TMath::TwoPi()*(i+1)/kN));//trace the next photon 
189     area+=(pos1-fTrkPos)*(pos2-fTrkPos);                               //add area of the triangle... 
190   }
191   return area;
192 }//FindRingArea()
193 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
194 Double_t AliHMPIDRecon::FindRingCkov(Int_t)
195 {
196 // 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
197 // collecting errors for all single Ckov candidates thetas. (Assuming they are independent)  
198 // Arguments: iNclus- total number of clusters in chamber for background estimation
199 //    Return: best estimation of track Theta ckov
200
201   Double_t wei = 0.;
202   Double_t weightThetaCerenkov = 0.;
203
204   Double_t ckovMin=9999.,ckovMax=0.;
205   Double_t sigma2 = 0;   //to collect error squared for this ring
206   
207   for(Int_t i=0;i<fPhotCnt;i++){//candidates loop
208     if(fPhotFlag[i] == 2){
209       if(fPhotCkov[i]<ckovMin) ckovMin=fPhotCkov[i];                         //find max and min Theta ckov from all candidates within probable window
210       if(fPhotCkov[i]>ckovMax) ckovMax=fPhotCkov[i]; 
211       weightThetaCerenkov += fPhotCkov[i]*fPhotWei[i];
212       wei += fPhotWei[i];                                                    //collect weight as sum of all candidate weghts   
213       
214       sigma2 += 1./Sigma2(fPhotCkov[i],fPhotPhi[i]);
215     }
216   }//candidates loop
217   
218   if(sigma2>0) fCkovSigma2=1./sigma2;
219   else         fCkovSigma2=1e10;  
220   
221   if(wei != 0.) weightThetaCerenkov /= wei; else weightThetaCerenkov = 0.;
222   return weightThetaCerenkov;
223 }//FindCkovRing()
224 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
225 Int_t AliHMPIDRecon::FlagPhot(Double_t ckov)
226 {
227 // Flag photon candidates if their individual ckov angle is inside the window around ckov angle returned by  HoughResponse()
228 // Arguments: ckov- value of most probable ckov angle for track as returned by HoughResponse()
229 //   Returns: number of photon candidates happened to be inside the window
230
231 // Photon Flag:  Flag = 0 initial set; 
232 //               Flag = 1 good candidate (charge compatible with photon); 
233 //               Flag = 2 photon used for the ring;
234   
235   Int_t steps = (Int_t)((ckov )/ fDTheta); //how many times we need to have fDTheta to fill the distance between 0  and thetaCkovHough
236
237   Double_t tmin = (Double_t)(steps - 1)*fDTheta;
238   Double_t tmax = (Double_t)(steps)*fDTheta;
239   Double_t tavg = 0.5*(tmin+tmax);
240
241   tmin = tavg - 0.5*fWindowWidth;  tmax = tavg + 0.5*fWindowWidth;
242
243   Int_t iInsideCnt = 0; //count photons which Theta ckov inside the window
244   for(Int_t i=0;i<fPhotCnt;i++){//photon candidates loop
245     if(fPhotCkov[i] >= tmin && fPhotCkov[i] <= tmax)    { 
246       fPhotFlag[i]=2;     
247       iInsideCnt++;
248     }
249   }
250   return iInsideCnt;
251 }//FlagPhot()
252 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
253 TVector2 AliHMPIDRecon::TracePhot(Double_t ckovThe,Double_t ckovPhi)const
254 {
255 // Trace a single Ckov photon from emission point somewhere in radiator up to photocathode taking into account ref indexes of materials it travereses
256 // Arguments: ckovThe,ckovPhi- photon ckov angles in DRS, [rad]    
257 //   Returns: distance between photon point on PC and track projection  
258   TRotation mtheta;   mtheta.RotateY(fTrkDir.Theta());
259   TRotation mphi;       mphi.RotateZ(fTrkDir.Phi());  
260   TRotation mrot=mphi*mtheta;
261   TVector3  dirCkov,dirCkovTors;   
262
263   dirCkovTors.SetMagThetaPhi(1,ckovThe,ckovPhi);                    //initially photon is directed according to requested ckov angle
264   dirCkov=mrot*dirCkovTors;                                         //now we know photon direction in LORS
265   return TraceForward(dirCkov);
266 }//TracePhot()
267 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
268 void AliHMPIDRecon::Propagate(const TVector3 dir,TVector3 &pos,Double_t z)const
269 {
270 // Finds an intersection point between a line and XY plane shifted along Z.
271 // Arguments:  dir,pos   - vector along the line and any point of the line
272 //             z         - z coordinate of plain 
273 //   Returns:  none
274 //   On exit:  pos is the position if this intesection if any
275   static TVector3 nrm(0,0,1); 
276          TVector3 pnt(0,0,z);
277   
278   TVector3 diff=pnt-pos;
279   Double_t sint=(nrm*diff)/(nrm*dir);
280   pos+=sint*dir;
281 }//Propagate()
282 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
283 void AliHMPIDRecon::Refract(TVector3 &dir,Double_t n1,Double_t n2)const
284 {
285 // Refract direction vector according to Snell law
286 // Arguments: 
287 //            n1 - ref idx of first substance
288 //            n2 - ref idx of second substance
289 //   Returns: none
290 //   On exit: dir is new direction
291   Double_t sinref=(n1/n2)*TMath::Sin(dir.Theta());
292   if(sinref>1.)    dir.SetXYZ(-999,-999,-999);
293   else             dir.SetTheta(TMath::ASin(sinref));
294 }//Refract()
295 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
296 Double_t AliHMPIDRecon::HoughResponse()
297 {
298 //
299 //
300 //       
301   Double_t kThetaMax=0.75;
302   Int_t nChannels = (Int_t)(kThetaMax/fDTheta+0.5);
303   TH1D *phots   = new TH1D("Rphot"  ,"phots"         ,nChannels,0,kThetaMax);
304   TH1D *photsw  = new TH1D("RphotWeighted" ,"photsw" ,nChannels,0,kThetaMax);
305   TH1D *resultw = new TH1D("resultw","resultw"       ,nChannels,0,kThetaMax);
306   Int_t nBin = (Int_t)(kThetaMax/fDTheta);
307   Int_t nCorrBand = (Int_t)(fWindowWidth/(2*fDTheta));
308   
309   for (Int_t i=0; i< fPhotCnt; i++){//photon cadidates loop
310     Double_t angle = fPhotCkov[i];  if(angle<0||angle>kThetaMax) continue;
311     phots->Fill(angle);
312     Int_t bin = (Int_t)(0.5+angle/(fDTheta));
313     Double_t weight=1.;
314     if(fIsWEIGHT){
315       Double_t lowerlimit = ((Double_t)bin)*fDTheta - 0.5*fDTheta;  Double_t upperlimit = ((Double_t)bin)*fDTheta + 0.5*fDTheta;   
316       Double_t diffArea = FindRingArea(upperlimit)-FindRingArea(lowerlimit);
317       if(diffArea>0) weight = 1./diffArea;
318     }
319     photsw->Fill(angle,weight);
320     fPhotWei[i]=weight;
321   }//photon candidates loop 
322    
323   for (Int_t i=1; i<=nBin;i++){
324     Int_t bin1= i-nCorrBand;
325     Int_t bin2= i+nCorrBand;
326     if(bin1<1) bin1=1;
327     if(bin2>nBin)bin2=nBin;
328     Double_t sumPhots=phots->Integral(bin1,bin2);
329     if(sumPhots<3) continue;                            // if less then 3 photons don't trust to this ring
330     Double_t sumPhotsw=photsw->Integral(bin1,bin2);
331     resultw->Fill((Double_t)((i+0.5)*fDTheta),sumPhotsw);
332   } 
333 // evaluate the "BEST" theta ckov as the maximum value of histogramm
334   Double_t *pVec = resultw->GetArray();
335   Int_t locMax = TMath::LocMax(nBin,pVec);
336   phots->Delete();photsw->Delete();resultw->Delete(); // Reset and delete objects
337   
338   return (Double_t)(locMax*fDTheta+0.5*fDTheta); //final most probable track theta ckov   
339 }//HoughResponse()
340 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
341 Double_t AliHMPIDRecon::Sigma2(Double_t ckovTh, Double_t ckovPh)const
342 {
343 // Analithical calculation of total error (as a sum of localization, geometrical and chromatic errors) on Cerenkov angle for a given Cerenkov photon 
344 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
345 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
346 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
347 //            MIP beta
348 //   Returns: absolute error on Cerenkov angle, [radians]    
349   
350   TVector3 v(-999,-999,-999);
351   Double_t trkBeta = 1./(TMath::Cos(ckovTh)*fRadNmean);
352
353   v.SetX(SigLoc (ckovTh,ckovPh,trkBeta));
354   v.SetY(SigGeom(ckovTh,ckovPh,trkBeta));
355   v.SetZ(SigCrom(ckovTh,ckovPh,trkBeta));
356
357   return v.Mag2();
358 }
359 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
360 Double_t AliHMPIDRecon::SigLoc(Double_t thetaC, Double_t phiC,Double_t betaM)const
361 {
362 // Analithical calculation of localization error (due to finite segmentation of PC) on Cerenkov angle for a given Cerenkov photon 
363 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
364 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
365 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
366 //            MIP beta
367 //   Returns: absolute error on Cerenkov angle, [radians]    
368   Double_t phiDelta = phiC - fTrkDir.Phi();
369
370   Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
371   Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
372   if (k<0) return 1e10;
373
374   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()));
375   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()));
376
377   Double_t kk = betaM*TMath::Sqrt(k)/(8*alpha);
378   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));
379   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));
380
381   return  TMath::Sqrt(0.2*0.2*dtdxc*dtdxc + 0.25*0.25*dtdyc*dtdyc);
382 }
383 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
384 Double_t AliHMPIDRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)const
385 {
386 // Analithical calculation of chromatic error (due to lack of knowledge of Cerenkov photon energy) on Cerenkov angle for a given Cerenkov photon 
387 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
388 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
389 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
390 //            MIP beta
391 //   Returns: absolute error on Cerenkov angle, [radians]    
392   Double_t phiDelta = phiC - fTrkDir.Phi();
393   Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
394
395   Double_t dtdn = TMath::Cos(fTrkDir.Theta())*fRadNmean*betaM*betaM/(alpha*TMath::Tan(thetaC));
396             
397   Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
398
399   return f*dtdn;
400 }//SigCrom()
401 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
402 Double_t AliHMPIDRecon::SigGeom(Double_t thetaC, Double_t phiC,Double_t betaM)const
403 {
404 // Analithical calculation of geometric error (due to lack of knowledge of creation point in radiator) on Cerenkov angle for a given Cerenkov photon 
405 // created by a given MIP. Formulae according to CERN-EP-2000-058 
406 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
407 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
408 //            MIP beta
409 //   Returns: absolute error on Cerenkov angle, [radians]    
410
411   Double_t phiDelta = phiC - fTrkDir.Phi();
412   Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
413
414   Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
415   if (k<0) return 1e10;
416
417   Double_t eTr = 0.5*1.5*betaM*TMath::Sqrt(k)/(8*alpha);
418   Double_t lambda = 1.-TMath::Sin(fTrkDir.Theta())*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiC)*TMath::Sin(phiC);
419
420   Double_t c = 1./(1.+ eTr*k/(alpha*alpha*TMath::Cos(thetaC)*TMath::Cos(thetaC)));
421   Double_t i = betaM*TMath::Tan(thetaC)*lambda*TMath::Power(k,1.5);
422   Double_t ii = 1.+eTr*betaM*i;
423
424   Double_t err = c * (i/(alpha*alpha*8) +  ii*(1.-lambda) / ( alpha*alpha*8*betaM*(1.+eTr)) );
425   Double_t trErr = 1.5/(TMath::Sqrt(12.)*TMath::Cos(fTrkDir.Theta()));
426
427   return trErr*err;
428 }//SigGeom()
429 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++