]> git.uio.no Git - u/mrichter/AliRoot.git/blame_incremental - HMPID/AliHMPIDTracker.cxx
Coverity fixed
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDTracker.cxx
... / ...
CommitLineData
1#include "AliHMPIDTracker.h" //class header
2#include "AliHMPIDtrack.h" //class header
3#include "AliHMPIDCluster.h" //GetTrackPoint(),PropagateBack()
4#include "AliHMPIDParam.h" //GetTrackPoint(),PropagateBack()
5#include "AliHMPIDPid.h" //Recon(),reconHTA()
6#include "AliHMPIDRecon.h" //Recon()
7#include "AliHMPIDRecoParamV1.h" //Recon()
8#include "AliHMPIDReconstructor.h"//Recon()
9#include "AliHMPIDReconHTA.h" //ReconHTA()
10#include <AliLog.h> //Recon()
11#include <AliESDEvent.h> //PropagateBack(),Recon()
12#include <AliESDtrack.h> //Intersect()
13#include <AliTracker.h>
14#include <AliRun.h> //GetTrackPoint(),PropagateBack()
15#include <AliTrackPointArray.h> //GetTrackPoint()
16#include <AliAlignObj.h> //GetTrackPoint()
17#include <AliCDBManager.h> //PropageteBack()
18#include <AliCDBEntry.h> //PropageteBack()
19#include "TTreeStream.h" // debug streamer
20//
21// HMPID base class fo tracking
22//.
23//.
24//.
25ClassImp(AliHMPIDTracker)
26//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27AliHMPIDTracker::AliHMPIDTracker():
28 AliTracker(),
29 fClu(new TObjArray(AliHMPIDParam::kMaxCh+1)),
30 fDebugStreamer(0)
31{
32// ctor. Create TObjArray of TClonesArray of AliHMPIDCluster
33//
34//
35 fClu->SetOwner(kTRUE);
36 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
37 fDebugStreamer = new TTreeSRedirector("HMPIDdebug.root");
38
39}//ctor
40
41
42AliHMPIDTracker::~AliHMPIDTracker(){
43 //
44 // destructor
45 //
46 delete fClu;
47 if (fDebugStreamer) delete fDebugStreamer;
48}
49
50//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
51Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
52{
53// Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
54// MIP cluster is reffered by index which is stored in AliESDtrack ???????
55// Arguments: idx- cluster index which is stored by HMPID in AliESDtrack
56// point- reference to the object where to store the point
57// Returns: status of operation if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track.
58 if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
59 Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
60 iClu = iClu%1000; //GetHMPIDcluIdx -> 1e+6*ch + 1e+3*clusize + cluIdx;
61 point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham));//layer and chamber number
62 TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
63 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
64 Float_t xyz[3];
65 pClu->GetGlobalXYZ(xyz);
66 Float_t cov[6];
67 pClu->GetGlobalCov(cov);
68 point.SetXYZ(xyz,cov);
69 return kTRUE;
70}//GetTrackPoint()
71//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
72Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
73{
74// Static method to find intersection in between given track and HMPID chambers
75// Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
76// Returns: intersected chamber ID or -1
77 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
78 for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){ //chambers loop
79 Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
80 if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;}
81 } //chambers loop
82 delete hmpTrk; hmpTrk=0x0;
83 return -1; //no intersection with HMPID chambers
84}//IntTrkCha()
85//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
86Int_t AliHMPIDTracker::IntTrkCha(Int_t ch,AliHMPIDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
87{
88// Static method to find intersection in between given track and HMPID chambers
89// Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm]
90// Returns: intersected chamber ID or -1
91 AliHMPIDParam *pParam=AliHMPIDParam::Instance();
92 Double_t p1[3],n1[3];
93 pParam->Norm(ch,n1);
94 pParam->Point(ch,p1,AliHMPIDParam::kRad); //point & norm for middle of radiator plane
95 Double_t p2[3],n2[3];
96 pParam->Norm(ch,n2);
97 pParam->Point(ch,p2,AliHMPIDParam::kPc); //point & norm for entrance to PC plane
98 if(pTrk->Intersect(p1,n1)==kFALSE) return -1; //try to intersect track with the middle of radiator
99 if(pTrk->Intersect(p2,n2)==kFALSE) return -1;
100 pParam->Mars2LorsVec(ch,n1,theta,phi); //track angles at RAD
101 pParam->Mars2Lors (ch,p1,xRa,yRa); //TRKxRAD position
102 pParam->Mars2Lors (ch,p2,xPc,yPc); //TRKxPC position
103 if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch; //return intersected chamber
104 return -1; //no intersection with HMPID chambers
105}//IntTrkCha()
106//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
107Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
108{
109// Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
110// Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
111// Returns: error code (currently ignored in AliReconstruction::RunTraking())
112 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
113 pCluTree->GetEntry(0);
114 return 0;
115}//LoadClusters()
116//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
117Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
118{
119// Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
120// Agruments: pEsd - pointer to ESD
121// Returns: error code
122 AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
123 AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
124 if(!pNmeanEnt) AliError("No Nmean C6F14 ");
125 if(!pQthreEnt) AliError("No Qthre");
126
127 return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());
128}//PropagateBack()
129//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
130Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
131{
132// Static method to reconstruct the Cherenkov angle for all valid tracks of a given event.
133// Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
134// Returns: error code, 0 if no errors
135//
136// Algortihm: Loop over tracks
137//
138// 1. Find the closest MIP cluster using fast tracks extrapolation method
139// 2. Propagate track to the MIP cluster using the STEER method
140// 3. Update the track information with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
141// 4. Propagate back the constrained track to the radiator radius ( not exact yet)
142// 5. Propagation in the last 10 cm with the fast method
143// 6. Set ESDtrack information
144// 7. Calculate the Cherenkov angle
145
146 const Double_t kMaxSnp=0.9; //maximal snp for prolongation
147 const Double_t kdRadiator=10; // distance between radiator and the plane
148 AliHMPIDRecon recon; //instance of reconstruction class, nothing important in ctor
149 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
150 Float_t xPc,yPc,xRa,yRa,theta,phi;
151
152 Double_t cluLORS[2]={0};
153 Int_t nMipClusTot=0;
154
155 Double_t qthre = 0; Double_t nmean=0; Int_t hvsec=0;
156 Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
157
158 Bool_t tsRight = kTRUE;
159
160 UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); //
161 UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); //
162 UInt_t ts = pEsd->GetTimeStamp();
163
164 if(ts<tsmin || ts>tsmax) {
165 AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
166 tsRight = kFALSE;
167 }
168
169 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
170
171 // Double_t bestChi2=99999;chi2=99999; //init. track matching params
172 Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1};
173
174 Bool_t isOkDcut=kFALSE;
175 Bool_t isOkQcut=kFALSE;
176 Bool_t isMatched=kFALSE;
177
178 AliHMPIDCluster *bestHmpCluster=0x0; //the best matching cluster
179 AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //get reconstructed track
180
181 if(!pTrk->IsOn(AliESDtrack::kTPCout)) continue;
182
183 if(pTrk->IsOn(AliESDtrack::kTPCrefit)) continue;
184 AliESDfriendTrack *ftrack= (AliESDfriendTrack *)pTrk->GetFriendTrack();
185 if (!ftrack) continue;
186 if (!ftrack->GetTPCOut()) continue;
187 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
188 AliHMPIDtrack *hmpTrkConstrained = 0; //create a hmpid track to be used for propagation and matching
189 hmpTrk->Set(ftrack->GetTPCOut()->GetX(), ftrack->GetTPCOut()->GetAlpha(),ftrack->GetTPCOut()->GetParameter(), ftrack->GetTPCOut()->GetCovariance());
190 //
191 bz=AliTracker::GetBz();
192
193 //initial flags for HMPID ESD infos
194 pTrk->SetHMPIDtrk(0,0,0,0); //no intersection found
195 pTrk->SetHMPIDmip(0,0,0,0); //store mip info in any case
196 pTrk->SetHMPIDcluIdx(99,99999); //chamber not found, mip not yet considered
197 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed); //ring reconstruction not yet performed
198
199 Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi); //find the intersected chamber for this track
200 if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;} //no intersection at all, go after next track
201
202 pTrk->SetHMPIDtrk(xPc,yPc,theta,phi); //store initial infos
203 pTrk->SetHMPIDcluIdx(ipCh,9999); //set chamber, index of cluster + cluster size
204
205 // track intersects the chamber ipCh: find the MIP
206
207 TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters
208 nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber
209 nClusCh[ipCh] = nMipClusTot;
210
211 if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;}
212
213 Int_t index=-1; //index of the "best" matching cluster
214 //
215 // 1. Find the closest MIP cluster using fast tracks extrapolation method
216 //
217 for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop
218
219 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster
220 // evaluate qThre
221 if(tsRight){
222 if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {
223 qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts); //
224 } else { // in the past just 1 qthre
225 hvsec = pParam->InHVSector(pClu->Y()); // per chamber
226 if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts); //
227 }
228 } else qthre = pParam->QCut();
229
230 //
231 if(pClu->Q()<qthre) continue; //charge compartible with MIP clusters
232 isOkQcut = kTRUE;
233 //
234 cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y(); //get the LORS coordinates of the cluster
235 Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
236
237 if(dist<dmin) {
238 dmin = dist;
239 index=iClu;
240 bestHmpCluster=pClu;
241 }
242 } // clusters loop
243
244 // moved down
245 /*if(!isOkQcut) {
246 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
247 delete hmpTrk;hmpTrk=0x0; continue;
248 }*/
249
250 //
251 // 2. Propagate track to the MIP cluster using the STEER method
252 //
253
254 if(!bestHmpCluster) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
255
256 TVector3 vG = pParam->Lors2Mars(ipCh,bestHmpCluster->X(),bestHmpCluster->Y());
257 Double_t gx = vG[0];
258 Double_t gy = vG[1];
259 Double_t gz = vG[2];
260 Double_t alpha=TMath::ATan2(gy,gx);
261 Double_t radiusH=TMath::Sqrt(gy*gy+gx*gx);
262 if (!(hmpTrk->Rotate(alpha,kTRUE))) continue;
263 if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrk,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
264 //
265 // 3. Update the track with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
266 //
267 AliExternalTrackParam trackC(*hmpTrk);
268 Double_t posC[2]={0,gz};
269 Double_t covC[3]={0.1*0.1, 0, 0.1*0.1};
270 trackC.Update(posC,covC);
271 //
272 // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
273 //
274 hmpTrkConstrained = new AliHMPIDtrack(*pTrk);
275 hmpTrkConstrained->Set(trackC.GetX(), trackC.GetAlpha(),trackC.GetParameter(), trackC.GetCovariance());
276 if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrkConstrained,radiusH-kdRadiator,pTrk->GetMass(),1,kFALSE,kMaxSnp,1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;}
277 //
278 // 5. Propagation in the last 10 cm with the fast method
279 //
280 Float_t xPc0=0, yPc0=0;
281
282 IntTrkCha(ipCh, hmpTrk, xPc0,yPc0,xRa,yRa,theta,phi);
283 IntTrkCha(ipCh, hmpTrkConstrained, xPc,yPc,xRa,yRa,theta,phi);
284 //
285 // 6. Set ESDtrack information
286 //
287 Int_t cluSiz = bestHmpCluster->Size();
288 pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0); //store mip info in any case
289 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set chamber, index of cluster + cluster size
290 pTrk->SetHMPIDtrk(xPc0,yPc0,theta,phi);
291 //
292 //
293 // Dump debug info if specified
294 //
295 if (AliHMPIDReconstructor::StreamLevel()>0) {
296 AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
297 AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);
298 if(!trackTPC->Rotate(alpha)) continue;
299 if(!trackCurrent->Rotate(alpha)) continue;
300 Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
301 Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
302 AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
303 if(!trackTPCNB->Rotate(alpha)) continue;
304 Bool_t statusTPCNB=kTRUE;
305 Double_t bfield[3]={0,0,0};
306 for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){
307 Double_t xyz[3];
308 trackTPCNB->GetXYZ(xyz);
309 GetBxByBz(xyz,bfield);
310 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield);
311 }
312 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield);
313 //
314 Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
315 Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());
316 //
317 AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
318 Double_t pos[2]={0,gz};
319 Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
320 Double_t chi2C = trackTPCConstrained->GetPredictedChi2(pos,cov);
321 trackTPCConstrained->Update(pos,cov);
322 (*fDebugStreamer)<<"track"<<
323 "rH="<<radiusH<< // radius of cluster
324 "angle="<<tanAlpha<< // tan of the local inlination angle
325 "dC="<<deltaC<< // delta of the curvature
326 "trackTPC.="<<trackTPC<< // TPC outer param extrapolated to the HMPID
327 "trackTPCNB.="<<trackTPCNB<< // TPC track prpagated with material budget correction
328 "chi2C="<<chi2C<<
329 "trackTPCC.="<<trackTPCConstrained<< // TPC outer param extrapolated to the HMPID constrained
330 "trackCurrent.="<<trackCurrent<< // current track extrapolated to the HMPID
331 "sTPC="<<statusTPC<< // status for the current TPC track
332 "sCurrent="<<statusCurrent<< // status for the current global track
333 "cl.="<<bestHmpCluster<< // HMPID cluster
334 //
335 "t.="<<pTrk<< // curent esd track
336 "ft.="<<ftrack<< // friend track
337 "hmpTrk.="<<hmpTrk<< // hmpid tracks as used in the following code
338 "hmpTrkC.="<<hmpTrkConstrained<< // constrained hmpid tracks as used in the following code
339 "gx="<<gx<< // global cluster position X
340 "gy="<<gy<< // Y
341 "gz="<<gz<< // Z
342 "\n";
343 }
344 //
345 //
346 //
347 if(!isOkQcut) {
348 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
349 delete hmpTrk;hmpTrk=0x0;
350 delete hmpTrkConstrained;hmpTrkConstrained=0x0;
351 continue;
352 }
353
354 if(AliHMPIDReconstructor::GetRecoParam()) //retrieve distance cut
355 {
356 if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE) //distance cut is fixed number
357 {
358 distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
359 }
360 else
361 {
362 for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar); //prevision: distance cut is function of momentum
363 distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
364 }
365 }
366 else {distCut=pParam->DistCut();}
367
368 //dmin recalculated
369
370 dmin = TMath::Sqrt((xPc0-bestHmpCluster->X())*(xPc0-bestHmpCluster->X())+(yPc0-bestHmpCluster->Y())*(yPc0-bestHmpCluster->Y()));
371
372 if(dmin < distCut) {
373 isOkDcut = kTRUE;
374 }
375 //isOkDcut = kTRUE; // switch OFF cut
376
377 if(!isOkDcut) {
378 pTrk->SetHMPIDsignal(pParam->kMipDistCut); //closest cluster with enough charge is still too far from intersection
379 }
380
381 if(isOkQcut*isOkDcut) isMatched = kTRUE; // MIP-Track matched !!
382
383 if(!isMatched) {delete hmpTrk;hmpTrk=0x0;delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;} // If matched continue...
384
385 Bool_t isOk = kTRUE;
386 if(!isOk) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
387 pTrk->SetOuterHmpParam(hmpTrkConstrained,AliESDtrack::kHMPIDout);
388
389 FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
390
391 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
392
393 //evaluate nMean
394 if(tsRight){
395 if(pNmean->GetEntries()==21) { //for backward compatibility
396 nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts); //C6F14 Nmean for this chamber
397 } else {
398 if(iRad < 0) {
399 nmean = -1;
400 } else {
401 Double_t tLow = ((TF1*)pNmean->At(6*ipCh+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber
402 Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber
403 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
404 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
405 }
406 if(nmean < 0){ //track didn' t pass through the radiator
407 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
408 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set index of cluster
409 delete hmpTrk;hmpTrk=0x0;
410 delete hmpTrkConstrained;hmpTrkConstrained=0x0;
411 continue;
412 }
413 }
414 } else nmean = pParam->MeanIdxRad();
415 //
416 // 7. Calculate the Cherenkov angle
417 //
418 recon.SetImpPC(xPc0,yPc0); //store track impact to PC
419 recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa); //search for Cerenkov angle of this track
420
421 Double_t thetaCkov = pTrk->GetHMPIDsignal();
422
423 if (AliHMPIDReconstructor::StreamLevel()>0) {
424 AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
425 AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);
426 if(!trackTPC->Rotate(alpha)) continue;
427 if(!trackCurrent->Rotate(alpha)) continue;
428 Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
429 Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
430 Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
431 Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());
432 //
433 AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
434 if(!trackTPCNB->Rotate(alpha)) continue;
435 Bool_t statusTPCNB=kTRUE;
436 Double_t bfield[3]={0,0,0};
437 for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){
438 Double_t xyz[3];
439 trackTPCNB->GetXYZ(xyz);
440 GetBxByBz(xyz,bfield);
441 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield);
442 }
443 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield);
444
445 AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
446 Double_t pos[2]={0,gz};
447 Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
448 Double_t chi2C = trackTPCConstrained->GetPredictedChi2(pos,cov);
449 trackTPCConstrained->Update(pos,cov);
450 (*fDebugStreamer)<<"track2"<<
451 "rH="<<radiusH<< // radius of cluster
452 "angle="<<tanAlpha<< // tan of the local inlination angle
453 "dC="<<deltaC<< // delta of the curvature
454 "trackTPC.="<<trackTPC<< // TPC outer param extrapolated to the HMPID
455 "trackTPCNB.="<<trackTPCNB<< // TPC outer param extrapolated to the HMPID
456 "chi2C="<<chi2C<<
457 "trackTPCC.="<<trackTPCConstrained<< // TPC outer param extrapolated to the HMPID constrained
458 "trackCurrent.="<<trackCurrent<< // current track extrapolated to the HMPID
459 "sTPC="<<statusTPC<< // status for the current TPC track
460 "sCurrent="<<statusCurrent<< // status for the current global track
461 "cl.="<<bestHmpCluster<< // HMPID cluster
462 //
463 "t.="<<pTrk<< // curent esd track
464 "ft.="<<ftrack<< // friend track
465 "hmpTrk.="<<hmpTrk<< // hmpid tracks as used in the following code
466 "hmpTrkC.="<<hmpTrkConstrained<< // constrained hmpid tracks as used in the following code
467 "gx="<<gx<< // global cluster position X
468 "gy="<<gy<< // Y
469 "gz="<<gz<< // Z
470 "thetaCkov="<<thetaCkov<< // Cherenkov angle
471 "\n";
472 }
473
474 if(pTrk->GetHMPIDsignal()<0) {
475 delete hmpTrk;hmpTrk=0x0;
476 delete hmpTrkConstrained; hmpTrkConstrained=0x0;
477 continue;}
478
479 AliHMPIDPid pID;
480 Double_t prob[5];
481 pID.FindPid(pTrk,5,prob);
482 pTrk->SetHMPIDpid(prob);
483// Printf(" Prob e- %6.2f mu %6.2f pi %6.2f k %6.2f p %6.2f",prob[0]*100,prob[1]*100,prob[2]*100,prob[3]*100,prob[4]*100);
484 delete hmpTrk; hmpTrk=0x0;
485 delete hmpTrkConstrained; hmpTrkConstrained=0x0;
486 }//iTrk
487
488 return 0; // error code: 0=no error;
489}//Recon()
490//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
491Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
492{
493// Static method to reconstruct Theta Ckov for all valid tracks of a given event.
494// Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time), pQthre - pointer to all function Qthre=f(time)
495// Returns: error code, 0 if no errors
496
497 AliHMPIDReconHTA reconHTA; //instance of reconstruction class, nothing important in ctor
498
499 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
500
501 Bool_t tsRight = kTRUE;
502
503 UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); //
504 UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); //
505 UInt_t ts = pEsd->GetTimeStamp();
506
507 if(ts<tsmin || ts>tsmax) {
508 AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
509 tsRight = kFALSE;
510 }
511
512 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
513
514 AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //here it is simulated or just empty track
515 Int_t ipCh;
516 ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000;
517 if(ipCh<0) continue;
518
519 TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters
520 Int_t nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber
521
522 Double_t qMip=-1;
523 Int_t chMip=-1;
524 Double_t xMip = 0;
525 Double_t yMip = 0;
526 Int_t indMip = 0;
527 Int_t cluMipSiz = 0;
528
529 for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop
530
531 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster
532 Double_t qClus = pClu->Q();
533 if(qClus>qMip) {
534 qMip = qClus;
535 chMip = pClu->Ch();
536 xMip = pClu->X();
537 yMip = pClu->Y();
538 indMip = iClu;
539 cluMipSiz = pClu->Size();
540 }
541 }//clus loop
542
543 if(chMip<0) return 1;
544
545 Int_t hvsec;
546 Double_t qthre=0;
547
548// evaluate qThre
549 if(tsRight){
550 if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) { // just for backward compatibility
551 qthre=((TF1*)pQthre->At(chMip))->Eval(ts); //
552 } else { // in the past just 1 qthre
553 hvsec = pParam->InHVSector(yMip); // per chamber
554 if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(ts); //
555 }
556 } else qthre = pParam->QCut();
557//
558 if(qMip<qthre) {
559 pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case
560 pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);
561 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
562 return 1; //charge compatible with MIP clusters
563 }
564
565 pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case
566 pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);
567
568 Double_t yRa = yMip; //just an approx...
569 Double_t nmean;
570
571 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
572
573 //evaluate nMean
574 if(tsRight){
575 if(pNmean->GetEntries()==21) { //for backward compatibility
576 nmean=((TF1*)pNmean->At(3*chMip))->Eval(ts); //C6F14 Nmean for this chamber
577 } else {
578 if(iRad < 0) {
579 nmean = -1;
580 } else {
581 Double_t tLow = ((TF1*)pNmean->At(6*chMip+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber
582 Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber
583 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
584 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
585 }
586 if(nmean < 0){ //track didn' t pass through the radiator
587 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
588 return 1;
589 }
590 }
591 } else nmean = pParam->MeanIdxRad();
592 //
593 if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) { //search for track parameters and Cerenkov angle of this track
594 AliHMPIDPid pID;
595 Double_t prob[5];
596 pID.FindPid(pTrk,5,prob);
597 pTrk->SetHMPIDpid(prob);
598 }
599// Printf(" Prob e- %6.2f mu %6.2f pi %6.2f k %6.2f p %6.2f",prob[0]*100,prob[1]*100,prob[2]*100,prob[3]*100,prob[4]*100);
600 }//iTrk
601
602 return 0; // error code: 0=no error;
603
604}//Recon()
605//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
606void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
607
608 // Publishes all pointers to clusters known to the tracker into the
609 // passed object array.
610 // The ownership is not transfered - the caller is not expected to delete
611 // the clusters
612
613 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
614 TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
615 for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
616 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);
617 array->AddLast(pClu);
618 }//cluster loop in iCh
619 pCluArr->Delete();
620 }//Ch loop
621
622 return;
623}
624//++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++