]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/MeanPtFluctuations/AliAnalysisTaskPtFluc.cxx
CMakeLists.txt are generated for the moment
[u/mrichter/AliRoot.git] / PWGCF / EBYE / MeanPtFluctuations / AliAnalysisTaskPtFluc.cxx
1 #include "TChain.h"
2 #include "TTree.h"
3 #include "TH1F.h"
4 #include "TH2F.h"
5 #include "TCanvas.h"
6 #include "TF1.h"
7 #include "TProfile.h"
8 #include "TList.h"
9 #include "iostream"
10
11 #include "AliAnalysisTaskSE.h"
12 #include "AliAnalysisManager.h"
13
14 #include "AliESD.h"
15 #include "AliESDEvent.h"
16 #include "AliESDInputHandler.h"
17 #include "AliESDtrackCuts.h"
18
19 #include "AliLog.h"
20
21 #include "AliAnalysisTaskPtFluc.h"
22
23
24 using namespace std;
25
26 // Analysis of Pt FLuctuations (pp)
27 // Author: Stefan Heckel
28 // Version of pp task:   8.0, 19.04.2011
29
30
31 ClassImp(AliAnalysisTaskPtFluc)
32
33 //________________________________________________________________________
34 AliAnalysisTaskPtFluc::AliAnalysisTaskPtFluc(const char *name)
35   :AliAnalysisTaskSE(name),
36   fESD(0),
37   fOutputList(0),
38   fPtSpec(0),
39   fMult(0),
40   fEta(0),
41   fEtaPhiPlus(0),
42   fEtaPhiMinus(0),
43   fVtxZ(0),
44   fVtxZCut(0),
45   fVtxZCont(0),
46   fVtxZTrackCuts(0),
47   fEventMeanPt(0),
48   fEventMeanPtSq(0),
49   fMultEventMeanPt(0),
50   fMultEventMeanPtSq(0),
51   fTwoPartCorrEv(0),
52   fTwoPartCorrEvSq(0),
53   fTwoPartCorrEvSample(0),
54   fTwoPartCorrEvSampleSq(0),
55   fESDTrackCuts(0),
56   fMaxVertexZ(0),
57   fNContributors(0),
58   fMC(0)
59 {
60   DefineOutput(1, TList::Class());
61 }
62
63 //________________________________________________________________________
64 AliAnalysisTaskPtFluc::~AliAnalysisTaskPtFluc()
65 {
66   if(fOutputList) delete fOutputList;  fOutputList =0;
67 }
68
69 //________________________________________________________________________
70 void AliAnalysisTaskPtFluc::UserCreateOutputObjects()
71 {
72   // Create histograms
73   // Called once
74
75   OpenFile(1, "RECREATE");
76   fOutputList = new TList();
77   fOutputList->SetOwner();
78
79
80   fPtSpec = new TH1F("fPtSpec","Pt spectrum",50,0,2.5);
81
82   fMult = new TH1F("fMult","Multiplicity distribution",120,-0.5,119.5);
83
84   fEta = new TH1F("fEta","Eta distribution",80,-2,2);
85
86   fEtaPhiPlus = new TH1F("fEtaPhiPlus","Phi distribution for positive eta",62,0,6.2);
87
88   fEtaPhiMinus = new TH1F("fEtaPhiMinus","Phi distribution for negative eta",62,0,6.2);
89
90   fVtxZ = new TH1F("fVtxZ","Vertex Z distribution before cuts",100,-20,20);
91
92   fVtxZCut = new TH1F("fVtxZCut","Vertex Z distribution after vtxZ cut",110,-11,11);
93
94   fVtxZCont = new TH1F("fVtxZCont","Vertex Z distribution after nCont cut",110,-11,11);
95
96   fVtxZTrackCuts = new TH1F("fVtxZTrackCuts","Vertex Z distribution after track cuts",110,-11,11);
97
98
99   fEventMeanPt = new TH1F("fEventMeanPt","Mean-Pt event by event",50,0,2.5);
100
101   fEventMeanPtSq = new TH1F("fEventMeanPtSq","Mean-Pt event by event squared",100,0,5);
102
103   fMultEventMeanPt = new TH1F("fMultEventMeanPt","Mean-Pt event by event vs. Multiplicity",100,-0.5,99.5);
104
105   fMultEventMeanPtSq = new TH1F("fMultEventMeanPtSq","Mean-Pt event by event squared vs. Multiplicity",100,-0.5,99.5);
106
107
108   fTwoPartCorrEv = new TH1F("fTwoPartCorrEv","Two-particle correlator vs. Multiplicity",100,-0.5,99.5);
109
110   fTwoPartCorrEvSq = new TH1F("fTwoPartCorrEvSq","Two-particle correlator squared vs. Multiplicity",100,-0.5,99.5);
111
112   fTwoPartCorrEvSample = new TH1F("fTwoPartCorrEvSample","Two-particle correlator vs. Multiplicity",100,-0.5,99.5);
113
114   fTwoPartCorrEvSampleSq = new TH1F("fTwoPartCorrEvSampleSq","Two-particle correlator squared vs. Multiplicity",100,-0.5,99.5);
115
116
117   // Add histograms to the output list
118   fOutputList->Add(fPtSpec);
119   fOutputList->Add(fMult);
120   fOutputList->Add(fEta);
121   fOutputList->Add(fEtaPhiPlus);
122   fOutputList->Add(fEtaPhiMinus);
123   fOutputList->Add(fVtxZ);
124   fOutputList->Add(fVtxZCut);
125   fOutputList->Add(fVtxZCont);
126   fOutputList->Add(fVtxZTrackCuts);
127   fOutputList->Add(fEventMeanPt);
128   fOutputList->Add(fEventMeanPtSq);
129   fOutputList->Add(fMultEventMeanPt);
130   fOutputList->Add(fMultEventMeanPtSq);
131   fOutputList->Add(fTwoPartCorrEv);
132   fOutputList->Add(fTwoPartCorrEvSq);
133   fOutputList->Add(fTwoPartCorrEvSample);
134   fOutputList->Add(fTwoPartCorrEvSampleSq);
135 }
136
137 //________________________________________________________________________
138 void AliAnalysisTaskPtFluc::UserExec(Option_t *)
139 {
140   // Main loop
141   // Called for each event
142
143
144   // --- ESD event handler ---
145
146   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
147
148   if (!esdH) {
149     Printf("ERROR: Could not get ESDInputHandler");
150   }
151   else {
152     fESD = esdH->GetEvent();
153   }
154
155   if (!fESD) {
156     Printf("ERROR: fESD not available");
157     return;
158   }
159
160
161   // --- End event handler ---
162
163
164   // --- Vertex cuts ---
165
166   // TPC+ITS vertex
167   const AliESDVertex* vtxESD = fESD->GetPrimaryVertexTracks();
168   Double_t vtxZ = vtxESD->GetZv();
169   Double_t vtxNCont = vtxESD->GetNContributors();
170
171 //   // TPConly vertex
172 //   const AliESDVertex* vtxESDTPC = fESD->GetPrimaryVertexTPC();
173 //   Double_t vtxZ = vtxESDTPC->GetZv();
174 //   Double_t vtxNCont = vtxESDTPC->GetNContributors();
175
176   fVtxZ->Fill(vtxZ); // VtxZ before cuts
177
178   // Event cut on the z-Position of the vertex
179   if(vtxZ > fMaxVertexZ || vtxZ < (-1.*fMaxVertexZ)) {
180 //     Printf("VertexZ out of range, Zv = %f",vtxZ);
181     return;
182   }
183
184   fVtxZCut->Fill(vtxZ); // VtxZ after cut on vtxZ
185
186   // Event cut on the number of contributors
187   if(vtxNCont < fNContributors) {
188 //     Printf("Vertex has no contributors");
189     return;
190   }
191 //   else {
192 //     Printf("Vertex nContributors = %i",vtxNCont);
193 //   }
194
195   fVtxZCont->Fill(vtxZ); // VtxZ after cut on nContributors
196
197   // --- End vertex cuts ---
198
199
200   Int_t nrTracks = 0;           // count for number of tracks which pass the track cuts
201   Int_t nrESDTracks = 0;        // number of tracks in the ESD file
202   Double_t sumPt = 0.;          // sum of track Pts
203   Double_t twoPartCorrEv = 0.;  // Two-particle correlator of one event for multiplicity bin analysis
204   Double_t twoPartCorrEvSample = 0.; // Two-part. corr. of one event for whole sample analysis
205
206   Double_t tracks[120] = {0.};  // array of track Pts, needed for the two-particle correlator
207
208   Double_t trackPt, trackEta, trackPhi;
209   Double_t eventMeanPt, eventMeanPtSq, evMptMult;
210   Double_t twoPartCorrPair, twoPartCorrEvSq;
211   Double_t twoPartCorrPairSample, twoPartCorrEvSampleSq;
212
213   Double_t *nbins = 0x0; // Mean pT values for multiplicity bin analysis
214   Double_t evMptSample;  // Mean pT value for whole sample analysis
215
216   Double_t energy; // beam energy
217
218   energy = fESD->GetBeamEnergy();
219
220
221 //   Printf("MC: %d,  %f",fMC,energy);
222
223
224 // --- Mean pT values ---
225
226 // !!!!! Have to be calculated in a first run - for each sample, which should be analysed, separately! !!!!!
227
228
229 // -- Mean pT values for the multiplicity bins (starting with N = 0, therefore the first value has to be 0.000) --
230
231 // MC 900 GeV 10e13 (Pythia6 - Perugia0) (Pt-Range: 0.15 - 2)
232 Double_t nbinsMC900[100] = {0.000, 0.456, 0.475, 0.490, 0.507, 0.523, 0.535, 0.546, 0.554, 0.561, 0.568, 0.572, 0.578, 0.582, 0.586, 0.590, 0.593, 0.596, 0.602, 0.601, 0.606, 0.605, 0.610, 0.618, 0.610, 0.620, 0.633, 0.615, 0.609, 0.620, 0.632, 0.634, 0.661, 0.647, 0.560, 0.671, 0.660, 0.678, 0.000, 0.612, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
233
234 // MC 7 TeV 10f6a (Pythia6 - Perugia0) (Pt-Range: 0.15 - 2)
235 Double_t nbinsMC7[100] = {0.000, 0.453, 0.475, 0.498, 0.523, 0.546, 0.563, 0.577, 0.587, 0.595, 0.603, 0.610, 0.616, 0.622, 0.627, 0.632, 0.636, 0.640, 0.644, 0.648, 0.651, 0.655, 0.658, 0.662, 0.664, 0.667, 0.669, 0.671, 0.674, 0.677, 0.679, 0.682, 0.683, 0.686, 0.688, 0.690, 0.692, 0.695, 0.696, 0.696, 0.699, 0.700, 0.701, 0.705, 0.709, 0.705, 0.701, 0.705, 0.714, 0.718, 0.727, 0.702, 0.725, 0.738, 0.724, 0.728, 0.710, 0.729, 0.711, 0.748, 0.702, 0.708, 0.732, 0.000, 0.683, 0.000, 0.000, 0.675, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
236
237
238 // // MC 900 GeV 10e12 (Phojet) (Pt-Range: 0.15 - 2)
239 // Double_t nbinsMC900[100] = {0.000, 0.478, 0.483, 0.488, 0.494, 0.498, 0.502, 0.505, 0.507, 0.509, 0.513, 0.516, 0.520, 0.522, 0.526, 0.528, 0.531, 0.535, 0.539, 0.540, 0.544, 0.540, 0.548, 0.548, 0.549, 0.556, 0.562, 0.549, 0.557, 0.550, 0.561, 0.559, 0.554, 0.547, 0.539, 0.564, 0.546, 0.603, 0.532, 0.583, 0.656, 0.531, 0.458, 0.547, 0.000, 0.633, 0.000, 0.603, 0.588, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
240
241 // // MC 7 TeV 10f6 (Phojet) (Pt-Range: 0.15 - 2)
242 // Double_t nbinsMC7[100] = {0.000, 0.488, 0.499, 0.512, 0.526, 0.536, 0.544, 0.549, 0.553, 0.557, 0.560, 0.563, 0.566, 0.569, 0.572, 0.574, 0.577, 0.579, 0.580, 0.582, 0.584, 0.585, 0.586, 0.587, 0.588, 0.589, 0.589, 0.590, 0.591, 0.591, 0.590, 0.591, 0.592, 0.591, 0.589, 0.591, 0.590, 0.588, 0.590, 0.589, 0.590, 0.588, 0.589, 0.589, 0.591, 0.588, 0.582, 0.580, 0.582, 0.586, 0.575, 0.575, 0.560, 0.570, 0.569, 0.604, 0.616, 0.549, 0.614, 0.611, 0.480, 0.628, 0.568, 0.000, 0.677, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
243
244
245 // Data 900 GeV 10c.pass3 (Pt-Range: 0.15 - 2)
246 Double_t nbinsData900[100] = {0.000, 0.490, 0.491, 0.495, 0.503, 0.512, 0.521, 0.529, 0.536, 0.542, 0.548, 0.551, 0.556, 0.561, 0.564, 0.567, 0.570, 0.573, 0.575, 0.578, 0.580, 0.582, 0.583, 0.590, 0.591, 0.591, 0.590, 0.598, 0.593, 0.608, 0.609, 0.609, 0.587, 0.615, 0.603, 0.613, 0.560, 0.583, 0.614, 0.649, 0.584, 0.584, 0.594, 0.000, 0.622, 0.000, 0.000, 0.703, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
247
248 // Data 2.76 TeV 11a.pass1 (Pt-Range: 0.15 - 2)
249 Double_t nbinsData276[100] = {0.000, 0.491, 0.500, 0.511, 0.524, 0.536, 0.548, 0.558, 0.567, 0.574, 0.581, 0.586, 0.592, 0.596, 0.601, 0.605, 0.609, 0.613, 0.616, 0.619, 0.622, 0.625, 0.628, 0.629, 0.632, 0.634, 0.637, 0.639, 0.641, 0.645, 0.645, 0.646, 0.649, 0.652, 0.655, 0.658, 0.656, 0.655, 0.659, 0.668, 0.660, 0.672, 0.672, 0.665, 0.669, 0.662, 0.666, 0.669, 0.677, 0.673, 0.595, 0.642, 0.700, 0.623, 0.670, 0.704, 0.694, 0.000, 0.000, 0.667, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
250
251 // Data 7 TeV 10d.pass2 + 10e.pass2 (Pt-Range: 0.15 - 2)
252 Double_t nbinsData7[100] = {0.000, 0.495, 0.498, 0.505, 0.516, 0.528, 0.539, 0.550, 0.559, 0.566, 0.573, 0.579, 0.584, 0.589, 0.593, 0.597, 0.601, 0.605, 0.608, 0.611, 0.614, 0.617, 0.619, 0.622, 0.624, 0.626, 0.629, 0.631, 0.633, 0.634, 0.636, 0.638, 0.640, 0.642, 0.643, 0.645, 0.646, 0.648, 0.649, 0.651, 0.652, 0.654, 0.654, 0.656, 0.657, 0.658, 0.659, 0.660, 0.662, 0.662, 0.663, 0.665, 0.664, 0.665, 0.668, 0.667, 0.670, 0.666, 0.670, 0.670, 0.677, 0.673, 0.672, 0.679, 0.680, 0.667, 0.682, 0.687, 0.675, 0.658, 0.677, 0.663, 0.686, 0.671, 0.682, 0.683, 0.667, 0.000, 0.000, 0.662, 0.720, 0.683, 0.720, 0.622, 0.000, 0.729, 0.000, 0.603, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000, 0.000};
253
254 // -- End mean pT values for the multiplicity bins --
255
256
257 // -- Selection of MC/Data and energy; whole sample values --
258
259 if (fMC) { // - MC -
260
261   if (energy < 1000.) {
262
263 //     Printf(" -- MC, 900 GeV -- ");
264
265     evMptSample = 0.5307; // MC 900 GeV 10e13 (Pythia6 - Perugia0) (Pt-Range: 0.15 - 2)
266 //     evMptSample = 0.5044; // MC 900 GeV 10e12 (Phojet) (Pt-Range: 0.15 - 2)
267
268     nbins = nbinsMC900;
269   }
270   else {
271
272 //     Printf(" -- MC, 7 TeV -- ");
273
274     evMptSample = 0.5847; // MC 7 TeV 10f6a (Pythia6 - Perugia0) (Pt-Range: 0.15 - 2)
275 //     evMptSample = 0.5517; // MC 7 TeV 10f6 (Phojet) (Pt-Range: 0.15 - 2)
276
277     nbins = nbinsMC7;
278   }
279
280 } // - End MC -
281 else { // - Data -
282
283   if (energy < 1000.) {
284
285 //     Printf(" -- Data, 900 GeV -- ");
286
287     evMptSample = 0.5267; // Data 900 GeV 10c.pass3 (Pt-Range: 0.15 - 2)
288     nbins = nbinsData900;
289   }
290   else if (energy > 1000. && energy < 3000.) {
291
292 //     Printf(" -- Data, 2.76 TeV -- ");
293
294     evMptSample = 0.5622; // Data 2.76 TeV 11a.pass1 (Pt-Range: 0.15 - 2)
295     nbins = nbinsData276;
296   }
297   else {
298
299 //     Printf(" -- Data, 7 TeV -- ");
300
301     evMptSample = 0.5743; // Data 7 TeV 10d.pass2 + 10e.pass2 (Pt-Range: 0.15 - 2)
302     nbins = nbinsData7;
303   }
304
305 } // - End data -
306
307 // -- End selection of MC/Data and energy; whole sample values --
308
309 // --- End mean pT values ---
310
311
312   nrESDTracks = fESD->GetNumberOfTracks();
313 //   if( nrESDTracks ) { printf("Found event with %i tracks.\n",nrESDTracks); }
314
315
316   // --- Loop over all tracks of one event ---
317
318   for (Int_t iTracks = 0; iTracks < nrESDTracks; iTracks++) {
319     AliESDtrack* track = fESD->GetTrack(iTracks);
320     if (!track) {
321       Printf("ERROR: Could not receive track %d\n", iTracks);
322       continue;
323     }
324     if(!fESDTrackCuts) {
325       Printf("ERROR: No esd track cut");
326       continue;
327     }
328     else {
329       if(!fESDTrackCuts->AcceptTrack(track))
330         continue;
331     }
332         trackPt = track->Pt();
333         fPtSpec->Fill(trackPt);
334         tracks[nrTracks] = trackPt;
335
336         trackEta = track->Eta();
337         fEta->Fill(trackEta);
338
339         trackPhi = track->Phi();
340
341         if (trackEta > 0.) {
342           fEtaPhiPlus->Fill(trackPhi);
343         }
344         else if (trackEta < 0.) {
345           fEtaPhiMinus->Fill(trackPhi);
346         }
347
348         sumPt = sumPt + trackPt;
349         nrTracks++;
350
351 //      printf("Track Pt = %.3f   Track Eta = %.3f   Track Phi = %3f\n",trackPt,trackEta,trackPhi);
352
353   } // --- End track loop ---
354
355
356   // --- Calculation of various values and filling of histograms (for remaining events with N_acc > 0) ---
357
358   if(nrTracks != 0) {
359
360     // VtxZ after all track cuts
361     fVtxZTrackCuts->Fill(vtxZ);
362
363     // Multiplicity distribution
364     fMult->Fill(nrTracks);
365
366     // Calculation of mean Pt and mean Pt Squared
367     eventMeanPt = sumPt / nrTracks;
368     eventMeanPtSq = eventMeanPt * eventMeanPt;
369
370     // Mean-Pt and Mean-Pt squared
371     fEventMeanPt->Fill(eventMeanPt);
372     fEventMeanPtSq->Fill(eventMeanPtSq);
373
374     // Mean-Pt and Mean-Pt squared depending on multiplicity
375     fMultEventMeanPt->Fill(nrTracks,eventMeanPt);
376     fMultEventMeanPtSq->Fill(nrTracks,eventMeanPtSq);
377
378 //     printf("nrTracks: %i   sumPt: %.8f   meanPt: %.8f   meanPtSq: %.8f\n",nrTracks,sumPt,eventMeanPt,eventMeanPtSq);
379
380
381     // --- Two-particle correlator ---
382
383     evMptMult = nbins[nrTracks];
384 //     printf("nrTracks = %3i   evMptMult = %.3f   evMptSample = %.3f \n\n",nrTracks,evMptMult,evMptSample);
385
386     for (int i=0; i<nrTracks; i++) {
387
388 //       printf("--- TrackPt = %.3f ",tracks[i]);
389
390       for (int j=i+1; j<nrTracks; j++) {
391
392         twoPartCorrPair = (tracks[i] - evMptMult) * (tracks[j] - evMptMult); // Multiplicity bins
393         twoPartCorrEv = twoPartCorrEv + twoPartCorrPair;
394
395         twoPartCorrPairSample = (tracks[i] - evMptSample) * (tracks[j] - evMptSample); // Whole sample
396         twoPartCorrEvSample = twoPartCorrEvSample + twoPartCorrPairSample;
397
398 //      printf("twoPartCorrPair = %.3f   twoPartCorrPairSample = %.3f \n",twoPartCorrPair,twoPartCorrPairSample);
399       }
400 //       printf("twoPartCorrEv = %.3f   twoPartCorrEvSample = %.3f \n",twoPartCorrEv,twoPartCorrEvSamples);
401     }
402 //     printf("Total Event: twoPartCorrEv = %.3f   twoPartCorrEvSample = %.3f \n\n",twoPartCorrEv,twoPartCorrEvSample);
403
404     // Multiplicity bins
405     twoPartCorrEvSq = twoPartCorrEv * twoPartCorrEv;
406     fTwoPartCorrEv->Fill(nrTracks,twoPartCorrEv);
407     fTwoPartCorrEvSq->Fill(nrTracks,twoPartCorrEvSq);
408
409     // Whole sample
410     twoPartCorrEvSampleSq = twoPartCorrEvSample * twoPartCorrEvSample;
411     fTwoPartCorrEvSample->Fill(nrTracks,twoPartCorrEvSample);
412     fTwoPartCorrEvSampleSq->Fill(nrTracks,twoPartCorrEvSampleSq);
413
414     // --- End two-particle correlator ---
415
416
417   } // --- End calculation of various values and filling of histograms ---
418
419
420   // Post output data
421   PostData(1, fOutputList);
422 }
423
424 void AliAnalysisTaskPtFluc::Terminate(Option_t *)
425 {
426   // Draw result to the screen
427   // Called once at the end of the query
428
429   fOutputList = dynamic_cast<TList*> (GetOutputData(1));
430   if (!fOutputList) {
431     Printf("ERROR: fOutputList not available\n");
432     return;
433   }
434
435 }