]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HMPID/AliHMPIDTracker.cxx
Make separate, specialized geometries for RPhi and RhoZ views.
[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():
21   AliTracker(),
22   fClu(new TObjArray(AliHMPIDParam::kMaxCh+1))  
23 {
24 // ctor. Create TObjArray of TClonesArray of AliHMPIDCluster  
25 // 
26 //  
27   fClu->SetOwner(kTRUE);
28   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) fClu->AddAt(new TClonesArray("AliHMPIDCluster"),i);
29 }//ctor
30 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++    
31 Bool_t AliHMPIDTracker::GetTrackPoint(Int_t idx, AliTrackPoint& point) const
32 {
33 // Interface callback methode invoked from AliReconstruction::WriteAlignmentData() to get position of MIP cluster in MARS associated to a current track.
34 // MIP cluster is reffered by index which is stored in AliESDtrack  ???????
35 // Arguments: idx- cluster index which is stored by HMPID in AliESDtrack
36 //            point- reference to the object where to store the point     
37 //   Returns: status of operation  if FALSE then AliReconstruction::WriteAlignmentData() do not store this point to array of points for current track. 
38   if(idx<0) return kFALSE; //no MIP cluster assigned to this track in PropagateBack()
39   Int_t iCham=idx/1000000; Int_t iClu=idx%1000000;
40   point.SetVolumeID(AliGeomManager::LayerToVolUID(AliGeomManager::kHMPID,iCham-1));//layer and chamber number
41   TClonesArray *pArr=(TClonesArray*)(*fClu)[iCham];
42   AliHMPIDCluster *pClu=(AliHMPIDCluster*)pArr->UncheckedAt(iClu);//get pointer to cluster
43   Double_t mars[3];
44   AliHMPIDParam::Instance()->Lors2Mars(iCham,pClu->X(),pClu->Y(),mars);
45   point.SetXYZ(mars[0],mars[1],mars[2]);
46   return kTRUE;
47 }//GetTrackPoint()
48 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
49 Int_t AliHMPIDTracker::IntTrkCha(AliESDtrack *pTrk,Float_t &xPc,Float_t &yPc,Float_t &xRa,Float_t &yRa,Float_t &theta,Float_t &phi)
50 {
51 // Static method to find intersection in between given track and HMPID chambers
52 // Arguments: pTrk- ESD track; xPc,yPc- track intersection with PC in LORS [cm]
53 //   Returns: intersected chamber ID or -1
54   AliHMPIDParam *pParam=AliHMPIDParam::Instance();
55   for(Int_t i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++){                              //chambers loop
56     Double_t p1[3],n1[3]; pParam->Norm(i,n1); pParam->Point(i,p1,AliHMPIDParam::kRad);          //point & norm  for middle of radiator plane
57     Double_t p2[3],n2[3]; pParam->Norm(i,n2); pParam->Point(i,p2,AliHMPIDParam::kPc);           //point & norm  for entrance to PC plane
58     if(pTrk->Intersect(p1,n1,-GetBz())==kFALSE) continue;                                       //try to intersect track with the middle of radiator
59     if(pTrk->Intersect(p2,n2,-GetBz())==kFALSE) continue;                                       //try to intersect track with PC
60     pParam->Mars2LorsVec(i,n1,theta,phi);                                                       //track angles at RAD
61     pParam->Mars2Lors   (i,p1,xRa,yRa);                                                         //TRKxRAD position
62     pParam->Mars2Lors   (i,p2,xPc,yPc);                                                         //TRKxPC position
63     if(AliHMPIDParam::IsInside(xPc,yPc,pParam->DistCut())==kTRUE) return i;                     //return intersected chamber  
64   }                                                                                             //chambers loop
65   return -1;                                                                                    //no intersection with HMPID chambers
66 }//IntTrkCha()
67 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
68 Int_t AliHMPIDTracker::LoadClusters(TTree *pCluTree)
69 {
70 // Interface callback methode invoked from AliReconstruction::RunTracking() to load HMPID clusters before PropagateBack() gets control. Done once per event.
71 // Arguments: pCluTree- pointer to clusters tree got by AliHMPIDLoader::LoadRecPoints("read") then AliHMPIDLoader::TreeR()
72 //   Returns: error code (currently ignored in AliReconstruction::RunTraking())    
73   for(int i=AliHMPIDParam::kMinCh;i<=AliHMPIDParam::kMaxCh;i++) pCluTree->SetBranchAddress(Form("HMPID%d",i),&((*fClu)[i]));
74   pCluTree->GetEntry(0);
75   return 0;  
76 }//LoadClusters()
77 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
78 Int_t AliHMPIDTracker::PropagateBack(AliESDEvent *pEsd)
79 {
80 // Interface pure virtual in AliTracker. Invoked from AliReconstruction::RunTracking() after invocation of AliTracker::LoadClusters() once per event
81 // Agruments: pEsd - pointer to ESD
82 //   Returns: error code    
83   AliCDBEntry *pNmeanEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Nmean"); //contains TObjArray of 42 TF1 + 1 EPhotMean
84   AliCDBEntry *pQthreEnt =AliCDBManager::Instance()->Get("HMPID/Calib/Qthre"); //contains TObjArray of 42 (7ch * 6sec) TF1
85   if(!pNmeanEnt) AliFatal("No Nmean C6F14 ");
86   if(!pQthreEnt) AliFatal("No Qthre");
87     
88   return Recon(pEsd,fClu,(TObjArray*)pNmeanEnt->GetObject(),(TObjArray*)pQthreEnt->GetObject());  
89 }//PropagateBack()
90 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
91 Int_t AliHMPIDTracker::Recon(AliESDEvent *pEsd,TObjArray *pClus,TObjArray *pNmean, TObjArray *pQthre)
92 {
93 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
94 // Arguments: pEsd- pointer ESD; pClu- pointer to clusters for all chambers; pNmean - pointer to all function Nmean=f(time)
95 //   Returns: error code, 0 if no errors   
96   AliHMPIDRecon recon;                                                                       //instance of reconstruction class, nothing important in ctor
97   Float_t xPc,yPc,xRa,yRa,theta,phi;
98   for(Int_t iTrk=0;iTrk<pEsd->GetNumberOfTracks();iTrk++){                                       //ESD tracks loop
99     AliESDtrack *pTrk = pEsd->GetTrack(iTrk);                                                    //get reconstructed track    
100     Int_t cham=IntTrkCha(pTrk,xPc,yPc,xRa,yRa,theta,phi);                                        //get chamber intersected by this track 
101     if(cham<0) {                                                                                 //no intersection at all, go after next track
102       pTrk->SetHMPIDtrk(0,0,0,0);                                                                //no intersection found
103       pTrk->SetHMPIDcluIdx   (99,99999);                                                         //chamber not found, mip not yet considered
104       pTrk->SetHMPIDsignal(AliHMPIDRecon::kNotPerformed);                                        //ring reconstruction not yet performed
105       continue;                                                                         
106     }
107     pTrk->SetHMPIDtrk(xRa,yRa,theta,phi);                                                        //store initial infos
108     Double_t nmean;
109     if(pNmean->GetEntries()==21) {                                                               //for backward compatibility
110        nmean=((TF1*)pNmean->At(3*cham))->Eval(pEsd->GetTimeStamp());                             //C6F14 Nmean for this chamber
111      } else {
112        Int_t iRad     = AliHMPIDParam::Radiator(yRa);                                            //evaluate the radiator involved
113        Double_t tLow  = ((TF1*)pNmean->At(6*cham+2*iRad  ))->Eval(pEsd->GetTimeStamp());         //C6F14 low  temp for this chamber
114        Double_t tHigh = ((TF1*)pNmean->At(6*cham+2*iRad+1))->Eval(pEsd->GetTimeStamp());         //C6F14 high temp for this chamber
115        Double_t tExp  = AliHMPIDParam::FindTemp(tLow,tHigh,yRa);                                 //estimated temp for that chamber at that y
116        nmean = AliHMPIDParam::NIdxRad(AliHMPIDParam::Instance()->GetEPhotMean(),tExp);           //mean ref idx @ a given temp
117      }
118     Double_t qthre = 0; 
119     if(pQthre->GetEntriesFast()==AliHMPIDParam::kMaxCh+1)                                        // just for backward compatibility
120       qthre=((TF1*)pQthre->At(cham))->Eval(pEsd->GetTimeStamp());                                //
121     else {                                                                                       // in the past just 1 qthre
122       Int_t hvsec = AliHMPIDParam::InHVSector(yPc);                                              //  per chamber
123       if (hvsec>=0)
124         qthre=((TF1*)pQthre->At(6*cham+hvsec))->Eval(pEsd->GetTimeStamp());                      //
125     }                                                                                            //
126     recon.SetImpPC(xPc,yPc);                                                                     //store track impact to PC
127     recon.CkovAngle(pTrk,(TClonesArray *)pClus->At(cham),nmean,qthre);                           //search for Cerenkov angle of this track
128 //    Printf("AliHMPIDTracker::Recon: nmean %f, qthre %f",nmean,qthre);
129   }                                                                                              //ESD tracks loop
130   return 0; // error code: 0=no error;
131 }//Recon()
132 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
133 Int_t AliHMPIDTracker::ReconHiddenTrk(Int_t iCh,Int_t iHVsec,AliESDtrack *pTrk,TClonesArray *pCluLst,TObjArray *pNmean,TObjArray *pQthre)
134 {
135 // Static method to reconstruct Theta Ckov for all valid tracks of a given event.
136 // 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)
137 //   Returns: error code, 0 if no errors
138   AliHMPIDReconHTA reconHTA;                                                                          //instance of reconstruction class, nothing important in ctor
139   Double_t nmean=((TF1*)pNmean->At(3*iCh))->Eval(0);                                            //C6F14 Nmean for this chamber
140   Double_t qthre = 0;
141   if(pQthre->GetEntriesFast()==AliHMPIDParam::kMaxCh+1)                                         //
142     qthre=((TF1*)pQthre->At(iCh))->Eval(0);                                                     //just for backward compatibi
143   else  qthre=((TF1*)pQthre->At(6*iCh+iHVsec))->Eval(0);                                        //
144   if(pCluLst->GetEntriesFast()<4) return 1;                                                     //min 4 clusters (3 + 1 mip) to find a ring! 
145   if(reconHTA.CkovHiddenTrk(pTrk,pCluLst,nmean,qthre)) return 0;                                   //search for track parameters and Cerenkov angle of this track
146   else return 1;                                                                                // error code: 0=no error,1=fit not performed;
147 }//Recon()
148 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++
149 void AliHMPIDTracker::FillClusterArray(TObjArray* array) const {
150   
151  // Publishes all pointers to clusters known to the tracker into the
152   // passed object array.
153   // The ownership is not transfered - the caller is not expected to delete
154   // the clusters
155  
156   for(Int_t iCh=AliHMPIDParam::kMinCh;iCh<=AliHMPIDParam::kMaxCh;iCh++){    
157     TClonesArray *pCluArr=(TClonesArray*)(*fClu)[iCh];
158     for (Int_t iClu=0; iClu<pCluArr->GetEntriesFast();iClu++){
159       AliHMPIDCluster *pClu=(AliHMPIDCluster*)pCluArr->UncheckedAt(iClu);    
160       array->AddLast(pClu);
161     }//cluster loop in iCh
162     pCluArr->Delete();
163   }//Ch loop
164     
165   return;
166 }
167 //++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++++