]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDTracker.cxx
Mip-track matching improved. Now Last point of HMPID stored and track re-fitted.
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDTracker.cxx
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 "AliHMPIDRecon.h"       //Recon()
6 #include "AliHMPIDReconHTA.h"    //ReconHTA()
7 #include <AliESDEvent.h>         //PropagateBack(),Recon()  
8 #include <AliESDtrack.h>         //Intersect()  
9 #include <AliRun.h>              //GetTrackPoint(),PropagateBack()  
10 #include <AliTrackPointArray.h>  //GetTrackPoint()
11 #include <AliAlignObj.h>         //GetTrackPoint()
12 #include <AliCDBManager.h>       //PropageteBack()
13 #include <AliCDBEntry.h>         //PropageteBack()
14 //.
15 // HMPID base class fo tracking
16 //.
17 //.
18 //.
19 ClassImp(AliHMPIDTracker)
20 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
21 AliHMPIDTracker::AliHMPIDTracker():
22   AliTracker(),
23   fClu(new TObjArray(AliHMPIDParam::kMaxCh+1))  
24 {
25 // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster  
26 // 
27 //  
28   fClu->SetOwner(kTRUE);
29   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
30 }//ctor
31 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
32 Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
33 {
34 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
35 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
36 // Arguments: idx- cluster index which is stored by HMPID in AliESDtrack
37 //            point- reference to the object where to store the point     
38 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
39   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
40   Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
41   point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham-1));//layer and chamber number
42   TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
43   AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
44   Double_t mars[3];
45   AliHMPIDParam::Instance()->Lors2Mars(iCham,pClu->X(),pClu->Y(),mars);
46   point.SetXYZ(mars[0],mars[1],mars[2]);
47   return kTRUE;
48 }//GetTrackPoint()
49 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
50 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
51 {
52 // Static method to find intersection in between given track and HMPID chambers
53 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
54 //   Returns: intersected chamber ID or -1
55   AliHMPIDParam *pParam=AliHMPIDParam::Instance();
56   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){                              //chambers loop
57     Double_t p1[3],n1[3]; pParam->Norm(i,n1); pParam->Point(i,p1,AliHMPIDParam::kRad);          //point & norm  for middle of radiator plane
58     Double_t p2[3],n2[3]; pParam->Norm(i,n2); pParam->Point(i,p2,AliHMPIDParam::kPc);           //point & norm  for entrance to PC plane
59     if(pTrk->Intersect(p1,n1,-GetBz())==kFALSE) continue;                                       //try to intersect track with the middle of radiator
60     if(pTrk->Intersect(p2,n2,-GetBz())==kFALSE) continue;                                       //try to intersect track with PC
61     pParam->Mars2LorsVec(i,n1,theta,phi);                                                       //track angles at RAD
62     pParam->Mars2Lors   (i,p1,xRa,yRa);                                                         //TRKxRAD position
63     pParam->Mars2Lors   (i,p2,xPc,yPc);                                                         //TRKxPC position
64     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return i;                     //return intersected chamber  
65   }                                                                                             //chambers loop
66   return -1;                                                                                    //no intersection with HMPID chambers
67 }//IntTrkCha()
68 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
69 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)
70 {
71 // Static method to find intersection in between given track and HMPID chambers
72 // Arguments: pTrk- HMPID track; xPc,yPc- track intersection with PC in LORS [cm]
73 //   Returns: intersected chamber ID or -1
74     AliHMPIDParam *pParam=AliHMPIDParam::Instance();
75     Double_t p1[3],n1[3]; 
76     pParam->Norm(ch,n1); 
77     pParam->Point(ch,p1,AliHMPIDParam::kRad);                                                    //point & norm  for middle of radiator plane
78     Double_t p2[3],n2[3]; 
79     pParam->Norm(ch,n2); 
80     pParam->Point(ch,p2,AliHMPIDParam::kPc);                                                     //point & norm  for entrance to PC plane
81     if(pTrk->Intersect(pTrk,p1,n1)==kFALSE) return -1;                                           //try to intersect track with the middle of radiator
82     if(pTrk->Intersect(pTrk,p2,n2)==kFALSE) return -1;   
83     pParam->Mars2LorsVec(ch,n1,theta,phi);                                                       //track angles at RAD
84     pParam->Mars2Lors   (ch,p1,xRa,yRa);                                                         //TRKxRAD position
85     pParam->Mars2Lors   (ch,p2,xPc,yPc);                                                         //TRKxPC position
86     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return ch;                     //return intersected chamber  
87   return -1;                                                                                     //no intersection with HMPID chambers
88 }//IntTrkCha()
89 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
90 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
91 {
92 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
93 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
94 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
95   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
96   pCluTree->GetEntry(0);
97   return 0;  
98 }//LoadClusters()
99 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
100 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
101 {
102 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
103 // Agruments: pEsd - pointer to ESD
104 //   Returns: error code    
105   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
106   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
107   if(!pNmeanEnt) AliFatal("No Nmean C6F14 ");
108   if(!pQthreEnt) AliFatal("No Qthre");
109     
110   return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());  
111 }//PropagateBack()
112 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
114 {
115 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
116 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
117 //   Returns: error code, 0 if no errors   
118   
119   AliHMPIDRecon recon;                                                                           //instance of reconstruction class, nothing important in ctor
120   Float_t xPc,yPc,xRa,yRa,theta,phi;
121   Double_t cluLORS[2]={0},cluMARS[3]={0},trkMARS[3]={0};
122 //  Double_t bestcluMARS[3]={0,0,0};
123   Double_t radClu,radInitTrk;   
124   Int_t nMipClusTot=0;
125   Double_t d3d=0,dmin=999999,bz=0;
126   Bool_t isMatched=kFALSE;
127   Int_t bestCluCh=-1;
128   Int_t cluSiz=0;
129   Double_t qthre = 0;   Double_t nmean=0; Int_t cham=0; Int_t hvsec=0;
130   Int_t index=0;                                                                                //index of the "best" matching cluster
131   Double_t bestChi2=-1;                                                                         //Chi2 of the "best" matching cluster
132   Double_t chi2=0;   
133   Int_t nClusCh[AliHMPIDParam::kMaxCh+1];
134   Bool_t isOkQcut=kFALSE;
135   Bool_t isOkDcut=kFALSE;
136   
137   AliHMPIDParam *pParam = AliHMPIDParam::Instance();                                             //Instance of AliHMPIDParam
138   
139   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                        //loop on the ESD tracks in the event
140     isMatched=kFALSE;bestCluCh=-1;dmin=999999;bestChi2=99999;chi2=99999;cluSiz=0;                 //init. track matching params
141     isOkQcut = kFALSE;
142     AliHMPIDCluster *bestHmpCluster=0x0;                                                          //the best matching cluster
143     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                     //get reconstructed track    
144     AliHMPIDtrack *hmpTrk = new AliHMPIDtrack(*pTrk);                                             //create a hmpid track to be used for propagation and matching 
145     bz=AliTracker::GetBz();  
146     
147     Int_t ipCh=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);
148     if(ipCh<0) {                                                                                 //no intersection at all, go after next track
149       pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
150       pTrk->SetHMPIDcluIdx   (99,99999);                                                         //chamber not found, mip not yet considered
151       pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
152       continue;                                                                         
153     }
154     
155 // track intersects the chamber ipCh: find the MIP          
156     
157     TClonesArray *pMipCluLst=(TClonesArray *)pClus->At(ipCh);                                   //get the list of clusters
158     nMipClusTot = pMipCluLst->GetEntries();                                                     //total number of clusters in the given chamber
159     nClusCh[ipCh] = nMipClusTot;
160     
161     for (Int_t iClu=0; iClu<nMipClusTot;iClu++) {                                               //clusters loop
162       
163       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pMipCluLst->UncheckedAt(iClu);                    //get the cluster
164 // evaluate qThre
165       if(pQthre->GetEntriesFast()==pParam->kMaxCh+1) {                                             // just for backward compatibility
166         qthre=((TF1*)pQthre->At(pClu->Ch()))->Eval(pEsd->GetTimeStamp());                          //
167       } else {                                                                                     // in the past just 1 qthre
168         hvsec = pParam->InHVSector(pClu->Y());                                              //  per chamber
169         if(hvsec>=0)
170           qthre=((TF1*)pQthre->At(6*cham+hvsec))->Eval(pEsd->GetTimeStamp());                      //
171       }                                                                                            //
172 //
173       if(pClu->Q()<qthre) continue;                                                                      //charge compartible with MIP clusters      
174       isOkQcut = kTRUE;
175
176       cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();                                            //get the LORS coordinates of the cluster
177       pParam->Lors2Mars(ipCh,cluLORS[0],cluLORS[1],cluMARS);              //convert cluster coors. from LORS to MARS
178       radClu=TMath::Sqrt(cluMARS[0]*cluMARS[0]+cluMARS[1]*cluMARS[1]);                       //radial distance of candidate cluster in MARS                                          
179       Double_t trkx0[3]; 
180       hmpTrk->GetXYZ(trkx0);                                                                 //get track position in MARS
181       radInitTrk=TMath::Sqrt(trkx0[0]*trkx0[0]+trkx0[1]*trkx0[1]);
182       hmpTrk->PropagateToR(radClu,10);
183       hmpTrk->GetXYZ(trkx0);                                                                   //get track position in MARS
184       hmpTrk->GetXYZAt(radClu,bz,trkMARS);                                                     //get the track coordinates at the rad distance after prop. 
185       d3d=TMath::Sqrt((cluMARS[0]-trkMARS[0])*(cluMARS[0]-trkMARS[0])+(cluMARS[1]-trkMARS[1])*(cluMARS[1]-trkMARS[1])+(cluMARS[2]-trkMARS[2])*(cluMARS[2]-trkMARS[2]));
186       chi2=hmpTrk->GetPredictedChi2(pClu);
187       if(dmin > d3d ) {                                                                         //to be saved for the moment...
188         cluSiz = pClu->Size();
189         dmin=d3d;
190         bestCluCh=ipCh;
191         bestHmpCluster=pClu;
192         index=iClu;
193         bestChi2=chi2;
194         cluLORS[0]=pClu->X(); cluLORS[1]=pClu->Y();
195 //        pParam->Lors2Mars(ipCh,cluLORS[0],cluLORS[1],bestcluMARS); 
196       }//global dmin cut 
197     }//clus loop
198
199     if(!isOkQcut) {
200       pTrk->SetHMPIDcluIdx(ipCh,9999);                                                          
201       pTrk->SetHMPIDsignal(pParam->kMipQdcCut);
202       continue;                                                                     
203     }
204     
205     if(dmin < pParam->DistCut()) {
206       isOkDcut = kTRUE;
207     }
208
209     if(!isOkDcut) {
210       pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                             //set chamber, index of cluster + cluster size
211       pTrk->SetHMPIDsignal(pParam->kMipDistCut);                                                //closest cluster with enough charge is still too far from intersection
212     }
213     
214     if(isOkQcut*isOkDcut) isMatched = kTRUE;                                                    // MIP-Track matched !!    
215     
216     if(!isMatched) continue;                                                                    // If matched continue...
217     
218     Int_t indexAll = 0;
219     for(Int_t iC=0;iC<bestCluCh;iC++) indexAll+=nClusCh[iC]; indexAll+=index;                    //to be verified...
220
221     Bool_t isOk = hmpTrk->Update(bestHmpCluster,bestChi2,indexAll);
222     if(!isOk) continue;
223     pTrk->SetOuterParam((AliExternalTrackParam*)&hmpTrk,AliESDtrack::kHMPIDout);                 
224
225 //    cham=IntTrkCha(bestCluCh,hmpTrk,xPc,yPc,xRa,yRa,theta,phi);
226     cham=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);
227     if(cham<0) {                                                                                  //no intersection at all, go after next track
228       pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
229       pTrk->SetHMPIDcluIdx   (99,99999);                                                         //chamber not found, mip not yet considered
230       pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
231       continue;                                                                         
232     }
233
234     pTrk->SetHMPIDtrk(xRa,yRa,theta,phi);                                                        //store initial infos
235     //evaluate nMean
236     if(pNmean->GetEntries()==21) {                                                              //for backward compatibility
237       nmean=((TF1*)pNmean->At(3*cham))->Eval(pEsd->GetTimeStamp());                             //C6F14 Nmean for this chamber
238     } else {
239       Int_t iRad     = pParam->Radiator(yRa);                                                   //evaluate the radiator involved
240       Double_t tLow  = ((TF1*)pNmean->At(6*cham+2*iRad  ))->Eval(pEsd->GetTimeStamp());         //C6F14 low  temp for this chamber
241       Double_t tHigh = ((TF1*)pNmean->At(6*cham+2*iRad+1))->Eval(pEsd->GetTimeStamp());         //C6F14 high temp for this chamber
242       Double_t tExp  = pParam->FindTemp(tLow,tHigh,yRa);                                        //estimated temp for that chamber at that y
243       nmean = pParam->NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);                  //mean ref idx @ a given temp
244       if(nmean < 0){                                                                            //track didn' t pass through the radiator
245          pTrk->SetHMPIDsignal(AliHMPIDRecon::kNoRad);                                           //set the appropriate flag
246          pTrk->SetHMPIDcluIdx(ipCh,index+1000*cluSiz);                                          //set index of cluster
247          continue;
248       }
249     }
250     //
251     recon.SetImpPC(xPc,yPc);                                                                     //store track impact to PC
252     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(cham),index,nmean);                           //search for Cerenkov angle of this track
253   }//iTrk
254
255   return 0; // error code: 0=no error;
256 }//Recon()
257 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
258 Int_t AliHMPIDTracker::ReconHiddenTrk(Int_t iCh,Int_t iHVsec,AliESDtrack *pTrk,TClonesArray *pCluLst,TObjArray *pNmean,TObjArray *pQthre)
259 {
260 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
261 // 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)
262 //   Returns: error code, 0 if no errors
263   AliHMPIDReconHTA reconHTA;                                                                          //instance of reconstruction class, nothing important in ctor
264   Double_t nmean=((TF1*)pNmean->At(3*iCh))->Eval(0);                                            //C6F14 Nmean for this chamber
265   Double_t qthre = 0;
266   if(pQthre->GetEntriesFast()==AliHMPIDParam::kMaxCh+1)                                         //
267     qthre=((TF1*)pQthre->At(iCh))->Eval(0);                                                     //just for backward compatibi
268   else  qthre=((TF1*)pQthre->At(6*iCh+iHVsec))->Eval(0);                                        //
269   if(pCluLst->GetEntriesFast()<4) return 1;                                                     //min 4 clusters (3 + 1 mip) to find a ring! 
270   if(reconHTA.CkovHiddenTrk(pTrk,pCluLst,nmean,qthre)) return 0;                                   //search for track parameters and Cerenkov angle of this track
271   else return 1;                                                                                // error code: 0=no error,1=fit not performed;
272 }//Recon()
273 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
274 void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
275   
276  // Publishes all pointers to clusters known to the tracker into the
277   // passed object array.
278   // The ownership is not transfered - the caller is not expected to delete
279   // the clusters
280  
281   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){    
282     TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
283     for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
284       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);    
285       array->AddLast(pClu);
286     }//cluster loop in iCh
287     pCluArr->Delete();
288   }//Ch loop
289     
290   return;
291 }
292 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++