]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWGUD/mcValidation/AliTrackletsTask.cxx
add few ebe dihadron and jet correlation to be compared to normal dihadron correlatio...
[u/mrichter/AliRoot.git] / PWGUD / mcValidation / AliTrackletsTask.cxx
CommitLineData
5e03d366 1/* $Id:$ */
2
3// -----------------------------------------------
4// Task to extract distributions
5// for traclets paramters
6// for a quick comparison
7// between MC and data
8// eta, phi, deltaPhi, deltaEta, vtxX, vtxY, vtxZ, SPD GFO
9// -----------------------------------------------
10
11
12#include "AliTrackletsTask.h"
13
14#include <TCanvas.h>
15#include <TFile.h>
16#include <TChain.h>
17#include <TMath.h>
18#include <TH1D.h>
19#include <TH2D.h>
20#include <TH3D.h>
21#include <TH1I.h>
22
23#include <AliLog.h>
24#include <AliESDEvent.h>
25#include <AliAnalysisManager.h>
26#include <AliESDInputHandler.h>
27#include <AliESDHeader.h>
28#include <AliESDVertex.h>
29#include <AliMultiplicity.h>
30
31ClassImp(AliTrackletsTask)
32
33AliTrackletsTask::AliTrackletsTask() :
34 AliAnalysisTask("AliTrackletsTask", ""),
35 fESD(0),
36 fOutput(0),
37 fNtracks(0x0),
38 fPhi(0x0),
39 fEtaPhi(0x0),
40 fDeltaPhi(0x0),
41 fDeltaTheta(0x0),
42 fVtxX(0x0),
43 fVtxY(0x0),
44 fVtxZ(0x0),
45 fVtx(0x0),
46 fVtxContributors(0x0)
47{
48 //
49 // Constructor. Initialization of pointers
50 //
51
52 // Define input and output slots here
53 DefineInput(0, TChain::Class());
54 DefineOutput(0, TList::Class());
55}
56
57//------------------------------------------------------------------
58
59AliTrackletsTask::~AliTrackletsTask()
60{
61 //
62 // Destructor
63 //
64
65 // histograms are in the output list and deleted when the output
66 // list is deleted by the TSelector dtor
67
68 if (fOutput) {
69 delete fOutput;
70 fOutput = 0;
71 }
72}
73
74//--------------------------------------------------------------------
75
76void AliTrackletsTask::ConnectInputData(Option_t *)
77{
78 //
79 // Connect ESD
80 //
81
82 Printf("AliTrackletsTask::ConnectInputData called");
83
84 AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
85
86 if (!esdH) {
87 Printf("ERROR: Could not get ESDInputHandler");
88 } else {
89 fESD = esdH->GetEvent();
90
91 TString branches("AliESDHeader Tracks AliMultiplicity");
92
93 // Enable only the needed branches
94 esdH->SetActiveBranches(branches);
95 }
96}
97
98//---------------------------------------------------------------------
99
100void AliTrackletsTask::CreateOutputObjects()
101{
102 //
103 // create result objects and add to output list
104 //
105
106 Printf("AliTrackletsTask::CreateOutputObjects");
107
108 fOutput = new TList;
109 fOutput->SetOwner();
110
111 fNtracks = new TH1I("fNtracks", "n. of ESD tracks", 1000, 0, 1000);
112 fOutput->Add(fNtracks);
113
114 fPhi = new TH1D("fPhi", "Phi (rad)", 720,0,2*TMath::Pi());
115 fOutput->Add(fPhi);
116
117 fEtaPhi = new TH2D("fEtaPhi", "Phi vs Eta;#eta;#phi [rad];count", 80, -4, 4, 18*5,0,2*TMath::Pi());
118 fOutput->Add(fEtaPhi);
119
120 fVtxX = new TH1D("fVtxX", "x SPD primary vertex; x [cm]; count", 100, -1, 1);
121 fOutput->Add(fVtxX);
122
123 fVtxY = new TH1D("fVtxY", "y SPD primary vertex; y [cm]; count", 100, -1, 1);
124 fOutput->Add(fVtxY);
125
126 fVtxZ = new TH1D("fVtxZ", "z SPD primary vertex; z [cm]; count", 100, -30, 30);
127 fOutput->Add(fVtxZ);
128
129 fVtx = new TH3D("fVtx", "SPD primary vertex; x [cm]; y [cm]; z [cm]; count", 100, -1, 1, 100,-1,1,100,-30,30);
130 fOutput->Add(fVtx);
131
132 fVtxContributors = new TH3D("fVtxContributors", "SPD primary vertex with N contributors > 0; x [cm]; y [cm]; z [cm]; count", 100, -1, 1, 100,-1,1,100,-30,30);
133 fOutput->Add(fVtxContributors);
134
135 fDeltaPhi = new TH1D("fDeltaPhi", "fDeltaPhi;#Delta #phi;Entries", 500, -1, 1);
136 fOutput->Add(fDeltaPhi);
137
138 fDeltaTheta = new TH1D("fDeltaTheta", "fDeltaTheta;#Delta #theta;Entries", 500, -0.2, 0.2);
139 fOutput->Add(fDeltaTheta);
140
141}
142
143//----------------------------------------------------------------------
144
145void AliTrackletsTask::Exec(Option_t*)
146{
147 //
148 // process the event
149 //
150
151 PostData(0, fOutput);
152
153 if (!fESD){
154 AliError("ESD branch not available");
155 return;
156 }
157
158 AliESDHeader* esdHeader = fESD->GetHeader();
159 if (!esdHeader){
160 Printf("ERROR: esdHeader could not be retrieved");
161 return;
162 }
163
164 Bool_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
165 if (isSelected == kFALSE) {
166 AliInfo("Event not selected by Physics analysis");
167 return;
168 }
169
170 Int_t ntrk = fESD->GetNumberOfTracks() ;
171 fNtracks->Fill(ntrk);
172
173 const AliMultiplicity* mult = fESD->GetMultiplicity();
174 if (!mult){
175 AliError("AliMultiplicity not found, returning");
176 return;
177 }
178
179 Int_t nTracklets = mult->GetNumberOfTracklets();
180
181 for (Int_t iTracklet = 0; iTracklet < nTracklets; iTracklet++){
182 Float_t eta = mult->GetEta(iTracklet);
183 Float_t phi = mult->GetPhi(iTracklet);
184 Float_t deltaPhi = mult->GetDeltaPhi(iTracklet);
185 Float_t deltaTheta = mult->GetDeltaTheta(iTracklet);
186 if (phi < 0) phi += TMath::Pi() * 2;
187 fPhi->Fill(phi);
188 fEtaPhi->Fill(eta, phi);
189 fDeltaPhi->Fill(deltaPhi);
190 fDeltaTheta->Fill(deltaTheta);
191 }
192
193 // vertex SPD
194 const AliESDVertex* vertexSPD = 0;
195 vertexSPD = fESD->GetPrimaryVertexSPD();
196 if (vertexSPD && (!(vertexSPD->IsFromVertexerZ() && vertexSPD->GetDispersion() > 0.02))) {
197 Double_t vtxX = vertexSPD->GetX();
198 Double_t vtxY = vertexSPD->GetY();
199 Double_t vtxZ = vertexSPD->GetZ();
200 fVtxX->Fill(vtxX);
201 fVtxY->Fill(vtxY);
202 fVtxZ->Fill(vtxZ);
203 fVtx->Fill(vtxX,vtxY,vtxZ);
204 if (vertexSPD->GetNContributors() > 0){
205 fVtxContributors->Fill(vtxX,vtxY,vtxZ);
206 }
207 }
208
209
210}
211
212//----------------------------------------------------------------------
213
214void AliTrackletsTask::Terminate(Option_t *)
215{
216 //
217 // Plotting distributions, and saving them in a file
218 //
219
220 fOutput = dynamic_cast<TList*> (GetOutputData(0));
221 if (!fOutput)
222 Printf("ERROR: fOutput not available");
223
224 fOutput->Print();
225
226 if (fOutput){
227 fNtracks = dynamic_cast<TH1I*>(fOutput->FindObject("fNtracks"));
228 fPhi = dynamic_cast<TH1D*>(fOutput->FindObject("fPhi"));
229 fEtaPhi = dynamic_cast<TH2D*>(fOutput->FindObject("fEtaPhi"));
230 fDeltaPhi = dynamic_cast<TH1D*>(fOutput->FindObject("fDeltaPhi"));
231 fDeltaTheta = dynamic_cast<TH1D*>(fOutput->FindObject("fDeltaTheta"));
232 fVtxX = dynamic_cast<TH1D*>(fOutput->FindObject("fVtxX"));
233 fVtxY = dynamic_cast<TH1D*>(fOutput->FindObject("fVtxY"));
234 fVtxZ = dynamic_cast<TH1D*>(fOutput->FindObject("fVtxZ"));
235 fVtx = dynamic_cast<TH3D*>(fOutput->FindObject("fVtx"));
236 fVtxContributors = dynamic_cast<TH3D*>(fOutput->FindObject("fVtxContributors"));
237
238 new TCanvas("Ntracks", " Ntracks ",50, 50, 550, 550) ;
239 fNtracks->Draw();
240 new TCanvas("phi", " phi ",50, 50, 550, 550) ;
241 fPhi->Draw();
242 new TCanvas("etaphi", " etaphi ",50, 50, 550, 550) ;
243 fEtaPhi->Draw();
244 new TCanvas("deltaPhi", " deltaPhi ",50, 50, 550, 550) ;
245 fDeltaPhi->Draw();
246 new TCanvas("deltaTheta", " deltaTheta ",50, 50, 550, 550) ;
247 fDeltaTheta->Draw();
248 new TCanvas("vtxX", " vtxX ",50, 50, 550, 550) ;
249 fVtxX->Draw();
250 new TCanvas("vtxY", " vtxY ",50, 50, 550, 550) ;
251 fVtxY->Draw();
252 new TCanvas("vtxZ", " vtxZ ",50, 50, 550, 550) ;
253 fVtxZ->Draw();
254
255 TFile* outputFile = new TFile("histograms.root", "RECREATE");
256 fNtracks->Write();
257 fPhi->Write();
258 fEtaPhi->Write();
259 fDeltaPhi->Write();
260 fDeltaTheta->Write();
261 fVtxX->Write();
262 fVtxY->Write();
263 fVtxZ->Write();
264 fVtx->Write();
265 fVtxContributors->Write();
266 outputFile->Write();
267 outputFile->Close();
268
269 }
270}
271
272
273
274
275