1 /**************************************************************************
2 * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
4 * Author: The ALICE Off-line Project. *
5 * Contributors are mentioned in the code where appropriate. *
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 **************************************************************************/
16 //////////////////////////////////////////////////////////////////////////
20 // HMPID class to perfom pattern recognition based on Hough transfrom //
21 // for single chamber //
22 //////////////////////////////////////////////////////////////////////////
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 #include <AliESDtrack.h> //CkovAngle()
31 const Double_t AliHMPIDRecon::fgkRadThick=1.5;
32 const Double_t AliHMPIDRecon::fgkWinThick=0.5;
33 const Double_t AliHMPIDRecon::fgkGapThick=8.0;
34 const Double_t AliHMPIDRecon::fgkWinIdx =1.5787;
35 const Double_t AliHMPIDRecon::fgkGapIdx =1.0005;
38 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
39 AliHMPIDRecon::AliHMPIDRecon():TTask("RichRec","RichPat"),
46 fTrkDir(TVector3(0,0,1)),fTrkPos(TVector2(30,40))
49 for (Int_t i=0; i<3000; i++) {
56 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
57 void AliHMPIDRecon::CkovAngle(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean)
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],
64 if(pCluLst->GetEntries()>200) fIsWEIGHT = kTRUE; // offset to take into account bkg in reconstruction
65 else fIsWEIGHT = kFALSE;
67 // Photon Flag: Flag = 0 initial set; Flag = 1 good candidate (charge compatible with photon); Flag = 2 photon used for the ring;
68 Float_t xPc,yPc,th,ph; pTrk->GetHMPIDtrk(xPc,yPc,th,ph); SetTrack(xPc,yPc,th,ph); //initialize this track
73 Float_t dMin=999,mipX=-1,mipY=-1;Int_t chId=-1,mipId=-1,mipQ=-1;
75 for (Int_t iClu=0; iClu<pCluLst->GetEntriesFast();iClu++){//clusters loop
76 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu); //get pointer to current cluster
78 if(pClu->Q()>100){ //charge compartible with MIP clusters
79 Float_t dX=xPc-pClu->X(),dY=yPc-pClu->Y(),d =TMath::Sqrt(dX*dX+dY*dY); //distance between current cluster and intersection point
80 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
81 }else{ //charge compartible with photon cluster
82 fPhotCkov[fPhotCnt]=FindPhotCkov(pClu->X(),pClu->Y()); //find ckov angle for this photon candidate
83 fPhotCnt++; //increment counter of photon candidates
86 Int_t iNacc=FlagPhot(HoughResponse()); //flag photons according to individual theta ckov with respect to most probable
88 pTrk->SetHMPIDmip (mipX,mipY,mipQ,iNacc); //store mip info
90 if(mipId==-1) {pTrk->SetHMPIDsignal (kMipQdcCut); return;} //no clusters with QDC more the threshold at all
91 if(dMin>1) {pTrk->SetHMPIDsignal (kMipDistCut); return;} //closest cluster with enough charge is still too far from intersection
92 pTrk->SetHMPIDcluIdx(chId,mipId);
93 if(iNacc<1) pTrk->SetHMPIDsignal(kNoPhotAccept); //no photon candidates is accepted
94 else pTrk->SetHMPIDsignal(FindRingCkov(pCluLst->GetEntries())); //find best Theta ckov for ring i.e. track
96 pTrk->SetHMPIDchi2(fCkovSigma2); //error squared
99 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
100 Double_t AliHMPIDRecon::FindPhotCkov(Double_t cluX,Double_t cluY)
102 // Finds Cerenkov angle for this photon candidate
103 // Arguments: cluX,cluY - position of cadidate's cluster
104 // Returns: Cerenkov angle
106 TVector2 pos(cluX,cluY); Double_t cluR=(pos-fTrkPos).Mod(); Double_t phi=FindPhotPhi(cluX,cluY);
107 Double_t ckov1=0,ckov2=0.75;
108 const Double_t kTol=0.05;
111 if(iIterCnt>=50) return -1;
112 Double_t ckov=0.5*(ckov1+ckov2);
113 Double_t dist=cluR-TracePhot(ckov,phi,pos); iIterCnt++; //get distance between trial point and cluster position
114 if (dist> kTol) ckov1=ckov; //cluster @ larger ckov
115 else if(dist<-kTol) ckov2=ckov; //cluster @ smaller ckov
116 else return ckov; //precision achived
119 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
120 Double_t AliHMPIDRecon::FindPhotPhi(Double_t cluX,Double_t cluY)
122 // Finds phi angle og photon candidate by considering the cluster's position of this candudate w.r.t track position
125 return fPhotPhi[fPhotCnt]=TMath::ATan2(cluY-fTrkPos.Y()-emiss*TMath::Tan(fTrkDir.Theta())*TMath::Sin(fTrkDir.Phi()),
126 cluX-fTrkPos.X()-emiss*TMath::Tan(fTrkDir.Theta())*TMath::Cos(fTrkDir.Phi()));
128 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
129 Double_t AliHMPIDRecon::FindRingArea(Double_t ckovAng)const
131 // Find area inside the cerenkov ring which lays inside PCs
132 // Arguments: ckovThe - cernkov
133 // Returns: area of the ring in cm^2 for given theta ckov
140 for(Int_t i=0;i<kN;i++){
141 TracePhot(ckovAng,Double_t(TMath::TwoPi()*i /kN),pos1);//trace this photon
142 TracePhot(ckovAng,Double_t(TMath::TwoPi()*(i+1)/kN),pos2);//trace this photon
143 area+=(pos1-fTrkPos)*(pos2-fTrkPos);
148 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
149 Double_t AliHMPIDRecon::FindRingCkov(Int_t)
151 // 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
152 // collecting errors for all single Ckov candidates thetas. (Assuming they are independent)
153 // Arguments: iNclus- total number of clusters in chamber for background estimation
154 // Return: best estimation of track Theta ckov
157 Double_t weightThetaCerenkov = 0.;
159 Double_t ckovMin=9999.,ckovMax=0.;
160 Double_t sigma2 = 0; //to collect error squared for this ring
162 for(Int_t i=0;i<fPhotCnt;i++){//candidates loop
163 if(fPhotFlag[i] == 2){
164 if(fPhotCkov[i]<=0) continue;//?????????????????Flag photos = 2 may imply CkovEta = 0??????????????
165 if(fPhotCkov[i]<ckovMin) ckovMin=fPhotCkov[i]; //find max and min Theta ckov from all candidates within probable window
166 if(fPhotCkov[i]>ckovMax) ckovMax=fPhotCkov[i];
167 weightThetaCerenkov += fPhotCkov[i]*fPhotWei[i]; wei += fPhotWei[i]; //collect weight as sum of all candidate weghts
169 sigma2 += 1./Sigma2(fPhotCkov[i],fPhotPhi[i]);
173 if(sigma2>0) fCkovSigma2=1./sigma2;
174 else fCkovSigma2=1e10;
176 if(wei != 0.) weightThetaCerenkov /= wei; else weightThetaCerenkov = 0.;
177 return weightThetaCerenkov;
179 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
180 Int_t AliHMPIDRecon::FlagPhot(Double_t ckov)
182 // Flag photon candidates if their individual ckov angle is inside the window around ckov angle returned by HoughResponse()
183 // Arguments: ckov- value of most probable ckov angle for track as returned by HoughResponse()
184 // Returns: number of photon candidates happened to be inside the window
187 Int_t steps = (Int_t)((ckov )/ fDTheta); //how many times we need to have fDTheta to fill the distance between 0 and thetaCkovHough
189 Double_t tmin = (Double_t)(steps - 1)*fDTheta;
190 Double_t tmax = (Double_t)(steps)*fDTheta;
191 Double_t tavg = 0.5*(tmin+tmax);
193 tmin = tavg - 0.5*fWindowWidth; tmax = tavg + 0.5*fWindowWidth;
195 Int_t iInsideCnt = 0; //count photons which Theta ckov inside the window
196 for(Int_t i=0;i<fPhotCnt;i++){//photon candidates loop
197 if(fPhotCkov[i] >= tmin && fPhotCkov[i] <= tmax) {
204 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
205 Double_t AliHMPIDRecon::TracePhot(Double_t ckovThe,Double_t ckovPhi,TVector2 &pos)const
207 // Trace a single Ckov photon from emission point somewhere in radiator up to photocathode taking into account ref indexes of materials it travereses
208 // Arguments: ckovThe,ckovPhi- photon ckov angles, [rad] (warning: not photon theta and phi)
209 // Returns: distance between photon point on PC and track projection
210 TRotation mtheta; mtheta.RotateY(fTrkDir.Theta());
211 TRotation mphi; mphi.RotateZ(fTrkDir.Phi());
212 TRotation mrot=mphi*mtheta;
214 TVector3 posCkov(fTrkPos.X(),fTrkPos.Y(),-0.5*fgkRadThick-fgkWinThick-fgkGapThick); //RAD: photon position is track position @ middle of RAD
215 TVector3 dirCkov; dirCkov.SetMagThetaPhi(1,ckovThe,ckovPhi); //initially photon is directed according to requested ckov angle
216 dirCkov=mrot*dirCkov; //now we know photon direction in LORS
217 dirCkov.SetPhi(ckovPhi);
218 if(dirCkov.Theta() > TMath::ASin(1./fRadNmean)) return -999;//total refraction on WIN-GAP boundary
220 Propagate(dirCkov,posCkov,-fgkWinThick-fgkGapThick); //go to RAD-WIN boundary remeber that z=0 is PC plane
221 Refract (dirCkov, fRadNmean,fgkWinIdx ); //RAD-WIN refraction
222 Propagate(dirCkov,posCkov,-fgkGapThick ); //go to WIN-GAP boundary
223 Refract (dirCkov, fgkWinIdx,fgkGapIdx ); //WIN-GAP refraction
224 Propagate(dirCkov,posCkov,0 ); //go to PC
226 pos.Set(posCkov.X(),posCkov.Y());
227 return (pos-fTrkPos).Mod();
229 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
230 void AliHMPIDRecon::Propagate(const TVector3 &dir,TVector3 &pos,Double_t z)const
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
236 // On exit: pos is the position if this intesection if any
237 static TVector3 nrm(0,0,1);
240 TVector3 diff=pnt-pos;
241 Double_t sint=(nrm*diff)/(nrm*dir);
244 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
245 void AliHMPIDRecon::Refract(TVector3 &dir,Double_t n1,Double_t n2)const
247 // Refract direction vector according to Snell law
249 // n1 - ref idx of first substance
250 // n2 - ref idx of second substance
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));
257 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
258 Double_t AliHMPIDRecon::HoughResponse()
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));
271 for (Int_t i=0; i< fPhotCnt; i++){//photon cadidates loop
272 Double_t angle = fPhotCkov[i]; if(angle<0||angle>kThetaMax) continue;
274 Int_t bin = (Int_t)(0.5+angle/(fDTheta));
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;
281 photsw->Fill(angle,weight);
283 }//photon candidates loop
285 for (Int_t i=1; i<=nBin;i++){
286 Int_t bin1= i-nCorrBand;
287 Int_t bin2= i+nCorrBand;
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);
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
300 return (Double_t)(locMax*fDTheta+0.5*fDTheta); //final most probable track theta ckov
302 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
303 Double_t AliHMPIDRecon::Sigma2(Double_t ckovTh, Double_t ckovPh)const
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]
310 // Returns: absolute error on Cerenkov angle, [radians]
312 TVector3 v(-999,-999,-999);
313 Double_t trkBeta = 1./(TMath::Cos(ckovTh)*fRadNmean);
315 v.SetX(SigLoc (ckovTh,ckovPh,trkBeta));
316 v.SetY(SigGeom(ckovTh,ckovPh,trkBeta));
317 v.SetZ(SigCrom(ckovTh,ckovPh,trkBeta));
321 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
322 Double_t AliHMPIDRecon::SigLoc(Double_t thetaC, Double_t phiC,Double_t betaM)const
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]
329 // Returns: absolute error on Cerenkov angle, [radians]
330 Double_t phiDelta = phiC - fTrkDir.Phi();
332 Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
333 Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
334 if (k<0) return 1e10;
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()));
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));
343 return TMath::Sqrt(0.2*0.2*dtdxc*dtdxc + 0.25*0.25*dtdyc*dtdyc);
345 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
346 Double_t AliHMPIDRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)const
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]
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());
357 Double_t dtdn = TMath::Cos(fTrkDir.Theta())*fRadNmean*betaM*betaM/(alpha*TMath::Tan(thetaC));
359 Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
363 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
364 Double_t AliHMPIDRecon::SigGeom(Double_t thetaC, Double_t phiC,Double_t betaM)const
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]
371 // Returns: absolute error on Cerenkov angle, [radians]
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());
376 Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
377 if (k<0) return 1e10;
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);
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;
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()));
391 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++