]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG/EMCAL/AliAnalysisTaskEmcalSample.cxx
Merge branch 'feature-movesplit'
[u/mrichter/AliRoot.git] / PWG / EMCAL / AliAnalysisTaskEmcalSample.cxx
1 // $Id$
2 //
3 // Emcal sample analysis task.
4 //
5 // Author: S.Aiola, M. Verweij
6
7 #include <TClonesArray.h>
8 #include <TH1F.h>
9 #include <TH2F.h>
10 #include <TH3F.h>
11 #include <TList.h>
12 #include <TLorentzVector.h>
13
14 #include "AliVCluster.h"
15 #include "AliAODCaloCluster.h"
16 #include "AliESDCaloCluster.h"
17 #include "AliVTrack.h"
18 #include "AliLog.h"
19 #include "AliParticleContainer.h"
20 #include "AliClusterContainer.h"
21 #include "AliPicoTrack.h"
22
23 #include "AliAnalysisTaskEmcalSample.h"
24
25 ClassImp(AliAnalysisTaskEmcalSample)
26
27 //________________________________________________________________________
28 AliAnalysisTaskEmcalSample::AliAnalysisTaskEmcalSample() : 
29   AliAnalysisTaskEmcal("AliAnalysisTaskEmcalSample", kTRUE),
30   fHistTracksPt(0),
31   fHistClustersPt(0),
32   fHistPtDEtaDPhiTrackClus(0),
33   fHistPtDEtaDPhiClusTrack(0),
34   fTracksCont(0),
35   fCaloClustersCont(0)
36 {
37   // Default constructor.
38
39   fHistTracksPt       = new TH1*[fNcentBins];
40   fHistClustersPt     = new TH1*[fNcentBins];
41
42   for (Int_t i = 0; i < fNcentBins; i++) {
43     fHistTracksPt[i] = 0;
44     fHistClustersPt[i] = 0;
45   }
46
47   SetMakeGeneralHistograms(kTRUE);
48 }
49
50 //________________________________________________________________________
51 AliAnalysisTaskEmcalSample::AliAnalysisTaskEmcalSample(const char *name) : 
52   AliAnalysisTaskEmcal(name, kTRUE),
53   fHistTracksPt(0),
54   fHistClustersPt(0),
55   fHistPtDEtaDPhiTrackClus(0),
56   fHistPtDEtaDPhiClusTrack(0),
57   fTracksCont(0),
58   fCaloClustersCont(0)
59 {
60   // Standard constructor.
61
62   fHistTracksPt       = new TH1*[fNcentBins];
63   fHistClustersPt     = new TH1*[fNcentBins];
64
65   for (Int_t i = 0; i < fNcentBins; i++) {
66     fHistTracksPt[i] = 0;
67     fHistClustersPt[i] = 0;
68   }
69
70   SetMakeGeneralHistograms(kTRUE);
71 }
72
73 //________________________________________________________________________
74 AliAnalysisTaskEmcalSample::~AliAnalysisTaskEmcalSample()
75 {
76   // Destructor.
77 }
78
79 //________________________________________________________________________
80 void AliAnalysisTaskEmcalSample::UserCreateOutputObjects()
81 {
82   // Create user output.
83
84   AliAnalysisTaskEmcal::UserCreateOutputObjects();
85
86   fTracksCont       = GetParticleContainer(0);
87   fCaloClustersCont = GetClusterContainer(0);
88   fTracksCont->SetClassName("AliVTrack");
89   fCaloClustersCont->SetClassName("AliVCluster");
90
91   TString histname;
92
93   for (Int_t i = 0; i < fNcentBins; i++) {
94     if (fParticleCollArray.GetEntriesFast()>0) {
95       histname = "fHistTracksPt_";
96       histname += i;
97       fHistTracksPt[i] = new TH1F(histname.Data(), histname.Data(), fNbins / 2, fMinBinPt, fMaxBinPt / 2);
98       fHistTracksPt[i]->GetXaxis()->SetTitle("p_{T,track} (GeV/c)");
99       fHistTracksPt[i]->GetYaxis()->SetTitle("counts");
100       fOutput->Add(fHistTracksPt[i]);
101     }
102
103     if (fClusterCollArray.GetEntriesFast()>0) {
104       histname = "fHistClustersPt_";
105       histname += i;
106       fHistClustersPt[i] = new TH1F(histname.Data(), histname.Data(), fNbins / 2, fMinBinPt, fMaxBinPt / 2);
107       fHistClustersPt[i]->GetXaxis()->SetTitle("p_{T,clus} (GeV/c)");
108       fHistClustersPt[i]->GetYaxis()->SetTitle("counts");
109       fOutput->Add(fHistClustersPt[i]);
110     }
111   }
112
113   histname = "fHistPtDEtaDPhiTrackClus";
114   fHistPtDEtaDPhiTrackClus = new TH3F(histname.Data(),Form("%s;#it{p}_{T}^{track};#Delta#eta;#Delta#varphi",histname.Data()),100,0.,100.,100,-0.1,0.1,100,-0.1,0.1);
115   fOutput->Add(fHistPtDEtaDPhiTrackClus);
116
117   histname = "fHistPtDEtaDPhiClusTrack";
118   fHistPtDEtaDPhiClusTrack = new TH3F(histname.Data(),Form("%s;#it{p}_{T}^{clus};#Delta#eta;#Delta#varphi",histname.Data()),100,0.,100.,100,-0.1,0.1,100,-0.1,0.1);
119   fOutput->Add(fHistPtDEtaDPhiClusTrack);
120
121   PostData(1, fOutput); // Post data for ALL output slots > 0 here.
122 }
123
124 //________________________________________________________________________
125 Bool_t AliAnalysisTaskEmcalSample::FillHistograms()
126 {
127   // Fill histograms.
128
129   if (fTracksCont) {
130     AliVTrack *track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle(0)); 
131     while(track) {
132       fHistTracksPt[fCentBin]->Fill(track->Pt()); 
133       track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle());
134     }
135   }
136   
137   if (fCaloClustersCont) {
138     AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
139     while(cluster) {
140       TLorentzVector nPart;
141       cluster->GetMomentum(nPart, fVertex);
142       fHistClustersPt[fCentBin]->Fill(nPart.Pt());
143       cluster = fCaloClustersCont->GetNextAcceptCluster();
144     }
145   }
146
147   CheckClusTrackMatching();
148
149   return kTRUE;
150 }
151
152 //________________________________________________________________________
153 void AliAnalysisTaskEmcalSample::CheckClusTrackMatching()
154 {
155   
156   if(!fTracksCont || !fCaloClustersCont)
157     return;
158
159   Double_t deta = 999;
160   Double_t dphi = 999;
161
162   //Get closest cluster to track
163   AliVTrack *track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle(0)); 
164   while(track) {
165     //Get matched cluster
166     Int_t emc1 = track->GetEMCALcluster();
167     if(fCaloClustersCont && emc1>=0) {
168       AliVCluster *clusMatch = fCaloClustersCont->GetCluster(emc1);
169       if(clusMatch) {
170         AliPicoTrack::GetEtaPhiDiff(track, clusMatch, dphi, deta);
171         fHistPtDEtaDPhiTrackClus->Fill(track->Pt(),deta,dphi);
172       }
173     }
174     track = static_cast<AliVTrack*>(fTracksCont->GetNextAcceptParticle());
175   }
176   
177   //Get closest track to cluster
178   AliVCluster *cluster = fCaloClustersCont->GetNextAcceptCluster(0); 
179   while(cluster) {
180     TLorentzVector nPart;
181     cluster->GetMomentum(nPart, fVertex);
182     fHistClustersPt[fCentBin]->Fill(nPart.Pt());
183     
184     //Get matched track
185     AliVTrack *mt = NULL;      
186     AliAODCaloCluster *acl = dynamic_cast<AliAODCaloCluster*>(cluster);
187     if(acl) {
188       if(acl->GetNTracksMatched()>1)
189         mt = static_cast<AliVTrack*>(acl->GetTrackMatched(0));
190     }
191     else {
192       AliESDCaloCluster *ecl = dynamic_cast<AliESDCaloCluster*>(cluster);
193       Int_t im = ecl->GetTrackMatchedIndex();
194       if(fTracksCont && im>=0) {
195         mt = static_cast<AliVTrack*>(fTracksCont->GetParticle(im));
196       }
197     }
198     if(mt) {
199       AliPicoTrack::GetEtaPhiDiff(mt, cluster, dphi, deta);
200       fHistPtDEtaDPhiClusTrack->Fill(nPart.Pt(),deta,dphi);
201       
202       /* //debugging
203          if(mt->IsEMCAL()) {
204          Int_t emc1 = mt->GetEMCALcluster();
205          Printf("current id: %d  emc1: %d",fCaloClustersCont->GetCurrentID(),emc1);
206          AliVCluster *clm = fCaloClustersCont->GetCluster(emc1);
207          AliPicoTrack::GetEtaPhiDiff(mt, clm, dphi, deta);
208          Printf("deta: %f dphi: %f",deta,dphi);
209          }
210       */
211     }
212     cluster = fCaloClustersCont->GetNextAcceptCluster();
213   }
214 }
215
216 //________________________________________________________________________
217 void AliAnalysisTaskEmcalSample::ExecOnce() {
218
219   AliAnalysisTaskEmcal::ExecOnce();
220
221   if (fTracksCont && fTracksCont->GetArray() == 0) fTracksCont = 0;
222   if (fCaloClustersCont && fCaloClustersCont->GetArray() == 0) fCaloClustersCont = 0;
223
224 }
225
226 //________________________________________________________________________
227 Bool_t AliAnalysisTaskEmcalSample::Run()
228 {
229   // Run analysis code here, if needed. It will be executed before FillHistograms().
230
231   return kTRUE;  // If return kFALSE FillHistogram() will NOT be executed.
232 }
233
234 //________________________________________________________________________
235 void AliAnalysisTaskEmcalSample::Terminate(Option_t *) 
236 {
237   // Called once at the end of the analysis.
238 }