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