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