]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDRecon.cxx
bugfix: corrected defines to use right default algorithms
[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 "AliHMPIDCluster.h" //CkovAngle()
26 #include <TRotation.h>       //TracePhot()
27 #include <TH1D.h>            //HoughResponse()
28 #include <TClonesArray.h>    //CkovAngle()
29 #include <AliESDtrack.h>     //CkovAngle()
30
31 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
32 AliHMPIDRecon::AliHMPIDRecon():TTask("RichRec","RichPat")
33 {
34 //..
35 //init of data members
36 //..
37   
38   fPhotCnt  = -1;
39   fPhotFlag = 0x0;
40   fPhotCkov = 0x0;
41   fPhotPhi  = 0x0;
42   fPhotWei  = 0x0;
43   fCkovSigma2 = 0;
44   fIsWEIGHT = kFALSE;
45   fDTheta   = 0.001;
46   fWindowWidth = 0.045;
47   fTrkDir = TVector3(0,0,1); // init just for test
48   fTrkPos = TVector2(30,40); // init just for test
49   
50   AliHMPIDParam *pParam=AliHMPIDParam::Instance();
51   fParam = pParam;
52   
53   fParam->SetRefIdx(fParam->MeanIdxRad()); // initialization of ref index to a default one
54 }
55 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
56 void AliHMPIDRecon::InitVars(Int_t n)
57 {
58 //..
59 //Init some variables
60 //..
61   if(n<0) return;
62   fPhotFlag = new Int_t[n];
63   fPhotCkov = new Double_t[n];
64   fPhotPhi  = new Double_t[n];
65   fPhotWei  = new Double_t[n];
66 //
67 }
68 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
69 void AliHMPIDRecon::DeleteVars()
70 {
71 //..
72 //Delete variables
73 //..
74   delete fPhotFlag;
75   delete fPhotCkov;
76   delete fPhotPhi;
77   delete fPhotWei;
78 }
79 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
80 void AliHMPIDRecon::CkovAngle(AliESDtrack *pTrk,TClonesArray *pCluLst,Double_t nmean,Double_t qthre)
81 {
82 // Pattern recognition method based on Hough transform
83 // Arguments:   pTrk     - track for which Ckov angle is to be found
84 //              pCluLst  - list of clusters for this chamber   
85 //   Returns:            - track ckov angle, [rad], 
86     
87   Int_t nClusTot = pCluLst->GetEntries();
88   if(nClusTot>fParam->MultCut()) fIsWEIGHT = kTRUE; // offset to take into account bkg in reconstruction
89   else                           fIsWEIGHT = kFALSE;
90
91   InitVars(nClusTot);
92   
93   Float_t xRa,yRa,th,ph;
94   pTrk->GetHMPIDtrk(xRa,yRa,th,ph);        //initialize this track: th and ph angles at middle of RAD 
95   SetTrack(xRa,yRa,th,ph);
96
97   fParam->SetRefIdx(nmean);
98
99   Float_t dMin=999,mipX=-1,mipY=-1;Int_t chId=-1,mipId=-1,mipQ=-1;                                                                           
100   fPhotCnt=0;                                                      
101   for (Int_t iClu=0; iClu<pCluLst->GetEntriesFast();iClu++){//clusters loop
102     AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluLst->UncheckedAt(iClu);                       //get pointer to current cluster    
103     chId=pClu->Ch();
104     if(pClu->Q()>qthre){                                                                      //charge compartible with MIP clusters      
105       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
106       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
107     }else{                                                                                    //charge compatible with photon cluster
108       Double_t thetaCer,phiCer;
109       if(FindPhotCkov(pClu->X(),pClu->Y(),thetaCer,phiCer)){                                  //find ckov angle for this  photon candidate
110         fPhotCkov[fPhotCnt]=thetaCer;                                                         //actual theta Cerenkov (in TRS)
111         fPhotPhi [fPhotCnt]=phiCer;                                                           //actual phi   Cerenkov (in TRS): -pi to come back to "unusual" ref system (X,Y,-Z)
112         //PH        Printf("photon n. %i reconstructed theta = %f",fPhotCnt,fPhotCkov[fPhotCnt]);
113         fPhotCnt++;                                                                           //increment counter of photon candidates
114       }
115     }
116   }//clusters loop
117   fMipPos.Set(mipX,mipY);
118   if(fPhotCnt<=3) pTrk->SetHMPIDsignal(kNoPhotAccept);                                        //no reconstruction with <=3 photon candidates
119   Int_t iNrec=FlagPhot(HoughResponse());                                                      //flag photons according to individual theta ckov with respect to most probable
120   pTrk->SetHMPIDmip(mipX,mipY,mipQ,iNrec);                                                    //store mip info 
121
122   if(mipId==-1)              {pTrk->SetHMPIDsignal(kMipQdcCut);  return;}                     //no clusters with QDC more the threshold at all
123   if(dMin>fParam->DistCut()) {pTrk->SetHMPIDsignal(kMipDistCut); return;}                     //closest cluster with enough charge is still too far from intersection
124   pTrk->SetHMPIDcluIdx(chId,mipId);                                                           //set index of cluster
125   if(iNrec<1){
126     pTrk->SetHMPIDsignal(kNoPhotAccept);                                                      //no photon candidates are accepted
127   }
128   else {
129     Double_t thetaC = FindRingCkov(pCluLst->GetEntries());                                    //find the best reconstructed theta Cherenkov
130 //    FindRingGeom(thetaC,2);
131     pTrk->SetHMPIDsignal(thetaC);                                                             //store theta Cherenkov
132     pTrk->SetHMPIDchi2(fCkovSigma2);                                                          //store errors squared
133   }
134
135   DeleteVars();
136 }//CkovAngle()
137 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
138 Bool_t AliHMPIDRecon::FindPhotCkov(Double_t cluX,Double_t cluY,Double_t &thetaCer,Double_t &phiCer)
139 {
140 // Finds Cerenkov angle  for this photon candidate
141 // Arguments: cluX,cluY - position of cadidate's cluster  
142 // Returns: Cerenkov angle 
143
144   TVector3 dirCkov;
145   
146   Double_t zRad= -0.5*fParam->RadThick()-0.5*fParam->WinThick();     //z position of middle of RAD
147   TVector3 rad(fTrkPos.X(),fTrkPos.Y(),zRad);                        //impact point at middle of RAD
148   TVector3  pc(cluX,cluY,0.5*fParam->WinThick()+fParam->GapIdx());   //mip at PC
149   Double_t cluR = TMath::Sqrt((cluX-fTrkPos.X())*(cluX-fTrkPos.X())+
150                               (cluY-fTrkPos.Y())*(cluY-fTrkPos.Y()));//ref. distance impact RAD-CLUSTER   
151   Double_t phi=(pc-rad).Phi();                                       //phi of photon
152     
153   Double_t ckov1=0;
154   Double_t ckov2=0.75+fTrkDir.Theta();                        //start to find theta cerenkov in DRS
155   const Double_t kTol=0.01;
156   Int_t iIterCnt = 0;
157   while(1){
158     if(iIterCnt>=50) return kFALSE;
159     Double_t ckov=0.5*(ckov1+ckov2);
160     dirCkov.SetMagThetaPhi(1,ckov,phi);
161     TVector2 posC=TraceForward(dirCkov);                      //trace photon with actual angles
162     Double_t dist=cluR-(posC-fTrkPos).Mod();                  //get distance between trial point and cluster position
163     if(posC.X()==-999) dist = - 999;                          //total reflection problem
164     iIterCnt++;                                               //counter step
165     if     (dist> kTol) ckov1=ckov;                           //cluster @ larger ckov
166     else if(dist<-kTol) ckov2=ckov;                           //cluster @ smaller ckov
167     else{                                                     //precision achived: ckov in DRS found
168       dirCkov.SetMagThetaPhi(1,ckov,phi);                     //
169       Lors2Trs(dirCkov,thetaCer,phiCer);                       //find ckov (in TRS:the effective Cherenkov angle!)
170       return kTRUE;
171     }
172   }
173 }//FindPhotTheta()
174 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
175 TVector2 AliHMPIDRecon::TraceForward(TVector3 dirCkov)const
176 {
177   //Trace forward a photon from (x,y) up to PC
178   // Arguments: dirCkov photon vector in LORS
179   //   Returns: pos of traced photon at PC
180   
181   TVector2 pos(-999,-999);
182   Double_t thetaCer = dirCkov.Theta();
183   if(thetaCer > TMath::ASin(1./fParam->GetRefIdx())) return pos;          //total refraction on WIN-GAP boundary
184   Double_t zRad= -0.5*fParam->RadThick()-0.5*fParam->WinThick();          //z position of middle of RAD
185   TVector3  posCkov(fTrkPos.X(),fTrkPos.Y(),zRad);                        //RAD: photon position is track position @ middle of RAD 
186   Propagate(dirCkov,posCkov,           -0.5*fParam->WinThick());          //go to RAD-WIN boundary  
187   Refract  (dirCkov,         fParam->GetRefIdx(),fParam->WinIdx());       //RAD-WIN refraction
188   Propagate(dirCkov,posCkov,            0.5*fParam->WinThick());          //go to WIN-GAP boundary
189   Refract  (dirCkov,         fParam->WinIdx(),fParam->GapIdx());          //WIN-GAP refraction
190   Propagate(dirCkov,posCkov,0.5*fParam->WinThick()+fParam->GapThick());   //go to PC
191   pos.Set(posCkov.X(),posCkov.Y());
192   return pos;
193 }//TraceForward()
194 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
195 void AliHMPIDRecon::Lors2Trs(TVector3 dirCkov,Double_t &thetaCer,Double_t &phiCer)const
196 {
197   //Theta Cerenkov reconstruction 
198   // Arguments: dirCkov photon vector in LORS
199   //   Returns: thetaCer of photon in TRS
200   //              phiCer of photon in TRS
201 //  TVector3 dirTrk;
202 //  dirTrk.SetMagThetaPhi(1,fTrkDir.Theta(),fTrkDir.Phi());
203 //  Double_t thetaCer = TMath::ACos(dirCkov*dirTrk);
204   TRotation mtheta;   mtheta.RotateY(-fTrkDir.Theta());
205   TRotation mphi;       mphi.RotateZ(-fTrkDir.Phi());
206   TRotation mrot=mtheta*mphi;
207   TVector3 dirCkovTRS;
208   dirCkovTRS=mrot*dirCkov;
209   phiCer  = dirCkovTRS.Phi();                                          //actual value of the phi of the photon
210   thetaCer= dirCkovTRS.Theta();                                        //actual value of thetaCerenkov of the photon
211 }
212 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
213 void AliHMPIDRecon::Trs2Lors(TVector3 dirCkov,Double_t &thetaCer,Double_t &phiCer)const
214 {
215   //Theta Cerenkov reconstruction 
216   // Arguments: dirCkov photon vector in TRS
217   //   Returns: thetaCer of photon in LORS
218   //              phiCer of photon in LORS
219   TRotation mtheta;   mtheta.RotateY(fTrkDir.Theta());
220   TRotation mphi;       mphi.RotateZ(fTrkDir.Phi());
221   TRotation mrot=mphi*mtheta;
222   TVector3 dirCkovLORS;
223   dirCkovLORS=mrot*dirCkov;
224   phiCer  = dirCkovLORS.Phi();                                          //actual value of the phi of the photon
225   thetaCer= dirCkovLORS.Theta();                                        //actual value of thetaCerenkov of the photon
226 }
227 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
228 void AliHMPIDRecon::FindRingGeom(Double_t ckovAng,Int_t level)
229 {
230 // Find area covered in the PC acceptance
231 // Arguments: ckovAng - cerenkov angle
232 //            level   - precision in finding area and portion of ring accepted (multiple of 50)
233 //   Returns: area of the ring in cm^2 for given theta ckov
234    
235   Int_t kN=50*level;
236   Int_t nPoints = 0;
237   Double_t area=0;
238   
239   Bool_t first=kFALSE;
240   TVector2 pos1;
241   
242   for(Int_t i=0;i<kN;i++){
243    if(!first) {
244       pos1=TracePhot(ckovAng,Double_t(TMath::TwoPi()*(i+1)/kN));                                    //find a good trace for the first photon
245      if(pos1.X()==-999) continue;                                                                   //no area: open ring                  
246      if(!fParam->IsInside(pos1.X(),pos1.Y(),0)) {
247        pos1 = IntWithEdge(fMipPos,pos1);                                                            // find the very first intersection...
248      } else {
249        if(!AliHMPIDParam::IsInDead(pos1.X(),pos1.Y())) nPoints++;                                   //photon is accepted if not in dead zone
250      }
251      first=kTRUE;
252      continue;
253    }
254    TVector2 pos2=TracePhot(ckovAng,Double_t(TMath::TwoPi()*(i+1)/kN));                              //trace the next photon
255    if(pos2.X()==-999) continue;                                                                     //no area: open ring             
256    if(!fParam->IsInside(pos2.X(),pos2.Y(),0)) {
257      pos2 = IntWithEdge(fMipPos,pos2);
258    } else {
259      if(!AliHMPIDParam::IsInDead(pos2.X(),pos2.Y())) nPoints++;                                     //photon is accepted if not in dead zone
260    }
261    area+=TMath::Abs((pos1-fMipPos).X()*(pos2-fMipPos).Y()-(pos1-fMipPos).Y()*(pos2-fMipPos).X());   //add area of the triangle...            
262    pos1 = pos2;
263   }
264 //---  find area and length of the ring;
265   fRingAcc = (Double_t)nPoints/(Double_t)kN;
266   area*=0.5;
267   fRingArea = area;
268 }//FindRingGeom()
269 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
270 TVector2 AliHMPIDRecon::IntWithEdge(TVector2 p1,TVector2 p2)const
271 {
272 // It finds the intersection of the line for 2 points traced as photons
273 // and the edge of a given PC
274 // Arguments: 2 points obtained tracing the photons
275 //   Returns: intersection point with detector (PC) edges
276
277   Double_t xmin = (p1.X()<p2.X())? p1.X():p2.X(); 
278   Double_t xmax = (p1.X()<p2.X())? p2.X():p1.X(); 
279   Double_t ymin = (p1.Y()<p2.Y())? p1.Y():p2.Y(); 
280   Double_t ymax = (p1.Y()<p2.Y())? p2.Y():p1.Y(); 
281   
282   Double_t m = TMath::Tan((p2-p1).Phi());
283   TVector2 pint;
284   //intersection with low  X
285   pint.Set((Double_t)(p1.X() + (0-p1.Y())/m),0.);
286   if(pint.X()>=0 && pint.X()<=fParam->SizeAllX() &&
287      pint.X()>=xmin && pint.X()<=xmax            &&
288      pint.Y()>=ymin && pint.Y()<=ymax) return pint;
289   //intersection with high X  
290   pint.Set((Double_t)(p1.X() + (fParam->SizeAllY()-p1.Y())/m),(Double_t)(fParam->SizeAllY()));
291   if(pint.X()>=0 && pint.X()<=fParam->SizeAllX() &&
292      pint.X()>=xmin && pint.X()<=xmax            &&
293      pint.Y()>=ymin && pint.Y()<=ymax) return pint;
294   //intersection with left Y  
295   pint.Set(0.,(Double_t)(p1.Y() + m*(0-p1.X())));
296   if(pint.Y()>=0 && pint.Y()<=fParam->SizeAllY() &&
297      pint.Y()>=ymin && pint.Y()<=ymax            &&
298      pint.X()>=xmin && pint.X()<=xmax) return pint;
299   //intersection with righ Y  
300   pint.Set((Double_t)(fParam->SizeAllX()),(Double_t)(p1.Y() + m*(fParam->SizeAllX()-p1.X())));
301   if(pint.Y()>=0 && pint.Y()<=fParam->SizeAllY() &&
302      pint.Y()>=ymin && pint.Y()<=ymax            &&
303      pint.X()>=xmin && pint.X()<=xmax) return pint;
304   return p1;
305 }//IntWithEdge()
306 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
307 Double_t AliHMPIDRecon::FindRingCkov(Int_t)
308 {
309 // 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
310 // collecting errors for all single Ckov candidates thetas. (Assuming they are independent)  
311 // Arguments: iNclus- total number of clusters in chamber for background estimation
312 //    Return: best estimation of track Theta ckov
313
314   Double_t wei = 0.;
315   Double_t weightThetaCerenkov = 0.;
316
317   Double_t ckovMin=9999.,ckovMax=0.;
318   Double_t sigma2 = 0;   //to collect error squared for this ring
319   
320   for(Int_t i=0;i<fPhotCnt;i++){//candidates loop
321     if(fPhotFlag[i] == 2){
322       if(fPhotCkov[i]<ckovMin) ckovMin=fPhotCkov[i];                         //find max and min Theta ckov from all candidates within probable window
323       if(fPhotCkov[i]>ckovMax) ckovMax=fPhotCkov[i]; 
324       weightThetaCerenkov += fPhotCkov[i]*fPhotWei[i];
325       wei += fPhotWei[i];                                                    //collect weight as sum of all candidate weghts   
326       
327       sigma2 += 1./Sigma2(fPhotCkov[i],fPhotPhi[i]);
328     }
329   }//candidates loop
330   
331   if(sigma2>0) fCkovSigma2=1./sigma2;
332   else         fCkovSigma2=1e10;  
333   
334   if(wei != 0.) weightThetaCerenkov /= wei; else weightThetaCerenkov = 0.;
335   return weightThetaCerenkov;
336 }//FindCkovRing()
337 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
338 Int_t AliHMPIDRecon::FlagPhot(Double_t ckov)
339 {
340 // Flag photon candidates if their individual ckov angle is inside the window around ckov angle returned by  HoughResponse()
341 // Arguments: ckov- value of most probable ckov angle for track as returned by HoughResponse()
342 //   Returns: number of photon candidates happened to be inside the window
343
344 // Photon Flag:  Flag = 0 initial set; 
345 //               Flag = 1 good candidate (charge compatible with photon); 
346 //               Flag = 2 photon used for the ring;
347   
348   Int_t steps = (Int_t)((ckov )/ fDTheta); //how many times we need to have fDTheta to fill the distance between 0  and thetaCkovHough
349
350   Double_t tmin = (Double_t)(steps - 1)*fDTheta;
351   Double_t tmax = (Double_t)(steps)*fDTheta;
352   Double_t tavg = 0.5*(tmin+tmax);
353
354   tmin = tavg - 0.5*fWindowWidth;  tmax = tavg + 0.5*fWindowWidth;
355
356   Int_t iInsideCnt = 0; //count photons which Theta ckov inside the window
357   for(Int_t i=0;i<fPhotCnt;i++){//photon candidates loop
358     fPhotFlag[i] = 0;
359     if(fPhotCkov[i] >= tmin && fPhotCkov[i] <= tmax)    { 
360       fPhotFlag[i]=2;     
361       iInsideCnt++;
362     }
363   }
364   return iInsideCnt;
365 }//FlagPhot()
366 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
367 TVector2 AliHMPIDRecon::TracePhot(Double_t ckovThe,Double_t ckovPhi)const
368 {
369 // Trace a single Ckov photon from emission point somewhere in radiator up to photocathode taking into account ref indexes of materials it travereses
370 // Arguments: ckovThe,ckovPhi- photon ckov angles in TRS, [rad]
371 //   Returns: distance between photon point on PC and track projection  
372   
373   Double_t theta,phi;
374   TVector3  dirTRS,dirLORS;   
375   dirTRS.SetMagThetaPhi(1,ckovThe,ckovPhi);                     //photon in TRS
376   Trs2Lors(dirTRS,theta,phi);
377   dirLORS.SetMagThetaPhi(1,theta,phi);                          //photon in LORS
378   return TraceForward(dirLORS);                                 //now foward tracing
379 }//TracePhot()
380 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
381 void AliHMPIDRecon::Propagate(const TVector3 dir,TVector3 &pos,Double_t z)const
382 {
383 // Finds an intersection point between a line and XY plane shifted along Z.
384 // Arguments:  dir,pos   - vector along the line and any point of the line
385 //             z         - z coordinate of plain 
386 //   Returns:  none
387 //   On exit:  pos is the position if this intesection if any
388   static TVector3 nrm(0,0,1); 
389          TVector3 pnt(0,0,z);
390   
391   TVector3 diff=pnt-pos;
392   Double_t sint=(nrm*diff)/(nrm*dir);
393   pos+=sint*dir;
394 }//Propagate()
395 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
396 void AliHMPIDRecon::Refract(TVector3 &dir,Double_t n1,Double_t n2)const
397 {
398 // Refract direction vector according to Snell law
399 // Arguments: 
400 //            n1 - ref idx of first substance
401 //            n2 - ref idx of second substance
402 //   Returns: none
403 //   On exit: dir is new direction
404   Double_t sinref=(n1/n2)*TMath::Sin(dir.Theta());
405   if(TMath::Abs(sinref)>1.) dir.SetXYZ(-999,-999,-999);
406   else             dir.SetTheta(TMath::ASin(sinref));
407 }//Refract()
408 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
409 Double_t AliHMPIDRecon::HoughResponse()
410 {
411 //
412 //    fIdxMip = mipId;
413
414 //       
415   Double_t kThetaMax=0.75;
416   Int_t nChannels = (Int_t)(kThetaMax/fDTheta+0.5);
417   TH1D *phots   = new TH1D("Rphot"  ,"phots"         ,nChannels,0,kThetaMax);
418   TH1D *photsw  = new TH1D("RphotWeighted" ,"photsw" ,nChannels,0,kThetaMax);
419   TH1D *resultw = new TH1D("resultw","resultw"       ,nChannels,0,kThetaMax);
420   Int_t nBin = (Int_t)(kThetaMax/fDTheta);
421   Int_t nCorrBand = (Int_t)(fWindowWidth/(2*fDTheta));
422   
423   for (Int_t i=0; i< fPhotCnt; i++){//photon cadidates loop
424     Double_t angle = fPhotCkov[i];  if(angle<0||angle>kThetaMax) continue;
425     phots->Fill(angle);
426     Int_t bin = (Int_t)(0.5+angle/(fDTheta));
427     Double_t weight=1.;
428     if(fIsWEIGHT){
429       Double_t lowerlimit = ((Double_t)bin)*fDTheta - 0.5*fDTheta;  Double_t upperlimit = ((Double_t)bin)*fDTheta + 0.5*fDTheta;
430       FindRingGeom(lowerlimit);
431       Double_t areaLow  = GetRingArea();
432       FindRingGeom(upperlimit);
433       Double_t areaHigh = GetRingArea();
434       Double_t diffArea = areaHigh - areaLow;
435       if(diffArea>0) weight = 1./diffArea;
436     }
437     photsw->Fill(angle,weight);
438     fPhotWei[i]=weight;
439   }//photon candidates loop 
440    
441   for (Int_t i=1; i<=nBin;i++){
442     Int_t bin1= i-nCorrBand;
443     Int_t bin2= i+nCorrBand;
444     if(bin1<1) bin1=1;
445     if(bin2>nBin)bin2=nBin;
446     Double_t sumPhots=phots->Integral(bin1,bin2);
447     if(sumPhots<3) continue;                            // if less then 3 photons don't trust to this ring
448     Double_t sumPhotsw=photsw->Integral(bin1,bin2);
449     resultw->Fill((Double_t)((i+0.5)*fDTheta),sumPhotsw);
450   } 
451 // evaluate the "BEST" theta ckov as the maximum value of histogramm
452   Double_t *pVec = resultw->GetArray();
453   Int_t locMax = TMath::LocMax(nBin,pVec);
454   delete phots;delete photsw;delete resultw; // Reset and delete objects
455   
456   return (Double_t)(locMax*fDTheta+0.5*fDTheta); //final most probable track theta ckov   
457 }//HoughResponse()
458 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
459 Double_t AliHMPIDRecon::Sigma2(Double_t ckovTh, Double_t ckovPh)const
460 {
461 // Analithical calculation of total error (as a sum of localization, geometrical and chromatic errors) on Cerenkov angle for a given Cerenkov photon 
462 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
463 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
464 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
465 //            MIP beta
466 //   Returns: absolute error on Cerenkov angle, [radians]    
467   
468   TVector3 v(-999,-999,-999);
469   Double_t trkBeta = 1./(TMath::Cos(ckovTh)*fParam->GetRefIdx());
470   
471   if(trkBeta > 1) trkBeta = 1;                 //protection against bad measured thetaCer  
472   if(trkBeta < 0) trkBeta = 0.0001;            //
473
474   v.SetX(SigLoc (ckovTh,ckovPh,trkBeta));
475   v.SetY(SigGeom(ckovTh,ckovPh,trkBeta));
476   v.SetZ(SigCrom(ckovTh,ckovPh,trkBeta));
477
478   return v.Mag2();
479 }
480 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
481 Double_t AliHMPIDRecon::SigLoc(Double_t thetaC, Double_t phiC,Double_t betaM)const
482 {
483 // Analithical calculation of localization error (due to finite segmentation of PC) on Cerenkov angle for a given Cerenkov photon 
484 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
485 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
486 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
487 //            MIP beta
488 //   Returns: absolute error on Cerenkov angle, [radians]    
489   
490   Double_t phiDelta = phiC - fTrkDir.Phi();
491
492   Double_t sint     = TMath::Sin(fTrkDir.Theta());
493   Double_t cost     = TMath::Cos(fTrkDir.Theta());
494   Double_t sinf     = TMath::Sin(fTrkDir.Phi());
495   Double_t cosf     = TMath::Cos(fTrkDir.Phi());
496   Double_t sinfd    = TMath::Sin(phiDelta);
497   Double_t cosfd    = TMath::Cos(phiDelta);
498   Double_t tantheta = TMath::Tan(thetaC);
499   
500   Double_t alpha =cost-tantheta*cosfd*sint;                                                 // formula (11)
501   Double_t k = 1.-fParam->GetRefIdx()*fParam->GetRefIdx()+alpha*alpha/(betaM*betaM);        // formula (after 8 in the text)
502   if (k<0) return 1e10;
503   Double_t mu =sint*sinf+tantheta*(cost*cosfd*sinf+sinfd*cosf);                             // formula (10)
504   Double_t e  =sint*cosf+tantheta*(cost*cosfd*cosf-sinfd*sinf);                             // formula (9)
505
506   Double_t kk = betaM*TMath::Sqrt(k)/(fParam->GapThick()*alpha);                            // formula (6) and (7)
507   Double_t dtdxc = kk*(k*(cosfd*cosf-cost*sinfd*sinf)-(alpha*mu/(betaM*betaM))*sint*sinfd); // formula (6)           
508   Double_t dtdyc = kk*(k*(cosfd*sinf+cost*sinfd*cosf)+(alpha* e/(betaM*betaM))*sint*sinfd); // formula (7)            pag.4
509
510   Double_t errX = 0.2,errY=0.25;                                                            //end of page 7
511   return  TMath::Sqrt(errX*errX*dtdxc*dtdxc + errY*errY*dtdyc*dtdyc);
512 }
513 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
514 Double_t AliHMPIDRecon::SigCrom(Double_t thetaC, Double_t phiC,Double_t betaM)const
515 {
516 // Analithical calculation of chromatic error (due to lack of knowledge of Cerenkov photon energy) on Cerenkov angle for a given Cerenkov photon 
517 // created by a given MIP. Fromulae according to CERN-EP-2000-058 
518 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
519 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
520 //            MIP beta
521 //   Returns: absolute error on Cerenkov angle, [radians]    
522   
523   Double_t phiDelta = phiC - fTrkDir.Phi();
524
525   Double_t sint     = TMath::Sin(fTrkDir.Theta());
526   Double_t cost     = TMath::Cos(fTrkDir.Theta());
527   Double_t cosfd    = TMath::Cos(phiDelta);
528   Double_t tantheta = TMath::Tan(thetaC);
529   
530   Double_t alpha =cost-tantheta*cosfd*sint;                                                 // formula (11)
531   Double_t dtdn = cost*fParam->GetRefIdx()*betaM*betaM/(alpha*tantheta);                    // formula (12)
532             
533 //  Double_t f = 0.00928*(7.75-5.635)/TMath::Sqrt(12.);
534   Double_t f = 0.0172*(7.75-5.635)/TMath::Sqrt(24.);
535
536   return f*dtdn;
537 }//SigCrom()
538 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
539 Double_t AliHMPIDRecon::SigGeom(Double_t thetaC, Double_t phiC,Double_t betaM)const
540 {
541 // Analithical calculation of geometric error (due to lack of knowledge of creation point in radiator) on Cerenkov angle for a given Cerenkov photon 
542 // created by a given MIP. Formulae according to CERN-EP-2000-058 
543 // Arguments: Cerenkov and azimuthal angles for Cerenkov photon, [radians]
544 //            dip and azimuthal angles for MIP taken at the entrance to radiator, [radians]        
545 //            MIP beta
546 //   Returns: absolute error on Cerenkov angle, [radians]    
547
548   Double_t phiDelta = phiC - fTrkDir.Phi();
549
550   Double_t sint     = TMath::Sin(fTrkDir.Theta());
551   Double_t cost     = TMath::Cos(fTrkDir.Theta());
552   Double_t sinf     = TMath::Sin(fTrkDir.Phi());
553   Double_t cosfd    = TMath::Cos(phiDelta);
554   Double_t costheta = TMath::Cos(thetaC);
555   Double_t tantheta = TMath::Tan(thetaC);
556   
557   Double_t alpha =cost-tantheta*cosfd*sint;                                                  // formula (11)
558   
559   Double_t k = 1.-fParam->GetRefIdx()*fParam->GetRefIdx()+alpha*alpha/(betaM*betaM);         // formula (after 8 in the text)
560   if (k<0) return 1e10;
561
562   Double_t eTr = 0.5*fParam->RadThick()*betaM*TMath::Sqrt(k)/(fParam->GapThick()*alpha);     // formula (14)
563   Double_t lambda = 1.-sint*sint*sinf*sinf;                                                  // formula (15)
564
565   Double_t c1 = 1./(1.+ eTr*k/(alpha*alpha*costheta*costheta));                              // formula (13.a)
566   Double_t c2 = betaM*TMath::Power(k,1.5)*tantheta*lambda/(fParam->GapThick()*alpha*alpha);  // formula (13.b)
567   Double_t c3 = (1.+eTr*k*betaM*betaM)/((1+eTr)*alpha*alpha);                                // formula (13.c)
568   Double_t c4 = TMath::Sqrt(k)*tantheta*(1-lambda)/(fParam->GapThick()*betaM);               // formula (13.d)
569   Double_t dtdT = c1 * (c2+c3*c4);
570   Double_t trErr = fParam->RadThick()/(TMath::Sqrt(12.)*cost);
571
572   return trErr*dtdT;
573 }//SigGeom()