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