]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGCF/EBYE/Fluctuations/AliHigherMomentsToyModel.cxx
new analysis task: multiplicity fluctuations (Maitreyee Mukherjee <maitreyee.mukherje...
[u/mrichter/AliRoot.git] / PWGCF / EBYE / Fluctuations / AliHigherMomentsToyModel.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: Nirbhay K. Behera                                              *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16
17 //=========================================================================//
18 //                                                                         //
19 //           Toy Model for Net-Charge Higher Moment Analysis               //
20 //              Author: Satyajit Jena || Nirbhay K. Behera                 //
21 //                      sjena@cern.ch || nbehera@cern.ch                   //
22 //                                                                         //
23 //=========================================================================//
24
25 #include "TChain.h"
26 #include "TList.h"
27 #include "TFile.h"
28 #include "TTree.h"
29 #include "TH1D.h"
30 #include "TH2F.h"
31 #include "THnSparse.h"
32 #include "TCanvas.h"
33 #include "TParticle.h"
34 #include <TDatabasePDG.h>
35 #include <TParticlePDG.h>
36
37 #include "AliAnalysisTaskSE.h"
38 #include "AliAnalysisManager.h"
39
40 #include "AliESDVertex.h"
41 #include "AliESDEvent.h"
42 #include "AliESDInputHandler.h"
43 #include "AliAODEvent.h"
44 #include "AliAODTrack.h"
45 #include "AliAODInputHandler.h"
46 #include "AliESD.h"
47 #include "AliESDEvent.h"
48 #include "AliAODEvent.h"
49 #include "AliStack.h"
50 #include "AliESDtrackCuts.h"
51 #include "AliAODMCHeader.h"
52 #include "AliAODMCParticle.h"
53 #include "AliMCEventHandler.h"
54 #include "AliMCEvent.h"
55 #include "AliMCParticle.h"
56 #include "AliStack.h"
57 #include "AliGenHijingEventHeader.h"
58 #include "AliGenEventHeader.h"
59 #include "AliPID.h"
60 #include "AliAODPid.h"                
61 #include "AliPIDResponse.h"
62 #include "AliAODpidUtil.h"        
63 #include "AliPIDCombined.h"
64
65 #include "AliHigherMomentsToyModel.h"
66
67 using std::cout;
68 using std::endl;
69
70
71 ClassImp(AliHigherMomentsToyModel)
72
73
74 //-----------------------------------------------------------------------
75 AliHigherMomentsToyModel::AliHigherMomentsToyModel( const char *name )
76 : AliAnalysisTaskSE( name ),
77   fListOfHistosQA(0),
78   fListOfHistos(0),
79   fAOD(0),
80   fArrayMC(0),
81   fPIDResponse(0),
82   fParticleSpecies(AliPID::kProton),
83   fAnalysisType("AOD"),
84   fCentralityEstimator("V0M"),
85   fCentrality(0),
86   fVxMax(3.),
87   fVyMax(3.),
88   fVzMax(10.),
89   fDCAxy(2.4),
90   fDCAz(3.),
91   fPtLowerLimit(0.3),
92   fPtHigherLimit(1.5),
93   fEtaLowerLimit(-0.8),
94   fEtaHigherLimit(0.8),
95   fRapidityCut(0.5),
96   fNSigmaCut(3.),
97   fTPCNClus(80),
98   fChi2perNDF(4.),
99   fAODtrackCutBit(1024),
100   fUsePid(kFALSE),
101   fEventCounter(0),
102   fHistDCA(0),
103   fTPCSig(0),
104   fTPCSigA(0),
105   fTHnCentNplusNminusCh(0),
106   fTHnCentNplusNminusChTruth(0),
107   fTHnCentNplusNminus(0)
108
109   
110   for ( Int_t i = 0; i < 13; i++) { 
111     fHistQA[i] = NULL;
112   }
113   
114   for ( Int_t i = 0; i < 5; i++ ){
115     fTHnCentNplusNminusPid[i] = NULL;
116     fTHnCentNplusNminusPidTruth[i] = NULL;
117   }
118
119   DefineOutput(1, TList::Class()); // Outpput....
120   DefineOutput(2, TList::Class()); 
121
122 }
123
124 AliHigherMomentsToyModel::~AliHigherMomentsToyModel() {
125   
126   if(fListOfHistosQA) delete fListOfHistosQA;
127   if(fListOfHistos) delete fListOfHistos;
128
129 }
130
131 //---------------------------------------------------------------------------------
132 void AliHigherMomentsToyModel::UserCreateOutputObjects() {
133   //For QA-Histograms
134   fListOfHistosQA = new TList();
135   fListOfHistosQA->SetOwner(kTRUE);
136   
137   fListOfHistos = new TList();
138   fListOfHistos->SetOwner(kTRUE);
139   
140   fEventCounter = new TH1D("fEventCounter","EventCounter", 10, 0.5,10.5);
141   fEventCounter->GetXaxis()->SetBinLabel(1,"Event Accesed");
142   fEventCounter->GetXaxis()->SetBinLabel(2,"Within 0-90% centrality");
143   fEventCounter->GetXaxis()->SetBinLabel(5,"Have a vertex");
144   fEventCounter->GetXaxis()->SetBinLabel(6,"After vertex Cut");
145   fEventCounter->GetXaxis()->SetBinLabel(7,"Event Analyzed");
146   fEventCounter->GetXaxis()->SetBinLabel(8,"Event Analysis finished");
147   fListOfHistosQA->Add(fEventCounter);
148   
149   fHistQA[0] = new TH1D("fHistQAvx", "Histo Vx After Cut", 400, -4., 4.);
150   fHistQA[1] = new TH1D("fHistQAvy", "Histo Vy After Cut", 400, -4., 4.);
151   fHistQA[2] = new TH1D("fHistQAvz", "Histo Vz After Cut", 500, -25., 25.);  
152   fHistQA[3] = new TH1D("fHistQAvxA", "Histo Vx All", 500, -5., 5.);
153   fHistQA[4] = new TH1D("fHistQAvyA", "Histo Vy All", 500, -5., 5.);
154   fHistQA[5] = new TH1D("fHistQAvzA", "Histo Vz All", 500, -25., 25.);
155   fHistQA[6] = new TH1D("fHistQADcaXyC", "Histo DCAxy after cut", 500, -5., 5.);
156   fHistQA[7] = new TH1D("fHistQADcaZC", "Histo DCAz after cut", 500, -5., 5.);   
157   fHistQA[8] = new TH1D("fHistQAPt","p_{T} distribution",600,0,6);
158   fHistQA[9] = new TH1D("fHistQAEta","#eta distribution",240,-1.2,1.2);
159   fHistQA[10] = new TH1D("fHistQAPhi","#phi distribution",340,0,6.8);
160   fHistQA[11] = new TH1D("fHistQANCls","Number of TPC cluster",200,0,200);
161   fHistQA[12] = new TH1D("fHistQAChi2","Chi2 per NDF",100,0,10);
162   
163   for(Int_t i = 0; i < 13; i++)
164     {
165       fListOfHistosQA->Add(fHistQA[i]);
166     }
167   
168   fHistDCA = new TH2D("fHistDCA","DCAxy Vs DCAz", 500, -5., 5., 500, -5., 5.);
169   fTPCSig = new TH2D("fTPCSig","TPC signal",200, 0.0, 10. ,1000,0.,1000);
170   fTPCSig->SetMarkerColor(kRed);
171   fTPCSigA = new TH2D("fTPCSigA","TPC signal all ",200, 0.0, 10. ,1000,0.,1000);
172   fListOfHistosQA->Add(fHistDCA);
173   fListOfHistosQA->Add(fTPCSig);
174   fListOfHistosQA->Add(fTPCSigA);
175   
176   
177   const Int_t nDim = 3;
178   const Int_t nPid = 5;
179   
180   TString Species[nPid] = {"Electron","Muon","Pion","Kaon","Proton"};
181   
182   Int_t fBins[nPid][nDim];
183   Double_t fMin[nPid][nDim];
184   Double_t fMax[nPid][nDim];
185   
186   for( Int_t i = 0; i < nPid; i++ ){
187     for( Int_t j = 0; j < nDim; j++ ){
188       fBins[i][j] = 0;
189       fMin[i][j] = 0.;
190       fMax[i][j] = 0.;
191     }
192   }
193   
194   Int_t fBinsCh[nDim] = {100, 1500, 1500};
195   Double_t fMinCh[nDim] = { -0.5, -0.5, -0.5 };
196   Double_t fMaxCh[nDim] = { 99.5, 1499.5, 1499.5};
197   
198   fTHnCentNplusNminusCh = new THnSparseD("fTHnCentNplusNminusCh","Cent-NplusChrg-NminusChrg", nDim, fBinsCh, fMinCh, fMaxCh); 
199   fTHnCentNplusNminusCh->GetAxis(0)->SetTitle("Centrality");
200   fTHnCentNplusNminusCh->GetAxis(1)->SetTitle("Nplus");
201   fTHnCentNplusNminusCh->GetAxis(2)->SetTitle("Nminus");
202   fListOfHistos->Add(fTHnCentNplusNminusCh);
203
204   if( fAnalysisType == "MCAOD" ){   
205     
206     fTHnCentNplusNminusChTruth = new THnSparseD("fTHnCentNplusNminusChTruth","Cent-NplusChrg-NminusChrg", nDim, fBinsCh, fMinCh, fMaxCh); 
207     fTHnCentNplusNminusChTruth->GetAxis(0)->SetTitle("Centrality");
208     fTHnCentNplusNminusChTruth->GetAxis(1)->SetTitle("Nplus");
209     fTHnCentNplusNminusChTruth->GetAxis(2)->SetTitle("Nminus");
210     fListOfHistos->Add(fTHnCentNplusNminusChTruth);
211   }//MCAOD---
212  
213   
214   TString hname1, hname11;
215   TString htitle1, axisTitle1, axisTitle2;
216   
217   
218   
219   if( fUsePid ){
220     
221     Int_t gPid = (Int_t)fParticleSpecies;
222     
223     if( gPid > 1 && gPid < 5 ){ 
224       //Pion----
225       fBins[2][0] = 100; fBins[2][1] = 1000; fBins[2][2] = 1000;
226       fMin[2][0] = -0.5; fMin[2][1] = -0.5; fMin[2][2] = -0.5;
227       fMax[2][0] = 99.5; fMax[2][1] = 999.5; fMax[2][2] = 999.5;
228       //Kaon------
229       fBins[3][0] = 100; fBins[3][1] = 600; fBins[3][2] = 600;
230       fMin[3][0] = -0.5; fMin[3][1] = -0.5; fMin[3][2] = -0.5;
231       fMax[3][0] = 99.5; fMax[3][1] = 599.5; fMax[3][2] = 599.5;
232       //Proton-----
233       fBins[4][0] = 100; fBins[4][1] = 400; fBins[4][2] = 400;
234       fMin[4][0] = -0.5; fMin[4][1] = -0.5; fMin[4][2] = -0.5;
235       fMax[4][0] = 99.5; fMax[4][1] = 399.5; fMax[4][2] = 399.5;
236       
237       hname1 = "fCentNplusNminusPid"; hname1 +=gPid;
238       htitle1 = Species[gPid]; htitle1 +=" And Neg-"; htitle1 +=Species[gPid];
239       axisTitle1 = Species[gPid];
240       axisTitle2 ="Neg-"; axisTitle2 += Species[gPid];
241       
242       fTHnCentNplusNminusPid[gPid] = new THnSparseD(hname1.Data(),htitle1.Data(), nDim, fBins[gPid], fMin[gPid],fMax[gPid]);
243       
244       fTHnCentNplusNminusPid[gPid]->GetAxis(0)->SetTitle("Centrality");
245       fTHnCentNplusNminusPid[gPid]->GetAxis(1)->SetTitle(axisTitle1.Data());
246       fTHnCentNplusNminusPid[gPid]->GetAxis(1)->SetTitle(axisTitle2.Data());
247       
248       fListOfHistos->Add(fTHnCentNplusNminusPid[gPid]);
249       
250       if( fAnalysisType == "MCAOD" ){
251         
252         hname11 = "fCentNplusNminusPidTruth"; hname11 +=gPid;
253         fTHnCentNplusNminusPidTruth[gPid] = new THnSparseD(hname11.Data(),htitle1.Data(), nDim, fBins[gPid], fMin[gPid],fMax[gPid]);
254         
255         fTHnCentNplusNminusPidTruth[gPid]->GetAxis(0)->SetTitle("Centrality");
256         fTHnCentNplusNminusPidTruth[gPid]->GetAxis(1)->SetTitle(axisTitle1.Data());
257         fTHnCentNplusNminusPidTruth[gPid]->GetAxis(1)->SetTitle(axisTitle2.Data());
258         
259         fListOfHistos->Add(fTHnCentNplusNminusPidTruth[gPid]);
260         
261       }//MCAOD-----
262     }//Pion, Koan and Proton-------
263     
264     else{
265       
266       Int_t fBinsX[nDim] = {100, 1500, 1500};
267       Double_t fMinX[nDim] = { -0.5, -0.5, -0.5 };
268       Double_t fMaxX[nDim] = { 99.5, 1499.5, 1499.5};
269       
270       fTHnCentNplusNminus = new THnSparseD("fTHnCentNplusNminus","Cent-NplusChrg-NminusChrg", nDim, fBinsX, fMinX, fMaxX); 
271       fTHnCentNplusNminus->GetAxis(0)->SetTitle("Centrality");
272       fTHnCentNplusNminus->GetAxis(1)->SetTitle("UnwantedPlus");
273       fTHnCentNplusNminus->GetAxis(2)->SetTitle("UnwantedMinus");
274       fListOfHistos->Add(fTHnCentNplusNminus);
275       
276     }//Unwanted particle -------
277     
278   }//fUsePid-------
279   
280   
281   
282   PostData(1, fListOfHistosQA);
283   PostData(2, fListOfHistos);
284   
285 }
286
287
288 //----------------------------------------------------------------------------------
289 void AliHigherMomentsToyModel::UserExec( Option_t * ){
290   
291   fEventCounter->Fill(1);
292
293   
294   if(fAnalysisType == "AOD") {
295     
296     doAODEvent();
297     
298   }//AOD--analysis-----
299
300   else if(fAnalysisType == "MCAOD") {
301   
302     doMCAODEvent();
303     
304   }
305   
306   else return;
307   
308   
309   
310   
311   fEventCounter->Fill(8);
312   
313   PostData(1, fListOfHistosQA);
314   PostData(2, fListOfHistos);
315   
316 }
317
318 //--------------------------------------------------------------------------------------
319 void AliHigherMomentsToyModel::doAODEvent(){
320  
321   //-------------------
322   Double_t nPlusCharge = 0.;
323   Double_t nMinusCharge = 0.;
324   Double_t nPartile = 0.;
325   Double_t nAntiParticle = 0.;
326   Int_t gPid = 0.;
327
328
329   //connect to the analysis mannager-----
330  
331   AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
332   if (!manager) {
333     cout<<"ERROR: Analysis manager not found."<<endl;
334     return;
335   }
336   //coneect to the inputHandler------------
337   AliAODInputHandler* inputHandler = dynamic_cast<AliAODInputHandler*> (manager->GetInputEventHandler());
338   if (!inputHandler) {
339     cout<<"ERROR: Input handler not found."<<endl;
340     return;
341   }
342   
343   fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
344   if (!fAOD)
345     {
346       cout<< "ERROR 01: AOD not found " <<endl;
347       return;
348     }
349   
350   fPIDResponse =inputHandler->GetPIDResponse();
351   
352   
353   if (!fPIDResponse){
354     AliFatal("This Task needs the PID response attached to the inputHandler");
355     return;
356   }  
357   
358   
359   AliAODHeader *aodHeader = fAOD->GetHeader();
360   
361   fCentrality = (Int_t)aodHeader->GetCentralityP()->GetCentralityPercentile(fCentralityEstimator.Data());
362   /* Int_t cent = -1;
363   cent =  aodHeader->GetCentralityP()->GetCentralityClass10(fCentralityEstimator.Data());
364   
365   if (cent == 0)
366     fCentrality = aodHeader->GetCentralityP()->GetCentralityClass5(fCentralityEstimator.Data());
367   else if (cent == 10 || cent == -1.)
368     fCentrality = -1;
369   else if (cent > 0 && cent < 9)
370     fCentrality = cent + 1;
371   */
372   if(fCentrality < 0 || fCentrality >= 91) return;
373   
374   fEventCounter->Fill(2);
375   
376   if(!ProperVertex()) return;   
377   
378   Int_t nTracks = fAOD->GetNumberOfTracks();
379   
380   TExMap *trackMap = new TExMap();//Mapping matrix----
381   
382   for(Int_t i = 0; i < nTracks; i++) {
383     
384     AliAODTrack* aodTrack = dynamic_cast<AliAODTrack *>(fAOD->GetTrack(i)); 
385     
386     if(!aodTrack) {
387       AliError(Form("ERROR: Could not retrieve AODtrack %d",i));
388       continue;
389     }
390     
391     Double_t tpcSignalAll = aodTrack->GetTPCsignal();
392     fTPCSigA->Fill(aodTrack->GetTPCmomentum(),tpcSignalAll);
393     
394     Int_t gID = aodTrack->GetID();
395     
396     if (aodTrack->TestFilterBit(1)) trackMap->Add(gID, i);//Global tracks-----
397   }//1st track loop----
398  
399  AliAODTrack* newAodTrack; 
400   
401   for( Int_t j = 0; j < nTracks; j++ )
402     {
403       
404       AliAODTrack* aodTrack1 = dynamic_cast<AliAODTrack *>(fAOD->GetTrack(j));
405       
406       if(!aodTrack1) {
407         AliError(Form("ERROR: Could not retrieve AODtrack %d",j));
408         continue;
409       }
410       
411       
412       if(!aodTrack1->TestFilterBit(fAODtrackCutBit)) continue;
413       
414       Int_t gID = aodTrack1->GetID();
415       
416       //if( aodTrack->GetID() != (-aodTrack1->GetID() -1) ) continue;
417       newAodTrack = gID >= 0 ? aodTrack1 : fAOD->GetTrack(trackMap->GetValue(-1-gID)); //Take those global track who corresponds to TPC only track
418       
419       Float_t dxy = 0., dz = 0.;
420       
421       dxy = aodTrack1->DCA();
422       dz  = aodTrack1->ZAtDCA();
423       
424       Double_t pt = aodTrack1->Pt();
425       Double_t eta = aodTrack1->Eta();
426       Double_t nclus = aodTrack1->GetTPCClusterInfo(2,1);
427       Double_t chi2ndf = aodTrack1->Chi2perNDF();
428       
429       /* 
430       if( fabs(dxy) > fDCAxy ) continue; 
431       if( fabs(dz) > fDCAz ) continue;
432       //Extra cut on DCA---( Similar to BF Task.. )           
433       if( fDCAxy !=0 && fDCAz !=0 ){
434         if( TMath::Sqrt((dxy*dxy)/(fDCAxy*fDCAxy)+(dz*dz)/(fDCAz*fDCAz)) > 1. ) continue;
435       }
436       */
437       if( pt < fPtLowerLimit || pt > fPtHigherLimit ) continue;
438       if( eta < fEtaLowerLimit || eta > fEtaHigherLimit ) continue;
439       if( nclus < fTPCNClus ) continue;
440       if( chi2ndf > fChi2perNDF ) continue;
441       
442       
443       fHistQA[6]->Fill(dxy);
444       fHistQA[7]->Fill(dz);
445       fHistQA[8]->Fill(pt);
446       fHistQA[9]->Fill(eta);
447       fHistQA[10]->Fill(aodTrack1->Phi());
448       fHistQA[11]->Fill(nclus);
449       fHistQA[12]->Fill(chi2ndf);
450       fHistDCA->Fill(dxy,dz);
451       
452       Short_t gCharge = aodTrack1->Charge();
453       
454       if(gCharge > 0) nPlusCharge += 1.;
455       if(gCharge < 0) nMinusCharge += 1.;
456       
457       if( fUsePid ) {
458         
459         gPid = (Int_t)fParticleSpecies;
460
461         Double_t rap =  newAodTrack->Y(AliPID::ParticleMass(fParticleSpecies));
462         Double_t tpcSignal = newAodTrack->GetTPCsignal();
463         //Double_t rap =  aodTrack1->Y(AliPID::ParticleMass(fParticleSpecies));
464         //Double_t tpcSignal = aodTrack1->GetTPCsignal();
465         
466         if( fabs(rap ) > fRapidityCut ) continue;//Rapidity cut
467         
468         fTPCSig->Fill(newAodTrack->GetTPCmomentum(),tpcSignal);
469         
470         Float_t nsigmaTPCPID = -999.;
471         Float_t nsigmaTOFPID = -999.;
472         //Float_t nsigmaTPCTOFPID = -999.;
473         
474         nsigmaTPCPID = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,fParticleSpecies));
475         nsigmaTOFPID = TMath::Abs(fPIDResponse->NumberOfSigmasTOF(newAodTrack,fParticleSpecies));
476         
477         if ( nsigmaTPCPID < fNSigmaCut  ){
478           
479           if (gCharge > 0) nPartile +=1.;
480           if( gCharge < 0 ) nAntiParticle +=1.;
481           
482         }
483       }//fUsepid-----
484       
485     }//--------- Track Loop to select with filterbit
486   
487   
488   
489   Double_t fContainerCh[3] = { fCentrality, nPlusCharge, nMinusCharge};
490   Double_t fContainerPid[3] = { fCentrality, nPartile, nAntiParticle};
491   
492   
493   fTHnCentNplusNminusCh->Fill(fContainerCh);
494   
495   if( fUsePid ){
496     gPid = (Int_t)fParticleSpecies;
497     fTHnCentNplusNminusPid[gPid]->Fill(fContainerPid);
498     
499     // cout << "nCentrality "<< fCentrality <<", nParticle="<< nPartile << ", nMinusParticle=" << nAntiParticle << endl;
500   }
501   
502   
503   //cout << "nCentrality "<< fCentrality <<", nPlusCharge="<< nPlusCharge << ", nMinusCharge=" << nMinusCharge << endl;
504   
505   fEventCounter->Fill(7);
506   
507   return;
508   
509 }
510 //--------------------------------------------------------------------------------------
511
512
513 //--------------------------------------------------------------------------------------
514 void AliHigherMomentsToyModel::doMCAODEvent(){
515   
516   
517   //---------
518   Double_t nPlusCharge = 0.;
519   Double_t nMinusCharge = 0.;
520
521   Double_t nPlusChargeTruth = 0.;
522   Double_t nMinusChargeTruth = 0.;
523
524   Double_t nPartile = 0.;
525   Double_t nAntiParticle = 0.;
526   Double_t nPartileTruth = 0.;
527   Double_t nAntiParticleTruth = 0.;
528
529   Int_t gPid = 0;
530   Int_t gPdgCode = AliPID::ParticleCode(fParticleSpecies);
531
532   //Connect to Anlaysis manager------
533   AliAnalysisManager* manager = AliAnalysisManager::GetAnalysisManager();
534   if (!manager) {
535     cout<<"ERROR: Analysis manager not found."<<endl;
536     return;
537   }
538   
539   AliAODInputHandler* inputHandler = dynamic_cast<AliAODInputHandler*> (manager->GetInputEventHandler());
540   if (!inputHandler) {
541     cout<<"ERROR: Input handler not found."<<endl;
542     return;
543   }
544   
545   //AOD event------
546   fAOD = dynamic_cast<AliAODEvent*>(InputEvent());
547   if (!fAOD)
548     {
549       cout<< "ERROR 01: AOD not found " <<endl;
550       return;
551     }
552   
553   fPIDResponse =inputHandler->GetPIDResponse();
554   
555   
556   if (!fPIDResponse){
557     AliFatal("This Task needs the PID response attached to the inputHandler");
558     return;
559   }  
560   // -- Get MC header
561   // ------------------------------------------------------------------
562   
563   fArrayMC = dynamic_cast<TClonesArray*>(fAOD->FindListObject(AliAODMCParticle::StdBranchName()));
564   if (!fArrayMC) {
565     AliFatal("Error: MC particles branch not found!\n");
566     return;
567   }
568   
569   AliAODMCHeader *mcHdr=NULL;
570   mcHdr=(AliAODMCHeader*)fAOD->GetList()->FindObject(AliAODMCHeader::StdBranchName());  
571   if(!mcHdr) {
572     printf("MC header branch not found!\n");
573     return;
574   }
575   
576   AliAODHeader *aodHeader = fAOD->GetHeader();
577   
578   fCentrality = (Int_t)aodHeader->GetCentralityP()->GetCentralityPercentile(fCentralityEstimator.Data());
579   
580   /*
581     Int_t cent = -1;
582     cent =  aodHeader->GetCentralityP()->GetCentralityClass10(fCentralityEstimator.Data());
583     
584     if (cent == 0)
585     fCentrality = aodHeader->GetCentralityP()->GetCentralityClass5(fCentralityEstimator.Data());
586     else if (cent == 10 || cent == -1.)
587     fCentrality = -1;
588     else if (cent > 0 && cent < 9)
589     fCentrality = cent + 1;
590   */
591   if( fCentrality < 0 || fCentrality >= 91) return;
592   
593   fEventCounter->Fill(2);
594   
595   
596   
597   if(!ProperVertex()) return;   
598   
599   Int_t nTracks = fAOD->GetNumberOfTracks();
600  
601   
602   TExMap *trackMap = new TExMap();//Mapping matrix----
603   
604   for(Int_t i = 0; i < nTracks; i++) {
605     
606     AliAODTrack* aodTrack = dynamic_cast<AliAODTrack *>(fAOD->GetTrack(i)); 
607     
608     if(!aodTrack) {
609       AliError(Form("ERROR: Could not retrieve AODtrack %d",i));
610       continue;
611     }
612     
613     Double_t tpcSignalAll = aodTrack->GetTPCsignal();
614     fTPCSigA->Fill(aodTrack->GetTPCmomentum(),tpcSignalAll);
615     
616     Int_t gID = aodTrack->GetID();
617     
618     if (aodTrack->TestFilterBit(1)) trackMap->Add(gID, i);//Global tracks-----
619     
620   }//1st track loop----
621   
622   AliAODTrack* newAodTrack; 
623   
624   for( Int_t j = 0; j < nTracks; j++ ){
625     
626     AliAODTrack* aodTrack1 = dynamic_cast<AliAODTrack *>(fAOD->GetTrack(j));
627     
628     if(!aodTrack1) {
629       AliError(Form("ERROR: Could not retrieve AODtrack %d",j));
630       continue;
631     }
632     
633    
634     if(!aodTrack1->TestFilterBit(fAODtrackCutBit)) continue;
635     
636     
637     Int_t gID = aodTrack1->GetID();
638     
639     //if( aodTrack->GetID() != (-aodTrack1->GetID() -1) ) continue;
640     
641     newAodTrack = gID >= 0 ? aodTrack1 : fAOD->GetTrack(trackMap->GetValue(-1-gID)); //Take those global track who corresponds to TPC only track
642     
643   
644     //cout << dxy << endl;
645     Double_t pt = aodTrack1->Pt();
646     Double_t eta = aodTrack1->Eta();
647     Double_t nclus = aodTrack1->GetTPCClusterInfo(2,1);
648     Double_t chi2ndf = aodTrack1->Chi2perNDF();
649   
650     Float_t dxy = 0., dz = 0.;
651     
652     if( fAODtrackCutBit == 128 ){
653      
654       dxy = aodTrack1->DCA();
655       dz  = aodTrack1->ZAtDCA();      
656       
657       if( fabs(dxy) > fDCAxy ) continue; 
658       if( fabs(dz) > fDCAz ) continue;
659       //Extra cut on DCA---( Similar to BF Task.. )           
660       if( fDCAxy !=0 && fDCAz !=0 ){
661         if( TMath::Sqrt((dxy*dxy)/(fDCAxy*fDCAxy)+(dz*dz)/(fDCAz*fDCAz)) > 1. ) continue;
662       }
663
664       fHistQA[6]->Fill(dxy);
665       fHistQA[7]->Fill(dz);
666       fHistDCA->Fill(dxy,dz);
667
668     }
669     
670
671     if( pt < fPtLowerLimit || pt > fPtHigherLimit ) continue;
672     if( eta < fEtaLowerLimit || eta > fEtaHigherLimit ) continue;
673     if( nclus < fTPCNClus ) continue;
674     if( chi2ndf > fChi2perNDF ) continue;
675     
676     
677     
678     fHistQA[8]->Fill(pt);
679     fHistQA[9]->Fill(eta);
680     fHistQA[10]->Fill(aodTrack1->Phi());
681     fHistQA[11]->Fill(nclus);
682     fHistQA[12]->Fill(chi2ndf);
683     
684     
685     Short_t gCharge = aodTrack1->Charge();
686     
687     if( gCharge == 0 ) continue;
688     
689     
690     if(gCharge > 0) nPlusCharge += 1.;
691     if(gCharge < 0) nMinusCharge += 1.;
692     
693     if( fUsePid ) {
694       
695       gPid = (Int_t)fParticleSpecies;;
696       Double_t rap =  newAodTrack->Y(AliPID::ParticleMass(fParticleSpecies));
697       Double_t tpcSignal = newAodTrack->GetTPCsignal();
698    
699       if( fabs(rap ) > fRapidityCut ) continue;//Rapidity cut
700       
701       fTPCSig->Fill(newAodTrack->GetTPCmomentum(),tpcSignal);
702       
703       Float_t nsigmaTPCPID = -999.;
704       Float_t nsigmaTOFPID = -999.;
705       //Float_t nsigmaTPCTOFPID = -999.;
706       
707       nsigmaTPCPID = TMath::Abs(fPIDResponse->NumberOfSigmasTPC(newAodTrack,fParticleSpecies));
708       nsigmaTOFPID = TMath::Abs(fPIDResponse->NumberOfSigmasTOF(newAodTrack,fParticleSpecies));
709       
710       if ( nsigmaTPCPID < fNSigmaCut  ){
711
712         
713         if (gCharge > 0) nPartile +=1.;
714         if( gCharge < 0 ) nAntiParticle +=1.;
715         
716       }//nSigmaCut-----
717     }//fUsepid-----
718    
719   }//--------- Track Loop to select with filterbit
720   
721   //cout << "Cent=" << fCentrality << " MC-PlusChrg=" << nPlusCharge << " MC-MinusChrg=" << nMinusCharge << endl;
722   
723   //cout << "nCentrality "<< fCentrality <<", nParticle="<< nPartile << ", nMinusParticle=" << nAntiParticle << endl; 
724   
725   Double_t fContainerCh[3] = { fCentrality, nPlusCharge, nMinusCharge};
726   Double_t fContainerPid[3] = { fCentrality, nPartile, nAntiParticle};
727   
728   fTHnCentNplusNminusCh->Fill(fContainerCh);  
729   
730   if( fUsePid ){
731
732     gPid = (Int_t)fParticleSpecies;
733     fTHnCentNplusNminusPid[gPid]->Fill(fContainerPid);
734     
735   }
736   
737   //===========================================
738   //--------MC Truth---------------------------
739   //===========================================
740
741   Int_t nMCTrack = fArrayMC->GetEntriesFast();
742   
743   for (Int_t iMC = 0; iMC < nMCTrack; iMC++) {
744     
745     AliAODMCParticle *partMC = (AliAODMCParticle*) fArrayMC->At(iMC);
746     
747     if(!partMC){
748       AliError(Form("ERROR: Could not retrieve AODMCtrack %d",iMC));
749       continue;
750     }
751     
752     if(partMC->Charge() == 0) continue;
753     if(!partMC->IsPhysicalPrimary()) continue;
754     
755     if (partMC->Eta() < fEtaLowerLimit || partMC->Eta() > fEtaHigherLimit) continue;
756     if (partMC->Pt() < fPtLowerLimit ||  partMC->Pt() > fPtHigherLimit) continue;
757     
758     Short_t gCharge = partMC->Charge();
759
760     if(gCharge > 0) nPlusChargeTruth += 1.;
761     if(gCharge < 0) nMinusChargeTruth += 1.;
762    
763     if(fUsePid){
764       
765       Double_t rap =  partMC->Y();
766       
767       if( fabs(rap ) > fRapidityCut ) continue;//Rapidity cut
768       
769       if(TMath::Abs(partMC->GetPdgCode()) != gPdgCode) continue;
770       
771       if(gCharge > 0) nPartileTruth += 1.;
772       if(gCharge < 0) nAntiParticleTruth += 1.;
773       
774     }//if(fUsePid) ----
775   
776   }//MC-Truth Track loop--
777   
778   Double_t fContainerChTruth[3] = { fCentrality, nPlusChargeTruth, nMinusChargeTruth };
779   Double_t fContainerPidTruth[3] = { fCentrality, nPartileTruth, nAntiParticleTruth };
780   
781   //cout << "Cent=" << fCentrality << " MC-PlusChrgT=" << nPlusChargeTruth << " MC-MinusChrgT=" << nMinusChargeTruth << endl;
782   
783   //cout << "nCentrality "<< fCentrality <<", nParticleT="<< nPartileTruth << ", nMinusParticleT=" << nAntiParticleTruth << endl; 
784   
785   //cout << "----------------------------" << endl;
786   fTHnCentNplusNminusChTruth->Fill(fContainerChTruth);  
787   
788   if( fUsePid ){
789     gPid = (Int_t)fParticleSpecies;
790     fTHnCentNplusNminusPidTruth[gPid]->Fill(fContainerPidTruth);
791     
792   }
793   
794   
795   fEventCounter->Fill(7);
796   
797   return;
798   
799 }
800 //---------------------------------------------------------------------------------------
801
802 //---------------------------------------------------------------------------------------
803 Bool_t AliHigherMomentsToyModel::ProperVertex(){
804   
805   Bool_t IsVtx = kFALSE;
806   
807   const AliAODVertex *vertex = fAOD->GetPrimaryVertex();
808   
809   if(vertex) {
810     Double32_t fCov[6];
811     vertex->GetCovarianceMatrix(fCov);
812     if(vertex->GetNContributors() > 0) {
813       if(fCov[5] != 0) {
814         
815         Double_t lvx = vertex->GetX();
816         Double_t lvy = vertex->GetY();
817         Double_t lvz = vertex->GetZ();
818         
819         fEventCounter->Fill(5);
820         
821         fHistQA[3]->Fill(lvx);fHistQA[4]->Fill(lvy);fHistQA[5]->Fill(lvz);
822         
823         if(TMath::Abs(lvx) < fVxMax) {
824           if(TMath::Abs(lvy) < fVyMax) {
825             if(TMath::Abs(lvz) < fVzMax) {
826               
827               fEventCounter->Fill(6);
828               fHistQA[0]->Fill(lvx);fHistQA[1]->Fill(lvy);fHistQA[2]->Fill(lvz);
829               
830               IsVtx = kTRUE; 
831               
832             }//Z-Vertex cut---
833           }//Y-vertex cut--
834         }//X-vertex cut---
835       }//Covariance------
836     }//Contributors check---
837   }//If vertex-----------
838   
839   return IsVtx;
840 }
841 //---------------------------------------------------------------------------------------
842
843
844
845 //------------------------------------------------------------------------
846 void AliHigherMomentsToyModel::Terminate( Option_t * ){
847   
848   Info("AliHigherMomentsToyModel"," Task Successfully finished");
849   
850 }
851
852 //------------------------------------------------------------------------