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
21 // HMPID base class fo tracking
25 ClassImp(AliHMPIDTracker)
26 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
27 AliHMPIDTracker::AliHMPIDTracker():
29 fClu(new TObjArray(AliHMPIDParam::kMaxCh+1)),
32 // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster
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");
42 AliHMPIDTracker::~AliHMPIDTracker(){
47 if (fDebugStreamer) delete fDebugStreamer;
50 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
51 Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
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 if(!pArr) return kFALSE;
64 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
65 if(!pClu) return kFALSE;
67 pClu->GetGlobalXYZ(xyz);
69 pClu->GetGlobalCov(cov);
70 point.SetXYZ(xyz,cov);
73 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
74 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
76 // Static method to find intersection in between given track and HMPID chambers
77 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
78 // Returns: intersected chamber ID or -1
79 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
80 for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){ //chambers loop
81 Int_t chInt = IntTrkCha(i,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
82 if(chInt>=0) {delete hmpTrk;hmpTrk=0x0;return chInt;}
84 delete hmpTrk; hmpTrk=0x0;
85 return -1; //no intersection with HMPID chambers
87 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
88 Int_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)
90 // Static method to find intersection in between given track and HMPID chambers
91 // Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm]
92 // Returns: intersected chamber ID or -1
93 AliHMPIDParam *pParam=AliHMPIDParam::Instance();
96 pParam->Point(ch,p1,AliHMPIDParam::kRad); //point & norm for middle of radiator plane
99 pParam->Point(ch,p2,AliHMPIDParam::kPc); //point & norm for entrance to PC plane
100 if(pTrk->Intersect(p1,n1)==kFALSE) return -1; //try to intersect track with the middle of radiator
101 if(pTrk->Intersect(p2,n2)==kFALSE) return -1;
102 pParam->Mars2LorsVec(ch,n1,theta,phi); //track angles at RAD
103 pParam->Mars2Lors (ch,p1,xRa,yRa); //TRKxRAD position
104 pParam->Mars2Lors (ch,p2,xPc,yPc); //TRKxPC position
105 if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch; //return intersected chamber
106 return -1; //no intersection with HMPID chambers
108 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
109 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
111 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
112 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
113 // Returns: error code (currently ignored in AliReconstruction::RunTraking())
114 for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
115 pCluTree->GetEntry(0);
118 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
119 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
121 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
122 // Agruments: pEsd - pointer to ESD
123 // Returns: error code
124 AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
125 AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
126 if(!pNmeanEnt) AliError("No Nmean C6F14 ");
127 if(!pQthreEnt) AliError("No Qthre");
129 return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());
131 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
132 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
134 // Static method to reconstruct the Cherenkov angle for all valid tracks of a given event.
135 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
136 // Returns: error code, 0 if no errors
138 // Algortihm: Loop over tracks
140 // 1. Find the closest MIP cluster using fast tracks extrapolation method
141 // 2. Propagate track to the MIP cluster using the STEER method
142 // 3. Update the track information with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
143 // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
144 // 5. Propagation in the last 10 cm with the fast method
145 // 6. Set ESDtrack information
146 // 7. Calculate the Cherenkov angle
148 const Double_t kMaxSnp=0.9; //maximal snp for prolongation
149 const Double_t kdRadiator=10; // distance between radiator and the plane
150 AliHMPIDRecon recon; //instance of reconstruction class, nothing important in ctor
151 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
152 Float_t xPc,yPc,xRa,yRa,theta,phi;
154 Double_t cluLORS[2]={0};
157 Double_t qthre = 0; Double_t nmean=0; Int_t hvsec=0;
158 // Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
160 Bool_t tsRight = kTRUE;
162 UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); //
163 UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); //
164 UInt_t ts = pEsd->GetTimeStamp();
166 if(ts<tsmin || ts>tsmax) {
167 AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
171 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
173 // Double_t bestChi2=99999;chi2=99999; //init. track matching params
174 // Double_t dmin=999999,bz=0,distCut=1,distParams[5]={1};
175 Double_t dmin=999999,distCut=1,distParams[5]={1};
177 Bool_t isOkDcut=kFALSE;
178 Bool_t isOkQcut=kFALSE;
179 Bool_t isMatched=kFALSE;
181 AliHMPIDCluster *bestHmpCluster=0x0; //the best matching cluster
182 AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //get reconstructed track
184 if(!pTrk->IsOn(AliESDtrack::kTPCout)) continue;
186 if(pTrk->IsOn(AliESDtrack::kTPCrefit)) continue;
187 AliESDfriendTrack *ftrack= (AliESDfriendTrack *)pTrk->GetFriendTrack();
188 if (!ftrack) continue;
189 if (!ftrack->GetTPCOut()) continue;
190 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
191 AliHMPIDtrack *hmpTrkConstrained = 0; //create a hmpid track to be used for propagation and matching
192 hmpTrk->Set(ftrack->GetTPCOut()->GetX(), ftrack->GetTPCOut()->GetAlpha(),ftrack->GetTPCOut()->GetParameter(), ftrack->GetTPCOut()->GetCovariance());
194 //bz=AliTracker::GetBz();
196 //initial flags for HMPID ESD infos
197 pTrk->SetHMPIDtrk(0,0,0,0); //no intersection found
198 pTrk->SetHMPIDmip(0,0,0,0); //store mip info in any case
199 pTrk->SetHMPIDcluIdx(99,99999); //chamber not found, mip not yet considered
200 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed); //ring reconstruction not yet performed
202 Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi); //find the intersected chamber for this track
203 if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;} //no intersection at all, go after next track
205 pTrk->SetHMPIDtrk(xPc,yPc,theta,phi); //store initial infos
206 pTrk->SetHMPIDcluIdx(ipCh,9999); //set chamber, index of cluster + cluster size
208 // track intersects the chamber ipCh: find the MIP
210 TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters
211 nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber
212 // nClusCh[ipCh] = nMipClusTot;
214 if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;}
216 Int_t index=-1; //index of the "best" matching cluster
218 // 1. Find the closest MIP cluster using fast tracks extrapolation method
220 for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop
222 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster
225 if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {
226 qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts); //
227 } else { // in the past just 1 qthre
228 hvsec = pParam->InHVSector(pClu->Y()); // per chamber
229 if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts); //
231 } else qthre = pParam->QCut();
234 if(pClu->Q()<qthre) continue; //charge compartible with MIP clusters
237 cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y(); //get the LORS coordinates of the cluster
238 Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
249 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
250 delete hmpTrk;hmpTrk=0x0; continue;
254 // 2. Propagate track to the MIP cluster using the STEER method
257 if(!bestHmpCluster) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
259 TVector3 vG = pParam->Lors2Mars(ipCh,bestHmpCluster->X(),bestHmpCluster->Y());
263 Double_t alpha=TMath::ATan2(gy,gx);
264 Double_t radiusH=TMath::Sqrt(gy*gy+gx*gx);
265 if (!(hmpTrk->Rotate(alpha,kTRUE))) continue;
266 if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrk,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
268 // 3. Update the track with MIP cluster (Improved angular and position resolution - to be used for Cherenkov angle calculation)
270 AliExternalTrackParam trackC(*hmpTrk);
271 Double_t posC[2]={0,gz};
272 Double_t covC[3]={0.1*0.1, 0, 0.1*0.1};
273 trackC.Update(posC,covC);
275 // 4. Propagate back the constrained track to the radiator radius ( not exact yet)
277 hmpTrkConstrained = new AliHMPIDtrack(*pTrk);
278 hmpTrkConstrained->Set(trackC.GetX(), trackC.GetAlpha(),trackC.GetParameter(), trackC.GetCovariance());
279 if(!AliTrackerBase::PropagateTrackToBxByBz(hmpTrkConstrained,radiusH-kdRadiator,pTrk->GetMass(),1,kFALSE,kMaxSnp,1)) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;}
281 // 5. Propagation in the last 10 cm with the fast method
283 Float_t xPc0=0, yPc0=0;
285 IntTrkCha(ipCh, hmpTrk, xPc0,yPc0,xRa,yRa,theta,phi);
286 IntTrkCha(ipCh, hmpTrkConstrained, xPc,yPc,xRa,yRa,theta,phi);
288 // 6. Set ESDtrack information
290 Int_t cluSiz = bestHmpCluster->Size();
291 pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0); //store mip info in any case
292 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set chamber, index of cluster + cluster size
293 pTrk->SetHMPIDtrk(xPc0,yPc0,theta,phi);
296 // Dump debug info if specified
298 if (AliHMPIDReconstructor::StreamLevel()>0) {
299 AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
300 AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);
301 if(!trackTPC->Rotate(alpha)) continue;
302 if(!trackCurrent->Rotate(alpha)) continue;
303 Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
304 Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
305 AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
306 if(!trackTPCNB->Rotate(alpha)) continue;
307 Bool_t statusTPCNB=kTRUE;
308 Double_t bfield[3]={0,0,0};
309 for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){
311 trackTPCNB->GetXYZ(xyz);
312 GetBxByBz(xyz,bfield);
313 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield);
315 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield);
317 Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
318 Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());
320 AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
321 Double_t pos[2]={0,gz};
322 Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
323 Double_t chi2C = trackTPCConstrained->GetPredictedChi2(pos,cov);
324 trackTPCConstrained->Update(pos,cov);
325 (*fDebugStreamer)<<"track"<<
326 "rH="<<radiusH<< // radius of cluster
327 "angle="<<tanAlpha<< // tan of the local inlination angle
328 "dC="<<deltaC<< // delta of the curvature
329 "trackTPC.="<<trackTPC<< // TPC outer param extrapolated to the HMPID
330 "trackTPCNB.="<<trackTPCNB<< // TPC track prpagated with material budget correction
332 "trackTPCC.="<<trackTPCConstrained<< // TPC outer param extrapolated to the HMPID constrained
333 "trackCurrent.="<<trackCurrent<< // current track extrapolated to the HMPID
334 "sTPC="<<statusTPC<< // status for the current TPC track
335 "sCurrent="<<statusCurrent<< // status for the current global track
336 "cl.="<<bestHmpCluster<< // HMPID cluster
338 "t.="<<pTrk<< // curent esd track
339 "ft.="<<ftrack<< // friend track
340 "hmpTrk.="<<hmpTrk<< // hmpid tracks as used in the following code
341 "hmpTrkC.="<<hmpTrkConstrained<< // constrained hmpid tracks as used in the following code
342 "gx="<<gx<< // global cluster position X
351 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
352 delete hmpTrk;hmpTrk=0x0;
353 delete hmpTrkConstrained;hmpTrkConstrained=0x0;
357 if(AliHMPIDReconstructor::GetRecoParam()) //retrieve distance cut
359 if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE) //distance cut is fixed number
361 distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
365 for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar); //prevision: distance cut is function of momentum
366 distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
369 else {distCut=pParam->DistCut();}
373 dmin = TMath::Sqrt((xPc0-bestHmpCluster->X())*(xPc0-bestHmpCluster->X())+(yPc0-bestHmpCluster->Y())*(yPc0-bestHmpCluster->Y()));
378 //isOkDcut = kTRUE; // switch OFF cut
381 pTrk->SetHMPIDsignal(pParam->kMipDistCut); //closest cluster with enough charge is still too far from intersection
384 if(isOkQcut*isOkDcut) isMatched = kTRUE; // MIP-Track matched !!
386 if(!isMatched) {delete hmpTrk;hmpTrk=0x0;delete hmpTrkConstrained;hmpTrkConstrained=0x0;continue;} // If matched continue...
389 if(!isOk) {delete hmpTrk;hmpTrk=0x0; delete hmpTrkConstrained;hmpTrkConstrained=0x0; continue;}
390 pTrk->SetOuterHmpParam(hmpTrkConstrained,AliESDtrack::kHMPIDout);
392 FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
394 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
398 if(pNmean->GetEntries()==21) { //for backward compatibility
399 nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts); //C6F14 Nmean for this chamber
404 Double_t tLow = ((TF1*)pNmean->At(6*ipCh+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber
405 Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber
406 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
407 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
409 if(nmean < 0){ //track didn' t pass through the radiator
410 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
411 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set index of cluster
412 delete hmpTrk;hmpTrk=0x0;
413 delete hmpTrkConstrained;hmpTrkConstrained=0x0;
417 } else nmean = pParam->MeanIdxRad();
419 // 7. Calculate the Cherenkov angle
421 recon.SetImpPC(xPc0,yPc0); //store track impact to PC
422 recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa); //search for Cerenkov angle of this track
424 Double_t thetaCkov = pTrk->GetHMPIDsignal();
426 if (AliHMPIDReconstructor::StreamLevel()>0) {
427 AliExternalTrackParam * trackTPC=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
428 AliExternalTrackParam * trackCurrent=new AliExternalTrackParam(*pTrk);
429 if(!trackTPC->Rotate(alpha)) continue;
430 if(!trackCurrent->Rotate(alpha)) continue;
431 Bool_t statusTPC= AliTracker::PropagateTrackToBxByBz(trackTPC,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
432 Bool_t statusCurrent=AliTracker::PropagateTrackToBxByBz(trackCurrent,radiusH,pTrk->GetMass(),1,kFALSE,kMaxSnp,-1);
433 Double_t tanAlpha=TMath::Tan(TMath::ASin(trackTPC->GetSnp()));
434 Double_t deltaC= trackTPC->GetC(AliTrackerBase::GetBz())-ftrack->GetTPCOut()->GetC(AliTrackerBase::GetBz());
436 AliExternalTrackParam * trackTPCNB=new AliExternalTrackParam(*(ftrack->GetTPCOut()));
437 if(!trackTPCNB->Rotate(alpha)) continue;
438 Bool_t statusTPCNB=kTRUE;
439 Double_t bfield[3]={0,0,0};
440 for (Double_t radius=trackTPCNB->GetX(); radius<radiusH; radius+=1){
442 trackTPCNB->GetXYZ(xyz);
443 GetBxByBz(xyz,bfield);
444 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radius,bfield);
446 statusTPCNB&=trackTPCNB->PropagateToBxByBz(radiusH,bfield);
448 AliExternalTrackParam * trackTPCConstrained= new AliExternalTrackParam(*trackTPC);
449 Double_t pos[2]={0,gz};
450 Double_t cov[3]={0.1*0.1, 0, 0.1*0.1};
451 Double_t chi2C = trackTPCConstrained->GetPredictedChi2(pos,cov);
452 trackTPCConstrained->Update(pos,cov);
453 (*fDebugStreamer)<<"track2"<<
454 "rH="<<radiusH<< // radius of cluster
455 "angle="<<tanAlpha<< // tan of the local inlination angle
456 "dC="<<deltaC<< // delta of the curvature
457 "trackTPC.="<<trackTPC<< // TPC outer param extrapolated to the HMPID
458 "trackTPCNB.="<<trackTPCNB<< // TPC outer param extrapolated to the HMPID
460 "trackTPCC.="<<trackTPCConstrained<< // TPC outer param extrapolated to the HMPID constrained
461 "trackCurrent.="<<trackCurrent<< // current track extrapolated to the HMPID
462 "sTPC="<<statusTPC<< // status for the current TPC track
463 "sCurrent="<<statusCurrent<< // status for the current global track
464 "cl.="<<bestHmpCluster<< // HMPID cluster
466 "t.="<<pTrk<< // curent esd track
467 "ft.="<<ftrack<< // friend track
468 "hmpTrk.="<<hmpTrk<< // hmpid tracks as used in the following code
469 "hmpTrkC.="<<hmpTrkConstrained<< // constrained hmpid tracks as used in the following code
470 "gx="<<gx<< // global cluster position X
473 "thetaCkov="<<thetaCkov<< // Cherenkov angle
477 if(pTrk->GetHMPIDsignal()<0) {
478 delete hmpTrk;hmpTrk=0x0;
479 delete hmpTrkConstrained; hmpTrkConstrained=0x0;
484 pID.FindPid(pTrk,nmean,5,prob);
485 pTrk->SetHMPIDpid(prob);
486 delete hmpTrk; hmpTrk=0x0;
487 delete hmpTrkConstrained; hmpTrkConstrained=0x0;
490 return 0; // error code: 0=no error;
492 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
493 Int_t AliHMPIDTracker::ReconFromKin(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
495 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
496 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
497 // Returns: error code, 0 if no errors
499 AliHMPIDRecon recon; //instance of reconstruction class, nothing important in ctor
500 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
501 Float_t xPc,yPc,xRa,yRa,theta,phi;
502 Double_t cluLORS[2]={0};
504 Double_t qthre = 0; Double_t nmean=0; Int_t hvsec=0;
505 //Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
507 Bool_t tsRight = kTRUE;
509 UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); //
510 UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); //
511 UInt_t ts = pEsd->GetTimeStamp();
513 if(ts<tsmin || ts>tsmax) {
514 AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
518 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
520 Double_t dmin=999999,distCut=1,distParams[5]={1};
522 Bool_t isOkDcut=kFALSE;
523 Bool_t isOkQcut=kFALSE;
524 Bool_t isMatched=kFALSE;
526 AliHMPIDCluster *bestHmpCluster=0x0; //the best matching cluster
527 AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //get reconstructed track
529 AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk); //create a hmpid track to be used for propagation and matching
530 // bz=AliTracker::GetBz();
531 //initial flags for HMPID ESD infos
532 pTrk->SetHMPIDtrk(0,0,0,0); //no intersection found
533 pTrk->SetHMPIDmip(0,0,0,0); //store mip info in any case
534 pTrk->SetHMPIDcluIdx(99,99999); //chamber not found, mip not yet considered
535 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed); //ring reconstruction not yet performed
537 Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi); //find the intersected chamber for this track
538 if(ipCh<0) {delete hmpTrk;hmpTrk=0x0;continue;} //no intersection at all, go after next track
540 pTrk->SetHMPIDtrk(xPc,yPc,theta,phi); //store initial infos
541 pTrk->SetHMPIDcluIdx(ipCh,9999); //set chamber, index of cluster + cluster size
543 // track intersects the chamber ipCh: find the MIP
545 TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters
546 nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber
547 // nClusCh[ipCh] = nMipClusTot;
549 if(nMipClusTot==0) {delete hmpTrk;hmpTrk=0x0;continue;}
551 Int_t index=-1; //index of the "best" matching cluster
553 for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop
555 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster
558 if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {
559 qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(ts); //
560 } else { // in the past just 1 qthre
561 hvsec = pParam->InHVSector(pClu->Y()); // per chamber
562 if(hvsec>=0) qthre=((TF1*)pQthre->At(6*ipCh+hvsec))->Eval(ts); //
564 } else qthre = pParam->QCut();
567 if(pClu->Q()<qthre) continue; //charge compartible with MIP clusters
570 cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y(); //get the LORS coordinates of the cluster
571 Double_t dist = TMath::Sqrt((xPc-cluLORS[0])*(xPc-cluLORS[0])+(yPc-cluLORS[1])*(yPc-cluLORS[1]));
581 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
582 delete hmpTrk;hmpTrk=0x0; continue;
585 Double_t radius = (pParam->Lors2Mars(ipCh,pParam->SizeAllX()/2,pParam->SizeAllY()/2)).Mag();
587 if(!AliTracker::PropagateTrackToBxByBz(hmpTrk,radius,pTrk->GetMass(),1,kFALSE)) {delete hmpTrk;hmpTrk=0x0;continue;}
589 if(!hmpTrk->PropagateTo(bestHmpCluster)) {delete hmpTrk;hmpTrk=0x0;continue;}
591 Int_t cluSiz = bestHmpCluster->Size();
592 pTrk->SetHMPIDmip(bestHmpCluster->X(),bestHmpCluster->Y(),(Int_t)bestHmpCluster->Q(),0); //store mip info in any case
593 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set chamber, index of cluster + cluster size
595 if(AliHMPIDReconstructor::GetRecoParam()) //retrieve distance cut
597 if(AliHMPIDReconstructor::GetRecoParam()->IsFixedDistCut()==kTRUE) //distance cut is fixed number
599 distCut=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDist();
603 for(Int_t ipar=0;ipar<5;ipar++) distParams[ipar]=AliHMPIDReconstructor::GetRecoParam()->GetHmpTrackMatchingDistParam(ipar); //prevision: distance cut is function of momentum
604 distCut=distParams[0]+distParams[1]*TMath::Power(distParams[2]*pTrk->GetP(),distParams[3]); //prevision: change functional form to be more precise
607 else {distCut=pParam->DistCut();}
614 pTrk->SetHMPIDsignal(pParam->kMipDistCut); //closest cluster with enough charge is still too far from intersection
617 if(isOkQcut*isOkDcut) isMatched = kTRUE; // MIP-Track matched !!
619 if(!isMatched) {delete hmpTrk;hmpTrk=0x0;continue;} // If matched continue...
621 Bool_t isOk = hmpTrk->Update(bestHmpCluster,0.1,0);
622 if(!isOk) {delete hmpTrk;hmpTrk=0x0;continue;}
623 pTrk->SetOuterHmpParam(hmpTrk,AliESDtrack::kHMPIDout);
625 FillResiduals(hmpTrk,bestHmpCluster,kFALSE);
627 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
631 if(pNmean->GetEntries()==21) { //for backward compatibility
632 nmean=((TF1*)pNmean->At(3*ipCh))->Eval(ts); //C6F14 Nmean for this chamber
637 Double_t tLow = ((TF1*)pNmean->At(6*ipCh+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber
638 Double_t tHigh = ((TF1*)pNmean->At(6*ipCh+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber
639 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
640 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
642 if(nmean < 0){ //track didn' t pass through the radiator
643 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
644 pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz); //set index of cluster
645 delete hmpTrk;hmpTrk=0x0;
649 } else nmean = pParam->MeanIdxRad();
652 recon.SetImpPC(xPc,yPc); //store track impact to PC
653 recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(ipCh),index,nmean,xRa,yRa); //search for Cerenkov angle of this track
655 if(pTrk->GetHMPIDsignal()<0) {delete hmpTrk;hmpTrk=0x0;continue;}
659 pID.FindPid(pTrk,nmean,5,prob);
660 // Printf("Tracker RefIndex = %f, chmaber = %i, prob1 = %f, prob2 = %f, prob3 = %f, prob4 = %f, prob5 = %f", nmean,ipCh,prob[0],prob[1],prob[2],prob[3],prob[4]);
661 pTrk->SetHMPIDpid(prob);
662 delete hmpTrk;hmpTrk=0x0;
665 return 0; // error code: 0=no error;
667 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
668 Int_t AliHMPIDTracker::ReconHiddenTrk(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
670 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
671 // 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)
672 // Returns: error code, 0 if no errors
674 AliHMPIDReconHTA reconHTA; //instance of reconstruction class, nothing important in ctor
676 AliHMPIDParam *pParam = AliHMPIDParam::Instance(); //Instance of AliHMPIDParam
678 Bool_t tsRight = kTRUE;
680 UInt_t tsmin = (UInt_t)((TF1*)pQthre->At(0))->GetXmin(); //
681 UInt_t tsmax = (UInt_t)((TF1*)pQthre->At(0))->GetXmax(); //
682 UInt_t ts = pEsd->GetTimeStamp();
684 if(ts<tsmin || ts>tsmax) {
685 AliWarning(Form(" in HMPID time stamp out of range!!! Please check!!! ts = %i",ts));
689 for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){ //loop on the ESD tracks in the event
691 AliESDtrack *pTrk = pEsd->GetTrack(iTrk); //here it is simulated or just empty track
693 ipCh = pTrk->GetHMPIDcluIdx();ipCh/=1000000;
696 TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh); //get the list of clusters
697 Int_t nMipClusTot = pMipCluLst->GetEntries(); //total number of clusters in the given chamber
706 for (Int_t iClu=0; iClu<nMipClusTot;iClu++) { //clusters loop
708 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu); //get the cluster
709 Double_t qClus = pClu->Q();
716 cluMipSiz = pClu->Size();
720 if(chMip<0) return 1;
727 if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) { // just for backward compatibility
728 qthre=((TF1*)pQthre->At(chMip))->Eval(ts); //
729 } else { // in the past just 1 qthre
730 hvsec = pParam->InHVSector(yMip); // per chamber
731 if(hvsec>=0) qthre=((TF1*)pQthre->At(6*chMip+hvsec))->Eval(ts); //
733 } else qthre = pParam->QCut();
736 pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case
737 pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);
738 pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
739 return 1; //charge compatible with MIP clusters
742 pTrk->SetHMPIDmip(xMip,yMip,(Int_t)qMip,0); //store mip info in any case
743 pTrk->SetHMPIDcluIdx(chMip,indMip+1000*cluMipSiz);
745 Double_t yRa = yMip; //just an approx...
748 Int_t iRad = pParam->Radiator(yRa); //evaluate the radiator involved
752 if(pNmean->GetEntries()==21) { //for backward compatibility
753 nmean=((TF1*)pNmean->At(3*chMip))->Eval(ts); //C6F14 Nmean for this chamber
758 Double_t tLow = ((TF1*)pNmean->At(6*chMip+2*iRad ))->Eval(ts); //C6F14 low temp for this chamber
759 Double_t tHigh = ((TF1*)pNmean->At(6*chMip+2*iRad+1))->Eval(ts); //C6F14 high temp for this chamber
760 Double_t tExp = pParam->FindTemp(tLow,tHigh,yRa); //estimated temp for that chamber at that y
761 nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp); //mean ref idx @ a given temp
763 if(nmean < 0){ //track didn' t pass through the radiator
764 pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad); //set the appropriate flag
768 } else nmean = pParam->MeanIdxRad();
770 if(!reconHTA.CkovHiddenTrk(pTrk,(TClonesArray *)pClus->At(ipCh),indMip,nmean)) { //search for track parameters and Cerenkov angle of this track
773 pID.FindPid(pTrk,nmean,5,prob);
774 pTrk->SetHMPIDpid(prob);
776 // 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);
779 return 0; // error code: 0=no error;
782 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
783 void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
785 // Publishes all pointers to clusters known to the tracker into the
786 // passed object array.
787 // The ownership is not transfered - the caller is not expected to delete
790 for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){
791 TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
792 for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
793 AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);
794 array->AddLast(pClu);
795 }//cluster loop in iCh
801 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++