]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/SPECTRA/AliProtonAnalysis.cxx
Adding some new QA plots (dca vs pT)
[u/mrichter/AliRoot.git] / PWG2 / SPECTRA / AliProtonAnalysis.cxx
1 /**************************************************************************
2  * Author: Panos Christakoglou.                                           *
3  * Contributors are mentioned in the code where appropriate.              *
4  *                                                                        *
5  * Permission to use, copy, modify and distribute this software and its   *
6  * documentation strictly for non-commercial purposes is hereby granted   *
7  * without fee, provided that the above copyright notice appears in all   *
8  * copies and that both the copyright notice and this permission notice   *
9  * appear in the supporting documentation. The authors make no claims     *
10  * about the suitability of this software for any purpose. It is          *
11  * provided "as is" without express or implied warranty.                  *
12  **************************************************************************/
13
14 /* $Id$ */
15
16 //-----------------------------------------------------------------
17 //                 AliProtonAnalysis class
18 //   This is the class to deal with the proton analysis
19 //   Origin: Panos Christakoglou | Panos.Christakoglou@cern.ch
20 //-----------------------------------------------------------------
21 #include <Riostream.h>
22 #include <TFile.h>
23 #include <TSystem.h>
24 #include <TF1.h>
25 #include <TH3F.h>
26 #include <TH2D.h>
27 #include <TH1D.h>
28 #include <TH1I.h>
29 #include <TParticle.h>
30 #include <TList.h>
31
32 #include <AliExternalTrackParam.h>
33 #include <AliAODEvent.h>
34 #include <AliESDEvent.h>
35 //#include <AliLog.h>
36 #include <AliPID.h>
37 #include <AliStack.h>
38 #include <AliCFContainer.h>
39 #include <AliCFEffGrid.h>
40 #include <AliCFDataGrid.h>
41 //#include <AliESDVertex.h>
42 class AliLog;
43 class AliESDVertex;
44
45 #include "AliProtonAnalysis.h"
46 #include "AliProtonAnalysisBase.h"
47
48 ClassImp(AliProtonAnalysis)
49
50 //____________________________________________________________________//
51 AliProtonAnalysis::AliProtonAnalysis() : 
52   TObject(), fProtonAnalysisBase(0),
53   fNBinsY(0), fMinY(0), fMaxY(0),
54   fNBinsPt(0), fMinPt(0), fMaxPt(0),
55   fProtonContainer(0), fAntiProtonContainer(0),
56   fHistEvents(0), fHistYPtProtons(0), fHistYPtAntiProtons(0), 
57   fHistEventStats(0),
58   fEffGridListProtons(0), fCorrectionListProtons2D(0), 
59   fEfficiencyListProtons1D(0), fCorrectionListProtons1D(0),
60   fEffGridListAntiProtons(0), fCorrectionListAntiProtons2D(0), 
61   fEfficiencyListAntiProtons1D(0), fCorrectionListAntiProtons1D(0),
62   fCorrectProtons(0), fCorrectAntiProtons(0),
63   fGlobalQAList(0), fQA2DList(0),
64   fQAProtonsAcceptedList(0), fQAProtonsRejectedList(0),
65   fQAAntiProtonsAcceptedList(0), fQAAntiProtonsRejectedList(0),
66   fInitQAFlag(kFALSE) {
67   //Default constructor
68 }
69
70 //____________________________________________________________________//
71 AliProtonAnalysis::AliProtonAnalysis(Int_t nbinsY, 
72                                      Float_t fLowY, Float_t fHighY,
73                                      Int_t nbinsPt, 
74                                      Float_t fLowPt, Float_t fHighPt) : 
75   TObject(), fProtonAnalysisBase(0),
76   fNBinsY(nbinsY), fMinY(fLowY), fMaxY(fHighY),
77   fNBinsPt(nbinsPt), fMinPt(fLowPt), fMaxPt(fHighPt),
78   fProtonContainer(0), fAntiProtonContainer(0),
79   fHistEvents(0), fHistYPtProtons(0), fHistYPtAntiProtons(0), fHistEventStats(0),
80   fEffGridListProtons(0), fCorrectionListProtons2D(0), 
81   fEfficiencyListProtons1D(0), fCorrectionListProtons1D(0),
82   fEffGridListAntiProtons(0), fCorrectionListAntiProtons2D(0), 
83   fEfficiencyListAntiProtons1D(0), fCorrectionListAntiProtons1D(0),
84   fCorrectProtons(0), fCorrectAntiProtons(0),
85   fGlobalQAList(0), fQA2DList(0),
86   fQAProtonsAcceptedList(0), fQAProtonsRejectedList(0),
87   fQAAntiProtonsAcceptedList(0), fQAAntiProtonsRejectedList(0),
88   fInitQAFlag(kFALSE) {
89   //Default constructor
90   fHistEvents = new TH1I("fHistEvents","Analyzed events",2,0.5,2.5);
91   fHistEvents->GetXaxis()->SetBinLabel(1,"Analyzed events");
92   fHistEvents->GetXaxis()->SetBinLabel(2,"Events with (anti)protons");
93
94   fHistYPtProtons = new TH2D("fHistYPtProtons","Protons",
95                              fNBinsY,fMinY,fMaxY,
96                              fNBinsPt,fMinPt,fMaxPt);
97   fHistYPtProtons->SetStats(kTRUE);
98   fHistYPtProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
99   if(fProtonAnalysisBase->GetEtaMode())
100     fHistYPtProtons->GetXaxis()->SetTitle("#eta");
101   else
102     fHistYPtProtons->GetXaxis()->SetTitle("y");
103   fHistYPtProtons->GetXaxis()->SetTitleColor(1);
104
105   fHistYPtAntiProtons = new TH2D("fHistYPtAntiProtons","Antiprotons",
106                                  fNBinsY,fMinY,fMaxY,
107                                  fNBinsPt,fMinPt,fMaxPt);
108   fHistYPtAntiProtons->SetStats(kTRUE);
109   fHistYPtAntiProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
110   if(fProtonAnalysisBase->GetEtaMode())
111     fHistYPtAntiProtons->GetXaxis()->SetTitle("#eta");
112   else
113     fHistYPtAntiProtons->GetXaxis()->SetTitle("y");
114   fHistYPtAntiProtons->GetXaxis()->SetTitleColor(1);
115
116   //setting up the containers
117   Int_t iBin[2];
118   iBin[0] = nbinsY;
119   iBin[1] = nbinsPt;
120   Double_t *binLimY = new Double_t[nbinsY+1];
121   Double_t *binLimPt = new Double_t[nbinsPt+1];
122   //values for bin lower bounds
123   for(Int_t i = 0; i <= nbinsY; i++) 
124     binLimY[i]=(Double_t)fLowY  + (fHighY - fLowY)  /nbinsY*(Double_t)i;
125   for(Int_t i = 0; i <= nbinsPt; i++) 
126     binLimPt[i]=(Double_t)fLowPt  + (fHighPt - fLowPt)  /nbinsPt*(Double_t)i;
127
128   fProtonContainer = new AliCFContainer("containerProtons",
129                                         "container for protons",
130                                         kNSteps,2,iBin);
131   fProtonContainer->SetBinLimits(0,binLimY); //rapidity or eta
132   fProtonContainer->SetBinLimits(1,binLimPt); //pT
133   fAntiProtonContainer = new AliCFContainer("containerAntiProtons",
134                                             "container for antiprotons",
135                                             kNSteps,2,iBin);
136   fAntiProtonContainer->SetBinLimits(0,binLimY); //rapidity or eta
137   fAntiProtonContainer->SetBinLimits(1,binLimPt); //pT
138
139   //Initialize the QA
140   if(!fInitQAFlag) InitQA();
141
142
143 //____________________________________________________________________//
144 AliProtonAnalysis::AliProtonAnalysis(Int_t nbinsY, Double_t *gY,
145                                      Int_t nbinsPt,Double_t *gPt) : 
146   TObject(), fProtonAnalysisBase(0),
147   fNBinsY(nbinsY), fMinY(gY[0]), fMaxY(gY[nbinsY]),
148   fNBinsPt(nbinsPt), fMinPt(gPt[0]), fMaxPt(gPt[nbinsPt]),
149   fProtonContainer(0), fAntiProtonContainer(0),
150   fHistEvents(0), fHistYPtProtons(0), fHistYPtAntiProtons(0), fHistEventStats(0),
151   fEffGridListProtons(0), fCorrectionListProtons2D(0), 
152   fEfficiencyListProtons1D(0), fCorrectionListProtons1D(0),
153   fEffGridListAntiProtons(0), fCorrectionListAntiProtons2D(0), 
154   fEfficiencyListAntiProtons1D(0), fCorrectionListAntiProtons1D(0),
155   fCorrectProtons(0), fCorrectAntiProtons(0),
156   fGlobalQAList(0), fQA2DList(0),
157   fQAProtonsAcceptedList(0), fQAProtonsRejectedList(0),
158   fQAAntiProtonsAcceptedList(0), fQAAntiProtonsRejectedList(0),
159   fInitQAFlag(kFALSE) {
160   //Default constructor
161   fHistEvents = new TH1I("fHistEvents","Analyzed events",2,0.5,2.5);
162   fHistEvents->GetXaxis()->SetBinLabel(1,"Analyzed events");
163   fHistEvents->GetXaxis()->SetBinLabel(2,"Events with (anti)protons");
164
165   fHistYPtProtons = new TH2D("fHistYPtProtons","Protons",
166                              fNBinsY,gY,fNBinsPt,gPt);
167   fHistYPtProtons->SetStats(kTRUE);
168   fHistYPtProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
169   if(fProtonAnalysisBase->GetEtaMode())
170     fHistYPtProtons->GetXaxis()->SetTitle("#eta");
171   else
172     fHistYPtProtons->GetXaxis()->SetTitle("y");
173   fHistYPtProtons->GetXaxis()->SetTitleColor(1);
174
175   fHistYPtAntiProtons = new TH2D("fHistYPtAntiProtons","Antiprotons",
176                                  fNBinsY,gY,fNBinsPt,gPt);
177   fHistYPtAntiProtons->SetStats(kTRUE);
178   fHistYPtAntiProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
179   if(fProtonAnalysisBase->GetEtaMode())
180     fHistYPtAntiProtons->GetXaxis()->SetTitle("#eta");
181   else
182     fHistYPtAntiProtons->GetXaxis()->SetTitle("y");
183   fHistYPtAntiProtons->GetXaxis()->SetTitleColor(1);
184
185   //setting up the containers
186   Int_t iBin[2];
187   iBin[0] = nbinsY;
188   iBin[1] = nbinsPt;
189   fProtonContainer = new AliCFContainer("containerProtons",
190                                         "container for protons",
191                                         kNSteps,2,iBin);
192   fProtonContainer->SetBinLimits(0,gY); //rapidity or eta
193   fProtonContainer->SetBinLimits(1,gPt); //pT
194   fAntiProtonContainer = new AliCFContainer("containerAntiProtons",
195                                             "container for antiprotons",
196                                             kNSteps,2,iBin);
197   fAntiProtonContainer->SetBinLimits(0,gY); //rapidity or eta
198   fAntiProtonContainer->SetBinLimits(1,gPt); //pT
199
200   //Initialize the QA
201   if(!fInitQAFlag) InitQA();
202
203
204 //____________________________________________________________________//
205 AliProtonAnalysis::~AliProtonAnalysis() {
206   //Default destructor
207   if(fProtonAnalysisBase) delete fProtonAnalysisBase;
208
209   if(fHistEvents) delete fHistEvents;
210   if(fHistYPtProtons) delete fHistYPtProtons;
211   if(fHistYPtAntiProtons) delete fHistYPtAntiProtons;
212   if(fHistEventStats) delete fHistEventStats;
213   if(fProtonContainer) delete fProtonContainer;
214   if(fAntiProtonContainer) delete fAntiProtonContainer;
215
216   if(fEffGridListProtons) delete fEffGridListProtons;
217   if(fCorrectionListProtons2D) delete fCorrectionListProtons2D;
218   if(fEfficiencyListProtons1D) delete fEfficiencyListProtons1D;
219   if(fCorrectionListProtons1D) delete fCorrectionListProtons1D;
220   if(fEffGridListAntiProtons) delete fEffGridListAntiProtons;
221   if(fCorrectionListAntiProtons2D) delete fCorrectionListAntiProtons2D;
222   if(fEfficiencyListAntiProtons1D) delete fEfficiencyListAntiProtons1D;
223   if(fCorrectionListAntiProtons1D) delete fCorrectionListAntiProtons1D;
224   if(fCorrectProtons) delete fCorrectProtons;
225   if(fCorrectAntiProtons) delete fCorrectAntiProtons;
226
227   //QA lists
228   if(fGlobalQAList) delete fGlobalQAList;
229   if(fQA2DList) delete fQA2DList;
230   if(fQAProtonsAcceptedList) delete fQAProtonsAcceptedList; 
231   if(fQAProtonsRejectedList) delete fQAProtonsRejectedList;
232   if(fQAAntiProtonsAcceptedList) delete fQAAntiProtonsAcceptedList; 
233   if(fQAAntiProtonsRejectedList) delete fQAAntiProtonsRejectedList;
234 }
235
236 //____________________________________________________________________//
237 void AliProtonAnalysis::InitAnalysisHistograms(Int_t nbinsY, 
238                                                Float_t fLowY, Float_t fHighY, 
239                                                Int_t nbinsPt, 
240                                                Float_t fLowPt, Float_t fHighPt) {
241   //Initializes the histograms
242   fNBinsY = nbinsY;
243   fMinY = fLowY;
244   fMaxY = fHighY;
245   fNBinsPt = nbinsPt;
246   fMinPt = fLowPt;
247   fMaxPt = fHighPt;
248
249   fHistEvents = new TH1I("fHistEvents","Analyzed events",2,0.5,2.5);
250   fHistEvents->GetXaxis()->SetBinLabel(1,"Analyzed events");
251   fHistEvents->GetXaxis()->SetBinLabel(2,"Events with (anti)protons");
252
253   fHistYPtProtons = new TH2D("fHistYPtProtons","Protons",
254                              fNBinsY,fMinY,fMaxY,
255                              fNBinsPt,fMinPt,fMaxPt);
256   fHistYPtProtons->SetStats(kTRUE);
257   fHistYPtProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
258   if(fProtonAnalysisBase->GetEtaMode())
259     fHistYPtProtons->GetXaxis()->SetTitle("#eta");
260   else
261     fHistYPtProtons->GetXaxis()->SetTitle("y");
262   fHistYPtProtons->GetXaxis()->SetTitleColor(1);
263
264   fHistYPtAntiProtons = new TH2D("fHistYPtAntiProtons","Antiprotons",
265                                  fNBinsY,fMinY,fMaxY,
266                                  fNBinsPt,fMinPt,fMaxPt);
267   fHistYPtAntiProtons->SetStats(kTRUE);
268   fHistYPtAntiProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
269   if(fProtonAnalysisBase->GetEtaMode())
270     fHistYPtAntiProtons->GetXaxis()->SetTitle("#eta");
271   else
272     fHistYPtAntiProtons->GetXaxis()->SetTitle("y");
273   fHistYPtAntiProtons->GetXaxis()->SetTitleColor(1);
274
275   //setting up the containers
276   Int_t iBin[2];
277   iBin[0] = nbinsY;
278   iBin[1] = nbinsPt;
279   Double_t *binLimY = new Double_t[nbinsY+1];
280   Double_t *binLimPt = new Double_t[nbinsPt+1];
281   //values for bin lower bounds
282   for(Int_t i = 0; i <= nbinsY; i++) 
283     binLimY[i]=(Double_t)fLowY  + (fHighY - fLowY)  /nbinsY*(Double_t)i;
284   for(Int_t i = 0; i <= nbinsPt; i++) 
285     binLimPt[i]=(Double_t)fLowPt  + (fHighPt - fLowPt)  /nbinsPt*(Double_t)i;
286
287   fProtonContainer = new AliCFContainer("containerProtons",
288                                         "container for protons",
289                                         kNSteps,2,iBin);
290   fProtonContainer->SetBinLimits(0,binLimY); //rapidity
291   fProtonContainer->SetBinLimits(1,binLimPt); //pT
292   fAntiProtonContainer = new AliCFContainer("containerAntiProtons",
293                                             "container for antiprotons",
294                                             kNSteps,2,iBin);
295   fAntiProtonContainer->SetBinLimits(0,binLimY); //rapidity
296   fAntiProtonContainer->SetBinLimits(1,binLimPt); //pT
297
298   //Initialize the QA
299   if(!fInitQAFlag) InitQA();
300 }
301
302 //____________________________________________________________________//
303 void AliProtonAnalysis::InitAnalysisHistograms(Int_t nbinsY, Double_t *gY, 
304                                                Int_t nbinsPt, Double_t *gPt) {
305   //Initializes the histograms using asymmetric values - global tracking
306   fNBinsY = nbinsY;
307   fMinY = gY[0];
308   fMaxY = gY[nbinsY];
309   fNBinsPt = nbinsPt;
310   fMinPt = gPt[0];
311   fMaxPt = gPt[nbinsPt];
312
313   fHistEvents = new TH1I("fHistEvents","Analyzed events",2,0.5,2.5);
314   fHistEvents->GetXaxis()->SetBinLabel(1,"Analyzed events");
315   fHistEvents->GetXaxis()->SetBinLabel(2,"Events with (anti)protons");
316
317   fHistYPtProtons = new TH2D("fHistYPtProtons","Protons",
318                              fNBinsY,gY,fNBinsPt,gPt);
319   fHistYPtProtons->SetStats(kTRUE);
320   fHistYPtProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
321   if(fProtonAnalysisBase->GetEtaMode())
322     fHistYPtProtons->GetXaxis()->SetTitle("#eta");
323   else
324     fHistYPtProtons->GetXaxis()->SetTitle("y");
325   fHistYPtProtons->GetXaxis()->SetTitleColor(1);
326
327   fHistYPtAntiProtons = new TH2D("fHistYPtAntiProtons","Antiprotons",
328                                  fNBinsY,gY,fNBinsPt,gPt);
329   fHistYPtAntiProtons->SetStats(kTRUE);
330   fHistYPtAntiProtons->GetYaxis()->SetTitle("P_{T} [GeV/c]");
331   if(fProtonAnalysisBase->GetEtaMode())
332     fHistYPtAntiProtons->GetXaxis()->SetTitle("#eta");
333   else
334     fHistYPtAntiProtons->GetXaxis()->SetTitle("y");
335   fHistYPtAntiProtons->GetXaxis()->SetTitleColor(1);
336
337   //setting up the containers
338   Int_t iBin[2];
339   iBin[0] = nbinsY;
340   iBin[1] = nbinsPt;
341
342   fProtonContainer = new AliCFContainer("containerProtons",
343                                         "container for protons",
344                                         kNSteps,2,iBin);
345   fProtonContainer->SetBinLimits(0,gY); //rapidity
346   fProtonContainer->SetBinLimits(1,gPt); //pT
347   fAntiProtonContainer = new AliCFContainer("containerAntiProtons",
348                                             "container for antiprotons",
349                                             kNSteps,2,iBin);
350   fAntiProtonContainer->SetBinLimits(0,gY); //rapidity
351   fAntiProtonContainer->SetBinLimits(1,gPt); //pT
352
353   //Initialize the QA
354   if(!fInitQAFlag) InitQA();
355 }
356
357 //____________________________________________________________________//
358 Bool_t AliProtonAnalysis::ReadFromFile(const char* filename) {
359   //Read the containers from the existing file
360   Bool_t status = kTRUE;
361
362   TFile *file = TFile::Open(filename);
363   if(!file) {
364     cout<<"Could not find the input file "<<filename<<endl;
365     status = kFALSE;
366   }
367
368   TList *list = (TList *)file->Get("outputList");
369   if(list) {
370     cout<<"Retrieving objects from the list "<<list->GetName()<<"..."<<endl; 
371     fHistYPtProtons = (TH2D *)list->At(0);
372     fHistYPtAntiProtons = (TH2D *)list->At(1);
373     fHistEvents = (TH1I *)list->At(2);
374     fProtonContainer = (AliCFContainer *)list->At(3);
375     fAntiProtonContainer = (AliCFContainer *)list->At(4);
376     fHistEventStats = (TH1F *)list->At(5);
377   }
378   else if(!list) {
379     cout<<"Retrieving objects from the file... "<<endl;
380     fHistYPtProtons = (TH2D *)file->Get("fHistYPtProtons");
381     fHistYPtAntiProtons = (TH2D *)file->Get("fHistYPtAntiProtons");
382     fHistEvents = (TH1I *)file->Get("fHistEvents");
383     fProtonContainer = (AliCFContainer *)file->Get("containerProtons");
384     fAntiProtonContainer = (AliCFContainer *)file->Get("containerAntiProtons");
385     fHistEventStats = (TH1F *)file->Get("fHistEventStats");
386   }
387   if((!fHistYPtProtons)||(!fHistYPtAntiProtons)||(!fHistEvents)
388      ||(!fProtonContainer)||(!fAntiProtonContainer)||(!fHistEventStats)) {
389     cout<<"Input containers were not found!!!"<<endl;
390     status = kFALSE;
391   }
392   else {
393     //fHistYPtProtons = fProtonContainer->ShowProjection(0,1,0);
394     //fHistYPtAntiProtons = fAntiProtonContainer->ShowProjection(0,1,0);
395     fHistYPtProtons->Sumw2();
396     fHistYPtAntiProtons->Sumw2();
397   }
398
399   return status;
400 }
401
402 //____________________________________________________________________//
403 TH1D *AliProtonAnalysis::GetProtonYHistogram() {
404   //Get the y histogram for protons
405   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
406
407   //TH1D *fYProtons = (TH1D *)fHistYPtProtons->ProjectionX("fYProtons",0,fHistYPtProtons->GetYaxis()->GetNbins(),"");
408   TH1D *fYProtons = fProtonContainer->ShowProjection(0,2); //variable-step
409    
410   fYProtons->SetStats(kFALSE);
411   fYProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dy)");
412   fYProtons->SetTitle("dN/dy protons");
413   fYProtons->SetMarkerStyle(kFullCircle);
414   fYProtons->SetMarkerColor(4);
415   if(nAnalyzedEvents > 0)
416   fYProtons->Scale(1./nAnalyzedEvents);
417   
418   return fYProtons;
419 }
420
421 //____________________________________________________________________//
422 TH1D *AliProtonAnalysis::GetAntiProtonYHistogram() {
423   //Get the y histogram for antiprotons
424   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
425   
426   //TH1D *fYAntiProtons = (TH1D *)fHistYPtAntiProtons->ProjectionX("fYAntiProtons",0,fHistYPtAntiProtons->GetYaxis()->GetNbins(),"");
427   TH1D *fYAntiProtons = fAntiProtonContainer->ShowProjection(0,2);//variable-step 
428  
429   fYAntiProtons->SetStats(kFALSE);
430   fYAntiProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dy)");
431   fYAntiProtons->SetTitle("dN/dy antiprotons");
432   fYAntiProtons->SetMarkerStyle(kFullCircle);
433   fYAntiProtons->SetMarkerColor(4);
434   if(nAnalyzedEvents > 0)
435     fYAntiProtons->Scale(1./nAnalyzedEvents);
436
437   return fYAntiProtons;
438 }
439
440 //____________________________________________________________________//
441 TH1D *AliProtonAnalysis::GetProtonPtHistogram() {
442   //Get the Pt histogram for protons
443   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
444   
445   //TH1D *fPtProtons = (TH1D *)fHistYPtProtons->ProjectionY("fPtProtons",0,fHistYPtProtons->GetXaxis()->GetNbins(),""); 
446   TH1D *fPtProtons = fProtonContainer->ShowProjection(1,2); //variable-step
447
448   fPtProtons->SetStats(kFALSE);
449   fPtProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dP_{T})");
450   fPtProtons->SetTitle("dN/dPt protons");
451   fPtProtons->SetMarkerStyle(kFullCircle);
452   fPtProtons->SetMarkerColor(4);
453   if(nAnalyzedEvents > 0)
454     fPtProtons->Scale(1./nAnalyzedEvents);
455
456   return fPtProtons;
457 }
458
459 //____________________________________________________________________//
460 TH1D *AliProtonAnalysis::GetAntiProtonPtHistogram() {
461   //Get the Pt histogram for antiprotons
462   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
463   
464   //TH1D *fPtAntiProtons = (TH1D *)fHistYPtAntiProtons->ProjectionY("fPtAntiProtons",0,fHistYPtProtons->GetXaxis()->GetNbins(),""); 
465   TH1D *fPtAntiProtons = fAntiProtonContainer->ShowProjection(1,2); //variable-step
466
467   fPtAntiProtons->SetStats(kFALSE);
468   fPtAntiProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dP_{T})");
469   fPtAntiProtons->SetTitle("dN/dPt antiprotons");
470   fPtAntiProtons->SetMarkerStyle(kFullCircle);
471   fPtAntiProtons->SetMarkerColor(4);
472   if(nAnalyzedEvents > 0)
473     fPtAntiProtons->Scale(1./nAnalyzedEvents);
474
475   return fPtAntiProtons;
476 }
477
478 //____________________________________________________________________//
479 TH1D *AliProtonAnalysis::GetProtonCorrectedYHistogram() {
480   //Get the corrected y histogram for protons
481   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
482
483   TH1D *fYProtons = fCorrectProtons->Project(0); //0: rapidity
484    
485   fYProtons->SetStats(kFALSE);
486   fYProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dy)");
487   fYProtons->GetXaxis()->SetTitle("y");
488   fYProtons->SetTitle("dN/dy protons");
489   fYProtons->SetMarkerStyle(kFullCircle);
490   fYProtons->SetMarkerColor(4);
491   if(nAnalyzedEvents > 0)
492     fYProtons->Scale(1./nAnalyzedEvents);
493   
494   return fYProtons;
495 }
496
497 //____________________________________________________________________//
498 TH1D *AliProtonAnalysis::GetAntiProtonCorrectedYHistogram() {
499   //Get the corrected y histogram for antiprotons
500   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
501
502   TH1D *fYAntiProtons = fCorrectAntiProtons->Project(0); //0: rapidity
503    
504   fYAntiProtons->SetStats(kFALSE);
505   fYAntiProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dy)");
506   fYAntiProtons->GetXaxis()->SetTitle("y");
507   fYAntiProtons->SetTitle("dN/dy protons");
508   fYAntiProtons->SetMarkerStyle(kFullCircle);
509   fYAntiProtons->SetMarkerColor(4);
510   if(nAnalyzedEvents > 0)
511     fYAntiProtons->Scale(1./nAnalyzedEvents);
512   
513   return fYAntiProtons;
514 }
515
516 //____________________________________________________________________//
517 TH1D *AliProtonAnalysis::GetProtonCorrectedPtHistogram() {
518   //Get the corrected Pt histogram for protons
519   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
520
521   TH1D *fPtProtons = fCorrectProtons->Project(0); //0: rapidity
522    
523   fPtProtons->SetStats(kFALSE);
524   fPtProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dP_{T})");
525   fPtProtons->GetXaxis()->SetTitle("P_{T} [GeV/c]");
526   fPtProtons->SetTitle("dN/dPt protons");
527   fPtProtons->SetMarkerStyle(kFullCircle);
528   fPtProtons->SetMarkerColor(4);
529   if(nAnalyzedEvents > 0)
530     fPtProtons->Scale(1./nAnalyzedEvents);
531   
532   return fPtProtons;
533 }
534
535 //____________________________________________________________________//
536 TH1D *AliProtonAnalysis::GetAntiProtonCorrectedPtHistogram() {
537   //Get the corrected Pt histogram for antiprotons
538   Int_t nAnalyzedEvents = GetNumberOfAnalyzedEvents();
539
540   TH1D *fPtAntiProtons = fCorrectAntiProtons->Project(0); //0: rapidity
541    
542   fPtAntiProtons->SetStats(kFALSE);
543   fPtAntiProtons->GetYaxis()->SetTitle("(1/N_{events})(dN/dP_{T})");
544   fPtAntiProtons->GetXaxis()->SetTitle("P_{T} [GeV/c]");
545   fPtAntiProtons->SetTitle("dN/dPt antiprotons");
546   fPtAntiProtons->SetMarkerStyle(kFullCircle);
547   fPtAntiProtons->SetMarkerColor(4);
548   if(nAnalyzedEvents > 0)
549     fPtAntiProtons->Scale(1./nAnalyzedEvents);
550   
551   return fPtAntiProtons;
552 }
553
554 //____________________________________________________________________//
555 TH1D *AliProtonAnalysis::GetYRatioHistogram() {
556   //Returns the rapidity dependence of the ratio (uncorrected)
557   TH1D *fYProtons = GetProtonYHistogram();
558   TH1D *fYAntiProtons = GetAntiProtonYHistogram();
559   
560   TH1D *hRatioY = new TH1D("hRatioY","",fYProtons->GetNbinsX(),fYProtons->GetXaxis()->GetXmin(),fYProtons->GetXaxis()->GetXmax());
561   hRatioY->Divide(fYAntiProtons,fYProtons,1.0,1.0);
562   hRatioY->SetMarkerStyle(kFullCircle);
563   hRatioY->SetMarkerColor(4);
564   hRatioY->GetYaxis()->SetTitle("#bar{p}/p");
565   hRatioY->GetYaxis()->SetTitleOffset(1.4);
566   hRatioY->GetXaxis()->SetTitle("y");
567   hRatioY->GetXaxis()->SetTitleColor(1);
568   hRatioY->SetStats(kFALSE);
569
570   return hRatioY;
571 }
572
573 //____________________________________________________________________//
574 TH1D *AliProtonAnalysis::GetYRatioCorrectedHistogram(TH2D *gCorrectionProtons, 
575                                                      TH2D *gCorrectionAntiProtons) {
576   //Returns the rapidity dependence of the ratio (corrected)
577   fHistYPtProtons->Multiply(gCorrectionProtons);
578   TH1D *fYProtons = GetProtonYHistogram();
579   fHistYPtAntiProtons->Multiply(gCorrectionAntiProtons);
580   TH1D *fYAntiProtons = GetAntiProtonYHistogram();
581   
582   TH1D *hRatioY = new TH1D("hRatioY","",fYProtons->GetNbinsX(),fYProtons->GetXaxis()->GetXmin(),fYProtons->GetXaxis()->GetXmax());
583   hRatioY->Divide(fYAntiProtons,fYProtons,1.0,1.0);
584   hRatioY->SetMarkerStyle(kFullCircle);
585   hRatioY->SetMarkerColor(4);
586   hRatioY->GetYaxis()->SetTitle("#bar{p}/p");
587   hRatioY->GetYaxis()->SetTitleOffset(1.4);
588   hRatioY->GetXaxis()->SetTitle("y");
589   hRatioY->GetXaxis()->SetTitleColor(1);
590   hRatioY->SetStats(kFALSE);
591
592   return hRatioY;
593 }
594
595 //____________________________________________________________________//
596 TH1D *AliProtonAnalysis::GetPtRatioHistogram() {
597   //Returns the pT dependence of the ratio (uncorrected)
598   TH1D *fPtProtons = GetProtonPtHistogram();
599   TH1D *fPtAntiProtons = GetAntiProtonPtHistogram();
600   
601   TH1D *hRatioPt = new TH1D("hRatioPt","",fPtProtons->GetNbinsX(),fPtProtons->GetXaxis()->GetXmin(),fPtProtons->GetXaxis()->GetXmax());
602   hRatioPt->Divide(fPtAntiProtons,fPtProtons,1.0,1.0);
603   hRatioPt->SetMarkerStyle(kFullCircle);
604   hRatioPt->SetMarkerColor(4);
605   hRatioPt->GetYaxis()->SetTitle("#bar{p}/p");
606   hRatioPt->GetYaxis()->SetTitleOffset(1.4);
607   hRatioPt->GetXaxis()->SetTitle("P_{T} [GeV/c]");
608   hRatioPt->GetXaxis()->SetTitleColor(1);
609   hRatioPt->SetStats(kFALSE);
610
611   return hRatioPt;
612 }
613
614 //____________________________________________________________________//
615 TH1D *AliProtonAnalysis::GetPtRatioCorrectedHistogram(TH2D *gCorrectionProtons, 
616                                                       TH2D *gCorrectionAntiProtons) {
617   //Returns the Pt dependence of the ratio (corrected)
618   fHistYPtProtons->Multiply(gCorrectionProtons);
619   TH1D *fPtProtons = GetProtonPtHistogram();
620   fHistYPtAntiProtons->Multiply(gCorrectionAntiProtons);
621   TH1D *fPtAntiProtons = GetAntiProtonPtHistogram();
622   
623   TH1D *hRatioPt = new TH1D("hRatioPt","",fPtProtons->GetNbinsX(),fPtProtons->GetXaxis()->GetXmin(),fPtProtons->GetXaxis()->GetXmax());
624   hRatioPt->Divide(fPtAntiProtons,fPtProtons,1.0,1.0);
625   hRatioPt->SetMarkerStyle(kFullCircle);
626   hRatioPt->SetMarkerColor(4);
627   hRatioPt->GetYaxis()->SetTitle("#bar{p}/p");
628   hRatioPt->GetYaxis()->SetTitleOffset(1.4);
629   hRatioPt->GetXaxis()->SetTitle("y");
630   hRatioPt->GetXaxis()->SetTitleColor(1);
631   hRatioPt->SetStats(kFALSE);
632
633   return hRatioPt;
634 }
635
636 //____________________________________________________________________//
637 TH1D *AliProtonAnalysis::GetYAsymmetryHistogram() {
638   //Returns the rapidity dependence of the asymmetry (uncorrected)
639   TH1D *fYProtons = GetProtonYHistogram();
640   TH1D *fYAntiProtons = GetAntiProtonYHistogram();
641   
642   TH1D *hsum = new TH1D("hsumY","",fYProtons->GetNbinsX(),fYProtons->GetXaxis()->GetXmin(),fYProtons->GetXaxis()->GetXmax());
643   hsum->Add(fYProtons,fYAntiProtons,1.0,1.0);
644
645   TH1D *hdiff = new TH1D("hdiffY","",fYProtons->GetNbinsX(),fYProtons->GetXaxis()->GetXmin(),fYProtons->GetXaxis()->GetXmax());
646   hdiff->Add(fYProtons,fYAntiProtons,1.0,-1.0);
647
648   TH1D *hAsymmetryY = new TH1D("hAsymmetryY","",fYProtons->GetNbinsX(),fYProtons->GetXaxis()->GetXmin(),fYProtons->GetXaxis()->GetXmax());
649   hAsymmetryY->Divide(hdiff,hsum,2.0,1.);
650   hAsymmetryY->SetMarkerStyle(kFullCircle);
651   hAsymmetryY->SetMarkerColor(4);
652   hAsymmetryY->GetYaxis()->SetTitle("A_{p}");
653   hAsymmetryY->GetYaxis()->SetTitleOffset(1.4);
654   hAsymmetryY->GetXaxis()->SetTitle("y");
655   hAsymmetryY->GetXaxis()->SetTitleColor(1);
656   hAsymmetryY->SetStats(kFALSE);
657
658   return hAsymmetryY;
659 }
660
661 //____________________________________________________________________//
662 TH1D *AliProtonAnalysis::GetPtAsymmetryHistogram() {
663   //Returns the pT dependence of the asymmetry (uncorrected)
664   TH1D *fPtProtons = GetProtonPtHistogram();
665   TH1D *fPtAntiProtons = GetAntiProtonPtHistogram();
666   
667   TH1D *hsum = new TH1D("hsumPt","",fPtProtons->GetNbinsX(),fPtProtons->GetXaxis()->GetXmin(),fPtProtons->GetXaxis()->GetXmax());
668   hsum->Add(fPtProtons,fPtAntiProtons,1.0,1.0);
669
670   TH1D *hdiff = new TH1D("hdiffPt","",fPtProtons->GetNbinsX(),fPtProtons->GetXaxis()->GetXmin(),fPtProtons->GetXaxis()->GetXmax());
671   hdiff->Add(fPtProtons,fPtAntiProtons,1.0,-1.0);
672
673   TH1D *hAsymmetryPt = new TH1D("hAsymmetryPt","",fPtProtons->GetNbinsX(),fPtProtons->GetXaxis()->GetXmin(),fPtProtons->GetXaxis()->GetXmax());
674   hAsymmetryPt->Divide(hdiff,hsum,2.0,1.);
675   hAsymmetryPt->SetMarkerStyle(kFullCircle);
676   hAsymmetryPt->SetMarkerColor(4);
677   hAsymmetryPt->GetYaxis()->SetTitle("A_{p}");
678   hAsymmetryPt->GetYaxis()->SetTitleOffset(1.4);
679   hAsymmetryPt->GetXaxis()->SetTitle("P_{T} [GeV/c]");
680   hAsymmetryPt->GetXaxis()->SetTitleColor(1);
681   hAsymmetryPt->SetStats(kFALSE);
682
683   return hAsymmetryPt;
684 }
685
686 //____________________________________________________________________//
687 void AliProtonAnalysis::Analyze(AliESDEvent* esd,
688                                 const AliESDVertex *vertex) {
689   //Main analysis part - ESD
690   Int_t nTracks = 0;
691   Int_t nIdentifiedProtons = 0, nIdentifiedAntiProtons = 0;
692   Int_t nSurvivedProtons = 0, nSurvivedAntiProtons = 0;
693
694   fHistEvents->Fill(1); //number of analyzed events
695   Double_t containerInput[2] ;
696   Double_t gPt = 0.0, gP = 0.0;
697   nTracks = esd->GetNumberOfTracks();
698   for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
699     AliESDtrack* track = esd->GetTrack(iTracks);
700     AliESDtrack trackTPC;
701
702     //in case it's a TPC only track relate it to the proper vertex
703     /*if(fProtonAnalysisBase->GetAnalysisMode()==AliProtonAnalysisBase::kTPC) {
704       Float_t p[2],cov[3];
705       track->GetImpactParametersTPC(p,cov);
706       if (p[0]==0 && p[1]==0)  
707         track->RelateToVertexTPC(((AliESDEvent*)esd)->GetPrimaryVertexTPC(),esd->GetMagneticField(),kVeryBig);
708       if (!track->FillTPCOnlyTrack(trackTPC)) {
709         continue;
710       }
711       track = &trackTPC ;
712       }*/
713
714     Int_t  fIdxInt[200];
715     Int_t nClustersTPC = track->GetTPCclusters(fIdxInt);
716     Int_t npointsTPCdEdx = track->GetTPCsignalN();
717     Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};  //The impact parameters and their covariance.
718
719     if((fProtonAnalysisBase->GetAnalysisMode()==AliProtonAnalysisBase::kTPC)||(fProtonAnalysisBase->GetAnalysisMode()==AliProtonAnalysisBase::kHybrid)) {
720       AliExternalTrackParam *tpcTrack = (AliExternalTrackParam *)track->GetTPCInnerParam();
721       if(!tpcTrack) continue;
722       gPt = tpcTrack->Pt();
723       gP = tpcTrack->P();
724     
725       tpcTrack->PropagateToDCA(vertex,
726                                esd->GetMagneticField(),
727                                100.,dca,cov);
728
729       if(fProtonAnalysisBase->IsPrimary(esd,vertex,track)) {
730         if(fProtonAnalysisBase->IsAccepted(esd,vertex,track)) {
731           ((TH2F *)(fQA2DList->At(0)))->Fill(gP,track->GetTPCsignal());
732           ((TH3F *)(fQA2DList->At(2)))->Fill(tpcTrack->Eta(),
733                                              tpcTrack->Phi()*180./TMath::Pi(),
734                                              npointsTPCdEdx);
735           ((TH3F *)(fQA2DList->At(4)))->Fill(tpcTrack->Eta(),
736                                              tpcTrack->Phi()*180./TMath::Pi(),
737                                              nClustersTPC);
738           ((TH3F *)(fQA2DList->At(6)))->Fill(gPt,
739                                              tpcTrack->Phi()*180./TMath::Pi(),
740                                              npointsTPCdEdx);
741           ((TH3F *)(fQA2DList->At(8)))->Fill(gPt,
742                                              tpcTrack->Phi()*180./TMath::Pi(),
743                                              nClustersTPC);     
744         }//quality cuts
745       }//primary cuts
746       
747       if(fProtonAnalysisBase->IsProton(track)) {
748         //Step: kStepIdentified
749         if(fProtonAnalysisBase->GetEtaMode())
750           containerInput[0] = tpcTrack->Eta();
751         else
752           containerInput[0] = fProtonAnalysisBase->Rapidity(tpcTrack->Px(),
753                                                             tpcTrack->Py(),
754                                                             tpcTrack->Pz());
755         containerInput[1] = gPt;
756         if(tpcTrack->Charge() > 0) {
757           nIdentifiedProtons += 1;
758           fProtonContainer->Fill(containerInput,kStepIdentified);
759         }//protons
760         else if(tpcTrack->Charge() < 0) {
761           nIdentifiedAntiProtons += 1;
762           fAntiProtonContainer->Fill(containerInput,kStepIdentified);
763         }//protons
764         
765         //Step: kStepIsPrimary
766         if(fProtonAnalysisBase->IsPrimary(esd,vertex,track)) {
767           if(tpcTrack->Charge() > 0)
768             fProtonContainer->Fill(containerInput,kStepIsPrimary);   
769           else if(tpcTrack->Charge() < 0) 
770             fAntiProtonContainer->Fill(containerInput,kStepIsPrimary);   
771           
772           FillQA(esd,vertex,track);       
773           
774           //Step: kStepSurvived
775           if(fProtonAnalysisBase->IsAccepted(esd,vertex,track)) {
776             ((TH2F *)(fQA2DList->At(1)))->Fill(gP,track->GetTPCsignal());
777             ((TH3F *)(fQA2DList->At(3)))->Fill(tpcTrack->Eta(),
778                                                tpcTrack->Phi()*180./TMath::Pi(),
779                                                npointsTPCdEdx);
780             ((TH3F *)(fQA2DList->At(5)))->Fill(tpcTrack->Eta(),
781                                                tpcTrack->Phi()*180./TMath::Pi(),
782                                                nClustersTPC);
783             ((TH3F *)(fQA2DList->At(7)))->Fill(gPt,
784                                                tpcTrack->Phi()*180./TMath::Pi(),
785                                                npointsTPCdEdx);
786             ((TH3F *)(fQA2DList->At(9)))->Fill(gPt,
787                                                tpcTrack->Phi()*180./TMath::Pi(),
788                                                nClustersTPC);
789             
790             if(tpcTrack->Charge() > 0) {
791               fProtonContainer->Fill(containerInput,kStepSurvived);   
792               ((TH2F *)(fQA2DList->At(10)))->Fill(tpcTrack->Eta(),
793                                                   tpcTrack->Phi()*180./TMath::Pi());
794               ((TH2F *)(fQA2DList->At(12)))->Fill(tpcTrack->Pt(),
795                                                   TMath::Abs(dca[0]));
796               ((TH2F *)(fQA2DList->At(13)))->Fill(tpcTrack->Pt(),
797                                                   TMath::Abs(dca[1]));
798             }//protons
799             else if(tpcTrack->Charge() < 0) {
800               fAntiProtonContainer->Fill(containerInput,kStepSurvived);   
801               ((TH2F *)(fQA2DList->At(11)))->Fill(tpcTrack->Eta(),
802                                                   tpcTrack->Phi()*180./TMath::Pi());
803               ((TH2F *)(fQA2DList->At(14)))->Fill(tpcTrack->Pt(),
804                                                   TMath::Abs(dca[0]));
805               ((TH2F *)(fQA2DList->At(15)))->Fill(tpcTrack->Pt(),
806                                                   TMath::Abs(dca[1]));
807             }//antiprotons
808             
809             //Step: kStepInPhaseSpace
810             if(fProtonAnalysisBase->IsInPhaseSpace(track)) {
811               if(tpcTrack->Charge() > 0) {
812                 nSurvivedProtons += 1;
813                 fHistYPtProtons->Fill(containerInput[0],
814                                       containerInput[1]);
815                 fProtonContainer->Fill(containerInput,kStepInPhaseSpace);   
816               }//protons
817               else if(tpcTrack->Charge() < 0) {
818                 nSurvivedAntiProtons += 1;
819                 fHistYPtAntiProtons->Fill(containerInput[0],
820                                           containerInput[1]);
821                 fAntiProtonContainer->Fill(containerInput,kStepInPhaseSpace);
822               }//antiprotons
823             }//Step: kStepInPhaseSpace
824           }//Step: kStepSurvived
825         }//Step: kStepIsPrimary
826       }//Step: kStepIdentified
827     }//TPC only tracks
828     else if(fProtonAnalysisBase->GetAnalysisMode() == AliProtonAnalysisBase::kGlobal) {
829       gPt = track->Pt();
830       gP = track->P();
831       track->PropagateToDCA(vertex,
832                             esd->GetMagneticField(),
833                             100.,dca,cov);
834       
835       if(fProtonAnalysisBase->IsPrimary(esd,vertex,track)) {
836         if(fProtonAnalysisBase->IsAccepted(esd,vertex,track)) {
837           ((TH2F *)(fQA2DList->At(0)))->Fill(gP,track->GetTPCsignal());
838           ((TH3F *)(fQA2DList->At(2)))->Fill(track->Eta(),
839                                              track->Phi()*180./TMath::Pi(),
840                                              npointsTPCdEdx);
841           ((TH3F *)(fQA2DList->At(4)))->Fill(track->Eta(),
842                                              track->Phi()*180./TMath::Pi(),
843                                              nClustersTPC);
844           ((TH3F *)(fQA2DList->At(6)))->Fill(gPt,
845                                              track->Phi()*180./TMath::Pi(),
846                                              npointsTPCdEdx);
847           ((TH3F *)(fQA2DList->At(8)))->Fill(gPt,
848                                              track->Phi()*180./TMath::Pi(),
849                                              nClustersTPC);     
850         }//quality cuts
851       }//primary cuts
852       
853       if(fProtonAnalysisBase->IsProton(track)) {
854         //Step: kStepIdentified
855         if(fProtonAnalysisBase->GetEtaMode())
856           containerInput[0] = track->Eta();
857         else
858           containerInput[0] = fProtonAnalysisBase->Rapidity(track->Px(),
859                                                             track->Py(),
860                                                             track->Pz());
861         containerInput[1] = gPt;
862         if(track->Charge() > 0) {
863           nIdentifiedProtons += 1;
864           fProtonContainer->Fill(containerInput,kStepIdentified);
865         }//protons
866         else if(track->Charge() < 0) {
867           nIdentifiedAntiProtons += 1;
868           fAntiProtonContainer->Fill(containerInput,kStepIdentified);
869         }//protons
870         
871         //Step: kStepIsPrimary
872         if(fProtonAnalysisBase->IsPrimary(esd,vertex,track)) {
873           if(track->Charge() > 0)
874             fProtonContainer->Fill(containerInput,kStepIsPrimary);   
875           else if(track->Charge() < 0) 
876             fAntiProtonContainer->Fill(containerInput,kStepIsPrimary);   
877           
878           FillQA(esd,vertex,track);       
879           
880           //Step: kStepSurvived
881           if(fProtonAnalysisBase->IsAccepted(esd,vertex,track)) {
882             ((TH2F *)(fQA2DList->At(1)))->Fill(gP,track->GetTPCsignal());
883             ((TH3F *)(fQA2DList->At(3)))->Fill(track->Eta(),
884                                                track->Phi()*180./TMath::Pi(),
885                                                npointsTPCdEdx);
886             ((TH3F *)(fQA2DList->At(5)))->Fill(track->Eta(),
887                                                track->Phi()*180./TMath::Pi(),
888                                                nClustersTPC);
889             ((TH3F *)(fQA2DList->At(7)))->Fill(gPt,
890                                                track->Phi()*180./TMath::Pi(),
891                                                npointsTPCdEdx);
892             ((TH3F *)(fQA2DList->At(9)))->Fill(gPt,
893                                                track->Phi()*180./TMath::Pi(),
894                                                nClustersTPC);
895             
896             if(track->Charge() > 0) {
897               fProtonContainer->Fill(containerInput,kStepSurvived);   
898               ((TH2F *)(fQA2DList->At(10)))->Fill(track->Eta(),
899                                                   track->Phi()*180./TMath::Pi());
900               ((TH2F *)(fQA2DList->At(12)))->Fill(track->Pt(),
901                                                   TMath::Abs(dca[0]));
902               ((TH2F *)(fQA2DList->At(13)))->Fill(track->Pt(),
903                                                   TMath::Abs(dca[1]));
904             }//protons
905             else if(track->Charge() < 0) {
906               fAntiProtonContainer->Fill(containerInput,kStepSurvived);   
907               ((TH2F *)(fQA2DList->At(11)))->Fill(track->Eta(),
908                                                   track->Phi()*180./TMath::Pi());
909               ((TH2F *)(fQA2DList->At(14)))->Fill(track->Pt(),
910                                                   TMath::Abs(dca[0]));
911               ((TH2F *)(fQA2DList->At(15)))->Fill(track->Pt(),
912                                                   TMath::Abs(dca[1]));
913             }//antiprotons
914             
915             //Step: kStepInPhaseSpace
916             if(fProtonAnalysisBase->IsInPhaseSpace(track)) {
917               if(track->Charge() > 0) {
918                 nSurvivedProtons += 1;
919                 fHistYPtProtons->Fill(containerInput[0],
920                                       containerInput[1]);
921                 fProtonContainer->Fill(containerInput,kStepInPhaseSpace);   
922               }//protons
923               else if(track->Charge() < 0) {
924                 nSurvivedAntiProtons += 1;
925                 fHistYPtAntiProtons->Fill(containerInput[0],
926                                           containerInput[1]);
927                 fAntiProtonContainer->Fill(containerInput,kStepInPhaseSpace);
928               }//antiprotons
929             }//Step: kStepInPhaseSpace
930           }//Step: kStepSurvived
931         }//Step: kStepIsPrimary
932       }//Step: kStepIdentified
933     }//Global tracking
934   }//track loop 
935   
936   if((nIdentifiedProtons > 0)||(nIdentifiedAntiProtons > 0))
937     fHistEvents->Fill(2); //number of analyzed events with at least one (anti)proton
938
939   if(fProtonAnalysisBase->GetDebugMode())
940     Printf("Initial number of tracks: %d | Identified (anti)protons: %d - %d | Survived (anti)protons: %d - %d",nTracks,nIdentifiedProtons,nIdentifiedAntiProtons,nSurvivedProtons,nSurvivedAntiProtons);
941 }
942
943 //____________________________________________________________________//
944 void AliProtonAnalysis::Analyze(AliAODEvent* const fAOD) {
945   //Main analysis part - AOD
946   fHistEvents->Fill(1); //number of analyzed events
947   Int_t nTracks = fAOD->GetNumberOfTracks();
948   for(Int_t iTracks = 0; iTracks < nTracks; iTracks++) {
949     AliAODTrack* track = fAOD->GetTrack(iTracks);
950     Double_t gPt = track->Pt();
951     Double_t gP = track->P();
952     
953     //pid
954     Double_t probability[10];
955     track->GetPID(probability);
956     Double_t rcc = 0.0;
957     for(Int_t i = 0; i < AliPID::kSPECIESN; i++) rcc += probability[i]*fProtonAnalysisBase->GetParticleFraction(i,gP);
958     if(rcc == 0.0) continue;
959     Double_t w[10];
960     for(Int_t i = 0; i < AliPID::kSPECIESN; i++) w[i] = probability[i]*fProtonAnalysisBase->GetParticleFraction(i,gP)/rcc;
961     Long64_t fParticleType = TMath::LocMax(AliPID::kSPECIESN,w);
962     if(fParticleType == 4) {
963       if(track->Charge() > 0) 
964         fHistYPtProtons->Fill(track->Y(fParticleType),gPt);
965       else if(track->Charge() < 0) 
966         fHistYPtAntiProtons->Fill(track->Y(fParticleType),gPt);
967     }//proton check
968   }//track loop 
969 }
970
971 //____________________________________________________________________//
972 void AliProtonAnalysis::Analyze(AliStack* const stack, 
973                                 Bool_t iInclusive) {
974   //Main analysis part - MC
975   fHistEvents->Fill(1); //number of analyzed events
976
977   Int_t nParticles = 0;
978   //inclusive protons - 
979   if(iInclusive) nParticles = stack->GetNtrack();
980   else nParticles = stack->GetNprimary();
981
982   for(Int_t i = 0; i < nParticles; i++) {
983     TParticle *particle = stack->Particle(i);
984     if(!particle) continue;
985
986     //in case of inclusive protons reject the secondaries from hadronic inter.
987     if(particle->GetUniqueID() == 13) continue;
988
989     if(TMath::Abs(particle->Eta()) > 1.0) continue;
990     if((particle->Pt() > fMaxPt)||(particle->Pt() < fMinPt)) continue;
991     if((fProtonAnalysisBase->Rapidity(particle->Px(),particle->Py(),particle->Pz()) > fMaxY)||(fProtonAnalysisBase->Rapidity(particle->Px(),particle->Py(),particle->Pz()) < fMinY)) continue;
992
993     Int_t pdgcode = particle->GetPdgCode();
994     if(pdgcode == 2212) fHistYPtProtons->Fill(fProtonAnalysisBase->Rapidity(particle->Px(),
995                                                                             particle->Py(),
996                                                                             particle->Pz()),
997                                               particle->Pt());
998     if(pdgcode == -2212) fHistYPtAntiProtons->Fill(fProtonAnalysisBase->Rapidity(particle->Px(),
999                                                                                  particle->Py(),
1000                                                                                  particle->Pz()),
1001                                                    particle->Pt());
1002   }//particle loop
1003 }
1004
1005 //____________________________________________________________________//
1006 Bool_t AliProtonAnalysis::PrintMean(TH1 *hist, Double_t edge) {
1007   //calculates the mean value of the ratio/asymmetry within \pm edge
1008   Double_t sum = 0.0;
1009   Int_t nentries = 0;
1010   //calculate the mean
1011   for(Int_t i = 0; i < hist->GetXaxis()->GetNbins(); i++) {
1012     Double_t x = hist->GetBinCenter(i+1);
1013     Double_t y = hist->GetBinContent(i+1);
1014     if(TMath::Abs(x) < edge) {
1015       sum += y;
1016       nentries += 1;
1017     }
1018   }
1019   Double_t mean = 0.0;
1020   if(nentries != 0)
1021     mean = sum/nentries;
1022
1023   //calculate the error
1024   for(Int_t i = 0; i < hist->GetXaxis()->GetNbins(); i++) {
1025     Double_t x = hist->GetBinCenter(i+1);
1026     Double_t y = hist->GetBinContent(i+1);
1027     if(TMath::Abs(x) < edge) {
1028       sum += TMath::Power((mean - y),2);
1029       nentries += 1;
1030     }
1031   }
1032
1033   Double_t error = 0.0;
1034   if(nentries != 0)
1035     error =  TMath::Sqrt(sum)/nentries;
1036
1037   cout<<"========================================="<<endl;
1038   cout<<"Input distribution: "<<hist->GetName()<<endl;
1039   cout<<"Interval used: -"<<edge<<" -> "<<edge<<endl;
1040   cout<<"Mean value :"<<mean<<endl;
1041   cout<<"Error: "<<error<<endl;
1042   cout<<"========================================="<<endl;
1043
1044   return 0;
1045 }
1046
1047 //____________________________________________________________________//
1048 Bool_t AliProtonAnalysis::PrintYields(TH1 *hist, Double_t edge) {
1049   //calculates the (anti)proton yields within the \pm edge
1050   Double_t sum = 0.0, sumerror = 0.0;
1051   Double_t error = 0.0;
1052   for(Int_t i = 0; i < hist->GetXaxis()->GetNbins(); i++) {
1053     Double_t x = hist->GetBinCenter(i+1);
1054     Double_t y = hist->GetBinContent(i+1);
1055     if(TMath::Abs(x) < edge) {
1056       sum += y;
1057       sumerror += TMath::Power(hist->GetBinError(i+1),2); 
1058     }
1059   }
1060
1061   error = TMath::Sqrt(sumerror);
1062
1063   cout<<"========================================="<<endl;
1064   cout<<"Input distribution: "<<hist->GetName()<<endl;
1065   cout<<"Interval used: -"<<edge<<" -> "<<edge<<endl;
1066   cout<<"Yields :"<<sum<<endl;
1067   cout<<"Error: "<<error<<endl;
1068   cout<<"========================================="<<endl;
1069
1070   return 0;
1071 }
1072
1073 //____________________________________________________________________//
1074 void AliProtonAnalysis::Correct(Int_t step) {
1075   //Applies the correction maps to the initial containers
1076   fCorrectProtons = new AliCFDataGrid("correctProtons",
1077                                       "corrected data",
1078                                       *fProtonContainer);
1079   fCorrectProtons->SetMeasured(0);
1080   fCorrectProtons->ApplyEffCorrection(*(AliCFEffGrid *)fEffGridListProtons->At(step));
1081
1082   fCorrectAntiProtons = new AliCFDataGrid("correctAntiProtons",
1083                                           "corrected data",
1084                                           *fAntiProtonContainer);
1085   fCorrectAntiProtons->SetMeasured(0);
1086   fCorrectAntiProtons->ApplyEffCorrection(*(AliCFEffGrid *)fEffGridListAntiProtons->At(step));
1087 }
1088
1089 //____________________________________________________________________//
1090 Bool_t AliProtonAnalysis::ReadCorrectionContainer(const char* filename) {
1091   // Reads the outout of the correction framework task
1092   // Creates the correction maps
1093   // Puts the results in the different TList objects
1094   Bool_t status = kTRUE;
1095
1096   TFile *file = TFile::Open(filename);
1097   if(!file) {
1098     cout<<"Could not find the input CORRFW file "<<filename<<endl;
1099     status = kFALSE;
1100   }
1101
1102   //________________________________________//
1103   //Protons
1104   fEffGridListProtons = new TList();
1105   fCorrectionListProtons2D = new TList(); 
1106   fEfficiencyListProtons1D = new TList(); 
1107   fCorrectionListProtons1D = new TList();
1108   
1109   AliCFContainer *corrfwContainerProtons = (AliCFContainer*) (file->Get("containerProtons"));
1110   if(!corrfwContainerProtons) {
1111     cout<<"CORRFW container for protons not found!"<<endl;
1112     status = kFALSE;
1113   }
1114   
1115   Int_t nSteps = corrfwContainerProtons->GetNStep();
1116   TH2D *gYPt[4];
1117   //currently the GRID is formed by the y-pT parameters
1118   //Add Vz as a next step
1119   Int_t iRap = 0, iPt = 1;
1120   AliCFEffGrid *effProtonsStep0Step1 = new AliCFEffGrid("eff10",
1121                                          "effProtonsStep0Step1",
1122                                          *corrfwContainerProtons);
1123   effProtonsStep0Step1->CalculateEfficiency(1,0); //eff= step1/step0
1124   fEffGridListProtons->Add(effProtonsStep0Step1);
1125   gYPt[0] = effProtonsStep0Step1->Project(iRap,iPt);
1126   fCorrectionListProtons2D->Add(gYPt[0]);
1127   
1128   AliCFEffGrid *effProtonsStep0Step2 = new AliCFEffGrid("eff20",
1129                                          "effProtonsStep0Step2",
1130                                          *corrfwContainerProtons);
1131   effProtonsStep0Step2->CalculateEfficiency(2,0); //eff= step2/step0
1132   fEffGridListProtons->Add(effProtonsStep0Step2);
1133   gYPt[1] = effProtonsStep0Step2->Project(iRap,iPt);
1134   fCorrectionListProtons2D->Add(gYPt[1]);
1135
1136   AliCFEffGrid *effProtonsStep0Step3 = new AliCFEffGrid("eff30",
1137                                          "effProtonsStep0Step3",
1138                                          *corrfwContainerProtons);
1139   effProtonsStep0Step3->CalculateEfficiency(3,0); //eff= step1/step0
1140   fEffGridListProtons->Add(effProtonsStep0Step3);
1141   gYPt[2] = effProtonsStep0Step3->Project(iRap,iPt);
1142   fCorrectionListProtons2D->Add(gYPt[2]);
1143
1144   TH1D *gEfficiency[2][3]; //efficiency as a function of pT and of y (raws-[2])
1145   TH1D *gCorrection[2][3]; //efficiency as a function of pT and of y (raws-[2])
1146   TString gTitle;
1147   //Get the projection of the efficiency maps
1148   for(Int_t iParameter = 0; iParameter < 2; iParameter++) {
1149     gEfficiency[iParameter][0] = effProtonsStep0Step1->Project(iParameter);
1150     gTitle = "ProtonsEfficiency_Parameter"; gTitle += iParameter+1;
1151     gTitle += "_Step0_Step1"; 
1152     gEfficiency[iParameter][0]->SetName(gTitle.Data());
1153     fEfficiencyListProtons1D->Add(gEfficiency[iParameter][0]);  
1154     gTitle = "ProtonsCorrection_Parameter"; gTitle += iParameter+1;
1155     gTitle += "_Step0_Step1"; 
1156     gCorrection[iParameter][0] = new TH1D(gTitle.Data(),
1157                                           gTitle.Data(),
1158                                           gEfficiency[iParameter][0]->GetNbinsX(),
1159                                           gEfficiency[iParameter][0]->GetXaxis()->GetXmin(),
1160                                           gEfficiency[iParameter][0]->GetXaxis()->GetXmax());
1161     //initialisation of the correction
1162     for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][0]->GetNbinsX(); iBin++)
1163       gCorrection[iParameter][0]->SetBinContent(iBin,1.0);
1164
1165     gEfficiency[iParameter][1] = effProtonsStep0Step2->Project(iParameter);
1166     gTitle = "ProtonsEfficiency_Parameter"; gTitle += iParameter+1;
1167     gTitle += "_Step0_Step2"; 
1168     gEfficiency[iParameter][1]->SetName(gTitle.Data());
1169     fEfficiencyListProtons1D->Add(gEfficiency[iParameter][1]);  
1170     gTitle = "ProtonsCorrection_Parameter"; gTitle += iParameter+1;
1171     gTitle += "_Step0_Step2"; 
1172     gCorrection[iParameter][1] = new TH1D(gTitle.Data(),
1173                                           gTitle.Data(),
1174                                           gEfficiency[iParameter][1]->GetNbinsX(),
1175                                           gEfficiency[iParameter][1]->GetXaxis()->GetXmin(),
1176                                           gEfficiency[iParameter][1]->GetXaxis()->GetXmax());
1177     //initialisation of the correction
1178     for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][1]->GetNbinsX(); iBin++)
1179       gCorrection[iParameter][1]->SetBinContent(iBin,1.0);
1180
1181     gEfficiency[iParameter][2] = effProtonsStep0Step3->Project(iParameter);
1182     gTitle = "ProtonsEfficiency_Parameter"; gTitle += iParameter+1;
1183     gTitle += "_Step0_Step3"; 
1184     gEfficiency[iParameter][2]->SetName(gTitle.Data());
1185     fEfficiencyListProtons1D->Add(gEfficiency[iParameter][2]);  
1186     gTitle = "ProtonsCorrection_Parameter"; gTitle += iParameter+1;
1187     gTitle += "_Step0_Step3"; 
1188     gCorrection[iParameter][2] = new TH1D(gTitle.Data(),
1189                                           gTitle.Data(),
1190                                           gEfficiency[iParameter][2]->GetNbinsX(),
1191                                           gEfficiency[iParameter][2]->GetXaxis()->GetXmin(),
1192                                           gEfficiency[iParameter][2]->GetXaxis()->GetXmax());
1193     //initialisation of the correction
1194     for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][2]->GetNbinsX(); iBin++)
1195       gCorrection[iParameter][2]->SetBinContent(iBin,1.0);
1196   }//parameter loop
1197   //Calculate the 1D correction parameters as a function of y and pT
1198   for(Int_t iParameter = 0; iParameter < 2; iParameter++) { 
1199     for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
1200       gCorrection[iParameter][iStep-1]->Divide(gEfficiency[iParameter][iStep-1]);
1201       fCorrectionListProtons1D->Add(gCorrection[iParameter][iStep-1]);  
1202     }
1203   }
1204
1205   //________________________________________//
1206   //AntiProtons
1207   fEffGridListAntiProtons = new TList();
1208   fCorrectionListAntiProtons2D = new TList(); 
1209   fEfficiencyListAntiProtons1D = new TList(); 
1210   fCorrectionListAntiProtons1D = new TList();
1211   
1212   AliCFContainer *corrfwContainerAntiProtons = (AliCFContainer*) (file->Get("containerAntiProtons"));
1213   if(!corrfwContainerAntiProtons) {
1214     cout<<"CORRFW container for antiprotons not found!"<<endl;
1215     status = kFALSE;
1216   }
1217   
1218   nSteps = corrfwContainerAntiProtons->GetNStep();
1219   //currently the GRID is formed by the y-pT parameters
1220   //Add Vz as a next step
1221   AliCFEffGrid *effAntiProtonsStep0Step1 = new AliCFEffGrid("eff10",
1222                                          "effAntiProtonsStep0Step1",
1223                                          *corrfwContainerAntiProtons);
1224   effAntiProtonsStep0Step1->CalculateEfficiency(1,0); //eff= step1/step0
1225   fEffGridListAntiProtons->Add(effAntiProtonsStep0Step1);
1226   gYPt[0] = effAntiProtonsStep0Step1->Project(iRap,iPt);
1227   fCorrectionListAntiProtons2D->Add(gYPt[0]);
1228   
1229   AliCFEffGrid *effAntiProtonsStep0Step2 = new AliCFEffGrid("eff20",
1230                                          "effAntiProtonsStep0Step2",
1231                                          *corrfwContainerAntiProtons);
1232   effAntiProtonsStep0Step2->CalculateEfficiency(2,0); //eff= step2/step0
1233   fEffGridListAntiProtons->Add(effAntiProtonsStep0Step2);
1234   gYPt[1] = effAntiProtonsStep0Step2->Project(iRap,iPt);
1235   fCorrectionListAntiProtons2D->Add(gYPt[1]);
1236
1237   AliCFEffGrid *effAntiProtonsStep0Step3 = new AliCFEffGrid("eff30",
1238                                          "effAntiProtonsStep0Step3",
1239                                          *corrfwContainerAntiProtons);
1240   effAntiProtonsStep0Step3->CalculateEfficiency(3,0); //eff= step1/step0
1241   fEffGridListAntiProtons->Add(effAntiProtonsStep0Step3);
1242   gYPt[2] = effAntiProtonsStep0Step3->Project(iRap,iPt);
1243   fCorrectionListAntiProtons2D->Add(gYPt[2]);
1244
1245   //Get the projection of the efficiency maps
1246   for(Int_t iParameter = 0; iParameter < 2; iParameter++) {
1247     gEfficiency[iParameter][0] = effAntiProtonsStep0Step1->Project(iParameter);
1248     gTitle = "AntiProtonsEfficiency_Parameter"; gTitle += iParameter+1;
1249     gTitle += "_Step0_Step1"; 
1250     gEfficiency[iParameter][0]->SetName(gTitle.Data());
1251     fEfficiencyListAntiProtons1D->Add(gEfficiency[iParameter][0]);  
1252     gTitle = "AntiProtonsCorrection_Parameter"; gTitle += iParameter+1;
1253     gTitle += "_Step0_Step1"; 
1254     gCorrection[iParameter][0] = new TH1D(gTitle.Data(),
1255                                           gTitle.Data(),
1256                                           gEfficiency[iParameter][0]->GetNbinsX(),
1257                                           gEfficiency[iParameter][0]->GetXaxis()->GetXmin(),
1258                                           gEfficiency[iParameter][0]->GetXaxis()->GetXmax());
1259     //initialisation of the correction
1260     for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][0]->GetNbinsX(); iBin++)
1261       gCorrection[iParameter][0]->SetBinContent(iBin,1.0);
1262
1263     gEfficiency[iParameter][1] = effAntiProtonsStep0Step2->Project(iParameter);
1264     gTitle = "AntiProtonsEfficiency_Parameter"; gTitle += iParameter+1;
1265     gTitle += "_Step0_Step2"; 
1266     gEfficiency[iParameter][1]->SetName(gTitle.Data());
1267     fEfficiencyListAntiProtons1D->Add(gEfficiency[iParameter][1]);  
1268     gTitle = "AntiProtonsCorrection_Parameter"; gTitle += iParameter+1;
1269     gTitle += "_Step0_Step2"; 
1270     gCorrection[iParameter][1] = new TH1D(gTitle.Data(),
1271                                           gTitle.Data(),
1272                                           gEfficiency[iParameter][1]->GetNbinsX(),
1273                                           gEfficiency[iParameter][1]->GetXaxis()->GetXmin(),
1274                                           gEfficiency[iParameter][1]->GetXaxis()->GetXmax());
1275     //initialisation of the correction
1276     for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][1]->GetNbinsX(); iBin++)
1277       gCorrection[iParameter][1]->SetBinContent(iBin,1.0);
1278
1279     gEfficiency[iParameter][2] = effAntiProtonsStep0Step3->Project(iParameter);
1280     gTitle = "AntiProtonsEfficiency_Parameter"; gTitle += iParameter+1;
1281     gTitle += "_Step0_Step3"; 
1282     gEfficiency[iParameter][2]->SetName(gTitle.Data());
1283     fEfficiencyListAntiProtons1D->Add(gEfficiency[iParameter][2]);  
1284     gTitle = "AntiProtonsCorrection_Parameter"; gTitle += iParameter+1;
1285     gTitle += "_Step0_Step3"; 
1286     gCorrection[iParameter][2] = new TH1D(gTitle.Data(),
1287                                           gTitle.Data(),
1288                                           gEfficiency[iParameter][2]->GetNbinsX(),
1289                                           gEfficiency[iParameter][2]->GetXaxis()->GetXmin(),
1290                                           gEfficiency[iParameter][2]->GetXaxis()->GetXmax());
1291     //initialisation of the correction
1292     for(Int_t iBin = 1; iBin <= gEfficiency[iParameter][2]->GetNbinsX(); iBin++)
1293       gCorrection[iParameter][2]->SetBinContent(iBin,1.0);
1294   }//parameter loop
1295   //Calculate the 1D correction parameters as a function of y and pT
1296   for(Int_t iParameter = 0; iParameter < 2; iParameter++) { 
1297     for(Int_t iStep = 1; iStep < nSteps; iStep++) { 
1298       gCorrection[iParameter][iStep-1]->Divide(gEfficiency[iParameter][iStep-1]);
1299       fCorrectionListAntiProtons1D->Add(gCorrection[iParameter][iStep-1]);  
1300     }
1301   }
1302
1303   return status;
1304 }
1305  
1306 //____________________________________________________________________//
1307 void AliProtonAnalysis::InitQA() {
1308   //Applies the correction maps to the initial containers
1309   fInitQAFlag = kTRUE;
1310   fGlobalQAList = new TList();
1311   fGlobalQAList->SetName("fGlobalQAList");
1312
1313   //========================================================//
1314   fQA2DList = new TList();
1315   fQA2DList->SetName("fQA2DList");
1316   fGlobalQAList->Add(fQA2DList);
1317
1318   //dEdx plots
1319   TH2F *gHistdEdxP = new TH2F("gHistdEdxP","dE/dx (TPC); P [GeV/c]; dE/dx [a.u]",1000,0.05,20.05,600,0,600);
1320   fQA2DList->Add(gHistdEdxP);
1321   TH2F *gHistProtonsdEdxP = new TH2F("gHistProtonsdEdxP","Accepted protons dE/dx (TPC); P [GeV/c]; dE/dx [a.u]",1000,0.05,20.05,600,0,600);
1322   fQA2DList->Add(gHistProtonsdEdxP);
1323
1324   //eta-phi-Npoints(dEdx)
1325   TH3F *gHistEtaPhiTPCdEdxNPoints = new TH3F("gHistEtaPhiTPCdEdxNPoints",
1326                                              ";#eta;#phi;N_{points}(TPC)",
1327                                              18,-0.9,0.9,
1328                                              180,0,360,
1329                                              100,0,200);
1330   gHistEtaPhiTPCdEdxNPoints->SetStats(kTRUE);
1331   fQA2DList->Add(gHistEtaPhiTPCdEdxNPoints);
1332   TH3F *gHistProtonsEtaPhiTPCdEdxNPoints = new TH3F("gHistProtonsEtaPhiTPCdEdxNPoints",
1333                                                     ";#eta;#phi;N_{points}(TPC)",
1334                                                     18,-0.9,0.9,
1335                                                     180,0,360,
1336                                                     100,0,200);
1337   gHistProtonsEtaPhiTPCdEdxNPoints->SetStats(kTRUE);
1338   fQA2DList->Add(gHistProtonsEtaPhiTPCdEdxNPoints);
1339
1340   //eta-phi-Npoints
1341   TH3F *gHistEtaPhiTPCNPoints = new TH3F("gHistEtaPhiTPCNPoints",
1342                                          ";#eta;#phi;N_{points}(TPC)",
1343                                          18,-0.9,0.9,
1344                                          180,0,360,
1345                                          100,0,200);
1346   gHistEtaPhiTPCNPoints->SetStats(kTRUE);
1347   fQA2DList->Add(gHistEtaPhiTPCNPoints);
1348   TH3F *gHistProtonsEtaPhiTPCNPoints = new TH3F("gHistProtonsEtaPhiTPCNPoints",
1349                                                 ";#eta;#phi;N_{points}(TPC)",
1350                                                 18,-0.9,0.9,
1351                                                 180,0,360,
1352                                                 100,0,200);
1353   gHistProtonsEtaPhiTPCNPoints->SetStats(kTRUE);
1354   fQA2DList->Add(gHistProtonsEtaPhiTPCNPoints);
1355
1356   //pt-phi-Npoints(dEdx)
1357   TH3F *gHistPtPhiTPCdEdxNPoints = new TH3F("gHistPtPhiTPCdEdxNPoints",
1358                                             ";P_{T} [GeV/c];#phi;N_{points}(TPC)",
1359                                             fNBinsPt,fMinPt,fMaxPt,
1360                                             180,0,360,
1361                                             100,0,200);
1362   gHistPtPhiTPCdEdxNPoints->SetStats(kTRUE);
1363   fQA2DList->Add(gHistPtPhiTPCdEdxNPoints);
1364   TH3F *gHistProtonsPtPhiTPCdEdxNPoints = new TH3F("gHistProtonsPtPhiTPCdEdxNPoints",
1365                                                     ";P_{T} [GeV/c];#phi;N_{points}(TPC)",
1366                                                     fNBinsPt,fMinPt,fMaxPt,
1367                                                     180,0,360,
1368                                                     100,0,200);
1369   gHistProtonsPtPhiTPCdEdxNPoints->SetStats(kTRUE);
1370   fQA2DList->Add(gHistProtonsPtPhiTPCdEdxNPoints);
1371
1372   //pt-phi-Npoints
1373   TH3F *gHistPtPhiTPCNPoints = new TH3F("gHistPtPhiTPCNPoints",
1374                                         ";P_{T} [GeV/c];#phi;N_{points}(TPC)",
1375                                         fNBinsPt,fMinPt,fMaxPt,
1376                                         180,0,360,
1377                                         100,0,200);
1378   gHistPtPhiTPCNPoints->SetStats(kTRUE);
1379   fQA2DList->Add(gHistPtPhiTPCNPoints);
1380   TH3F *gHistProtonsPtPhiTPCNPoints = new TH3F("gHistProtonsPtPhiTPCNPoints",
1381                                                ";P_{T} [GeV/c];#phi;N_{points}(TPC)",
1382                                                fNBinsPt,fMinPt,fMaxPt,
1383                                                180,0,360,
1384                                                100,0,200);
1385   gHistProtonsPtPhiTPCNPoints->SetStats(kTRUE);
1386   fQA2DList->Add(gHistProtonsPtPhiTPCNPoints);
1387
1388   //eta-phi for protons & antiprotons
1389   TH2F *gHistProtonsEtaPhi = new TH2F("gHistProtonsEtaPhi",
1390                                       ";#eta;#phi",
1391                                       18,-0.9,0.9,
1392                                       180,0,360);
1393   gHistProtonsEtaPhi->SetStats(kTRUE);
1394   fQA2DList->Add(gHistProtonsEtaPhi);
1395   TH2F *gHistAntiProtonsEtaPhi = new TH2F("gHistAntiProtonsEtaPhi",
1396                                           ";#eta;#phi",
1397                                           18,-0.9,0.9,
1398                                           180,0,360);
1399   gHistAntiProtonsEtaPhi->SetStats(kTRUE);
1400   fQA2DList->Add(gHistAntiProtonsEtaPhi);
1401
1402   //dca vs pT for protons & antiprotons
1403   TH2F *gHistProtonsDCAxyPt = new TH2F("gHistProtonsDCAxyPt",
1404                                        ";P_{T} [GeV/c];dca_{xy} [cm]",
1405                                        16,0.3,1.1,
1406                                        1000,0,10);
1407   gHistProtonsDCAxyPt->SetStats(kTRUE);
1408   fQA2DList->Add(gHistProtonsDCAxyPt);
1409   TH2F *gHistProtonsDCAzPt = new TH2F("gHistProtonsDCAzPt",
1410                                       ";P_{T} [GeV/c];dca_{z} [cm]",
1411                                       16,0.3,1.1,
1412                                       1000,0,10);
1413   gHistProtonsDCAzPt->SetStats(kTRUE);
1414   fQA2DList->Add(gHistProtonsDCAzPt);
1415   TH2F *gHistAntiProtonsDCAxyPt = new TH2F("gHistAntiProtonsDCAxyPt",
1416                                            ";P_{T} [GeV/c];dca_{xy} [cm]",
1417                                            16,0.3,1.1,
1418                                            1000,0,10);
1419   gHistAntiProtonsDCAxyPt->SetStats(kTRUE);
1420   fQA2DList->Add(gHistAntiProtonsDCAxyPt);
1421   TH2F *gHistAntiProtonsDCAzPt = new TH2F("gHistAntiProtonsDCAzPt",
1422                                           ";P_{T} [GeV/c];dca_{z} [cm]",
1423                                           16,0.3,1.1,
1424                                           1000,0,10);
1425   gHistAntiProtonsDCAzPt->SetStats(kTRUE);
1426   fQA2DList->Add(gHistAntiProtonsDCAzPt);
1427
1428   //========================================================//
1429   fQAProtonsAcceptedList = new TList();
1430   fQAProtonsAcceptedList->SetName("fQAProtonsAcceptedList");
1431   fGlobalQAList->Add(fQAProtonsAcceptedList);
1432   //Accepted protons
1433   TH1F *gProtonsITSClustersPass = new TH1F("gProtonsITSClustersPass",
1434                                             ";N_{clusters} (ITS);Entries",
1435                                             7,0,7);
1436   fQAProtonsAcceptedList->Add(gProtonsITSClustersPass);
1437   TH1F *gProtonsChi2PerClusterITSPass = new TH1F("gProtonsChi2PerClusterITSPass",
1438                                                  ";x^{2}/N_{clusters} (ITS);Entries",
1439                                                  100,0,4);
1440   fQAProtonsAcceptedList->Add(gProtonsChi2PerClusterITSPass);
1441   TH1F *gProtonsTPCClustersPass = new TH1F("gProtonsTPCClustersPass",
1442                                            ";N_{clusters} (TPC);Entries",
1443                                            100,0,200);
1444   fQAProtonsAcceptedList->Add(gProtonsTPCClustersPass);
1445   TH1F *gProtonsChi2PerClusterTPCPass = new TH1F("gProtonsChi2PerClusterTPCPass",
1446                                                  ";x^{2}/N_{clusters} (TPC);Entries",
1447                                                  100,0,4);
1448   fQAProtonsAcceptedList->Add(gProtonsChi2PerClusterTPCPass);
1449   TH1F *gProtonsExtCov11Pass = new TH1F("gProtonsExtCov11Pass",
1450                                         ";#sigma_{y} [cm];Entries",
1451                                         100,0,4);
1452   fQAProtonsAcceptedList->Add(gProtonsExtCov11Pass);
1453   TH1F *gProtonsExtCov22Pass = new TH1F("gProtonsExtCov22Pass",
1454                                         ";#sigma_{z} [cm];Entries",
1455                                         100,0,4);
1456   fQAProtonsAcceptedList->Add(gProtonsExtCov22Pass);
1457   TH1F *gProtonsExtCov33Pass = new TH1F("gProtonsExtCov33Pass",
1458                                         ";#sigma_{sin(#phi)};Entries",
1459                                         100,0,4);
1460   fQAProtonsAcceptedList->Add(gProtonsExtCov33Pass);
1461   TH1F *gProtonsExtCov44Pass = new TH1F("gProtonsExtCov44Pass",
1462                                         ";#sigma_{tan(#lambda)};Entries",
1463                                         100,0,4);
1464   fQAProtonsAcceptedList->Add(gProtonsExtCov44Pass);
1465   TH1F *gProtonsExtCov55Pass = new TH1F("gProtonsExtCov55Pass",
1466                                         ";#sigma_{1/P_{T}} [GeV/c]^{-1};Entries",
1467                                         100,0,4);
1468   fQAProtonsAcceptedList->Add(gProtonsExtCov55Pass);
1469   TH1F *gProtonsSigmaToVertexPass = new TH1F("gProtonsSigmaToVertexPass",
1470                                              ";#sigma_{Vertex};Entries",
1471                                              100,0,10);
1472   fQAProtonsAcceptedList->Add(gProtonsSigmaToVertexPass);
1473   TH1F *gProtonsSigmaToVertexTPCPass = new TH1F("gProtonsSigmaToVertexTPCPass",
1474                                                 ";#sigma_{Vertex};Entries",
1475                                                 100,0,10);
1476   fQAProtonsAcceptedList->Add(gProtonsSigmaToVertexTPCPass);
1477   TH1F *gProtonsDCAXYPass = new TH1F("gProtonsDCAXYPass",
1478                                      ";DCA_{xy} [cm];Entries",
1479                                      100,0,20);
1480   fQAProtonsAcceptedList->Add(gProtonsDCAXYPass);
1481   TH1F *gProtonsDCAXYTPCPass = new TH1F("gProtonsDCAXYTPCPass",
1482                                         ";DCA_{xy} [cm];Entries",
1483                                         100,0,20);
1484   fQAProtonsAcceptedList->Add(gProtonsDCAXYTPCPass);
1485   TH1F *gProtonsDCAZPass = new TH1F("gProtonsDCAZPass",
1486                                     ";DCA_{z} [cm];Entries",
1487                                     100,0,20);
1488   fQAProtonsAcceptedList->Add(gProtonsDCAZPass);
1489   TH1F *gProtonsDCAZTPCPass = new TH1F("gProtonsDCAZTPCPass",
1490                                        ";DCA_{z} [cm];Entries",
1491                                        100,0,20);
1492   fQAProtonsAcceptedList->Add(gProtonsDCAZTPCPass);
1493   TH1F *gProtonsConstrainChi2Pass = new TH1F("gProtonsConstrainChi2Pass",
1494                                              ";Log_{10}(#chi^{2});Entries",
1495                                              100,-10,10);
1496   fQAProtonsAcceptedList->Add(gProtonsConstrainChi2Pass);
1497   TH1F *gProtonsITSRefitPass = new TH1F("gProtonsITSRefitPass",
1498                                         "",10,-1,1);
1499   fQAProtonsAcceptedList->Add(gProtonsITSRefitPass);
1500   TH1F *gProtonsTPCRefitPass = new TH1F("gProtonsTPCRefitPass",
1501                                         "",10,-1,1);
1502   fQAProtonsAcceptedList->Add(gProtonsTPCRefitPass);
1503   TH1F *gProtonsESDpidPass = new TH1F("gProtonsESDpidPass",
1504                                       "",10,-1,1);
1505   fQAProtonsAcceptedList->Add(gProtonsESDpidPass);
1506   TH1F *gProtonsTPCpidPass = new TH1F("gProtonsTPCpidPass",
1507                                       "",10,-1,1);
1508   fQAProtonsAcceptedList->Add(gProtonsTPCpidPass);
1509   TH1F *gProtonsPointOnITSLayer1Pass = new TH1F("gProtonsPointOnITSLayer1Pass",
1510                                                 "",10,-1,1);
1511   fQAProtonsAcceptedList->Add(gProtonsPointOnITSLayer1Pass);
1512   TH1F *gProtonsPointOnITSLayer2Pass = new TH1F("gProtonsPointOnITSLayer2Pass",
1513                                                 "",10,-1,1);
1514   fQAProtonsAcceptedList->Add(gProtonsPointOnITSLayer2Pass);
1515   TH1F *gProtonsPointOnITSLayer3Pass = new TH1F("gProtonsPointOnITSLayer3Pass",
1516                                                 "",10,-1,1);
1517   fQAProtonsAcceptedList->Add(gProtonsPointOnITSLayer3Pass);
1518   TH1F *gProtonsPointOnITSLayer4Pass = new TH1F("gProtonsPointOnITSLayer4Pass",
1519                                                 "",10,-1,1);
1520   fQAProtonsAcceptedList->Add(gProtonsPointOnITSLayer4Pass);
1521   TH1F *gProtonsPointOnITSLayer5Pass = new TH1F("gProtonsPointOnITSLayer5Pass",
1522                                                 "",10,-1,1);
1523   fQAProtonsAcceptedList->Add(gProtonsPointOnITSLayer5Pass);
1524   TH1F *gProtonsPointOnITSLayer6Pass = new TH1F("gProtonsPointOnITSLayer6Pass",
1525                                                 "",10,-1,1);
1526   fQAProtonsAcceptedList->Add(gProtonsPointOnITSLayer6Pass);
1527   TH1F *gProtonsNumberOfTPCdEdxPointsPass = new TH1F("gProtonsNumberOfTPCdEdxPointsPass","",100,0,200);
1528   fQAProtonsAcceptedList->Add(gProtonsNumberOfTPCdEdxPointsPass);
1529
1530   //========================================================//  
1531   fQAProtonsRejectedList = new TList();
1532   fQAProtonsRejectedList->SetName("fQAProtonsRejectedList");
1533   fGlobalQAList->Add(fQAProtonsRejectedList);
1534   //Rejected protons
1535   TH1F *gProtonsITSClustersReject = new TH1F("gProtonsITSClustersReject",
1536                                              ";N_{clusters} (ITS);Entries",
1537                                              7,0,7);
1538   gProtonsITSClustersReject->SetFillColor(kRed-2);
1539   fQAProtonsRejectedList->Add(gProtonsITSClustersReject);
1540   TH1F *gProtonsChi2PerClusterITSReject = new TH1F("gProtonsChi2PerClusterITSReject",
1541                                                    ";x^{2}/N_{clusters} (ITS);Entries",
1542                                                    100,0,4);
1543   gProtonsChi2PerClusterITSReject->SetFillColor(kRed-2);
1544   fQAProtonsRejectedList->Add(gProtonsChi2PerClusterITSReject);
1545   TH1F *gProtonsTPCClustersReject = new TH1F("gProtonsTPCClustersReject",
1546                                              ";N_{clusters} (TPC);Entries",
1547                                              100,0,200);
1548   gProtonsTPCClustersReject->SetFillColor(kRed-2);
1549   fQAProtonsRejectedList->Add(gProtonsTPCClustersReject);
1550   TH1F *gProtonsChi2PerClusterTPCReject = new TH1F("gProtonsChi2PerClusterTPCReject",
1551                                                    ";x^{2}/N_{clusters} (TPC);Entries",
1552                                                    100,0,4);
1553   gProtonsChi2PerClusterTPCReject->SetFillColor(kRed-2);
1554   fQAProtonsRejectedList->Add(gProtonsChi2PerClusterTPCReject);
1555   TH1F *gProtonsExtCov11Reject = new TH1F("gProtonsExtCov11Reject",
1556                                           ";#sigma_{y} [cm];Entries",
1557                                           100,0,4);
1558   gProtonsExtCov11Reject->SetFillColor(kRed-2);
1559   fQAProtonsRejectedList->Add(gProtonsExtCov11Reject);
1560   TH1F *gProtonsExtCov22Reject = new TH1F("gProtonsExtCov22Reject",
1561                                           ";#sigma_{z} [cm];Entries",
1562                                           100,0,4);
1563   gProtonsExtCov22Reject->SetFillColor(kRed-2);
1564   fQAProtonsRejectedList->Add(gProtonsExtCov22Reject);
1565   TH1F *gProtonsExtCov33Reject = new TH1F("gProtonsExtCov33Reject",
1566                                           ";#sigma_{sin(#phi)};Entries",
1567                                           100,0,4);
1568   gProtonsExtCov33Reject->SetFillColor(kRed-2);
1569   fQAProtonsRejectedList->Add(gProtonsExtCov33Reject);
1570   TH1F *gProtonsExtCov44Reject = new TH1F("gProtonsExtCov44Reject",
1571                                           ";#sigma_{tan(#lambda)};Entries",
1572                                           100,0,4);
1573   gProtonsExtCov44Reject->SetFillColor(kRed-2);
1574   fQAProtonsRejectedList->Add(gProtonsExtCov44Reject);
1575   TH1F *gProtonsExtCov55Reject = new TH1F("gProtonsExtCov55Reject",
1576                                           ";#sigma_{1/P_{T}} [GeV/c]^{-1};Entries",
1577                                           100,0,4);
1578   gProtonsExtCov55Reject->SetFillColor(kRed-2);
1579   fQAProtonsRejectedList->Add(gProtonsExtCov55Reject);
1580   TH1F *gProtonsSigmaToVertexReject = new TH1F("gProtonsSigmaToVertexReject",
1581                                                ";#sigma_{Vertex};Entries",
1582                                                100,0,10);
1583   gProtonsSigmaToVertexReject->SetFillColor(kRed-2);
1584   fQAProtonsRejectedList->Add(gProtonsSigmaToVertexReject);
1585   TH1F *gProtonsSigmaToVertexTPCReject = new TH1F("gProtonsSigmaToVertexTPCReject",
1586                                                   ";#sigma_{Vertex};Entries",
1587                                                   100,0,10);
1588   gProtonsSigmaToVertexTPCReject->SetFillColor(kRed-2);
1589   fQAProtonsRejectedList->Add(gProtonsSigmaToVertexTPCReject);
1590   TH1F *gProtonsDCAXYReject = new TH1F("gProtonsDCAXYReject",
1591                                        ";DCA_{xy} [cm];Entries",
1592                                        100,0,20);
1593   gProtonsDCAXYReject->SetFillColor(kRed-2);
1594   fQAProtonsRejectedList->Add(gProtonsDCAXYReject);
1595   TH1F *gProtonsDCAXYTPCReject = new TH1F("gProtonsDCAXYTPCReject",
1596                                           ";DCA_{xy} [cm];Entries",
1597                                           100,0,20);
1598   gProtonsDCAXYTPCReject->SetFillColor(kRed-2);
1599   fQAProtonsRejectedList->Add(gProtonsDCAXYTPCReject);
1600   TH1F *gProtonsDCAZReject = new TH1F("gProtonsDCAZReject",
1601                                       ";DCA_{z} [cm];Entries",
1602                                       100,0,20);
1603   gProtonsDCAZReject->SetFillColor(kRed-2);
1604   fQAProtonsRejectedList->Add(gProtonsDCAZReject);
1605   TH1F *gProtonsDCAZTPCReject = new TH1F("gProtonsDCAZTPCReject",
1606                                          ";DCA_{z} [cm];Entries",
1607                                          100,0,20);
1608   gProtonsDCAZTPCReject->SetFillColor(kRed-2);
1609   fQAProtonsRejectedList->Add(gProtonsDCAZTPCReject);
1610   TH1F *gProtonsConstrainChi2Reject = new TH1F("gProtonsConstrainChi2Reject",
1611                                                ";Log_{10}(#chi^{2});Entries",
1612                                                100,-10,10);
1613   gProtonsConstrainChi2Reject->SetFillColor(kRed-2);
1614   fQAProtonsRejectedList->Add(gProtonsConstrainChi2Reject);
1615   TH1F *gProtonsITSRefitReject = new TH1F("gProtonsITSRefitReject",
1616                                           "",10,-1,1);
1617   gProtonsITSRefitReject->SetFillColor(kRed-2);
1618   fQAProtonsRejectedList->Add(gProtonsITSRefitReject);
1619   TH1F *gProtonsTPCRefitReject = new TH1F("gProtonsTPCRefitReject",
1620                                           "",10,-1,1);
1621   gProtonsTPCRefitReject->SetFillColor(kRed-2);
1622   fQAProtonsRejectedList->Add(gProtonsTPCRefitReject);
1623   TH1F *gProtonsESDpidReject = new TH1F("gProtonsESDpidReject",
1624                                         "",10,-1,1);
1625   gProtonsESDpidReject->SetFillColor(kRed-2);
1626   fQAProtonsRejectedList->Add(gProtonsESDpidReject);
1627   TH1F *gProtonsTPCpidReject = new TH1F("gProtonsTPCpidReject",
1628                                         "",10,-1,1);
1629   gProtonsTPCpidReject->SetFillColor(kRed-2);
1630   fQAProtonsRejectedList->Add(gProtonsTPCpidReject);
1631   TH1F *gProtonsPointOnITSLayer1Reject = new TH1F("gProtonsPointOnITSLayer1Reject",
1632                                                   "",10,-1,1);
1633   gProtonsPointOnITSLayer1Reject->SetFillColor(kRed-2);
1634   fQAProtonsRejectedList->Add(gProtonsPointOnITSLayer1Reject);
1635   TH1F *gProtonsPointOnITSLayer2Reject = new TH1F("gProtonsPointOnITSLayer2Reject",
1636                                                   "",10,-1,1);
1637   gProtonsPointOnITSLayer2Reject->SetFillColor(kRed-2);
1638   fQAProtonsRejectedList->Add(gProtonsPointOnITSLayer2Reject);
1639   TH1F *gProtonsPointOnITSLayer3Reject = new TH1F("gProtonsPointOnITSLayer3Reject",
1640                                                   "",10,-1,1);
1641   gProtonsPointOnITSLayer3Reject->SetFillColor(kRed-2);
1642   fQAProtonsRejectedList->Add(gProtonsPointOnITSLayer3Reject);
1643   TH1F *gProtonsPointOnITSLayer4Reject = new TH1F("gProtonsPointOnITSLayer4Reject",
1644                                                   "",10,-1,1);
1645   gProtonsPointOnITSLayer4Reject->SetFillColor(kRed-2);
1646   fQAProtonsRejectedList->Add(gProtonsPointOnITSLayer4Reject);
1647   TH1F *gProtonsPointOnITSLayer5Reject = new TH1F("gProtonsPointOnITSLayer5Reject",
1648                                                   "",10,-1,1);
1649   gProtonsPointOnITSLayer5Reject->SetFillColor(kRed-2);
1650   fQAProtonsRejectedList->Add(gProtonsPointOnITSLayer5Reject);
1651   TH1F *gProtonsPointOnITSLayer6Reject = new TH1F("gProtonsPointOnITSLayer6Reject",
1652                                                   "",10,-1,1);
1653   gProtonsPointOnITSLayer6Reject->SetFillColor(kRed-2);
1654   fQAProtonsRejectedList->Add(gProtonsPointOnITSLayer6Reject);
1655   TH1F *gProtonsNumberOfTPCdEdxPointsReject = new TH1F("gProtonsNumberOfTPCdEdxPointsReject","",100,0,200);
1656   gProtonsNumberOfTPCdEdxPointsReject->SetFillColor(kRed-2);
1657   fQAProtonsRejectedList->Add(gProtonsNumberOfTPCdEdxPointsReject);
1658     
1659   //========================================================//
1660   fQAAntiProtonsAcceptedList = new TList();
1661   fQAAntiProtonsAcceptedList->SetName("fQAAntiProtonsAcceptedList");
1662   fGlobalQAList->Add(fQAAntiProtonsAcceptedList);
1663   //Accepted antiprotons
1664   TH1F *gAntiProtonsITSClustersPass = new TH1F("gAntiProtonsITSClustersPass",
1665                                                ";N_{clusters} (ITS);Entries",
1666                                                7,0,7);
1667   fQAAntiProtonsAcceptedList->Add(gAntiProtonsITSClustersPass);
1668   TH1F *gAntiProtonsChi2PerClusterITSPass = new TH1F("gAntiProtonsChi2PerClusterITSPass",
1669                                                      ";x^{2}/N_{clusters} (ITS);Entries",
1670                                                      100,0,4);
1671   fQAAntiProtonsAcceptedList->Add(gAntiProtonsChi2PerClusterITSPass);
1672   TH1F *gAntiProtonsTPCClustersPass = new TH1F("gAntiProtonsTPCClustersPass",
1673                                                ";N_{clusters} (TPC);Entries",
1674                                                100,0,200);
1675   fQAAntiProtonsAcceptedList->Add(gAntiProtonsTPCClustersPass);
1676   TH1F *gAntiProtonsChi2PerClusterTPCPass = new TH1F("gAntiProtonsChi2PerClusterTPCPass",
1677                                                      ";x^{2}/N_{clusters} (TPC);Entries",
1678                                                      100,0,4);
1679   fQAAntiProtonsAcceptedList->Add(gAntiProtonsChi2PerClusterTPCPass);
1680   TH1F *gAntiProtonsExtCov11Pass = new TH1F("gAntiProtonsExtCov11Pass",
1681                                             ";#sigma_{y} [cm];Entries",
1682                                             100,0,4);
1683   fQAAntiProtonsAcceptedList->Add(gAntiProtonsExtCov11Pass);
1684   TH1F *gAntiProtonsExtCov22Pass = new TH1F("gAntiProtonsExtCov22Pass",
1685                                             ";#sigma_{z} [cm];Entries",
1686                                             100,0,4);
1687   fQAAntiProtonsAcceptedList->Add(gAntiProtonsExtCov22Pass);
1688   TH1F *gAntiProtonsExtCov33Pass = new TH1F("gAntiProtonsExtCov33Pass",
1689                                             ";#sigma_{sin(#phi)};Entries",
1690                                             100,0,4);
1691   fQAAntiProtonsAcceptedList->Add(gAntiProtonsExtCov33Pass);
1692   TH1F *gAntiProtonsExtCov44Pass = new TH1F("gAntiProtonsExtCov44Pass",
1693                                             ";#sigma_{tan(#lambda)};Entries",
1694                                             100,0,4);
1695   fQAAntiProtonsAcceptedList->Add(gAntiProtonsExtCov44Pass);
1696   TH1F *gAntiProtonsExtCov55Pass = new TH1F("gAntiProtonsExtCov55Pass",
1697                                             ";#sigma_{1/P_{T}} [GeV/c]^{-1};Entries",
1698                                             100,0,4);
1699   fQAAntiProtonsAcceptedList->Add(gAntiProtonsExtCov55Pass);
1700   TH1F *gAntiProtonsSigmaToVertexPass = new TH1F("gAntiProtonsSigmaToVertexPass",
1701                                                  ";#sigma_{Vertex};Entries",
1702                                                  100,0,10);
1703   fQAAntiProtonsAcceptedList->Add(gAntiProtonsSigmaToVertexPass);
1704   TH1F *gAntiProtonsSigmaToVertexTPCPass = new TH1F("gAntiProtonsSigmaToVertexTPCPass",
1705                                                     ";#sigma_{Vertex};Entries",
1706                                                     100,0,10);
1707   fQAAntiProtonsAcceptedList->Add(gAntiProtonsSigmaToVertexTPCPass);
1708   TH1F *gAntiProtonsDCAXYPass = new TH1F("gAntiProtonsDCAXYPass",
1709                                          ";DCA_{xy} [cm];Entries",
1710                                          100,0,20);
1711   fQAAntiProtonsAcceptedList->Add(gAntiProtonsDCAXYPass);
1712   TH1F *gAntiProtonsDCAXYTPCPass = new TH1F("gAntiProtonsDCAXYTPCPass",
1713                                             ";DCA_{xy} [cm];Entries",
1714                                             100,0,20);
1715   fQAAntiProtonsAcceptedList->Add(gAntiProtonsDCAXYTPCPass);
1716   TH1F *gAntiProtonsDCAZPass = new TH1F("gAntiProtonsDCAZPass",
1717                                         ";DCA_{z} [cm];Entries",
1718                                         100,0,20);
1719   fQAAntiProtonsAcceptedList->Add(gAntiProtonsDCAZPass);
1720   TH1F *gAntiProtonsDCAZTPCPass = new TH1F("gAntiProtonsDCAZTPCPass",
1721                                            ";DCA_{z} [cm];Entries",
1722                                            100,0,20);
1723   fQAAntiProtonsAcceptedList->Add(gAntiProtonsDCAZTPCPass);
1724   TH1F *gAntiProtonsConstrainChi2Pass = new TH1F("gAntiProtonsConstrainChi2Pass",
1725                                                  ";Log_{10}(#chi^{2});Entries",
1726                                                  100,-10,10);
1727   fQAAntiProtonsAcceptedList->Add(gAntiProtonsConstrainChi2Pass);
1728   TH1F *gAntiProtonsITSRefitPass = new TH1F("gAntiProtonsITSRefitPass",
1729                                             "",10,-1,1);
1730   fQAAntiProtonsAcceptedList->Add(gAntiProtonsITSRefitPass);
1731   TH1F *gAntiProtonsTPCRefitPass = new TH1F("gAntiProtonsTPCRefitPass",
1732                                             "",10,-1,1);
1733   fQAAntiProtonsAcceptedList->Add(gAntiProtonsTPCRefitPass);
1734   TH1F *gAntiProtonsESDpidPass = new TH1F("gAntiProtonsESDpidPass",
1735                                           "",10,-1,1);
1736   fQAAntiProtonsAcceptedList->Add(gAntiProtonsESDpidPass);
1737   TH1F *gAntiProtonsTPCpidPass = new TH1F("gAntiProtonsTPCpidPass",
1738                                           "",10,-1,1);
1739   fQAAntiProtonsAcceptedList->Add(gAntiProtonsTPCpidPass);
1740   TH1F *gAntiProtonsPointOnITSLayer1Pass = new TH1F("gAntiProtonsPointOnITSLayer1Pass",
1741                                                     "",10,-1,1);
1742   fQAAntiProtonsAcceptedList->Add(gAntiProtonsPointOnITSLayer1Pass);
1743   TH1F *gAntiProtonsPointOnITSLayer2Pass = new TH1F("gAntiProtonsPointOnITSLayer2Pass",
1744                                                     "",10,-1,1);
1745   fQAAntiProtonsAcceptedList->Add(gAntiProtonsPointOnITSLayer2Pass);
1746   TH1F *gAntiProtonsPointOnITSLayer3Pass = new TH1F("gAntiProtonsPointOnITSLayer3Pass",
1747                                                     "",10,-1,1);
1748   fQAAntiProtonsAcceptedList->Add(gAntiProtonsPointOnITSLayer3Pass);
1749   TH1F *gAntiProtonsPointOnITSLayer4Pass = new TH1F("gAntiProtonsPointOnITSLayer4Pass",
1750                                                     "",10,-1,1);
1751   fQAAntiProtonsAcceptedList->Add(gAntiProtonsPointOnITSLayer4Pass);
1752   TH1F *gAntiProtonsPointOnITSLayer5Pass = new TH1F("gAntiProtonsPointOnITSLayer5Pass",
1753                                                     "",10,-1,1);
1754   fQAAntiProtonsAcceptedList->Add(gAntiProtonsPointOnITSLayer5Pass);
1755   TH1F *gAntiProtonsPointOnITSLayer6Pass = new TH1F("gAntiProtonsPointOnITSLayer6Pass",
1756                                                     "",10,-1,1);
1757   fQAAntiProtonsAcceptedList->Add(gAntiProtonsPointOnITSLayer6Pass);
1758   TH1F *gAntiProtonsNumberOfTPCdEdxPointsPass = new TH1F("gAntiProtonsNumberOfTPCdEdxPointsPass","",100,0,200);
1759   fQAAntiProtonsAcceptedList->Add(gAntiProtonsNumberOfTPCdEdxPointsPass);
1760   
1761   //========================================================//
1762   fQAAntiProtonsRejectedList = new TList();
1763   fQAAntiProtonsRejectedList->SetName("fQAAntiProtonsRejectedList");
1764   fGlobalQAList->Add(fQAAntiProtonsRejectedList);
1765   //Rejected antiprotons
1766   TH1F *gAntiProtonsITSClustersReject = new TH1F("gAntiProtonsITSClustersReject",
1767                                                  ";N_{clusters} (ITS);Entries",
1768                                                  7,0,7);
1769   gAntiProtonsITSClustersReject->SetFillColor(kRed-2);
1770   fQAAntiProtonsRejectedList->Add(gAntiProtonsITSClustersReject);
1771   TH1F *gAntiProtonsChi2PerClusterITSReject = new TH1F("gAntiProtonsChi2PerClusterITSReject",
1772                                                        ";x^{2}/N_{clusters} (ITS);Entries",
1773                                                        100,0,4);
1774   gAntiProtonsChi2PerClusterITSReject->SetFillColor(kRed-2);
1775   fQAAntiProtonsRejectedList->Add(gAntiProtonsChi2PerClusterITSReject);
1776   TH1F *gAntiProtonsTPCClustersReject = new TH1F("gAntiProtonsTPCClustersReject",
1777                                                  ";N_{clusters} (TPC);Entries",
1778                                                  100,0,200);
1779   gAntiProtonsTPCClustersReject->SetFillColor(kRed-2);
1780   fQAAntiProtonsRejectedList->Add(gAntiProtonsTPCClustersReject);
1781   TH1F *gAntiProtonsChi2PerClusterTPCReject = new TH1F("gAntiProtonsChi2PerClusterTPCReject",
1782                                                        ";x^{2}/N_{clusters} (TPC);Entries",
1783                                                        100,0,4);
1784   gAntiProtonsChi2PerClusterTPCReject->SetFillColor(kRed-2);
1785   fQAAntiProtonsRejectedList->Add(gAntiProtonsChi2PerClusterTPCReject);
1786   TH1F *gAntiProtonsExtCov11Reject = new TH1F("gAntiProtonsExtCov11Reject",
1787                                               ";#sigma_{y} [cm];Entries",
1788                                               100,0,4);
1789   gAntiProtonsExtCov11Reject->SetFillColor(kRed-2);
1790   fQAAntiProtonsRejectedList->Add(gAntiProtonsExtCov11Reject);
1791   TH1F *gAntiProtonsExtCov22Reject = new TH1F("gAntiProtonsExtCov22Reject",
1792                                               ";#sigma_{z} [cm];Entries",
1793                                               100,0,4);
1794   gAntiProtonsExtCov22Reject->SetFillColor(kRed-2);
1795   fQAAntiProtonsRejectedList->Add(gAntiProtonsExtCov22Reject);
1796   TH1F *gAntiProtonsExtCov33Reject = new TH1F("gAntiProtonsExtCov33Reject",
1797                                               ";#sigma_{sin(#phi)};Entries",
1798                                               100,0,4);
1799   gAntiProtonsExtCov33Reject->SetFillColor(kRed-2);
1800   fQAAntiProtonsRejectedList->Add(gAntiProtonsExtCov33Reject);
1801   TH1F *gAntiProtonsExtCov44Reject = new TH1F("gAntiProtonsExtCov44Reject",
1802                                               ";#sigma_{tan(#lambda)};Entries",
1803                                               100,0,4);
1804   gAntiProtonsExtCov44Reject->SetFillColor(kRed-2);
1805   fQAAntiProtonsRejectedList->Add(gAntiProtonsExtCov44Reject);
1806   TH1F *gAntiProtonsExtCov55Reject = new TH1F("gAntiProtonsExtCov55Reject",
1807                                               ";#sigma_{1/P_{T}} [GeV/c]^{-1};Entries",
1808                                               100,0,4);
1809   gAntiProtonsExtCov55Reject->SetFillColor(kRed-2);
1810   fQAAntiProtonsRejectedList->Add(gAntiProtonsExtCov55Reject);
1811   TH1F *gAntiProtonsSigmaToVertexReject = new TH1F("gAntiProtonsSigmaToVertexReject",
1812                                                    ";#sigma_{Vertex};Entries",
1813                                                    100,0,10);
1814   gAntiProtonsSigmaToVertexReject->SetFillColor(kRed-2);
1815   fQAAntiProtonsRejectedList->Add(gAntiProtonsSigmaToVertexReject);
1816   TH1F *gAntiProtonsSigmaToVertexTPCReject = new TH1F("gAntiProtonsSigmaToVertexTPCReject",
1817                                                       ";#sigma_{Vertex};Entries",
1818                                                       100,0,10);
1819   gAntiProtonsSigmaToVertexTPCReject->SetFillColor(kRed-2);
1820   fQAAntiProtonsRejectedList->Add(gAntiProtonsSigmaToVertexTPCReject);
1821   TH1F *gAntiProtonsDCAXYReject = new TH1F("gAntiProtonsDCAXYReject",
1822                                            ";DCA_{xy} [cm];Entries",
1823                                            100,0,20);
1824   gAntiProtonsDCAXYReject->SetFillColor(kRed-2);
1825   fQAAntiProtonsRejectedList->Add(gAntiProtonsDCAXYReject);
1826   TH1F *gAntiProtonsDCAXYTPCReject = new TH1F("gAntiProtonsDCAXYTPCReject",
1827                                               ";DCA_{xy} [cm];Entries",
1828                                               100,0,20);
1829   gAntiProtonsDCAXYTPCReject->SetFillColor(kRed-2);
1830   fQAAntiProtonsRejectedList->Add(gAntiProtonsDCAXYTPCReject);
1831   TH1F *gAntiProtonsDCAZReject = new TH1F("gAntiProtonsDCAZReject",
1832                                           ";DCA_{z} [cm];Entries",
1833                                           100,0,20);
1834   gAntiProtonsDCAZReject->SetFillColor(kRed-2);
1835   fQAAntiProtonsRejectedList->Add(gAntiProtonsDCAZReject);
1836   TH1F *gAntiProtonsDCAZTPCReject = new TH1F("gAntiProtonsDCAZTPCReject",
1837                                              ";DCA_{z} [cm];Entries",
1838                                              100,0,20);
1839   gAntiProtonsDCAZTPCReject->SetFillColor(kRed-2);
1840   fQAAntiProtonsRejectedList->Add(gAntiProtonsDCAZTPCReject);
1841   TH1F *gAntiProtonsConstrainChi2Reject = new TH1F("gAntiProtonsConstrainChi2Reject",
1842                                                    ";Log_{10}(#chi^{2});Entries",
1843                                                    100,-10,10);
1844   gAntiProtonsConstrainChi2Reject->SetFillColor(kRed-2);
1845   fQAAntiProtonsRejectedList->Add(gAntiProtonsConstrainChi2Reject);
1846   TH1F *gAntiProtonsITSRefitReject = new TH1F("gAntiProtonsITSRefitReject",
1847                                               "",10,-1,1);
1848   gAntiProtonsITSRefitReject->SetFillColor(kRed-2);
1849   fQAAntiProtonsRejectedList->Add(gAntiProtonsITSRefitReject);
1850   TH1F *gAntiProtonsTPCRefitReject = new TH1F("gAntiProtonsTPCRefitReject",
1851                                               "",10,-1,1);
1852   gAntiProtonsTPCRefitReject->SetFillColor(kRed-2);
1853   fQAAntiProtonsRejectedList->Add(gAntiProtonsTPCRefitReject);
1854   TH1F *gAntiProtonsESDpidReject = new TH1F("gAntiProtonsESDpidReject",
1855                                             "",10,-1,1);
1856   gAntiProtonsESDpidReject->SetFillColor(kRed-2);
1857   fQAAntiProtonsRejectedList->Add(gAntiProtonsESDpidReject);
1858   TH1F *gAntiProtonsTPCpidReject = new TH1F("gAntiProtonsTPCpidReject",
1859                                             "",10,-1,1);
1860   gAntiProtonsTPCpidReject->SetFillColor(kRed-2);
1861   fQAAntiProtonsRejectedList->Add(gAntiProtonsTPCpidReject);
1862   TH1F *gAntiProtonsPointOnITSLayer1Reject = new TH1F("gAntiProtonsPointOnITSLayer1Reject",
1863                                                       "",10,-1,1);
1864   gAntiProtonsPointOnITSLayer1Reject->SetFillColor(kRed-2);
1865   fQAAntiProtonsRejectedList->Add(gAntiProtonsPointOnITSLayer1Reject);
1866   TH1F *gAntiProtonsPointOnITSLayer2Reject = new TH1F("gAntiProtonsPointOnITSLayer2Reject",
1867                                                       "",10,-1,1);
1868   gAntiProtonsPointOnITSLayer2Reject->SetFillColor(kRed-2);
1869   fQAAntiProtonsRejectedList->Add(gAntiProtonsPointOnITSLayer2Reject);
1870   TH1F *gAntiProtonsPointOnITSLayer3Reject = new TH1F("gAntiProtonsPointOnITSLayer3Reject",
1871                                                       "",10,-1,1);
1872   gAntiProtonsPointOnITSLayer3Reject->SetFillColor(kRed-2);
1873   fQAAntiProtonsRejectedList->Add(gAntiProtonsPointOnITSLayer3Reject);
1874   TH1F *gAntiProtonsPointOnITSLayer4Reject = new TH1F("gAntiProtonsPointOnITSLayer4Reject",
1875                                                       "",10,-1,1);
1876   gAntiProtonsPointOnITSLayer4Reject->SetFillColor(kRed-2);
1877   fQAAntiProtonsRejectedList->Add(gAntiProtonsPointOnITSLayer4Reject);
1878   TH1F *gAntiProtonsPointOnITSLayer5Reject = new TH1F("gAntiProtonsPointOnITSLayer5Reject",
1879                                                       "",10,-1,1);
1880   gAntiProtonsPointOnITSLayer5Reject->SetFillColor(kRed-2);
1881   fQAAntiProtonsRejectedList->Add(gAntiProtonsPointOnITSLayer5Reject);
1882   TH1F *gAntiProtonsPointOnITSLayer6Reject = new TH1F("gAntiProtonsPointOnITSLayer6Reject",
1883                                                       "",10,-1,1);
1884   gAntiProtonsPointOnITSLayer6Reject->SetFillColor(kRed-2);
1885   fQAAntiProtonsRejectedList->Add(gAntiProtonsPointOnITSLayer6Reject);
1886   TH1F *gAntiProtonsNumberOfTPCdEdxPointsReject = new TH1F("gAntiProtonsNumberOfTPCdEdxPointsReject","",100,0,200);
1887   gAntiProtonsNumberOfTPCdEdxPointsReject->SetFillColor(kRed-2);
1888   fQAAntiProtonsRejectedList->Add(gAntiProtonsNumberOfTPCdEdxPointsReject);
1889 }
1890
1891 //____________________________________________________________________//
1892 void AliProtonAnalysis::FillQA(AliESDEvent *esd,
1893                                const AliESDVertex *vertex, 
1894                                AliESDtrack* track) {
1895   //Fills the QA histograms
1896   Double_t gPt = 0.0, gPx = 0.0, gPy = 0.0, gPz = 0.0;
1897   Double_t dca[2] = {0.0,0.0}, cov[3] = {0.0,0.0,0.0};  //The impact parameters and their covariance.
1898
1899   if((fProtonAnalysisBase->GetAnalysisMode()==AliProtonAnalysisBase::kTPC)||(fProtonAnalysisBase->GetAnalysisMode()==AliProtonAnalysisBase::kHybrid)) {
1900     AliExternalTrackParam *tpcTrack = (AliExternalTrackParam *)track->GetTPCInnerParam();
1901     if(!tpcTrack) {
1902       gPt = 0.0; gPx = 0.0; gPy = 0.0; gPz = 0.0;
1903       dca[0] = -100.; dca[1] = -100.;
1904       cov[0] = -100.; cov[1] = -100.; cov[2] = -100.;
1905     }
1906     else {
1907       gPt = tpcTrack->Pt();
1908       gPx = tpcTrack->Px();
1909       gPy = tpcTrack->Py();
1910       gPz = tpcTrack->Pz();
1911       tpcTrack->PropagateToDCA(vertex,
1912                                esd->GetMagneticField(),
1913                                100.,dca,cov);
1914     }
1915   }
1916   else{
1917     gPt = track->Pt();
1918     gPx = track->Px();
1919     gPy = track->Py();
1920     gPz = track->Pz();
1921     track->PropagateToDCA(vertex,
1922                           esd->GetMagneticField(),
1923                           100.,dca,cov);
1924   }
1925
1926   Int_t  fIdxInt[200];
1927   Int_t nClustersITS = track->GetITSclusters(fIdxInt);
1928   Int_t nClustersTPC = track->GetTPCclusters(fIdxInt);
1929
1930   Float_t chi2PerClusterITS = -1;
1931   if (nClustersITS!=0)
1932     chi2PerClusterITS = track->GetITSchi2()/Float_t(nClustersITS);
1933   Float_t chi2PerClusterTPC = -1;
1934   if (nClustersTPC!=0)
1935     chi2PerClusterTPC = track->GetTPCchi2()/Float_t(nClustersTPC);
1936
1937   Double_t extCov[15];
1938   track->GetExternalCovariance(extCov);
1939   
1940   //protons
1941   if(track->Charge() > 0) {
1942     if(fProtonAnalysisBase->IsUsedMinITSClusters()) {
1943       if(nClustersITS < fProtonAnalysisBase->GetMinITSClusters()) {
1944         ((TH1F *)(fQAProtonsRejectedList->At(0)))->Fill(nClustersITS);
1945       }
1946       else if(nClustersITS >= fProtonAnalysisBase->GetMinITSClusters()) 
1947         ((TH1F *)(fQAProtonsAcceptedList->At(0)))->Fill(nClustersITS);
1948     }//ITS clusters
1949     if(fProtonAnalysisBase->IsUsedMaxChi2PerITSCluster()) {
1950       if(chi2PerClusterITS > fProtonAnalysisBase->GetMaxChi2PerITSCluster()) {
1951         ((TH1F *)(fQAProtonsRejectedList->At(1)))->Fill(chi2PerClusterITS);
1952       }
1953       else if(chi2PerClusterITS <= fProtonAnalysisBase->GetMaxChi2PerITSCluster())
1954         ((TH1F *)(fQAProtonsAcceptedList->At(1)))->Fill(chi2PerClusterITS);
1955     }//chi2 per ITS cluster
1956     if(fProtonAnalysisBase->IsUsedMinTPCClusters()) {
1957       if(nClustersTPC < fProtonAnalysisBase->GetMinTPCClusters()) {
1958         ((TH1F *)(fQAProtonsRejectedList->At(2)))->Fill(nClustersTPC);
1959       }
1960       else if(nClustersTPC >= fProtonAnalysisBase->GetMinTPCClusters()) {
1961         ((TH1F *)(fQAProtonsAcceptedList->At(2)))->Fill(nClustersTPC);
1962       }
1963     }//TPC clusters
1964     if(fProtonAnalysisBase->IsUsedMaxChi2PerTPCCluster()) {
1965       if(chi2PerClusterTPC > fProtonAnalysisBase->GetMaxChi2PerTPCCluster()) {
1966         ((TH1F *)(fQAProtonsRejectedList->At(3)))->Fill(chi2PerClusterTPC);
1967       }
1968       else if(chi2PerClusterTPC <= fProtonAnalysisBase->GetMaxChi2PerTPCCluster())
1969         ((TH1F *)(fQAProtonsAcceptedList->At(3)))->Fill(chi2PerClusterTPC);
1970     }//chi2 per TPC cluster
1971     if(fProtonAnalysisBase->IsUsedMaxCov11()) {
1972       if(extCov[0] > fProtonAnalysisBase->GetMaxCov11()) {
1973         ((TH1F *)(fQAProtonsRejectedList->At(4)))->Fill(extCov[0]);
1974       }
1975       else if(extCov[0] <= fProtonAnalysisBase->GetMaxCov11())
1976         ((TH1F *)(fQAProtonsAcceptedList->At(4)))->Fill(extCov[0]);
1977     }//cov11
1978     if(fProtonAnalysisBase->IsUsedMaxCov22()) {
1979       if(extCov[2] > fProtonAnalysisBase->GetMaxCov22()) {
1980         ((TH1F *)(fQAProtonsRejectedList->At(5)))->Fill(extCov[2]);
1981       }
1982       else if(extCov[2] <= fProtonAnalysisBase->GetMaxCov22())
1983         ((TH1F *)(fQAProtonsAcceptedList->At(5)))->Fill(extCov[2]);
1984     }//cov11
1985     if(fProtonAnalysisBase->IsUsedMaxCov33()) {
1986       if(extCov[5] > fProtonAnalysisBase->GetMaxCov33()) {
1987         ((TH1F *)(fQAProtonsRejectedList->At(6)))->Fill(extCov[5]);
1988       }
1989       else if(extCov[5] <= fProtonAnalysisBase->GetMaxCov33())
1990         ((TH1F *)(fQAProtonsAcceptedList->At(6)))->Fill(extCov[5]);
1991     }//cov11
1992     if(fProtonAnalysisBase->IsUsedMaxCov44()) {
1993       if(extCov[9] > fProtonAnalysisBase->GetMaxCov44()) {
1994         ((TH1F *)(fQAProtonsRejectedList->At(7)))->Fill(extCov[9]);
1995       }
1996       else if(extCov[9] <= fProtonAnalysisBase->GetMaxCov44())
1997         ((TH1F *)(fQAProtonsAcceptedList->At(7)))->Fill(extCov[9]);
1998     }//cov11
1999     if(fProtonAnalysisBase->IsUsedMaxCov55()) {
2000       if(extCov[14] > fProtonAnalysisBase->GetMaxCov55()) {
2001         ((TH1F *)(fQAProtonsRejectedList->At(8)))->Fill(extCov[14]);
2002       }
2003       else if(extCov[14] <= fProtonAnalysisBase->GetMaxCov55())
2004         ((TH1F *)(fQAProtonsAcceptedList->At(8)))->Fill(extCov[14]);
2005     }//cov55
2006     if(fProtonAnalysisBase->IsUsedMaxSigmaToVertex()) {
2007       if(fProtonAnalysisBase->GetSigmaToVertex(track) > fProtonAnalysisBase->GetMaxSigmaToVertex()) {
2008         ((TH1F *)(fQAProtonsRejectedList->At(9)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2009       }
2010       else if(fProtonAnalysisBase->GetSigmaToVertex(track) <= fProtonAnalysisBase->GetMaxSigmaToVertex())
2011         ((TH1F *)(fQAProtonsAcceptedList->At(9)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2012     }//sigma to vertex
2013     if(fProtonAnalysisBase->IsUsedMaxSigmaToVertexTPC()) {
2014       if(fProtonAnalysisBase->GetSigmaToVertex(track) > fProtonAnalysisBase->GetMaxSigmaToVertexTPC()) {
2015         ((TH1F *)(fQAProtonsRejectedList->At(10)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2016       }
2017       else if(fProtonAnalysisBase->GetSigmaToVertex(track) <= fProtonAnalysisBase->GetMaxSigmaToVertexTPC())
2018         ((TH1F *)(fQAProtonsAcceptedList->At(10)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2019     }//sigma to vertex TPC
2020     if(fProtonAnalysisBase->IsUsedMaxDCAXY()) {
2021       if(TMath::Abs(dca[0]) > fProtonAnalysisBase->GetMaxDCAXY()) {
2022         ((TH1F *)(fQAProtonsRejectedList->At(11)))->Fill(TMath::Abs(dca[0]));
2023       }
2024       else if(TMath::Abs(dca[0]) <= fProtonAnalysisBase->GetMaxDCAXY())
2025         ((TH1F *)(fQAProtonsAcceptedList->At(11)))->Fill(TMath::Abs(dca[0]));
2026     }//DCA xy global tracking
2027     if(fProtonAnalysisBase->IsUsedMaxDCAXYTPC()) {
2028       if(TMath::Abs(dca[0]) > fProtonAnalysisBase->GetMaxDCAXYTPC()) {
2029         ((TH1F *)(fQAProtonsRejectedList->At(12)))->Fill(TMath::Abs(dca[0]));
2030       }
2031       else if(TMath::Abs(dca[0]) <= fProtonAnalysisBase->GetMaxDCAXYTPC())
2032         ((TH1F *)(fQAProtonsAcceptedList->At(12)))->Fill(TMath::Abs(dca[0]));
2033     }//DCA xy TPC tracking
2034     if(fProtonAnalysisBase->IsUsedMaxDCAZ()) {
2035       if(TMath::Abs(dca[1]) > fProtonAnalysisBase->GetMaxDCAZ()) {
2036         ((TH1F *)(fQAProtonsRejectedList->At(13)))->Fill(TMath::Abs(dca[1]));
2037       }
2038       else if(TMath::Abs(dca[1]) <= fProtonAnalysisBase->GetMaxDCAZ())
2039         ((TH1F *)(fQAProtonsAcceptedList->At(13)))->Fill(TMath::Abs(dca[1]));
2040     }//DCA z global tracking
2041     if(fProtonAnalysisBase->IsUsedMaxDCAZTPC()) {
2042       if(TMath::Abs(dca[1]) > fProtonAnalysisBase->GetMaxDCAZTPC()) {
2043         ((TH1F *)(fQAProtonsRejectedList->At(14)))->Fill(TMath::Abs(dca[1]));
2044       }
2045       else if(TMath::Abs(dca[1]) <= fProtonAnalysisBase->GetMaxDCAZTPC())
2046         ((TH1F *)(fQAProtonsAcceptedList->At(14)))->Fill(TMath::Abs(dca[1]));
2047     }//DCA z TPC tracking
2048     if(fProtonAnalysisBase->IsUsedMaxConstrainChi2()) {
2049       if(track->GetConstrainedChi2() > 0) {
2050         if(TMath::Log(track->GetConstrainedChi2()) > fProtonAnalysisBase->GetMaxConstrainChi2()) {
2051           ((TH1F *)(fQAProtonsRejectedList->At(15)))->Fill(TMath::Log(track->GetConstrainedChi2()));
2052         }
2053         else if(TMath::Log(track->GetConstrainedChi2()) <= fProtonAnalysisBase->GetMaxConstrainChi2())
2054           ((TH1F *)(fQAProtonsAcceptedList->At(15)))->Fill(TMath::Log(track->GetConstrainedChi2()));
2055       }
2056     }//constrain chi2 - vertex
2057     if(fProtonAnalysisBase->IsUsedITSRefit()) {
2058       if ((track->GetStatus() & AliESDtrack::kITSrefit) == 0) {
2059         ((TH1F *)(fQAProtonsRejectedList->At(16)))->Fill(0);
2060       }
2061       else if((track->GetStatus() & AliESDtrack::kITSrefit) != 0)
2062         ((TH1F *)(fQAProtonsAcceptedList->At(16)))->Fill(0);
2063     }//ITS refit
2064     if(fProtonAnalysisBase->IsUsedTPCRefit()) {
2065       if ((track->GetStatus() & AliESDtrack::kTPCrefit) == 0) {
2066         ((TH1F *)(fQAProtonsRejectedList->At(17)))->Fill(0);
2067       }
2068       else if((track->GetStatus() & AliESDtrack::kTPCrefit) != 0)
2069         ((TH1F *)(fQAProtonsAcceptedList->At(17)))->Fill(0);
2070     }//TPC refit
2071     if(fProtonAnalysisBase->IsUsedESDpid()) {
2072       if ((track->GetStatus() & AliESDtrack::kESDpid) == 0) {
2073         ((TH1F *)(fQAProtonsRejectedList->At(18)))->Fill(0);
2074       }
2075       else if((track->GetStatus() & AliESDtrack::kESDpid) != 0)
2076         ((TH1F *)(fQAProtonsAcceptedList->At(18)))->Fill(0);
2077     }//ESD pid
2078     if(fProtonAnalysisBase->IsUsedTPCpid()) {
2079       if ((track->GetStatus() & AliESDtrack::kTPCpid) == 0) {
2080         ((TH1F *)(fQAProtonsRejectedList->At(19)))->Fill(0);
2081       }
2082       else if((track->GetStatus() & AliESDtrack::kTPCpid) != 0)
2083         ((TH1F *)(fQAProtonsAcceptedList->At(19)))->Fill(0);
2084     }//TPC pid
2085     if(fProtonAnalysisBase->IsUsedPointOnITSLayer1()) {
2086       if(!track->HasPointOnITSLayer(0)) {
2087         ((TH1F *)(fQAProtonsRejectedList->At(20)))->Fill(0);
2088       }
2089       else if(track->HasPointOnITSLayer(0))
2090         ((TH1F *)(fQAProtonsAcceptedList->At(20)))->Fill(0);
2091     }//point on SPD1
2092     if(fProtonAnalysisBase->IsUsedPointOnITSLayer2()) {
2093       if(!track->HasPointOnITSLayer(1)) {
2094         ((TH1F *)(fQAProtonsRejectedList->At(21)))->Fill(0);
2095       }
2096       else if(track->HasPointOnITSLayer(1))
2097         ((TH1F *)(fQAProtonsAcceptedList->At(21)))->Fill(0);
2098     }//point on SPD2
2099     if(fProtonAnalysisBase->IsUsedPointOnITSLayer3()) {
2100       if(!track->HasPointOnITSLayer(2)) {
2101         ((TH1F *)(fQAProtonsRejectedList->At(22)))->Fill(0);
2102       }
2103       else if(track->HasPointOnITSLayer(2))
2104         ((TH1F *)(fQAProtonsAcceptedList->At(22)))->Fill(0);
2105     }//point on SDD1
2106     if(fProtonAnalysisBase->IsUsedPointOnITSLayer4()) {
2107       if(!track->HasPointOnITSLayer(3)) {
2108         ((TH1F *)(fQAProtonsRejectedList->At(23)))->Fill(0);
2109       }
2110       else if(track->HasPointOnITSLayer(3))
2111         ((TH1F *)(fQAProtonsAcceptedList->At(23)))->Fill(0);
2112     }//point on SDD2
2113     if(fProtonAnalysisBase->IsUsedPointOnITSLayer5()) {
2114       if(!track->HasPointOnITSLayer(4)) {
2115         ((TH1F *)(fQAProtonsRejectedList->At(24)))->Fill(0);
2116       }
2117       else if(track->HasPointOnITSLayer(4))
2118         ((TH1F *)(fQAProtonsAcceptedList->At(24)))->Fill(0);
2119     }//point on SSD1
2120     if(fProtonAnalysisBase->IsUsedPointOnITSLayer6()) {
2121       if(!track->HasPointOnITSLayer(5)) {
2122         ((TH1F *)(fQAProtonsRejectedList->At(25)))->Fill(0);
2123       }
2124       else if(track->HasPointOnITSLayer(5))
2125         ((TH1F *)(fQAProtonsAcceptedList->At(25)))->Fill(0);
2126     }//point on SSD2
2127     if(fProtonAnalysisBase->IsUsedMinTPCdEdxPoints()) {
2128       if(track->GetTPCsignalN() < fProtonAnalysisBase->GetMinTPCdEdxPoints()) {
2129         ((TH1F *)(fQAProtonsRejectedList->At(26)))->Fill(track->GetTPCsignalN());
2130       }
2131       if(track->GetTPCsignalN() >= fProtonAnalysisBase->GetMinTPCdEdxPoints())
2132         ((TH1F *)(fQAProtonsAcceptedList->At(26)))->Fill(track->GetTPCsignalN());
2133     }//number of TPC points for the dE/dx
2134   }//protons
2135
2136   //antiprotons
2137   if(track->Charge() < 0) {
2138     if(fProtonAnalysisBase->IsUsedMinITSClusters()) {
2139       if(nClustersITS < fProtonAnalysisBase->GetMinITSClusters()) {
2140         ((TH1F *)(fQAAntiProtonsRejectedList->At(0)))->Fill(nClustersITS);
2141       }
2142       else if(nClustersITS >= fProtonAnalysisBase->GetMinITSClusters()) 
2143         ((TH1F *)(fQAAntiProtonsAcceptedList->At(0)))->Fill(nClustersITS);
2144     }//ITS clusters
2145     if(fProtonAnalysisBase->IsUsedMaxChi2PerITSCluster()) {
2146       if(chi2PerClusterITS > fProtonAnalysisBase->GetMaxChi2PerITSCluster()) {
2147         ((TH1F *)(fQAAntiProtonsRejectedList->At(1)))->Fill(chi2PerClusterITS);
2148       }
2149       else if(chi2PerClusterITS <= fProtonAnalysisBase->GetMaxChi2PerITSCluster())
2150         ((TH1F *)(fQAAntiProtonsAcceptedList->At(1)))->Fill(chi2PerClusterITS);
2151     }//chi2 per ITS cluster
2152     if(fProtonAnalysisBase->IsUsedMinTPCClusters()) {
2153       if(nClustersTPC < fProtonAnalysisBase->GetMinTPCClusters()) {
2154         ((TH1F *)(fQAAntiProtonsRejectedList->At(2)))->Fill(nClustersTPC);
2155       }
2156       else if(nClustersTPC >= fProtonAnalysisBase->GetMinTPCClusters()) {
2157         ((TH1F *)(fQAAntiProtonsAcceptedList->At(2)))->Fill(nClustersTPC);
2158       }
2159     }//TPC clusters
2160     if(fProtonAnalysisBase->IsUsedMaxChi2PerTPCCluster()) {
2161       if(chi2PerClusterTPC > fProtonAnalysisBase->GetMaxChi2PerTPCCluster()) {
2162         ((TH1F *)(fQAAntiProtonsRejectedList->At(3)))->Fill(chi2PerClusterTPC);
2163       }
2164       else if(chi2PerClusterTPC <= fProtonAnalysisBase->GetMaxChi2PerTPCCluster())
2165         ((TH1F *)(fQAAntiProtonsAcceptedList->At(3)))->Fill(chi2PerClusterTPC);
2166     }//chi2 per TPC cluster
2167     if(fProtonAnalysisBase->IsUsedMaxCov11()) {
2168       if(extCov[0] > fProtonAnalysisBase->GetMaxCov11()) {
2169         ((TH1F *)(fQAAntiProtonsRejectedList->At(4)))->Fill(extCov[0]);
2170       }
2171       else if(extCov[0] <= fProtonAnalysisBase->GetMaxCov11())
2172         ((TH1F *)(fQAAntiProtonsAcceptedList->At(4)))->Fill(extCov[0]);
2173     }//cov11
2174     if(fProtonAnalysisBase->IsUsedMaxCov22()) {
2175       if(extCov[2] > fProtonAnalysisBase->GetMaxCov22()) {
2176         ((TH1F *)(fQAAntiProtonsRejectedList->At(5)))->Fill(extCov[2]);
2177       }
2178       else if(extCov[2] <= fProtonAnalysisBase->GetMaxCov22())
2179         ((TH1F *)(fQAAntiProtonsAcceptedList->At(5)))->Fill(extCov[2]);
2180     }//cov11
2181     if(fProtonAnalysisBase->IsUsedMaxCov33()) {
2182       if(extCov[5] > fProtonAnalysisBase->GetMaxCov33()) {
2183         ((TH1F *)(fQAAntiProtonsRejectedList->At(6)))->Fill(extCov[5]);
2184       }
2185       else if(extCov[5] <= fProtonAnalysisBase->GetMaxCov33())
2186         ((TH1F *)(fQAAntiProtonsAcceptedList->At(6)))->Fill(extCov[5]);
2187     }//cov11
2188     if(fProtonAnalysisBase->IsUsedMaxCov44()) {
2189       if(extCov[9] > fProtonAnalysisBase->GetMaxCov44()) {
2190         ((TH1F *)(fQAAntiProtonsRejectedList->At(7)))->Fill(extCov[9]);
2191       }
2192       else if(extCov[9] <= fProtonAnalysisBase->GetMaxCov44())
2193         ((TH1F *)(fQAAntiProtonsAcceptedList->At(7)))->Fill(extCov[9]);
2194     }//cov11
2195     if(fProtonAnalysisBase->IsUsedMaxCov55()) {
2196       if(extCov[14] > fProtonAnalysisBase->GetMaxCov55()) {
2197         ((TH1F *)(fQAAntiProtonsRejectedList->At(8)))->Fill(extCov[14]);
2198       }
2199       else if(extCov[14] <= fProtonAnalysisBase->GetMaxCov55())
2200         ((TH1F *)(fQAAntiProtonsAcceptedList->At(8)))->Fill(extCov[14]);
2201     }//cov55
2202     if(fProtonAnalysisBase->IsUsedMaxSigmaToVertex()) {
2203       if(fProtonAnalysisBase->GetSigmaToVertex(track) > fProtonAnalysisBase->GetMaxSigmaToVertex()) {
2204         ((TH1F *)(fQAAntiProtonsRejectedList->At(9)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2205       }
2206       else if(fProtonAnalysisBase->GetSigmaToVertex(track) <= fProtonAnalysisBase->GetMaxSigmaToVertex())
2207         ((TH1F *)(fQAAntiProtonsAcceptedList->At(9)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2208     }//sigma to vertex
2209     if(fProtonAnalysisBase->IsUsedMaxSigmaToVertexTPC()) {
2210       if(fProtonAnalysisBase->GetSigmaToVertex(track) > fProtonAnalysisBase->GetMaxSigmaToVertexTPC()) {
2211         ((TH1F *)(fQAAntiProtonsRejectedList->At(10)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2212       }
2213       else if(fProtonAnalysisBase->GetSigmaToVertex(track) <= fProtonAnalysisBase->GetMaxSigmaToVertexTPC())
2214         ((TH1F *)(fQAAntiProtonsAcceptedList->At(10)))->Fill(fProtonAnalysisBase->GetSigmaToVertex(track));
2215     }//sigma to vertex TPC
2216     if(fProtonAnalysisBase->IsUsedMaxDCAXY()) {
2217       if(TMath::Abs(dca[0]) > fProtonAnalysisBase->GetMaxDCAXY()) {
2218         ((TH1F *)(fQAAntiProtonsRejectedList->At(11)))->Fill(TMath::Abs(dca[0]));
2219       }
2220       else if(TMath::Abs(dca[0]) <= fProtonAnalysisBase->GetMaxDCAXY())
2221         ((TH1F *)(fQAAntiProtonsAcceptedList->At(11)))->Fill(TMath::Abs(dca[0]));
2222     }//DCA xy global tracking
2223     if(fProtonAnalysisBase->IsUsedMaxDCAXYTPC()) {
2224       if(TMath::Abs(dca[0]) > fProtonAnalysisBase->GetMaxDCAXYTPC()) {
2225         ((TH1F *)(fQAAntiProtonsRejectedList->At(12)))->Fill(TMath::Abs(dca[0]));
2226       }
2227       else if(TMath::Abs(dca[0]) <= fProtonAnalysisBase->GetMaxDCAXYTPC())
2228         ((TH1F *)(fQAAntiProtonsAcceptedList->At(12)))->Fill(TMath::Abs(dca[0]));
2229     }//DCA xy TPC tracking
2230     if(fProtonAnalysisBase->IsUsedMaxDCAZ()) {
2231       if(TMath::Abs(dca[1]) > fProtonAnalysisBase->GetMaxDCAZ()) {
2232         ((TH1F *)(fQAAntiProtonsRejectedList->At(13)))->Fill(TMath::Abs(dca[1]));
2233       }
2234       else if(TMath::Abs(dca[1]) <= fProtonAnalysisBase->GetMaxDCAZ())
2235         ((TH1F *)(fQAAntiProtonsAcceptedList->At(13)))->Fill(TMath::Abs(dca[1]));
2236     }//DCA z global tracking
2237     if(fProtonAnalysisBase->IsUsedMaxDCAZTPC()) {
2238       if(TMath::Abs(dca[1]) > fProtonAnalysisBase->GetMaxDCAZTPC()) {
2239         ((TH1F *)(fQAAntiProtonsRejectedList->At(14)))->Fill(TMath::Abs(dca[1]));
2240       }
2241       else if(TMath::Abs(dca[1]) <= fProtonAnalysisBase->GetMaxDCAZTPC())
2242         ((TH1F *)(fQAAntiProtonsAcceptedList->At(14)))->Fill(TMath::Abs(dca[1]));
2243     }//DCA z TPC tracking
2244     if(fProtonAnalysisBase->IsUsedMaxConstrainChi2()) {
2245       if(track->GetConstrainedChi2() > 0) {
2246         if(TMath::Log(track->GetConstrainedChi2()) > fProtonAnalysisBase->GetMaxConstrainChi2()) {
2247           ((TH1F *)(fQAAntiProtonsRejectedList->At(15)))->Fill(TMath::Log(track->GetConstrainedChi2()));
2248         }
2249         else if(TMath::Log(track->GetConstrainedChi2()) <= fProtonAnalysisBase->GetMaxConstrainChi2())
2250           ((TH1F *)(fQAAntiProtonsAcceptedList->At(15)))->Fill(TMath::Log(track->GetConstrainedChi2()));
2251       }
2252     }//constrain chi2 - vertex
2253     if(fProtonAnalysisBase->IsUsedITSRefit()) {
2254       if ((track->GetStatus() & AliESDtrack::kITSrefit) == 0) {
2255         ((TH1F *)(fQAAntiProtonsRejectedList->At(16)))->Fill(0);
2256       }
2257       else if((track->GetStatus() & AliESDtrack::kITSrefit) != 0)
2258         ((TH1F *)(fQAAntiProtonsAcceptedList->At(16)))->Fill(0);
2259     }//ITS refit
2260     if(fProtonAnalysisBase->IsUsedTPCRefit()) {
2261       if ((track->GetStatus() & AliESDtrack::kTPCrefit) == 0) {
2262         ((TH1F *)(fQAAntiProtonsRejectedList->At(17)))->Fill(0);
2263       }
2264       else if((track->GetStatus() & AliESDtrack::kTPCrefit) != 0)
2265         ((TH1F *)(fQAAntiProtonsAcceptedList->At(17)))->Fill(0);
2266     }//TPC refit
2267     if(fProtonAnalysisBase->IsUsedESDpid()) {
2268       if ((track->GetStatus() & AliESDtrack::kESDpid) == 0) {
2269         ((TH1F *)(fQAAntiProtonsRejectedList->At(18)))->Fill(0);
2270       }
2271       else if((track->GetStatus() & AliESDtrack::kESDpid) != 0)
2272         ((TH1F *)(fQAAntiProtonsAcceptedList->At(18)))->Fill(0);
2273     }//ESD pid
2274     if(fProtonAnalysisBase->IsUsedTPCpid()) {
2275       if ((track->GetStatus() & AliESDtrack::kTPCpid) == 0) {
2276         ((TH1F *)(fQAAntiProtonsRejectedList->At(19)))->Fill(0);
2277       }
2278       else if((track->GetStatus() & AliESDtrack::kTPCpid) != 0)
2279         ((TH1F *)(fQAAntiProtonsAcceptedList->At(19)))->Fill(0);
2280     }//TPC pid
2281     if(fProtonAnalysisBase->IsUsedPointOnITSLayer1()) {
2282       if(!track->HasPointOnITSLayer(0)) {
2283         ((TH1F *)(fQAAntiProtonsRejectedList->At(20)))->Fill(0);
2284       }
2285       else if(track->HasPointOnITSLayer(0))
2286         ((TH1F *)(fQAAntiProtonsAcceptedList->At(20)))->Fill(0);
2287     }//point on SPD1
2288     if(fProtonAnalysisBase->IsUsedPointOnITSLayer2()) {
2289       if(!track->HasPointOnITSLayer(1)) {
2290         ((TH1F *)(fQAAntiProtonsRejectedList->At(21)))->Fill(0);
2291       }
2292       else if(track->HasPointOnITSLayer(1))
2293         ((TH1F *)(fQAAntiProtonsAcceptedList->At(21)))->Fill(0);
2294     }//point on SPD2
2295     if(fProtonAnalysisBase->IsUsedPointOnITSLayer3()) {
2296       if(!track->HasPointOnITSLayer(2)) {
2297         ((TH1F *)(fQAAntiProtonsRejectedList->At(22)))->Fill(0);
2298       }
2299       else if(track->HasPointOnITSLayer(2))
2300         ((TH1F *)(fQAAntiProtonsAcceptedList->At(22)))->Fill(0);
2301     }//point on SDD1
2302     if(fProtonAnalysisBase->IsUsedPointOnITSLayer4()) {
2303       if(!track->HasPointOnITSLayer(3)) {
2304         ((TH1F *)(fQAAntiProtonsRejectedList->At(23)))->Fill(0);
2305       }
2306       else if(track->HasPointOnITSLayer(3))
2307         ((TH1F *)(fQAAntiProtonsAcceptedList->At(23)))->Fill(0);
2308     }//point on SDD2
2309     if(fProtonAnalysisBase->IsUsedPointOnITSLayer5()) {
2310       if(!track->HasPointOnITSLayer(4)) {
2311         ((TH1F *)(fQAAntiProtonsRejectedList->At(24)))->Fill(0);
2312       }
2313       else if(track->HasPointOnITSLayer(4))
2314         ((TH1F *)(fQAAntiProtonsAcceptedList->At(24)))->Fill(0);
2315     }//point on SSD1
2316     if(fProtonAnalysisBase->IsUsedPointOnITSLayer6()) {
2317       if(!track->HasPointOnITSLayer(5)) {
2318         ((TH1F *)(fQAAntiProtonsRejectedList->At(25)))->Fill(0);
2319       }
2320       else if(track->HasPointOnITSLayer(5))
2321         ((TH1F *)(fQAAntiProtonsAcceptedList->At(25)))->Fill(0);
2322     }//point on SSD2
2323     if(fProtonAnalysisBase->IsUsedMinTPCdEdxPoints()) {
2324       if(track->GetTPCsignalN() < fProtonAnalysisBase->GetMinTPCdEdxPoints()) {
2325         ((TH1F *)(fQAAntiProtonsRejectedList->At(26)))->Fill(track->GetTPCsignalN());
2326       }
2327       if(track->GetTPCsignalN() >= fProtonAnalysisBase->GetMinTPCdEdxPoints())
2328         ((TH1F *)(fQAAntiProtonsAcceptedList->At(26)))->Fill(track->GetTPCsignalN());
2329     }//number of TPC points for the dE/dx
2330   }//antiprotons
2331 }