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