]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/EMCAL/AliEmcalClusTrackMatcherTask.cxx
Fix when not being able to open the input file
[u/mrichter/AliRoot.git] / PWG / EMCAL / AliEmcalClusTrackMatcherTask.cxx
1 // $Id$
2 //
3 // Track/cluster matcher
4 // 
5 // Author: C.Loizides, S.Aiola
6
7 #include "AliEmcalClusTrackMatcherTask.h"
8
9 #include <TClonesArray.h>
10
11 #include "AliEmcalParticle.h"
12 #include "AliLog.h"
13 #include "AliPicoTrack.h"
14 #include "AliVCluster.h"
15 #include "AliVTrack.h"
16
17 ClassImp(AliEmcalClusTrackMatcherTask)
18
19 //________________________________________________________________________
20 AliEmcalClusTrackMatcherTask::AliEmcalClusTrackMatcherTask() : 
21   AliAnalysisTaskEmcal("AliEmcalClusTrackMatcherTask",kFALSE),
22   fMaxDistance(0.06)
23 {
24   // Constructor.
25
26   for(Int_t icent=0; icent<8; ++icent) {
27     for(Int_t ipt=0; ipt<9; ++ipt) {
28       for(Int_t ieta=0; ieta<2; ++ieta) {
29         fHistMatchEta[icent][ipt][ieta] = 0;
30         fHistMatchPhi[icent][ipt][ieta] = 0;
31       }
32     }
33   }
34 }
35
36 //________________________________________________________________________
37 AliEmcalClusTrackMatcherTask::AliEmcalClusTrackMatcherTask(const char *name, Bool_t histo) : 
38   AliAnalysisTaskEmcal(name,histo),
39   fMaxDistance(0.06)
40 {
41   // Standard constructor.
42
43   for(Int_t icent=0; icent<8; ++icent) {
44     for(Int_t ipt=0; ipt<9; ++ipt) {
45       for(Int_t ieta=0; ieta<2; ++ieta) {
46         fHistMatchEta[icent][ipt][ieta] = 0;
47         fHistMatchPhi[icent][ipt][ieta] = 0;
48       }
49     }
50   }
51 }
52
53 //________________________________________________________________________
54 AliEmcalClusTrackMatcherTask::~AliEmcalClusTrackMatcherTask()
55 {
56   // Destructor.
57 }
58
59 //________________________________________________________________________
60 void AliEmcalClusTrackMatcherTask::UserCreateOutputObjects()
61 {
62   // Create my user objects.
63
64   if (!fCreateHisto)
65     return;
66
67   AliAnalysisTaskEmcal::UserCreateOutputObjects();
68
69   for(Int_t icent=0; icent<8; ++icent) {
70     for(Int_t ipt=0; ipt<9; ++ipt) {
71       for(Int_t ieta=0; ieta<2; ++ieta) {
72         TString nameEta(Form("fHistMatchEta_%i_%i_%i",icent,ipt,ieta));
73         fHistMatchEta[icent][ipt][ieta] = new TH1F(nameEta, nameEta, 400, -0.2, 0.2);
74         TString namePhi(Form("fHistMatchPhi_%i_%i_%i",icent,ipt,ieta));
75         fHistMatchPhi[icent][ipt][ieta] = new TH1F(namePhi, namePhi, 400, -0.2, 0.2);
76         fOutput->Add(fHistMatchEta[icent][ipt][ieta]);
77         fOutput->Add(fHistMatchPhi[icent][ipt][ieta]);
78       }
79     }
80   }
81
82   PostData(1, fOutput);
83 }
84
85 //________________________________________________________________________
86 Int_t AliEmcalClusTrackMatcherTask::GetMomBin(Double_t p) const
87 {
88   // Get momenum bin.
89
90   Int_t pbin=-1;
91   if (p<0.5) 
92     pbin=0;
93   else if (p>=0.5 && p<1.0) 
94     pbin=1;
95   else if (p>=1.0 && p<1.5) 
96     pbin=2;
97   else if (p>=1.5 && p<2.) 
98     pbin=3;
99   else if (p>=2. && p<3.) 
100     pbin=4;
101   else if (p>=3. && p<4.) 
102     pbin=5;
103   else if (p>=4. && p<5.) 
104     pbin=6;
105   else if (p>=5. && p<8.) 
106     pbin=7;
107   else if (p>=8.) 
108     pbin=8;
109
110   return pbin;
111 }
112
113 //________________________________________________________________________
114 Bool_t AliEmcalClusTrackMatcherTask::Run() 
115 {
116   // Run the matching for the selected options.
117   
118   const Double_t maxd2 = fMaxDistance*fMaxDistance;
119
120   const Int_t nC = fCaloClusters->GetEntries();
121   const Int_t nT = fTracks->GetEntries();
122
123   // set the links between tracks and clusters
124   for (Int_t c = 0; c < nC; ++c) {
125     AliEmcalParticle *partC = static_cast<AliEmcalParticle*>(fCaloClusters->At(c));
126     if (!partC)
127       continue;
128     if (!AcceptEmcalPart(partC))
129       continue;
130     for (Int_t t = 0; t < nT; ++t) {
131       AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(t));
132       if (!partT)
133         continue;
134       if (!AcceptEmcalPart(partT))
135         continue;
136       AliVCluster *clust = partC->GetCluster();
137       AliVTrack   *track = partT->GetTrack()  ;
138       Double_t deta = 999;
139       Double_t dphi = 999;
140       AliPicoTrack::GetEtaPhiDiff(track, clust, dphi, deta);
141       Double_t d2 = deta * deta + dphi * dphi;
142       if (d2 > maxd2)
143         continue;
144       Double_t d = TMath::Sqrt(d2);
145       partC->AddMatchedObj(t, d);
146       partT->AddMatchedObj(c, d);
147
148       if (fCreateHisto) {
149         Int_t mombin = GetMomBin(track->P());
150         Int_t centbinch = fCentBin;
151         if (track->Charge()<0) 
152           centbinch += 4;
153         Int_t etabin = 0;
154         if(track->Eta() > 0) 
155           etabin = 1;
156             
157         fHistMatchEta[centbinch][mombin][etabin]->Fill(deta);
158         fHistMatchPhi[centbinch][mombin][etabin]->Fill(dphi);
159       }
160     }
161   }
162
163   for (Int_t c = 0; c < nC; ++c) {
164     AliEmcalParticle *partC = static_cast<AliEmcalParticle*>(fCaloClusters->At(c));
165     if (!partC)
166       continue;
167     if (partC->GetNumberOfMatchedObj() <= 0)
168       continue;
169     const UInt_t matchedId = partC->GetMatchedObjId();
170     AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(matchedId));
171     AliVCluster *clust = partC->GetCluster();
172     AliVTrack   *track = partT->GetTrack()  ;
173     Double_t deta = 999;
174     Double_t dphi = 999;
175     AliPicoTrack::GetEtaPhiDiff(track, clust, dphi, deta);
176     clust->SetEmcCpvDistance(matchedId);
177     clust->SetTrackDistance(deta, dphi);
178   }
179
180   for (Int_t t = 0; t < nT; ++t) {
181     AliEmcalParticle *partT = static_cast<AliEmcalParticle*>(fTracks->At(t));
182     if (!partT)
183       continue;
184     if (partT->GetNumberOfMatchedObj() <= 0)
185       continue;
186     AliVTrack *track = partT->GetTrack();
187     track->SetEMCALcluster(partT->GetMatchedObjId());
188   }
189
190   return kTRUE;
191 }