]> git.uio.no Git - u/mrichter/AliRoot.git/blame - HMPID/AliHMPIDRecon.cxx
Algorithm to reconstruct rings without tracking infos included
[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
a591e55f 24#include "AliHMPIDRecon.h" //class header
25#include "AliHMPIDParam.h" //CkovAngle()
d3da6dc4 26#include "AliHMPIDCluster.h" //CkovAngle()
43400d2d 27#include <TMinuit.h> //FitEllipse()
a591e55f 28#include <TRotation.h> //TracePhot()
29#include <TH1D.h> //HoughResponse()
30#include <TClonesArray.h> //CkovAngle()
31#include <AliESDtrack.h> //CkovAngle()
d3da6dc4 32
33const Double_t AliHMPIDRecon::fgkRadThick=1.5;
34const Double_t AliHMPIDRecon::fgkWinThick=0.5;
35const Double_t AliHMPIDRecon::fgkGapThick=8.0;
d3da6dc4 36const Double_t AliHMPIDRecon::fgkWinIdx =1.5787;
37const Double_t AliHMPIDRecon::fgkGapIdx =1.0005;
38
d3da6dc4 39//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
40AliHMPIDRecon::AliHMPIDRecon():TTask("RichRec","RichPat"),
abb5f786 41 fRadNmean(1.292),
d3da6dc4 42 fPhotCnt(-1),
43 fCkovSigma2(0),
44 fIsWEIGHT(kFALSE),
45 fDTheta(0.001),
46 fWindowWidth(0.045),
47 fTrkDir(TVector3(0,0,1)),fTrkPos(TVector2(30,40))
48{
49// main ctor
50 for (Int_t i=0; i<3000; i++) {
51 fPhotFlag[i] = 0;
52 fPhotCkov[i] = -1;
53 fPhotPhi [i] = -1;
54 fPhotWei [i] = 0;
55 }
56}
57//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
a591e55f 58void AliHMPIDRecon::CkovAngle(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean)
d3da6dc4 59{
60// Pattern recognition method based on Hough transform
59280a5a 61// Arguments: pTrk - track for which Ckov angle is to be found
62// pCluLst - list of clusters for this chamber
63// Returns: - track ckov angle, [rad],
a591e55f 64
65 AliHMPIDParam *pParam=AliHMPIDParam::Instance();
d3da6dc4 66
a591e55f 67 if(pCluLst->GetEntries()>pParam->MultCut()) fIsWEIGHT = kTRUE; // offset to take into account bkg in reconstruction
68 else fIsWEIGHT = kFALSE;
d3da6dc4 69
a591e55f 70 Float_t xRa,yRa,th,ph;
71 pTrk->GetHMPIDtrk(xRa,yRa,th,ph); //initialize this track: th and ph angles at middle of RAD
a591e55f 72 SetTrack(xRa,yRa,th,ph);
73
abb5f786 74 fRadNmean=nmean;
d3da6dc4 75
59280a5a 76 Float_t dMin=999,mipX=-1,mipY=-1;Int_t chId=-1,mipId=-1,mipQ=-1;
d3da6dc4 77 fPhotCnt=0;
78 for (Int_t iClu=0; iClu<pCluLst->GetEntriesFast();iClu++){//clusters loop
79 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu); //get pointer to current cluster
59280a5a 80 chId=pClu->Ch();
a591e55f 81 if(pClu->Q()>pParam->QCut()){ //charge compartible with MIP clusters
82 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
83 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
84 }else{ //charge compatible with photon cluster
85 Double_t thetaCer,phiCer;
86 if(FindPhotCkov(pClu->X(),pClu->Y(),thetaCer,phiCer)){ //find ckov angle for this photon candidate
87 fPhotCkov[fPhotCnt]=thetaCer; //actual theta Cerenkov (in TRS)
b4ad85e9 88 fPhotPhi [fPhotCnt]=phiCer; //actual phi Cerenkov (in TRS): -pi to come back to "unusual" ref system (X,Y,-Z)
a591e55f 89 fPhotCnt++; //increment counter of photon candidates
90 }
59280a5a 91 }
d3da6dc4 92 }//clusters loop
a591e55f 93 Int_t iNacc=FlagPhot(HoughResponse()); //flag photons according to individual theta ckov with respect to most probable
94 pTrk->SetHMPIDmip(mipX,mipY,mipQ,iNacc); //store mip info
59280a5a 95
a591e55f 96 if(mipId==-1) {pTrk->SetHMPIDsignal(kMipQdcCut); return;} //no clusters with QDC more the threshold at all
97 if(dMin>pParam->DistCut()) {pTrk->SetHMPIDsignal(kMipDistCut); return;} //closest cluster with enough charge is still too far from intersection
98 pTrk->SetHMPIDcluIdx(chId,mipId); //set index of cluster
99 if(iNacc<1) pTrk->SetHMPIDsignal(kNoPhotAccept); //no photon candidates is accepted
100 else pTrk->SetHMPIDsignal(FindRingCkov(pCluLst->GetEntries())); //find best Theta ckov for ring i.e. track
59280a5a 101
a591e55f 102 pTrk->SetHMPIDchi2(fCkovSigma2); //errors squared
d3da6dc4 103
43400d2d 104}//CkovAngle()
d3da6dc4 105//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
a591e55f 106Bool_t AliHMPIDRecon::FindPhotCkov(Double_t cluX,Double_t cluY,Double_t &thetaCer,Double_t &phiCer)
d3da6dc4 107{
108// Finds Cerenkov angle for this photon candidate
109// Arguments: cluX,cluY - position of cadidate's cluster
a591e55f 110// Returns: Cerenkov angle
d3da6dc4 111
a591e55f 112 TVector3 dirCkov;
113
67a1c24c 114 Double_t zRad= -0.5*fgkRadThick-0.5*fgkWinThick; //z position of middle of RAD
115 TVector3 rad(fTrkPos.X(),fTrkPos.Y(),zRad); //impact point at middle of RAD
116 TVector3 pc(cluX,cluY,0.5*fgkWinThick+fgkGapIdx); //mip at PC
a591e55f 117 Double_t cluR = TMath::Sqrt((cluX-fTrkPos.X())*(cluX-fTrkPos.X())+
118 (cluY-fTrkPos.Y())*(cluY-fTrkPos.Y()));//ref. distance impact RAD-CLUSTER
67a1c24c 119 Double_t phi=(pc-rad).Phi(); //phi of photon
a591e55f 120
b4ad85e9 121 Double_t ckov1=0;
67a1c24c 122 Double_t ckov2=0.75+fTrkDir.Theta(); //start to find theta cerenkov in DRS
b4ad85e9 123 const Double_t kTol=0.01;
d3da6dc4 124 Int_t iIterCnt = 0;
125 while(1){
a591e55f 126 if(iIterCnt>=50) return kFALSE;
d3da6dc4 127 Double_t ckov=0.5*(ckov1+ckov2);
67a1c24c 128 dirCkov.SetMagThetaPhi(1,ckov,phi);
a591e55f 129 TVector2 posC=TraceForward(dirCkov); //trace photon with actual angles
130 Double_t dist=cluR-(posC-fTrkPos).Mod(); //get distance between trial point and cluster position
131 if(posC.X()==-999) dist = - 999; //total reflection problem
132 iIterCnt++; //counter step
b4ad85e9 133 if (dist> kTol) ckov1=ckov; //cluster @ larger ckov
d3da6dc4 134 else if(dist<-kTol) ckov2=ckov; //cluster @ smaller ckov
a591e55f 135 else{ //precision achived: ckov in DRS found
136 dirCkov.SetMagThetaPhi(1,ckov,phi); //
137 RecPhot(dirCkov,thetaCer,phiCer); //find ckov (in TRS:the effective Cherenkov angle!)
138 return kTRUE;
139 }
d3da6dc4 140 }
141}//FindPhotTheta()
142//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
a591e55f 143TVector2 AliHMPIDRecon::TraceForward(TVector3 dirCkov)const
d3da6dc4 144{
a591e55f 145 //Trace forward a photon from (x,y) up to PC
146 // Arguments: dirCkov photon vector in LORS
147 // Returns: pos of traced photon at PC
148 TVector2 pos(-999,-999);
67a1c24c 149 Double_t thetaCer = dirCkov.Theta();
150 if(thetaCer > TMath::ASin(1./fRadNmean)) return pos; //total refraction on WIN-GAP boundary
151 Double_t zRad= -0.5*fgkRadThick-0.5*fgkWinThick; //z position of middle of RAD
152 TVector3 posCkov(fTrkPos.X(),fTrkPos.Y(),zRad); //RAD: photon position is track position @ middle of RAD
153 Propagate(dirCkov,posCkov, -0.5*fgkWinThick); //go to RAD-WIN boundary
154 Refract (dirCkov, fRadNmean,fgkWinIdx); //RAD-WIN refraction
155 Propagate(dirCkov,posCkov, 0.5*fgkWinThick); //go to WIN-GAP boundary
156 Refract (dirCkov, fgkWinIdx,fgkGapIdx); //WIN-GAP refraction
157 Propagate(dirCkov,posCkov,0.5*fgkWinThick+fgkGapThick); //go to PC
a591e55f 158 pos.Set(posCkov.X(),posCkov.Y());
159 return pos;
160}//TraceForward()
161//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
162void AliHMPIDRecon::RecPhot(TVector3 dirCkov,Double_t &thetaCer,Double_t &phiCer)
163{
164 //Theta Cerenkov reconstruction
165 // Arguments: (x,y) of initial point in LORS, dirCkov photon vector in LORS
166 // Returns: thetaCer theta cerenkov reconstructed
167// TVector3 dirTrk;
168// dirTrk.SetMagThetaPhi(1,fTrkDir.Theta(),fTrkDir.Phi());
169// Double_t thetaCer = TMath::ACos(dirCkov*dirTrk);
170 TRotation mtheta; mtheta.RotateY(- fTrkDir.Theta());
171 TRotation mphi; mphi.RotateZ(- fTrkDir.Phi());
172 TRotation mrot=mtheta*mphi;
173 TVector3 dirCkovTRS;
174 dirCkovTRS=mrot*dirCkov;
175 phiCer = dirCkovTRS.Phi(); //actual value of the phi of the photon
176 thetaCer= dirCkovTRS.Theta(); //actual value of thetaCerenkov of the photon
d3da6dc4 177}
178//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
179Double_t AliHMPIDRecon::FindRingArea(Double_t ckovAng)const
180{
181// Find area inside the cerenkov ring which lays inside PCs
a591e55f 182// Arguments: ckovAng - cerenkov angle
d3da6dc4 183// Returns: area of the ring in cm^2 for given theta ckov
184
d3da6dc4 185 const Int_t kN=100;
186 Double_t area=0;
187 for(Int_t i=0;i<kN;i++){
a591e55f 188 TVector2 pos1=TracePhot(ckovAng,Double_t(TMath::TwoPi()*i /kN));//trace this photon
189 TVector2 pos2=TracePhot(ckovAng,Double_t(TMath::TwoPi()*(i+1)/kN));//trace the next photon
190 area+=(pos1-fTrkPos)*(pos2-fTrkPos); //add area of the triangle...
d3da6dc4 191 }
192 return area;
193}//FindRingArea()
194//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
195Double_t AliHMPIDRecon::FindRingCkov(Int_t)
196{
197// 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
198// collecting errors for all single Ckov candidates thetas. (Assuming they are independent)
199// Arguments: iNclus- total number of clusters in chamber for background estimation
200// Return: best estimation of track Theta ckov
201
202 Double_t wei = 0.;
203 Double_t weightThetaCerenkov = 0.;
204
205 Double_t ckovMin=9999.,ckovMax=0.;
206 Double_t sigma2 = 0; //to collect error squared for this ring
207
208 for(Int_t i=0;i<fPhotCnt;i++){//candidates loop
209 if(fPhotFlag[i] == 2){
a591e55f 210 if(fPhotCkov[i]<ckovMin) ckovMin=fPhotCkov[i]; //find max and min Theta ckov from all candidates within probable window
d3da6dc4 211 if(fPhotCkov[i]>ckovMax) ckovMax=fPhotCkov[i];
a591e55f 212 weightThetaCerenkov += fPhotCkov[i]*fPhotWei[i];
213 wei += fPhotWei[i]; //collect weight as sum of all candidate weghts
d3da6dc4 214
d3da6dc4 215 sigma2 += 1./Sigma2(fPhotCkov[i],fPhotPhi[i]);
216 }
217 }//candidates loop
218
219 if(sigma2>0) fCkovSigma2=1./sigma2;
220 else fCkovSigma2=1e10;
221
b4ad85e9 222 if(wei != 0.) weightThetaCerenkov /= wei; else weightThetaCerenkov = 0.;
d3da6dc4 223 return weightThetaCerenkov;
224}//FindCkovRing()
225//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
226Int_t AliHMPIDRecon::FlagPhot(Double_t ckov)
227{
228// Flag photon candidates if their individual ckov angle is inside the window around ckov angle returned by HoughResponse()
229// Arguments: ckov- value of most probable ckov angle for track as returned by HoughResponse()
230// Returns: number of photon candidates happened to be inside the window
231
a591e55f 232// Photon Flag: Flag = 0 initial set;
233// Flag = 1 good candidate (charge compatible with photon);
234// Flag = 2 photon used for the ring;
d3da6dc4 235
236 Int_t steps = (Int_t)((ckov )/ fDTheta); //how many times we need to have fDTheta to fill the distance between 0 and thetaCkovHough
237
238 Double_t tmin = (Double_t)(steps - 1)*fDTheta;
239 Double_t tmax = (Double_t)(steps)*fDTheta;
240 Double_t tavg = 0.5*(tmin+tmax);
241
242 tmin = tavg - 0.5*fWindowWidth; tmax = tavg + 0.5*fWindowWidth;
243
244 Int_t iInsideCnt = 0; //count photons which Theta ckov inside the window
245 for(Int_t i=0;i<fPhotCnt;i++){//photon candidates loop
246 if(fPhotCkov[i] >= tmin && fPhotCkov[i] <= tmax) {
247 fPhotFlag[i]=2;
248 iInsideCnt++;
249 }
250 }
251 return iInsideCnt;
252}//FlagPhot()
253//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
a591e55f 254TVector2 AliHMPIDRecon::TracePhot(Double_t ckovThe,Double_t ckovPhi)const
d3da6dc4 255{
256// Trace a single Ckov photon from emission point somewhere in radiator up to photocathode taking into account ref indexes of materials it travereses
a591e55f 257// Arguments: ckovThe,ckovPhi- photon ckov angles in DRS, [rad]
d3da6dc4 258// Returns: distance between photon point on PC and track projection
259 TRotation mtheta; mtheta.RotateY(fTrkDir.Theta());
260 TRotation mphi; mphi.RotateZ(fTrkDir.Phi());
261 TRotation mrot=mphi*mtheta;
a591e55f 262 TVector3 dirCkov,dirCkovTors;
263
264 dirCkovTors.SetMagThetaPhi(1,ckovThe,ckovPhi); //initially photon is directed according to requested ckov angle
265 dirCkov=mrot*dirCkovTors; //now we know photon direction in LORS
266 return TraceForward(dirCkov);
267}//TracePhot()
d3da6dc4 268//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
a591e55f 269void AliHMPIDRecon::Propagate(const TVector3 dir,TVector3 &pos,Double_t z)const
d3da6dc4 270{
271// Finds an intersection point between a line and XY plane shifted along Z.
272// Arguments: dir,pos - vector along the line and any point of the line
273// z - z coordinate of plain
274// Returns: none
275// On exit: pos is the position if this intesection if any
276 static TVector3 nrm(0,0,1);
277 TVector3 pnt(0,0,z);
278
279 TVector3 diff=pnt-pos;
280 Double_t sint=(nrm*diff)/(nrm*dir);
281 pos+=sint*dir;
282}//Propagate()
283//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
284void AliHMPIDRecon::Refract(TVector3 &dir,Double_t n1,Double_t n2)const
285{
286// Refract direction vector according to Snell law
287// Arguments:
288// n1 - ref idx of first substance
289// n2 - ref idx of second substance
290// Returns: none
291// On exit: dir is new direction
67a1c24c 292 Double_t sinref=(n1/n2)*TMath::Sin(dir.Theta());
d3da6dc4 293 if(sinref>1.) dir.SetXYZ(-999,-999,-999);
67a1c24c 294 else dir.SetTheta(TMath::ASin(sinref));
d3da6dc4 295}//Refract()
296//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
297Double_t AliHMPIDRecon::HoughResponse()
298{
299//
300//
301//
302 Double_t kThetaMax=0.75;
303 Int_t nChannels = (Int_t)(kThetaMax/fDTheta+0.5);
304 TH1D *phots = new TH1D("Rphot" ,"phots" ,nChannels,0,kThetaMax);
305 TH1D *photsw = new TH1D("RphotWeighted" ,"photsw" ,nChannels,0,kThetaMax);
306 TH1D *resultw = new TH1D("resultw","resultw" ,nChannels,0,kThetaMax);
307 Int_t nBin = (Int_t)(kThetaMax/fDTheta);
308 Int_t nCorrBand = (Int_t)(fWindowWidth/(2*fDTheta));
309
310 for (Int_t i=0; i< fPhotCnt; i++){//photon cadidates loop
311 Double_t angle = fPhotCkov[i]; if(angle<0||angle>kThetaMax) continue;
312 phots->Fill(angle);
313 Int_t bin = (Int_t)(0.5+angle/(fDTheta));
314 Double_t weight=1.;
315 if(fIsWEIGHT){
316 Double_t lowerlimit = ((Double_t)bin)*fDTheta - 0.5*fDTheta; Double_t upperlimit = ((Double_t)bin)*fDTheta + 0.5*fDTheta;
317 Double_t diffArea = FindRingArea(upperlimit)-FindRingArea(lowerlimit);
318 if(diffArea>0) weight = 1./diffArea;
319 }
320 photsw->Fill(angle,weight);
321 fPhotWei[i]=weight;
322 }//photon candidates loop
323
324 for (Int_t i=1; i<=nBin;i++){
325 Int_t bin1= i-nCorrBand;
326 Int_t bin2= i+nCorrBand;
327 if(bin1<1) bin1=1;
328 if(bin2>nBin)bin2=nBin;
329 Double_t sumPhots=phots->Integral(bin1,bin2);
330 if(sumPhots<3) continue; // if less then 3 photons don't trust to this ring
331 Double_t sumPhotsw=photsw->Integral(bin1,bin2);
332 resultw->Fill((Double_t)((i+0.5)*fDTheta),sumPhotsw);
333 }
334// evaluate the "BEST" theta ckov as the maximum value of histogramm
335 Double_t *pVec = resultw->GetArray();
336 Int_t locMax = TMath::LocMax(nBin,pVec);
337 phots->Delete();photsw->Delete();resultw->Delete(); // Reset and delete objects
338
339 return (Double_t)(locMax*fDTheta+0.5*fDTheta); //final most probable track theta ckov
340}//HoughResponse()
341//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
342Double_t AliHMPIDRecon::Sigma2(Double_t ckovTh, Double_t ckovPh)const
343{
344// Analithical calculation of total error (as a sum of localization, geometrical and chromatic errors) on Cerenkov angle for a given Cerenkov photon
345// created by a given MIP. Fromulae according to CERN-EP-2000-058
346// Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
347// dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]
348// MIP beta
349// Returns: absolute error on Cerenkov angle, [radians]
350
351 TVector3 v(-999,-999,-999);
abb5f786 352 Double_t trkBeta = 1./(TMath::Cos(ckovTh)*fRadNmean);
d3da6dc4 353
354 v.SetX(SigLoc (ckovTh,ckovPh,trkBeta));
355 v.SetY(SigGeom(ckovTh,ckovPh,trkBeta));
356 v.SetZ(SigCrom(ckovTh,ckovPh,trkBeta));
357
358 return v.Mag2();
359}
360//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
361Double_t AliHMPIDRecon::SigLoc(Double_t thetaC, Double_t phiC,Double_t betaM)const
362{
363// Analithical calculation of localization error (due to finite segmentation of PC) on Cerenkov angle for a given Cerenkov photon
364// created by a given MIP. Fromulae according to CERN-EP-2000-058
365// Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
366// dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]
367// MIP beta
368// Returns: absolute error on Cerenkov angle, [radians]
369 Double_t phiDelta = phiC - fTrkDir.Phi();
370
371 Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
abb5f786 372 Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
d3da6dc4 373 if (k<0) return 1e10;
374
375 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()));
376 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()));
377
378 Double_t kk = betaM*TMath::Sqrt(k)/(8*alpha);
379 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));
380 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));
381
382 return TMath::Sqrt(0.2*0.2*dtdxc*dtdxc + 0.25*0.25*dtdyc*dtdyc);
383}
384//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
385Double_t AliHMPIDRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)const
386{
387// Analithical calculation of chromatic error (due to lack of knowledge of Cerenkov photon energy) on Cerenkov angle for a given Cerenkov photon
388// created by a given MIP. Fromulae according to CERN-EP-2000-058
389// Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
390// dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]
391// MIP beta
392// Returns: absolute error on Cerenkov angle, [radians]
393 Double_t phiDelta = phiC - fTrkDir.Phi();
394 Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
395
abb5f786 396 Double_t dtdn = TMath::Cos(fTrkDir.Theta())*fRadNmean*betaM*betaM/(alpha*TMath::Tan(thetaC));
d3da6dc4 397
398 Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
399
400 return f*dtdn;
401}//SigCrom()
402//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
403Double_t AliHMPIDRecon::SigGeom(Double_t thetaC, Double_t phiC,Double_t betaM)const
404{
405// Analithical calculation of geometric error (due to lack of knowledge of creation point in radiator) on Cerenkov angle for a given Cerenkov photon
406// created by a given MIP. Formulae according to CERN-EP-2000-058
407// Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
408// dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]
409// MIP beta
410// Returns: absolute error on Cerenkov angle, [radians]
411
412 Double_t phiDelta = phiC - fTrkDir.Phi();
413 Double_t alpha =TMath::Cos(fTrkDir.Theta())-TMath::Tan(thetaC)*TMath::Cos(phiDelta)*TMath::Sin(fTrkDir.Theta());
414
abb5f786 415 Double_t k = 1.-fRadNmean*fRadNmean+alpha*alpha/(betaM*betaM);
d3da6dc4 416 if (k<0) return 1e10;
417
418 Double_t eTr = 0.5*1.5*betaM*TMath::Sqrt(k)/(8*alpha);
419 Double_t lambda = 1.-TMath::Sin(fTrkDir.Theta())*TMath::Sin(fTrkDir.Theta())*TMath::Sin(phiC)*TMath::Sin(phiC);
420
421 Double_t c = 1./(1.+ eTr*k/(alpha*alpha*TMath::Cos(thetaC)*TMath::Cos(thetaC)));
422 Double_t i = betaM*TMath::Tan(thetaC)*lambda*TMath::Power(k,1.5);
423 Double_t ii = 1.+eTr*betaM*i;
424
425 Double_t err = c * (i/(alpha*alpha*8) + ii*(1.-lambda) / ( alpha*alpha*8*betaM*(1.+eTr)) );
426 Double_t trErr = 1.5/(TMath::Sqrt(12.)*TMath::Cos(fTrkDir.Theta()));
427
428 return trErr*err;
429}//SigGeom()
430//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
43400d2d 431//
432//
433//
434//
435//
436//
437//
438//
439//
440//
441//
442// From here Hidden track algorithm....
443//
444//
445//
446//
447//
448//
449//
450//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
451void AliHMPIDRecon::CkovHiddenTrk(TClonesArray *pCluLst,Double_t nmean)
452{
453// Pattern recognition method without any infos from tracking...
454// The method finds in the chmber the cluster with the highest charge
455// compatibile with a MIP, then the strategy is applied
456// Arguments: pCluLst - list of clusters for this chamber
457// Returns: - track ckov angle, [rad],
458
459 fRadNmean=nmean;
460
461 Float_t mipX=-1,mipY=-1;Int_t mipId=-1,mipQ=-1;
462 fPhotCnt=0;
463 Double_t qRef = 0;
464 fNClu = pCluLst->GetEntriesFast();
465 for (Int_t iClu=0;iClu<fNClu;iClu++){ //clusters loop
466 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu); //get pointer to current cluster
467 fXClu[iClu] = pClu->X();fYClu[iClu] = pClu->Y(); //store x,y for fitting procedure
468 if(pClu->Q()>qRef){ //searching the highest charge to select a MIP
469 qRef = pClu->Q();
470 mipId=iClu; mipX=pClu->X();mipY=pClu->Y();mipQ=(Int_t)pClu->Q();
471 }
472 }//clusters loop
473
474 fIdxMip = mipId;
475 fMipX = mipX; fMipY=mipY;
476 DoRecHiddenTrk();
477
478}//CkovHiddenTrk()
479//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
480void AliHMPIDRecon::DoRecHiddenTrk()
481{
482// Pattern recognition method without any infos from tracking...
483// Before a preclustering filter to avoid part of the noise
484// Then only ellipsed-rings are fitted (no possibility,
485// for the moment, to reconstruct very inclined tracks
486// Finally a fitting with (th,ph) free starting by very close values
487// previously evaluated.
488// Arguments: none
489// Returns: none
490 Double_t phiRec;
491 CluPreFilter();
492 if(!FitEllipse(phiRec)) {Printf("Not an ellipse, bye!");return;}
493 Printf("--->now it starts the free fit with phi = %f",phiRec*TMath::RadToDeg());
494 return FitFree(phiRec);
495}
496//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
497void AliHMPIDRecon::CluPreFilter()
498{
499// Filter of bkg clusters
500// based on elliptical-shapes...
501//
502 ;
503}
504//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
505Bool_t AliHMPIDRecon::FitEllipse(Double_t &phiRec)
506{
507//Fit a set of clusters with an analitical conical section function:
508 //
509 // Ax^2 + B*y^2 + 2Hxy + 2Gx + 2Fy + 1 = 0 ---> conical section
510 //
511 // H*H - A*B > 0 hyperbola
512 // < 0 ellipse
513 // = 0 parabola
514 //
515 // tan 2alfa = 2H/(A-B) alfa=angle of rotation
516 //
517 // coordinate of the centre of the conical section:
518 // x = x' + a
519 // y = y' + b
520 //
521 // HF - BG
522 // a = ---------
523 // AB - H^2
524 //
525 // HG - AF
526 // b = --------
527 // AB - H^2
528
529 Double_t cA,cB,cF,cG,cH;
530 Double_t aArg=-1,parStep,parLow,parHigh; Int_t iErrFlg; //tmp vars for TMinuit
531
532 if(!gMinuit) gMinuit = new TMinuit(100); //init MINUIT with this number of parameters (5 params)
533 gMinuit->mncler(); // reset Minuit list of paramters
534 gMinuit->SetObjectFit((TObject*)this); gMinuit->SetFCN(AliHMPIDRecon::FunMinEl); //set fit function
535 gMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg); //suspend all printout from TMinuit
536 gMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg); //suspend all warning printout from TMinuit
537
538 Double_t d1,d2,d3;
539 TString sName;
540
541 gMinuit->mnparm(0," A ",1,parStep=0.01,parLow=0,parHigh=0,iErrFlg);
542 gMinuit->mnparm(1," B ",1,parStep=0.01,parLow=0,parHigh=0,iErrFlg);
543 gMinuit->mnparm(2," H ",1,parStep=0.01,parLow=0,parHigh=0,iErrFlg);
544 gMinuit->mnparm(3," G ",1,parStep=0.01,parLow=0,parHigh=0,iErrFlg);
545 gMinuit->mnparm(4," F ",1,parStep=0.01,parLow=0,parHigh=0,iErrFlg);
546
547 gMinuit->mnexcm("SIMPLEX" ,&aArg,0,iErrFlg);
548 gMinuit->mnexcm("MIGRAD" ,&aArg,0,iErrFlg);
549 gMinuit->mnpout(0,sName,cA,d1,d2,d3,iErrFlg);
550 gMinuit->mnpout(1,sName,cB,d1,d2,d3,iErrFlg);
551 gMinuit->mnpout(2,sName,cH,d1,d2,d3,iErrFlg);
552 gMinuit->mnpout(3,sName,cG,d1,d2,d3,iErrFlg);
553 gMinuit->mnpout(4,sName,cF,d1,d2,d3,iErrFlg);
554 delete gMinuit;
555
556 Double_t i2 = cA*cB-cH*cH; //quartic invariant : i2 > 0 ellipse, i2 < 0 hyperbola
557 Double_t aX = (cH*cF-cB*cG)/i2; //x centre of the canonical section
558 Double_t bY = (cH*cG-cA*cF)/i2; //y centre of the canonical section
559 Double_t alfa1 = TMath::ATan(2*cH/(cA-cB)); //alpha = angle of rotation of the conical section
560 if(alfa1<0) alfa1+=TMath::Pi();
561 alfa1*=0.5;
562 Double_t alfa2 = alfa1+TMath::Pi();
563 Double_t phiref = TMath::ATan2(bY-fMipY,aX-fMipX); //evaluate in a unique way the angle of rotation comapring it
564 if(phiref<0) phiref+=TMath::TwoPi(); //with the vector that poinst to the centre from the mip
565 if(i2<0) phiref+=TMath::Pi();
566 if(phiref>TMath::TwoPi()) phiref-=TMath::TwoPi();
567
568// Printf(" alfa1 %f",alfa1*TMath::RadToDeg());
569// Printf(" alfa2 %f",alfa2*TMath::RadToDeg());
570// Printf(" firef %f",phiref*TMath::RadToDeg());
571 if(TMath::Abs(alfa1-phiref)<TMath::Abs(alfa2-phiref)) phiRec = alfa1; else phiRec = alfa2;
572
573// cout << " phi reconstructed " << phiRec*TMath::RadToDeg() << endl;
574 return (i2>0);
575//
576}
577//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
578void AliHMPIDRecon::FitFree(Double_t phiRec)
579{
580// Fit performed by minimizing RMS/sqrt(n) of the
581// photons reconstructed. First phi is fixed and theta
582// is fouond, then (th,ph) of the track
583// as free parameters
584// Arguments: PhiRec phi of the track
585// Returns: none
586 Double_t aArg=-1,parStep,parLow,parHigh; Int_t iErrFlg; //tmp vars for TMinuit
587 if(!gMinuit) gMinuit = new TMinuit(100); //init MINUIT with this number of parameters (5 params)
588 gMinuit->mncler(); // reset Minuit list of paramters
589 gMinuit->SetObjectFit((TObject*)this); gMinuit->SetFCN(AliHMPIDRecon::FunMinPhot); //set fit function
590 gMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg); //suspend all printout from TMinuit
591 gMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg); //suspend all warning printout from TMinuit
592
593 Double_t d1,d2,d3;
594 TString sName;
595 Double_t th,ph;
596
597 gMinuit->mnexcm("SET PRI",&aArg,1,iErrFlg); //suspend all printout from TMinuit
598 gMinuit->mnexcm("SET NOW",&aArg,0,iErrFlg);
599
600 gMinuit->mnparm(0," theta ",0.01,parStep=0.01,parLow=0,parHigh=TMath::PiOver2(),iErrFlg);
601 gMinuit->mnparm(1," phi ",phiRec,parStep=0.01,parLow=0,parHigh=TMath::TwoPi(),iErrFlg);
602
603 gMinuit->FixParameter(1);
604 gMinuit->mnexcm("SIMPLEX" ,&aArg,0,iErrFlg);
605 gMinuit->mnexcm("MIGRAD" ,&aArg,0,iErrFlg);
606 gMinuit->Release(1);
607 gMinuit->mnexcm("MIGRAD" ,&aArg,0,iErrFlg);
608
609 gMinuit->mnpout(0,sName,th,d1,d2,d3,iErrFlg);
610 gMinuit->mnpout(1,sName,ph,d1,d2,d3,iErrFlg);
611
612 Printf(" reconstr. theta %f phi %f",th*TMath::RadToDeg(),ph*TMath::RadToDeg());
613}
614//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
615Double_t AliHMPIDRecon::FunConSect(Double_t *c,Double_t x,Double_t y)
616{
617 return c[0]*x*x+c[1]*y*y+2*c[2]*x*y+2*c[3]*x+2*c[4]*y+1;
618}
619//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
620void AliHMPIDRecon::FunMinEl(Int_t &/* */,Double_t* /* */,Double_t &f,Double_t *par,Int_t /* */)
621{
622 AliHMPIDRecon *pRec=(AliHMPIDRecon*)gMinuit->GetObjectFit();
623 Double_t minFun = 0;
624 Int_t np = pRec->NClu();
625 for(Int_t i=0;i<np;i++) {
626 if(i==pRec->IdxMip()) continue;
627 Double_t el = pRec->FunConSect(par,pRec->XClu(i),pRec->YClu(i));
628 minFun +=el*el;
629 }
630 f = minFun;
631}
632//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
633void AliHMPIDRecon::FunMinPhot(Int_t &/* */,Double_t* /* */,Double_t &f,Double_t *par,Int_t /* */)
634{
635 AliHMPIDRecon *pRec=(AliHMPIDRecon*)gMinuit->GetObjectFit();
636
637 Double_t sizeCh = fgkRadThick+fgkWinThick+fgkGapThick;
638 Double_t thTrk = par[0];
639 Double_t phTrk = par[1];
640 Double_t xrad = pRec->MipX() - sizeCh*TMath::Tan(thTrk)*TMath::Cos(phTrk);
641 Double_t yrad = pRec->MipY() - sizeCh*TMath::Tan(thTrk)*TMath::Sin(phTrk);
642
643 pRec->SetTrack(xrad,yrad,thTrk,phTrk);
644
645 Double_t meanCkov=0;
646 Double_t meanCkov2=0;
647 Double_t thetaCer,phiCer;
648 for(Int_t i=0;i<pRec->NClu();i++) {
649 pRec->FindPhotCkov(pRec->XClu(i),pRec->YClu(i),thetaCer,phiCer);
650
651 meanCkov += thetaCer;
652 meanCkov2 += thetaCer*thetaCer;
653 }
654 meanCkov/=pRec->NClu();
655 Double_t rms = TMath::Sqrt(meanCkov2/pRec->NClu() - meanCkov*meanCkov);
656 f = rms/TMath::Sqrt(pRec->NClu());
657 pRec->SetCkovFit(meanCkov);
658}//FunMinPhot()