]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDTracker.cxx
Raw data correctly written
[u/mrichter/AliRoot.git] / HMPID / AliHMPIDTracker.cxx
1 #include "AliHMPIDTracker.h"     //class header
2 #include "AliHMPIDCluster.h"     //GetTrackPoint(),PropagateBack() 
3 #include "AliHMPIDParam.h"       //GetTrackPoint(),PropagateBack()
4 #include "AliHMPIDRecon.h"       //Recon()
5 #include "AliHMPIDReconHTA.h"    //ReconHTA()
6 #include <AliESDEvent.h>         //PropagateBack(),Recon()  
7 #include <AliESDtrack.h>         //Intersect()  
8 #include <AliRun.h>              //GetTrackPoint(),PropagateBack()  
9 #include <AliTrackPointArray.h>  //GetTrackPoint()
10 #include <AliAlignObj.h>         //GetTrackPoint()
11 #include <AliCDBManager.h>       //PropageteBack()
12 #include <AliCDBEntry.h>         //PropageteBack()
13 //.
14 // HMPID base class fo tracking
15 //.
16 //.
17 //.
18 ClassImp(AliHMPIDTracker)
19 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
20 AliHMPIDTracker::AliHMPIDTracker():AliTracker()
21 {
22 // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster  
23 // 
24 //  
25   fClu=new TObjArray(AliHMPIDParam::kMaxCh+1);  fClu->SetOwner(kTRUE);
26   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
27 }//ctor
28 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
29 Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
30 {
31 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
32 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
33 // Arguments: idx- cluster index which is stored by HMPID in AliESDtrack
34 //            point- reference to the object where to store the point     
35 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
36   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
37   Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
38   point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham-1));//layer and chamber number
39   TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
40   AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
41   Double_t mars[3];
42   AliHMPIDParam::Instance()->Lors2Mars(iCham,pClu->X(),pClu->Y(),mars);
43   point.SetXYZ(mars[0],mars[1],mars[2]);
44   return kTRUE;
45 }//GetTrackPoint()
46 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
47 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc)
48 {
49 // Static method to find intersection in between given track and HMPID chambers
50 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
51 //   Returns: intersected chamber ID or -1
52   AliHMPIDParam *pParam=AliHMPIDParam::Instance();
53   Float_t xRa=0,yRa=0,theta=0,phi=0;                                                            //track intersection at PC and angles at RAD, LORS  
54   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){                              //chambers loop
55     Double_t p1[3],n1[3]; pParam->Norm(i,n1); pParam->Point(i,p1,AliHMPIDParam::kRad);          //point & norm  for middle of radiator plane
56     Double_t p2[3],n2[3]; pParam->Norm(i,n2); pParam->Point(i,p2,AliHMPIDParam::kPc);           //point & norm  for entrance to PC plane
57     if(pTrk->Intersect(p1,n1,-GetBz())==kFALSE) continue;                                       //try to intersect track with the middle of radiator
58     if(pTrk->Intersect(p2,n2,-GetBz())==kFALSE) continue;                                       //try to intersect track with PC
59     pParam->Mars2LorsVec(i,n1,theta,phi);                                                       //track angles at RAD
60     pParam->Mars2Lors   (i,p1,xRa,yRa);                                                         //TRKxRAD position
61     pParam->Mars2Lors   (i,p2,xPc,yPc);                                                         //TRKxPC position
62     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kFALSE) continue;                    //not in active area  
63     pTrk->SetHMPIDtrk      (xRa,yRa,theta,phi);                                                 //store track intersection info
64     pTrk->SetHMPIDcluIdx   (i,0);
65     return i;
66   }                                                                                             //chambers loop
67   return -1; //no intersection with HMPID chambers
68 }//IntTrkCha()
69 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
70 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
71 {
72 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
73 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
74 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
75   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
76   pCluTree->GetEntry(0);
77   return 0;  
78 }//LoadClusters()
79 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
80 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
81 {
82 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
83 // Agruments: pEsd - pointer to ESD
84 //   Returns: error code    
85   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 21 TF1
86   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 7 TF1
87   if(!pNmeanEnt) AliFatal("No Nmean C6F14 ");
88   if(!pQthreEnt) AliFatal("No Qthre");
89     
90   return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());  
91 }//PropagateBack()
92 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
93 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
94 {
95 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
96 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
97 //   Returns: error code, 0 if no errors   
98   AliHMPIDRecon recon;                                                                       //instance of reconstruction class, nothing important in ctor
99   Float_t xPc,yPc;
100   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                       //ESD tracks loop
101     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                    //get reconstructed track    
102     Int_t cham=IntTrkCha(pTrk,xPc,yPc);                                                          //get chamber intersected by this track 
103     if(cham<0) continue;                                                                         //no intersection at all, go after next track
104     Double_t nmean=((TF1*)pNmean->At(3*cham))->Eval(pEsd->GetTimeStamp());                       //C6F14 Nmean for this chamber
105     Double_t qthre=((TF1*)pQthre->At(cham))  ->Eval(pEsd->GetTimeStamp());                       //Qthre for this chamber
106     recon.SetImpPC(xPc,yPc);                                                                     //store track impact to PC
107     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(cham),nmean,qthre);                           //search for Cerenkov angle of this track
108 //    Printf("AliHMPIDTracker::Recon: nmean %f, qthre %f",nmean,qthre);
109   }                                                                                              //ESD tracks loop
110   return 0; // error code: 0=no error;
111 }//Recon()
112 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
113 Int_t AliHMPIDTracker::ReconHiddenTrk(Int_t iCh,AliESDtrack *pTrk,TClonesArray *pCluLst,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), pQthre - pointer to all function Qthre=f(time)
117 //   Returns: error code, 0 if no errors
118   AliHMPIDReconHTA reconHTA;                                                                          //instance of reconstruction class, nothing important in ctor
119   Double_t nmean=((TF1*)pNmean->At(3*iCh))->Eval(0);                                            //C6F14 Nmean for this chamber
120   Double_t qthre=((TF1*)pQthre->At(iCh))  ->Eval(0);                                            //C6F14 Nmean for this chamber
121   if(pCluLst->GetEntriesFast()<4) return 1;                                                     //min 4 clusters (3 + 1 mip) to find a ring! 
122   if(reconHTA.CkovHiddenTrk(pTrk,pCluLst,nmean,qthre)) return 0;                                   //search for track parameters and Cerenkov angle of this track
123   else return 1;                                                                                // error code: 0=no error,1=fit not performed;
124 }//Recon()
125 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++