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