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