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