]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/SPECTRA/PiKaPr/TPCTOF/AliAnalysisCombinedHadronSpectra.cxx
a5670aec07cb4a25654650cb2585bfae32c1eb8b
[u/mrichter/AliRoot.git] / PWGLF / SPECTRA / PiKaPr / TPCTOF / AliAnalysisCombinedHadronSpectra.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
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, proviyaded 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 purapose. It is         *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 ///////////////////////////////////////////////////////////////////////////
17 //                                                                       //
18 //                                                                       //
19 // Analysis for identified charged hadron spectra.                       //
20 //                                                                       //
21 //                                                                       //
22 ///////////////////////////////////////////////////////////////////////////
23
24 #include "Riostream.h"
25 #include "TChain.h"
26 #include "TTree.h"
27 #include "TH1F.h"
28 #include "TH2D.h"
29 #include "TH3F.h"
30 #include "TList.h"
31 #include "TMath.h"
32 #include "TCanvas.h"
33 #include "TObjArray.h"
34 #include "TF1.h"
35 #include "TFile.h"
36
37 #include "AliAnalysisTaskSE.h"
38 #include "AliAnalysisManager.h"
39
40 #include "AliHeader.h"
41 #include "AliGenPythiaEventHeader.h"
42 #include "AliGenCocktailEventHeader.h"
43
44 #include "AliPID.h"
45 #include "AliESDtrackCuts.h"
46 #include "AliESDVertex.h"
47 #include "AliESDEvent.h"
48 #include "AliESDInputHandler.h"
49 #include "AliESDtrack.h"
50 #include "AliESDpid.h"
51 #include "AliCentrality.h"
52 #include "AliESDUtils.h"
53 #include "AliMultiplicity.h"
54
55 #include "AliMCEventHandler.h"
56 #include "AliMCEvent.h"
57 #include "AliStack.h"
58
59 #include "AliLog.h"
60
61 #include "AliAnalysisCombinedHadronSpectra.h"
62
63
64 ClassImp(AliAnalysisCombinedHadronSpectra)
65
66 //________________________________________________________________________
67 AliAnalysisCombinedHadronSpectra::AliAnalysisCombinedHadronSpectra() 
68   : AliAnalysisTaskSE("TaskChargedHadron"), fESD(0), fListHist(0), fESDtrackCuts(0),fESDTrackCutsMult(0),fESDpid(0),
69     fMCtrue(0),
70     fOnlyQA(0),
71     fUseHBTmultiplicity(0),
72     fUseTPConlyTracks(0),
73     fSaveMotherPDG(0),
74     fAlephParameters(),
75     fHistRealTracks(0),
76     fHistMCparticles(0),
77     fHistPidQA(0),
78     fHistMult(0),
79     fHistCentrality(0)
80 {
81   // default Constructor
82   /* fast compilation test
83      gSystem->Load("libANALYSIS");
84      gSystem->Load("libANALYSISalice");
85      .L /d/alice09/akalweit/train/trunk/akalweit_hadronSpectra/AliAnalysisCombinedHadronSpectra.cxx++
86    */
87 }
88
89
90 //________________________________________________________________________
91 AliAnalysisCombinedHadronSpectra::AliAnalysisCombinedHadronSpectra(const char *name) 
92   : AliAnalysisTaskSE(name), fESD(0), fListHist(0), fESDtrackCuts(0),fESDTrackCutsMult(0),fESDpid(0),
93     fMCtrue(0),
94     fOnlyQA(0),
95     fUseHBTmultiplicity(0),
96     fUseTPConlyTracks(0),
97     fSaveMotherPDG(0),
98     fAlephParameters(),
99     fHistRealTracks(0),
100     fHistMCparticles(0),
101     fHistPidQA(0),
102     fHistMult(0),
103     fHistCentrality(0)
104 {
105   //
106   // standard constructur which should be used
107   //
108   Printf("*** CONSTRUCTOR CALLED ****");
109   //
110   fMCtrue = kTRUE; 
111   fOnlyQA = kFALSE;
112   fUseHBTmultiplicity = kTRUE;
113   fUseTPConlyTracks = kFALSE;
114   /* real */
115   fAlephParameters[0] = 0.0283086;
116   fAlephParameters[1] = 2.63394e+01;
117   fAlephParameters[2] = 5.04114e-11;
118   fAlephParameters[3] = 2.12543e+00;
119   fAlephParameters[4] = 4.88663e+00;
120   //
121   // initialize PID object
122   //
123   //fESDpid = new AliESDpid();
124   //
125   // create track cuts
126   //
127   fESDtrackCuts = new AliESDtrackCuts("AliESDtrackCuts","AliESDtrackCuts");
128   //
129   //Initialize();
130   // Output slot #0 writes into a TList container
131   DefineOutput(1, TList::Class());
132
133 }
134
135
136 //________________________________________________________________________
137 void AliAnalysisCombinedHadronSpectra::Initialize()
138 {
139   //
140   // updating parameters in case of changes
141   //
142   // 1. pid parameters
143   //
144   //fESDpid->GetTPCResponse().SetBetheBlochParameters(fAlephParameters[0],fAlephParameters[1],fAlephParameters[2],fAlephParameters[3],fAlephParameters[4]);
145   //
146   // 2. track cuts
147   //
148   /*
149   fESDtrackCuts->SetMaxCovDiagonalElements(2, 2, 0.5, 0.5, 2);  // BEWARE STANDARD VALUES ARE: 2, 2, 0.5, 0.5, 2
150   fESDtrackCuts->SetMaxNsigmaToVertex(3);
151   fESDtrackCuts->SetRequireSigmaToVertex(kTRUE);
152   fESDtrackCuts->SetAcceptKinkDaughters(kFALSE);
153   fESDtrackCuts->SetMinNClustersTPC(70);
154   fESDtrackCuts->SetMaxChi2PerClusterTPC(4);
155   fESDtrackCuts->SetMaxDCAToVertexXY(3);
156   fESDtrackCuts->SetMaxDCAToVertexZ(3);
157   fESDtrackCuts->SetRequireTPCRefit(kTRUE);
158   fESDtrackCuts->SetRequireITSRefit(kTRUE);
159   fESDtrackCuts->SetClusterRequirementITS(AliESDtrackCuts::kSPD, AliESDtrackCuts::kAny); //TEMPORARY <-> REMOVE
160   fESDtrackCuts->SetMinNClustersITS(3);
161   */
162   //fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kTRUE); // kTRUE = sel. primaries --> patch for the moment, do TFractionFitter later
163
164
165   if (!fUseTPConlyTracks) {
166           fESDtrackCuts = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kFALSE);
167           fESDtrackCuts->SetMaxDCAToVertexXY(3);
168           fESDtrackCuts->SetMaxDCAToVertexZ(2);
169           fESDtrackCuts->SetEtaRange(-0.9,0.9);
170   }
171   else {
172           //fESDtrackCuts = AliESDtrackCuts::GetStandardTPCOnlyTrackCuts();
173           fESDtrackCuts->SetMinNClustersTPC(70);
174           fESDtrackCuts->SetMaxChi2PerClusterTPC(4);
175           fESDtrackCuts->SetAcceptKinkDaughters(kFALSE);
176           fESDtrackCuts->SetRequireTPCRefit(kFALSE);
177
178           fESDtrackCuts->SetMaxDCAToVertexXY(15);
179           fESDtrackCuts->SetMaxDCAToVertexZ(6);
180           fESDtrackCuts->SetDCAToVertex2D(kFALSE);
181           fESDtrackCuts->SetRequireSigmaToVertex(kFALSE);
182
183           fESDtrackCuts->SetEtaRange(-0.9,0.9);
184   }
185   //
186   //
187   //
188   //
189   fESDTrackCutsMult = AliESDtrackCuts::GetStandardITSTPCTrackCuts2010(kTRUE);
190   fESDTrackCutsMult->SetEtaRange(-1.2,+1.2);
191   fESDTrackCutsMult->SetPtRange(0.15,1e10);
192
193 }
194
195
196 //________________________________________________________________________
197 void AliAnalysisCombinedHadronSpectra::UserCreateOutputObjects() 
198 {
199   // Create histograms
200   // Called once
201   fListHist = new TList();
202   fListHist->SetOwner(kTRUE);
203   //
204   const Int_t kPtBins = 35;
205   const Int_t kMultBins = 11;
206   const Int_t kDcaBins = 76;
207   const Float_t kDcaBinsTPConlyFactor = 5; //need to change binning of DCA plot for tpconly
208   // sort pT-bins ..
209   Double_t binsPt[kPtBins+1] = {0., 0.05, 0.1, 0.15, 0.2, 0.25, 0.3, 0.35, 0.4, 0.45, 0.5, 0.55, 0.6, 0.65, 0.7, 0.75, 0.8, 0.85, 0.9, 0.95, 1.0, 1.1, 1.2, 1.3, 1.4, 1.5, 1.6, 1.7, 1.8, 1.9, 2.0, 2.2, 2.4, 2.6, 2.8, 3.0};
210
211   Double_t binsDca[kDcaBins+1] =  {-3,-2.85,-2.7,-2.55,-2.4,-2.25,-2.1,-1.95,-1.8,-1.65,-1.5,-1.35,-1.2,-1.05,-0.9,-0.75,-0.6,-0.45,-0.3,-0.285,-0.27,-0.255,-0.24,-0.225,-0.21,-0.195,-0.18,-0.165,-0.15,-0.135,-0.12,-0.105,-0.09,-0.075,-0.06,-0.045,-0.03,-0.015,0,0.015,0.03,0.045,0.06,0.075,0.09,0.105,0.12,0.135,0.15,0.165,0.18,0.195,0.21,0.225,0.24,0.255,0.27,0.285,0.3,0.45,0.6,0.75,0.9,1.05,1.2,1.35,1.5,1.65,1.8,1.95,2.1,2.25,2.4,2.55,2.7,2.85,3};
212
213   // DCA bins borders get multiplied by constant factor for TPConlyTracks
214   Double_t binsDcaTPConly[kDcaBins+1];
215   for (Int_t i = 0; i< kDcaBins+1; i++) {
216         binsDcaTPConly[i] = kDcaBinsTPConlyFactor * binsDca[i];
217   }
218
219   //
220   // create the histograms with all necessary information --> it is filled 4x for each particle assumption
221   //
222   // (0.) assumed particle: 0. pion, 1. kaon, 2. proton, 3. deuteron
223   // (1.) multiplicity or centrality -- number of accepted ESD tracks per events (deprecated), but now classes from 1 to 10, 0: Min. Bias
224   // (2.) pT
225   // (3.) sign
226   // (4.) rapidity --> filled 4x
227   // (5.)  pull TPC dEx --> filled 4x
228   // (6.) has valid TOF pid signal
229   // (7.) nsigma TOF --> filled 4x
230   // (8..) dca_xy
231   // (9.) CODE -- only MC 0-generated, 1-true rec. primaries, 2-misident prim, 3-second weak, 4-second material, 5-misident sec, 6-sec. K0, 7-sec. lambda, 8-sec sigma+
232   //
233   //                              0,           1,           2,  3,       4,   5,    6,   7,       8
234   Int_t    binsHistReal[9] = {   3,   kMultBins,     kPtBins,  2,      10,   50,    2,  80, kDcaBins};
235   Double_t xminHistReal[9] = {-0.5,        -0.5,           0, -2,    -0.5,   -5,- 0.5,  -8,       -3};
236   Double_t xmaxHistReal[9] = { 2.5,        10.5,           3,  2,     0.5,    5,  1.5,   8,        3};
237   fHistRealTracks = new THnSparseF("fHistRealTracks","real tracks",9,binsHistReal,xminHistReal,xmaxHistReal);
238   //
239   fHistRealTracks->GetAxis(2)->Set(kPtBins, binsPt);
240
241   //different DCAxy binning for TPConlyTracks
242   if (!fUseTPConlyTracks) fHistRealTracks->GetAxis(8)->Set(kDcaBins, binsDca);
243   else  fHistRealTracks->GetAxis(8)->Set(kDcaBins, binsDcaTPConly);
244   fListHist->Add(fHistRealTracks);
245   //
246   //                      0.ptot,1.tpcSig,2.hasTOF, 3. assumed part., 4. nclDedx, 5. nSigmaTPC (4x), 6. nSigmaTOF (4x), 7. centrality
247   fHistPidQA = new TH3D("fHistPidQA","PID QA",500,0.1,10,1000,0,1000,2,-2,2);
248   BinLogAxis(fHistPidQA);
249   fListHist->Add(fHistPidQA);
250   
251   //                            0,            1,           2,  3,      4,   5,    6,   7,        8,    9
252   Int_t    binsHistMC[10] = {   3,    kMultBins,     kPtBins,  2,     10,  50,    2,  80, kDcaBins,    6};
253   Double_t xminHistMC[10] = {-0.5,         -0.5,           0, -2,   -0.5,  -5,- 0.5,  -8,       -3, -0.5};
254   Double_t xmaxHistMC[10] = { 2.5,         10.5,           3,  2,    0.5,   5,  1.5,   8,        3,  5.5};
255
256   //different binning for CODE axis, if we want to save motherPDG
257   if (fSaveMotherPDG) {
258     binsHistMC[9] = 9;
259     xmaxHistMC[9] = 8.5;
260   }
261
262   fHistMCparticles = new THnSparseF("fHistMCparticles","MC histogram",10,binsHistMC,xminHistMC,xmaxHistMC);
263   fHistMCparticles->GetAxis(2)->Set(kPtBins, binsPt);
264
265   //different DCAxy binning for TPConlyTracks
266   if (!fUseTPConlyTracks) fHistMCparticles->GetAxis(8)->Set(kDcaBins, binsDca);
267   else fHistMCparticles->GetAxis(8)->Set(kDcaBins, binsDcaTPConly);
268
269   fListHist->Add(fHistMCparticles);
270   //
271   fHistMult = new TH2D("fHistMult", "control histogram to count number of events", 502, -2.5, 499.5,4,-0.5,3.5);
272   fHistCentrality = new TH1D("fHistCentrality", "control histogram to count number of events", 22, -1.5, 20.5);
273   fListHist->Add(fHistMult);
274   fListHist->Add(fHistCentrality);
275   
276 }
277
278 //________________________________________________________________________
279 void AliAnalysisCombinedHadronSpectra::UserExec(Option_t *) 
280 {
281   //
282   // main event loop
283   //
284   if (!fESDpid) fESDpid = ((AliESDInputHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->GetESDpid();
285   if (!fESDpid) {
286     fESDpid = new AliESDpid(); // HACK FOR MC PBPB --> PLEASE REMOVE AS SOON AS POSSIBLE
287     fESDpid->GetTPCResponse().SetBetheBlochParameters(1.28778e+00/50., 3.13539e+01, TMath::Exp(-3.16327e+01), 1.87901e+00, 6.41583e+00);
288   }
289   //AliLog::SetGlobalLogLevel(AliLog::kError);
290   //
291   // Check Monte Carlo information and other access first:
292   //
293   AliMCEventHandler* eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
294   if (!eventHandler) {
295     //Printf("ERROR: Could not retrieve MC event handler");
296     fMCtrue = kFALSE;
297   }
298   //
299   AliMCEvent* mcEvent = 0x0;
300   AliStack* stack = 0x0;
301   if (eventHandler) mcEvent = eventHandler->MCEvent();
302   if (!mcEvent) {
303     //Printf("ERROR: Could not retrieve MC event");
304     if (fMCtrue) return;
305   }
306   if (fMCtrue) {
307     stack = mcEvent->Stack();
308     if (!stack) return;
309   }
310   //
311   fESD = dynamic_cast<AliESDEvent*>( InputEvent() );
312   if (!fESD) {
313     //Printf("ERROR: fESD not available");
314     return;
315   }
316   
317   if (!fESDtrackCuts) {
318     Printf("ERROR: fESDtrackCuts not available");
319     return;
320   }
321   //
322   // check if event is selected by physics selection class
323   //
324   Bool_t isSelected = kTRUE; // for reasons of backward compatibility --> check is now in AddTask macro
325   /*
326   Bool_t isSelected = kFALSE;
327   isSelected = (((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected() & AliVEvent::kMB) == AliVEvent::kMB;
328   */
329   //
330   // monitor vertex position
331   //
332   const AliESDVertex *vertex = fESD->GetPrimaryVertexTracks();
333   if(vertex->GetNContributors()<1) {
334     // SPD vertex
335     vertex = fESD->GetPrimaryVertexSPD();
336     if(vertex->GetNContributors()<1) vertex = 0x0;
337   }  
338   //
339   // small track loop to determine trigger Pt, multiplicity or centrality
340   //
341   Double_t triggerPt = 0;
342   for (Int_t i=0;i<fESD->GetNumberOfTracks();++i) {
343     AliESDtrack *track =fESD->GetTrack(i);
344     if (!fESDTrackCutsMult->AcceptTrack(track)) continue;
345     if (track->Pt() > triggerPt) triggerPt = track->Pt();
346   }
347   Int_t trackCounter = fESDTrackCutsMult->CountAcceptedTracks(fESD);
348   //
349   // 2nd multiplicity estimator SPD hits; replaces multiplicity from HBT
350   //
351   Float_t spdCorr = -1;
352   const AliMultiplicity *mult = fESD->GetMultiplicity();
353   Float_t nClusters[6]={0.0,0.0,0.0,0.0,0.0,0.0};
354   for(Int_t ilay=0; ilay<6; ilay++) nClusters[ilay] = (Float_t)mult->GetNumberOfITSClusters(ilay);
355   if (vertex) spdCorr = AliESDUtils::GetCorrSPD2(nClusters[1],vertex->GetZ());
356   //
357   Float_t centrality = -1;
358   //
359   // IMPORTANT CENTRALITY DEFINITION FOR pp
360   //
361   if (!fUseHBTmultiplicity) {
362     // test binning to check if multiplicity-dependence is due to diffractive events being s
363     if (trackCounter >= 0  && trackCounter <= 0)  centrality = 0;
364     if (trackCounter >= 1  && trackCounter <= 1)  centrality = 1;
365     if (trackCounter >= 2  && trackCounter <= 2)  centrality = 2;
366     if (trackCounter >= 3  && trackCounter <= 3)  centrality = 3;
367     if (trackCounter >= 4  && trackCounter <= 4)  centrality = 4;
368     if (trackCounter >= 5  && trackCounter <= 5)  centrality = 5;
369     // this was all the original bin 1 being [0..5]
370     if (trackCounter >= 6  && trackCounter <= 9)  centrality = 6;
371     if (trackCounter >= 10 && trackCounter <= 14) centrality = 7;
372     if (trackCounter >= 15 && trackCounter <= 22) centrality = 8;
373     if (trackCounter >= 23 && trackCounter <= 32) centrality = 9;
374     if (trackCounter >= 33 && trackCounter <= 42) centrality = 10;
375     /*
376     if (trackCounter >= 43 && trackCounter <= 52) centrality = 7
377     if (trackCounter >= 53 && trackCounter <= 62) centrality = 8;
378     if (trackCounter >= 63 && trackCounter <= 72) centrality = 9;
379     if (trackCounter >= 73 && trackCounter <= 82) centrality = 10;
380     */
381   } else {
382     if (spdCorr >= 0  && spdCorr <=  2)  centrality  = 0;
383     if (spdCorr >= 3  && spdCorr <=  5)  centrality  = 1;
384     if (spdCorr >= 6  && spdCorr <=  8)  centrality  = 2;
385     if (spdCorr >= 9  && spdCorr <= 11)  centrality  = 3;
386     if (spdCorr >= 12 && spdCorr <= 14)  centrality  = 4;
387     if (spdCorr >= 15 && spdCorr <= 16)  centrality  = 5;
388     // this was all the original bin 1 being [0..16]
389     if (spdCorr >= 17 && spdCorr <= 30)  centrality =  6;
390     if (spdCorr >= 31 && spdCorr <= 45)  centrality =  7;
391     if (spdCorr >= 46 && spdCorr <= 68)  centrality =  8;
392     if (spdCorr >= 69 && spdCorr <= 97)  centrality =  9;
393     if (spdCorr >= 98)                   centrality = 10;
394     /*
395     if (spdCorr >= 17 && spdCorr <= 30)  centrality = 2;
396     if (spdCorr >= 31 && spdCorr <= 45)  centrality = 3;
397     if (spdCorr >= 46 && spdCorr <= 68)  centrality = 4;
398     if (spdCorr >= 69 && spdCorr <= 97)  centrality = 5;
399     if (spdCorr >= 98)                   centrality = 6;
400     */
401   }
402   //
403   Int_t rootS = fESD->GetBeamEnergy() < 1000 ? 0 : 1;
404   if (fESD->GetEventSpecie() == 4) { // PbPb
405     rootS = 2;
406     AliCentrality *esdCentrality = fESD->GetCentrality();
407     centrality = esdCentrality->GetCentralityClass10("V0M") + 1; // centrality percentile determined with V0
408     if (TMath::Abs(centrality - 1) < 1e-5) {
409       centrality = esdCentrality->GetCentralityClass5("V0M");
410     }
411   }
412   Int_t nContributors = 0;
413   if (fESD->GetPrimaryVertexTPC()) nContributors = fESD->GetPrimaryVertexTPC()->GetNContributors();
414   //
415   Int_t processtype = 0;
416   Int_t processCode = 0;
417   //
418   // important change: fill generated only after vertex cut in case of heavy-ions
419   //
420   if (!vertex && fESD->GetEventSpecie() == 4) {
421     fHistMult->Fill(-1, processCode);
422     PostData(1, fListHist);
423     return;
424   } else {
425     if (TMath::Abs(vertex->GetZv()) > 10 && fESD->GetEventSpecie() == 4) {
426       fHistMult->Fill(-1, processCode);
427       PostData(1, fListHist);
428       return;
429     }
430   }
431   //
432   if (fMCtrue) {
433     //
434     //
435     //
436     AliHeader * header = mcEvent->Header();
437     processtype = GetPythiaEventProcessType(header);
438     // non diffractive
439     if (processtype !=92 && processtype !=93 && processtype != 94) processCode = 1;
440     // single diffractive
441     if ((processtype == 92 || processtype == 93)) processCode = 2;
442     // double diffractive
443     if (processtype == 94) processCode = 3;
444     //
445     for(Int_t i = 0; i < stack->GetNtrack(); i++) {
446       TParticle * trackMC = stack->Particle(i);
447       Int_t pdg = trackMC->GetPdgCode();
448       //
449       Double_t xv = trackMC->Vx();
450       Double_t yv = trackMC->Vy();
451       Double_t zv = trackMC->Vz();
452       Double_t dxy = 0;
453       dxy = TMath::Sqrt(xv*xv + yv*yv); // so stupid to avoid warnings
454       Double_t dz = 0;
455       dz = TMath::Abs(zv); // so stupid to avoid warnings
456       //
457       // vertex cut - selection of primaries
458       //
459       //if (dxy > 3 || dz > 10) continue; // fixed cut at 3cm in r
460       //
461       if (!stack->IsPhysicalPrimary(i)) continue;
462       //
463       // fill MC histograms here...
464       // 
465       Double_t rap = trackMC->Y();
466       Double_t pT  = trackMC->Pt();
467       Int_t sign = pdg < 0 ? -1 : 1; // only works for charged pi,K,p !!
468 //      Double_t transMass = TMath::Sqrt(trackMC->Pt()*trackMC->Pt() + trackMC->GetMass()*trackMC->GetMass()) - trackMC->GetMass();
469       //
470       Int_t iPart = -1;
471       if (TMath::Abs(pdg) == 211)  iPart = 0; // select Pi+/Pi- only
472       if (TMath::Abs(pdg) == 321)  iPart = 1; // select K+/K- only
473       if (TMath::Abs(pdg) == 2212) iPart = 2; // select p+/p- only
474       if (iPart == -1) continue;
475       //
476       Double_t vecHistMC[10] = {iPart, centrality,  pT, sign, rap, 0, 1, 0, dxy, 0};
477       if (!fOnlyQA) fHistMCparticles->Fill(vecHistMC);
478     }
479   }
480   //
481   if (!isSelected && !fOnlyQA) {
482     PostData(1, fListHist);
483     return;
484   }
485   //
486   if (!vertex) {
487     fHistMult->Fill(-1, processCode);
488     PostData(1, fListHist);
489     return;
490   } else {
491     if (TMath::Abs(vertex->GetZv()) > 10) {
492       fHistMult->Fill(-1, processCode);
493       PostData(1, fListHist);
494       return;
495     }
496   }
497   //
498   // count events after physics selection and after vertex selection
499   //
500   //cout << "MULTIPLICITY " << trackCounter << " " << fESD->GetEventNumberInFile() <<endl;
501   fHistMult->Fill(trackCounter, processCode);
502   fHistCentrality->Fill(centrality);
503
504
505   //***************************************************
506   // track loop
507   //***************************************************
508   //const Float_t kNsigmaCut = 3;
509   //const Float_t k2sigmaCorr = 1/(0.5*(TMath::Erf(kNsigmaCut/sqrt(2))-TMath::Erf(-kNsigmaCut/sqrt(2))))/*1/0.9545*/;
510   //
511   Float_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z for the vertex cut
512   //
513   for (Int_t i=0;i<fESD->GetNumberOfTracks();++i) {
514         
515         AliESDtrack *track = 0;
516     AliESDtrack *trackForTOF = 0; //normal track for all TOF information needed when using tpconly-tracks
517         
518         //normal tracks, if tpconly flag is set, use tpconlytracks
519         if (!fUseTPConlyTracks){
520         track =fESD->GetTrack(i); 
521         }
522         else {
523                 track = fESDtrackCuts->GetTPCOnlyTrack(fESD,i);
524                 if (!track) continue;
525                 trackForTOF = fESD->GetTrack(i);
526         }
527     //
528     if (!track->GetInnerParam()) {
529                 if (fUseTPConlyTracks) {delete track; track = 0;} //need to delete tpconlytrack
530                 continue;
531         }
532     Double_t ptot = track->GetInnerParam()->GetP(); // momentum for dEdx determination
533     Double_t pT = track->Pt();
534     track->GetImpactParameters(dca, cov);
535     //
536     //
537     // cut for dead regions in the detector
538     // if (track->Eta() > 0.1 && (track->Eta() < 0.2 && track->Phi() > 0.1 && track->Phi() < 0.1) continue;
539     //
540     // 2.a) apply some standard track cuts according to general recommendations
541     //
542     if (!fESDtrackCuts->AcceptTrack(track)) {
543                 if (fUseTPConlyTracks) {delete track; track = 0;} //need to delete tpconlytrack
544                 continue;
545         }
546
547         UInt_t status = 0;
548     if (!fUseTPConlyTracks) status = track->GetStatus();
549         else status = trackForTOF->GetStatus();
550     Bool_t hasTOFout  = status&AliESDtrack::kTOFout; 
551     Bool_t hasTOFtime = status&AliESDtrack::kTIME;
552     Bool_t hasTOFpid  = status&AliESDtrack::kTOFpid;
553     Bool_t hasTOF     = kFALSE;
554     if (hasTOFout && hasTOFtime && hasTOFpid) hasTOF = kTRUE;
555         Float_t length = 0.;
556     if (!fUseTPConlyTracks) length = track->GetIntegratedLength(); 
557     else length = trackForTOF->GetIntegratedLength();
558
559     if (length < 350.) hasTOF = kFALSE;
560     //
561     // calculate rapidities and kinematics
562     // 
563     //
564     Double_t pvec[3];
565     track->GetPxPyPz(pvec);
566     Double_t energyPion = TMath::Sqrt(track->GetP()*track->GetP() + AliPID::ParticleMass(AliPID::kPion)*AliPID::ParticleMass(AliPID::kPion));
567     Double_t energyKaon = TMath::Sqrt(track->GetP()*track->GetP() + AliPID::ParticleMass(AliPID::kKaon)*AliPID::ParticleMass(AliPID::kKaon));
568     Double_t energyProton = TMath::Sqrt(track->GetP()*track->GetP() + AliPID::ParticleMass(AliPID::kProton)*AliPID::ParticleMass(AliPID::kProton));
569     Double_t energyDeuteron = TMath::Sqrt(track->GetP()*track->GetP() + 4*AliPID::ParticleMass(AliPID::kProton)*AliPID::ParticleMass(AliPID::kProton));
570     //
571     Double_t rapPion = 0.5*TMath::Log((energyPion + pvec[2])/(energyPion - pvec[2]));
572     Double_t rapKaon = 0.5*TMath::Log((energyKaon + pvec[2])/(energyKaon - pvec[2]));
573     Double_t rapProton = 0.5*TMath::Log((energyProton + pvec[2])/(energyProton - pvec[2]));
574     Double_t rapDeuteron = 0.5*TMath::Log((energyDeuteron + pvec[2])/(energyDeuteron - pvec[2]));
575     //
576 //    Double_t transMassPion = TMath::Sqrt(track->Pt()*track->Pt() + AliPID::ParticleMass(AliPID::kPion)*AliPID::ParticleMass(AliPID::kPion))      -  AliPID::ParticleMass(AliPID::kPion);
577 //    Double_t transMassKaon = TMath::Sqrt(track->Pt()*track->Pt() + AliPID::ParticleMass(AliPID::kKaon)*AliPID::ParticleMass(AliPID::kKaon))     -  AliPID::ParticleMass(AliPID::kKaon);
578  //   Double_t transMassProton = TMath::Sqrt(track->Pt()*track->Pt() + AliPID::ParticleMass(AliPID::kProton)*AliPID::ParticleMass(AliPID::kProton))    -  AliPID::ParticleMass(AliPID::kProton);
579 //    Double_t transMassDeuteron = TMath::Sqrt(track->Pt()*track->Pt() + 4*AliPID::ParticleMass(AliPID::kProton)*AliPID::ParticleMass(AliPID::kProton))    -  2*AliPID::ParticleMass(AliPID::kProton);
580     //
581     // 3. make the PID
582     //
583     Double_t sign = track->GetSign();   
584     Double_t tpcSignal = track->GetTPCsignal();
585     //
586     //
587     // 3.a. calculate expected signals in nsigma
588     //
589     //  
590     // (0.) assumed particle: 0. pion, 1. kaon, 2. proton, 3. deuteron
591     // (1.) multiplicity or centrality -- number of accepted ESD tracks per events (deprecated), but now classes from 1 to 10, 0: Min. Bias
592     // (2.) pT
593     // (3.) sign
594     // (4.) rapidity --> filled 4x
595     // (5.)  pull TPC dEx --> filled 4x
596     // (6.) has valid TOF pid signal
597     // (7.) nsigma TOF --> filled 4x
598     // (8..) dca_xy
599     // (9.) CODE -- only MC 0-generated, 1-true rec. primaries, 2-misident, 3-second weak, 4-second material
600     //
601 //    Double_t transMass[4] = {transMassPion,transMassKaon,transMassProton,transMassDeuteron};
602     Double_t rap[4] = {rapPion,rapKaon,rapProton,rapDeuteron};
603     Double_t pullsTPC[4] = {fESDpid->NumberOfSigmasTPC(track,AliPID::kPion),
604                             fESDpid->NumberOfSigmasTPC(track,AliPID::kKaon),
605                             fESDpid->NumberOfSigmasTPC(track,AliPID::kProton),
606                             0}; // ASK FOR PUTTING THE DEUTERON TO AliPID !!!!!!!!!!!!!!
607     Float_t time0 = fESDpid->GetTOFResponse().GetTimeZero();
608     //fESDpid->GetTOFResponse().SetTimeResolution(130.);
609     Double_t pullsTOF[4] ={0.,0.,0.,0.};
610     if (!fUseTPConlyTracks) {
611                          pullsTOF[0] = fESDpid->NumberOfSigmasTOF(track,AliPID::kPion, time0);
612                              pullsTOF[1] = fESDpid->NumberOfSigmasTOF(track,AliPID::kKaon, time0);
613                              pullsTOF[2] = fESDpid->NumberOfSigmasTOF(track,AliPID::kProton, time0);
614                              pullsTOF[3] = 0; // ASK FOR PUTTING THE DEUTERON TO AliPID !!!!!!!!!!!!!!;
615         }
616         else {
617                          pullsTOF[0] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kPion, time0);
618                              pullsTOF[1] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kKaon, time0);
619                              pullsTOF[2] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kProton, time0);
620                              pullsTOF[3] = 0; // ASK FOR PUTTING THE DEUTERON TO AliPID !!!!!!!!!!!!!!;
621         }
622
623     //
624 //    Double_t tpcQA[4] = {fESDpid->NumberOfSigmasTPC(track,AliPID::kElectron),
625 //                       fESDpid->NumberOfSigmasTPC(track,AliPID::kPion),
626 //                       fESDpid->NumberOfSigmasTPC(track,AliPID::kKaon),
627 //                       fESDpid->NumberOfSigmasTPC(track,AliPID::kProton)};
628
629     Double_t tofQA[4] = {0.,0.,0.,0.}; 
630     if (!fUseTPConlyTracks) {
631                  tofQA[0] = fESDpid->NumberOfSigmasTOF(track,AliPID::kElectron, time0);
632                          tofQA[0] = fESDpid->NumberOfSigmasTOF(track,AliPID::kPion, time0);
633                  tofQA[0] = fESDpid->NumberOfSigmasTOF(track,AliPID::kKaon, time0);
634                  tofQA[0] = fESDpid->NumberOfSigmasTOF(track,AliPID::kProton, time0);
635         }
636         else{
637                  tofQA[0] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kElectron, time0);
638                          tofQA[0] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kPion, time0);
639                  tofQA[0] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kKaon, time0);
640                  tofQA[0] = fESDpid->NumberOfSigmasTOF(trackForTOF,AliPID::kProton, time0);
641         }
642
643     //
644     for(Int_t iPart = 0; iPart < 3; iPart++) { // loop over assumed particle type
645       //                              0,           1,    2,    3,           4,               5,      6,              7,     8
646       Double_t vecHistReal[9]  = {iPart,  centrality,   pT, sign,  rap[iPart], pullsTPC[iPart], hasTOF, pullsTOF[iPart], dca[0]};
647       if (!fOnlyQA) fHistRealTracks->Fill(vecHistReal);
648       //
649       // using MC truth for precise efficiencies...
650       //
651       if (fMCtrue && !fOnlyQA) {
652         Int_t code = 9; // code: 0-generated, 1-true rec. primaries, 2-misident, 3-second weak, 4-second material
653         Int_t assumedPdg = 0;//2212(proton); 321(Kaon); 211(pion);
654         Int_t motherCode = -1;
655         if (iPart == 0) assumedPdg = 211;
656         if (iPart == 1) assumedPdg = 321;
657         if (iPart == 2) assumedPdg = 2212;
658         //
659         //
660         TParticle *trackMC = stack->Particle(TMath::Abs(track->GetLabel()));
661         Int_t pdg = TMath::Abs(trackMC->GetPdgCode());
662         //
663         if (pdg != assumedPdg && stack->IsPhysicalPrimary(TMath::Abs(track->GetLabel()))) code = 2;
664         if (pdg != assumedPdg && stack->IsSecondaryFromWeakDecay(TMath::Abs(track->GetLabel()))) code = 5;
665         if (pdg == assumedPdg && stack->IsPhysicalPrimary(TMath::Abs(track->GetLabel()))) code = 1;
666         if (pdg == assumedPdg && stack->IsSecondaryFromWeakDecay(TMath::Abs(track->GetLabel()))) {
667           code = 3;
668           if (fSaveMotherPDG){
669             TParticle *trackMother =  stack->Particle(TMath::Abs(trackMC->GetFirstMother()));
670             if (trackMother->GetPdgCode() == 310) motherCode = 6; //K0
671             if (trackMother->GetPdgCode() == 3122) motherCode = 7; //Lambda
672             if (trackMother->GetPdgCode() == 3222) motherCode = 8; //Sigma+
673           }
674         }
675         if (pdg == assumedPdg && stack->IsSecondaryFromMaterial(TMath::Abs(track->GetLabel()))) code = 4;
676         
677         //
678         // muons need special treatment, because they are indistinguishable from pions
679         //
680         if (iPart == 0 && pdg == 13  && stack->IsPhysicalPrimary(TMath::Abs(track->GetLabel()))) code = 1;
681         if (iPart == 0 && pdg == 13  && stack->IsSecondaryFromWeakDecay(TMath::Abs(track->GetLabel()))) code = 3;
682         //
683         // check TOF mismatch on MC basis with TOF label
684         //
685         Int_t tofLabel[3];
686         if (!fUseTPConlyTracks) track->GetTOFLabel(tofLabel);
687         else trackForTOF->GetTOFLabel(tofLabel);
688         if (TMath::Abs(track->GetLabel()) != TMath::Abs(tofLabel[0])) hasTOF = kFALSE;
689         //
690         // IMPORTANT BIG PROBLEM HERE THE PROBABLILITY TO HAVE A PID SIGNAL MUST BE IN !!!!!!!!!!!!
691         //
692         //                              0,           1,   2,    3,           4,               5,      6,               7,      8,   9
693         Double_t vectorHistMC[10] = {iPart,  centrality,  pT, sign,  rap[iPart], pullsTPC[iPart], hasTOF, pullsTOF[iPart], dca[0], code};
694         if (!fOnlyQA) { 
695           fHistMCparticles->Fill(vectorHistMC);
696           if (motherCode != -1 && fSaveMotherPDG) { //if mother of weak decay is K0, lambda or sigma+ add track again with this information
697             Double_t vectorHistMCmother[10] = {iPart,  centrality,  pT, sign,  rap[iPart], pullsTPC[iPart], hasTOF, pullsTOF[iPart], dca[0], motherCode};
698             fHistMCparticles->Fill(vectorHistMCmother);
699           }
700         }
701       }
702       //
703       //
704       Int_t tpcShared = track->GetTPCnclsS();
705       if (TMath::Abs(track->Eta()) < 0.8 && iPart == 0 && tpcShared < 4) fHistPidQA->Fill(ptot,tpcSignal,sign);
706     } // end loop over assumed particle type
707
708           //need to delete tpconlytrack
709           if (fUseTPConlyTracks){
710                 delete track; 
711                 track = 0;     
712           }
713
714   } // end of track loop
715   
716   // Post output data  
717   PostData(1, fListHist);
718   
719 }      
720
721
722 //________________________________________________________________________
723 void AliAnalysisCombinedHadronSpectra::Terminate(Option_t *) 
724 {
725   // Draw result to the screen
726   // Called once at the end of the query
727   Printf("*** CONSTRUCTOR CALLED ****");
728
729 }
730
731
732 //________________________________________________________________________
733 Bool_t AliAnalysisCombinedHadronSpectra::SelectOnImpPar(AliESDtrack* t) {
734   //
735   // cut on transverse impact parameter // DEPRECATED
736   //
737   Float_t d0z0[2],covd0z0[3];
738   t->GetImpactParameters(d0z0,covd0z0);
739   Float_t sigma= 0.0050+0.0060/TMath::Power(t->Pt(),0.9);
740   Float_t d0max = 7.*sigma;
741   //
742   Float_t sigmaZ = 0.0146+0.0070/TMath::Power(t->Pt(),1.114758);
743   if (t->Pt() > 1) sigmaZ = 0.0216;
744   Float_t d0maxZ = 5.*sigmaZ;
745   //
746   if(TMath::Abs(d0z0[0]) < d0max && TMath::Abs(d0z0[1]) < d0maxZ) return kTRUE;
747   return kFALSE;
748 }
749
750
751 //________________________________________________________________________
752 void AliAnalysisCombinedHadronSpectra::BinLogAxis(const TH1 *h) {
753   //
754   // Method for the correct logarithmic binning of histograms
755   //
756   TAxis *axis = h->GetXaxis();
757   int bins = axis->GetNbins();
758
759   Double_t from = axis->GetXmin();
760   Double_t to = axis->GetXmax();
761   Double_t *newBins = new Double_t[bins + 1];
762    
763   newBins[0] = from;
764   Double_t factor = pow(to/from, 1./bins);
765   
766   for (int i = 1; i <= bins; i++) {
767    newBins[i] = factor * newBins[i-1];
768   }
769   axis->Set(bins, newBins);
770   delete [] newBins;
771   
772 }
773
774
775 //________________________________________________________________________
776 Int_t AliAnalysisCombinedHadronSpectra::GetPythiaEventProcessType(const AliHeader* aHeader, const Bool_t adebug) const {
777   //
778   // get the process type of the event.
779   //
780
781   // can only read pythia headers, either directly or from cocktalil header
782   AliGenPythiaEventHeader* pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(aHeader->GenEventHeader());
783
784   if (!pythiaGenHeader) {
785
786     AliGenCocktailEventHeader* genCocktailHeader = dynamic_cast<AliGenCocktailEventHeader*>(aHeader->GenEventHeader());
787     if (!genCocktailHeader) {
788       //printf("AliAnalysisCombinedHadronSpectra::GetProcessType : Unknown header type (not Pythia or Cocktail). \n");
789       return -1;
790     }
791
792     TList* headerList = genCocktailHeader->GetHeaders();
793     if (!headerList) {
794       return -1;
795     }
796
797     for (Int_t i=0; i<headerList->GetEntries(); i++) {
798       pythiaGenHeader = dynamic_cast<AliGenPythiaEventHeader*>(headerList->At(i));
799       if (pythiaGenHeader)
800         break;
801     }
802
803     if (!pythiaGenHeader) {
804       //printf("AliAnalysisCombinedHadronSpectra::GetProcessType : Could not find Pythia header. \n");
805       return -1;
806     }
807   }
808
809   if (adebug) {
810     //printf("AliAnalysisCombinedHadronSpectra::GetProcessType : Pythia process type found: %d \n",pythiaGenHeader->ProcessType());
811   }
812
813   return pythiaGenHeader->ProcessType();
814 }