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