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