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