]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliPWG4HighPtQATPConly.cxx
Missing packages produce a Fatal
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliPWG4HighPtQATPConly.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: The ALICE Off-line Project.                                    *
5  * Contributors are mentioned in the code where appropriate.              *
6  *                                                                        *
7  * Permission to use, copy, modify and distribute this software and its   *
8  * documentation strictly for non-commercial purposes is hereby granted   *
9  * without fee, provided that the above copyright notice appears in all   *
10  * copies and that both the copyright notice and this permission notice   *
11  * appear in the supporting documentation. The authors make no claims     *
12  * about the suitability of this software for any purpose. It is          *
13  * provided "as is" without express or implied warranty.                  *
14  **************************************************************************/
15
16 //-----------------------------------------------------------------------
17 // This class compares the global reconstruction with the TPConly 
18 // reconstruction
19 // Momentum resolution is stored as function of track cuts and pt.
20 // Output: Histograms for different set of cuts
21 //-----------------------------------------------------------------------
22 // Author : Marta Verweij - UU
23 //-----------------------------------------------------------------------
24
25 #ifndef ALIPWG4HIGHPTQATPCONLY_CXX
26 #define ALIPWG4HIGHPTQATPCONLY_CXX
27
28 #include "AliPWG4HighPtQATPConly.h"
29
30 #include "TVector3.h"
31 #include <iostream>
32 #include "TH1.h"
33 #include "TH2.h"
34 #include "TH3.h"
35 #include "TList.h"
36 #include "TChain.h"
37 #include "TH3F.h"
38 #include <Bytes.h>
39 #include <TTree.h>
40
41 #include "AliAnalysisManager.h"
42 #include "AliESDInputHandler.h"
43 #include "AliESDtrack.h"
44 #include "AliESDfriend.h"
45 #include "AliESDfriendTrack.h"
46 #include "AliESDtrackCuts.h"
47 #include "AliExternalTrackParam.h"
48 #include "AliLog.h"
49 //#include "AliAnalysisHelperJetTasks.h"
50
51 #include "AliStack.h"
52 #include "TParticle.h"
53 #include "TH1I.h"
54 #include "AliMCEvent.h"
55 #include "AliMCEventHandler.h"
56
57 using namespace std; //required for resolving the 'cout' symbol
58
59 ClassImp(AliPWG4HighPtQATPConly)
60
61 AliPWG4HighPtQATPConly::AliPWG4HighPtQATPConly(): AliAnalysisTask("AliPWG4HighPtQATPConly", ""), 
62   fESD(0), 
63   fESDfriend(0), 
64   fMC(0),
65   fTrackCuts(0), 
66   fTrackCutsITS(0),
67   fMaxCosmicAngle(0.002),
68   fNEventAll(0),
69   fNEventSel(0),
70   fPtAll(0),
71   fPtSel(0),
72   fPtAllminPtTPCvsPtAll(0),
73   fPtAllminPtTPCvsPtAllNPointTPC(0),
74   fPtAllminPtTPCvsPtAllNPointTPCS(0),
75   fPtAllminPtTPCvsPtAllDCAR(0),
76   fPtAllminPtTPCvsPtAllDCAZ(0),
77   fPtAllminPtTPCvsPtAllPhi(0),
78   fPtAllminPtTPCvsPtAllNPointITS(0),
79   fPtAllminPtTPCvsPtAllNSigmaToVertex(0),
80   fPtAllminPtTPCvsPtAllChi2C(0),
81   fPtAllminPtTPCvsPtAllRel1PtUncertainty(0),
82   fPtAllminPtTPCvsPtAllChi2PerNClusTPC(0),
83   fPtAllminPtTPCvsPtAllChi2PerNClusITS(0),
84   fEtaPhiOutliers(0),
85   fPtSelITSouter(0),
86   fPtITSouterminPtTPCvsPtAll(0),
87   fPtITSouterminPtTPCvsPtAllNPointTPC(0),
88   fPtITSouterminPtTPCvsPtAllNPointTPCS(0),
89   fPtITSouterminPtTPCvsPtAllDCAR(0),
90   fPtITSouterminPtTPCvsPtAllDCAZ(0),
91   fPtITSouterminPtTPCvsPtAllPhi(0),
92   fPtITSouterminPtTPCvsPtAllNPointITS(0),
93   fPtITSouterminPtTPCvsPtAllNSigmaToVertex(0),
94   fPtITSouterminPtTPCvsPtAllChi2C(0),
95   fPtITSouterminPtTPCvsPtAllRel1PtUncertainty(0),
96   fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC(0),
97   fPtITSouterminPtTPCvsPtAllChi2PerNClusITS(0),
98   fPtITSouterminPtTPCvsPtAllITSLayer0(0),
99   fPtITSouterminPtTPCvsPtAllITSLayer1(0),
100   fPtITSouterminPtTPCvsPtAllITSLayer2(0),
101   fPtITSouterminPtTPCvsPtAllITSLayer3(0),
102   fPtITSouterminPtTPCvsPtAllITSLayer4(0),
103   fPtITSouterminPtTPCvsPtAllITSLayer5(0),
104   fPtITSouterminPtTPCvsPtAllNoSPD(0),
105   fPtITSouterminPtTPCvsPtAllNoSDD(0),
106   fPtITSouterminPtTPCvsPtAllNoSSD(0),
107   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0(0),
108   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1(0),
109   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2(0),
110   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3(0),
111   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4(0),
112   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5(0),
113   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD(0),
114   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD(0),
115   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD(0),
116   fHistList(0),
117   fPtAllTPC(0),
118   fPtSelTPC(0),
119   fPtSelTPCITS(0),
120   fHistListTPC(0),
121   fPtSelITS(0),
122   fPtITSminPtTPCvsPtITS(0),
123   fPtITSminPtTPCvsPtITSNPointTPC(0),
124   fPtITSminPtTPCvsPtITSNPointTPCS(0),
125   fPtITSminPtTPCvsPtITSDCAR(0),
126   fPtITSminPtTPCvsPtITSDCAZ(0),
127   fPtITSminPtTPCvsPtITSPhi(0),
128   fPtITSminPtTPCvsPtITSNPointITS(0),
129   fPtITSminPtTPCvsPtITSNSigmaToVertex(0),
130   fPtITSminPtTPCvsPtITSChi2C(0),
131   fPtITSminPtTPCvsPtITSRel1PtUncertainty(0),
132   fPtITSminPtTPCvsPtITSChi2PerNClusTPC(0),
133   fPtITSminPtTPCvsPtITSChi2PerNClusITS(0),
134   fPtRel1PtUncertaintyChi2PerClusTPC(0),
135   fPtNPointTPCSChi2PerClusTPC(0),
136   fPtNPointTPCSRel1PtUncertainty(0),
137   fPtRel1PtUncertaintyChi2PerClusTPCSharedSel(0),
138   fHistListITS(0),
139   fPtCosmicCandidates(0),
140   fDeltaPtCosmicCandidates(0),
141   fDeltaPhi(0),
142   fDeltaEta(0),
143   fHistListCosmics(0)
144 {
145
146 }
147 //________________________________________________________________________
148 AliPWG4HighPtQATPConly::AliPWG4HighPtQATPConly(const char *name): 
149   AliAnalysisTask(name, ""), 
150   fESD(0),
151   fESDfriend(0), 
152   fMC(0),
153   fTrackCuts(),
154   fTrackCutsITS(),
155   fMaxCosmicAngle(0.002),
156   fNEventAll(0),
157   fNEventSel(0),
158   fPtAll(0),
159   fPtSel(0),
160   fPtAllminPtTPCvsPtAll(0),
161   fPtAllminPtTPCvsPtAllNPointTPC(0),
162   fPtAllminPtTPCvsPtAllNPointTPCS(0),
163   fPtAllminPtTPCvsPtAllDCAR(0),
164   fPtAllminPtTPCvsPtAllDCAZ(0),
165   fPtAllminPtTPCvsPtAllPhi(0),
166   fPtAllminPtTPCvsPtAllNPointITS(0),
167   fPtAllminPtTPCvsPtAllNSigmaToVertex(0),
168   fPtAllminPtTPCvsPtAllChi2C(0),
169   fPtAllminPtTPCvsPtAllRel1PtUncertainty(0),
170   fPtAllminPtTPCvsPtAllChi2PerNClusTPC(0),
171   fPtAllminPtTPCvsPtAllChi2PerNClusITS(0),
172   fEtaPhiOutliers(0),
173   fPtSelITSouter(0),
174   fPtITSouterminPtTPCvsPtAll(0),
175   fPtITSouterminPtTPCvsPtAllNPointTPC(0),
176   fPtITSouterminPtTPCvsPtAllNPointTPCS(0),
177   fPtITSouterminPtTPCvsPtAllDCAR(0),
178   fPtITSouterminPtTPCvsPtAllDCAZ(0),
179   fPtITSouterminPtTPCvsPtAllPhi(0),
180   fPtITSouterminPtTPCvsPtAllNPointITS(0),
181   fPtITSouterminPtTPCvsPtAllNSigmaToVertex(0),
182   fPtITSouterminPtTPCvsPtAllChi2C(0),
183   fPtITSouterminPtTPCvsPtAllRel1PtUncertainty(0),
184   fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC(0),
185   fPtITSouterminPtTPCvsPtAllChi2PerNClusITS(0),
186   fPtITSouterminPtTPCvsPtAllITSLayer0(0),
187   fPtITSouterminPtTPCvsPtAllITSLayer1(0),
188   fPtITSouterminPtTPCvsPtAllITSLayer2(0),
189   fPtITSouterminPtTPCvsPtAllITSLayer3(0),
190   fPtITSouterminPtTPCvsPtAllITSLayer4(0),
191   fPtITSouterminPtTPCvsPtAllITSLayer5(0),
192   fPtITSouterminPtTPCvsPtAllNoSPD(0),
193   fPtITSouterminPtTPCvsPtAllNoSDD(0),
194   fPtITSouterminPtTPCvsPtAllNoSSD(0),
195   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0(0),
196   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1(0),
197   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2(0),
198   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3(0),
199   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4(0),
200   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5(0),
201   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD(0),
202   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD(0),
203   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD(0),
204   fHistList(0),
205   fPtAllTPC(0),
206   fPtSelTPC(0),
207   fPtSelTPCITS(0),
208   fHistListTPC(0),
209   fPtSelITS(0),
210   fPtITSminPtTPCvsPtITS(0),
211   fPtITSminPtTPCvsPtITSNPointTPC(0),
212   fPtITSminPtTPCvsPtITSNPointTPCS(0),
213   fPtITSminPtTPCvsPtITSDCAR(0),
214   fPtITSminPtTPCvsPtITSDCAZ(0),
215   fPtITSminPtTPCvsPtITSPhi(0),
216   fPtITSminPtTPCvsPtITSNPointITS(0),
217   fPtITSminPtTPCvsPtITSNSigmaToVertex(0),
218   fPtITSminPtTPCvsPtITSChi2C(0),
219   fPtITSminPtTPCvsPtITSRel1PtUncertainty(0),
220   fPtITSminPtTPCvsPtITSChi2PerNClusTPC(0),
221   fPtITSminPtTPCvsPtITSChi2PerNClusITS(0),
222   fPtRel1PtUncertaintyChi2PerClusTPC(0),
223   fPtNPointTPCSChi2PerClusTPC(0),
224   fPtNPointTPCSRel1PtUncertainty(0),
225   fPtRel1PtUncertaintyChi2PerClusTPCSharedSel(0),
226   fHistListITS(0),
227   fPtCosmicCandidates(0),
228   fDeltaPtCosmicCandidates(0),
229   fDeltaPhi(0),
230   fDeltaEta(0),
231   fHistListCosmics(0)
232 {
233   //
234   // Constructor. Initialization of Inputs and Outputs
235   //
236   Info("AliPWG4HighPtQATPConly","Calling Constructor");
237   // Input slot #0 works with a TChain ESD
238   DefineInput(0, TChain::Class());
239   // Output slot #0 writes into a TList
240   DefineOutput(0, TList::Class());
241   // Output slot #1 writes into a TList
242   DefineOutput(1, TList::Class());
243   // Output slot #2 writes into a TList
244   DefineOutput(2, TList::Class());
245   // Output slot #3 writes into a TList
246   DefineOutput(3, TList::Class());
247 }
248
249 //________________________________________________________________________
250 void AliPWG4HighPtQATPConly::ConnectInputData(Option_t *) 
251 {
252   // Connect ESD here
253   // Called once
254   AliDebug(2,Form(">> AliPWG4HighPtSpectra::ConnectInputData \n"));
255   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
256   if (!tree) {
257     AliDebug(2,Form( "ERROR: Could not read chain from input slot 0 \n"));
258     return;
259   } 
260   
261   AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
262   
263   if (!esdH) {
264     AliDebug(2,Form("ERROR: Could not get ESDInputHandler")); 
265     return;
266   } else
267     fESD = esdH->GetEvent();
268   
269  AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
270  // AliMCEventHandler* mcH = dynamic_cast<AliMCEventHandler*>
271  //                        (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
272   if (!eventHandler) {
273     AliDebug(2,Form( "ERROR: Could not retrieve MC event handler \n"));
274   }
275   else
276     fMC = eventHandler->MCEvent();
277
278   //attach the ESD friend
279    //  tree->SetBranchStatus("*", kTRUE);
280 //   tree->SetBranchStatus("Tracks*", kTRUE);
281 //   tree->SetBranchStatus("ESDfriend*", kTRUE);
282   //  fESD->ReadFromTree(tree);
283
284   fESDfriend = (AliESDfriend*)fESD->FindListObject("AliESDfriend");
285   if (!fESDfriend)
286     {
287       // works for both, we just want to avoid setting the branch adress twice
288       // in case of the new ESD
289       tree->SetBranchAddress("ESDfriend.",&fESDfriend);
290     }
291   
292 }
293
294 //________________________________________________________________________
295 void AliPWG4HighPtQATPConly::CreateOutputObjects() {
296   //Create output objects
297   AliDebug(2,Form(">> AliPWG4HighPtQATPConly::CreateOutputObjects \n")); 
298
299   Bool_t oldStatus = TH1::AddDirectoryStatus();
300   TH1::AddDirectory(kFALSE); 
301
302   OpenFile(0);
303   fHistList = new TList();
304   OpenFile(1);
305   fHistListTPC = new TList();
306   OpenFile(2);
307   fHistListITS = new TList();
308   OpenFile(3);
309   fHistListCosmics = new TList();
310
311
312   Int_t fgkNPhiBins=18;
313   Float_t kMinPhi = 0.;
314   Float_t kMaxPhi = 2.*TMath::Pi();
315   
316   Float_t fgkPtMin=0.;
317   Float_t fgkPtMax=100.;
318   Int_t fgkNPtBins=(int)(fgkPtMax-fgkPtMin);
319
320   Float_t fgkChi2PerClusMin = 0.;
321   Float_t fgkChi2PerClusMax = 3.5;
322   Int_t fgkChi2PerClusBins = (int)(fgkChi2PerClusMax*10.);
323   
324
325   Int_t fgkResPtBins=80;
326
327   fNEventAll = new TH1F("fNEventAll","NEventAll",1,-0.5,0.5);
328   fHistList->Add(fNEventAll);
329   fNEventSel = new TH1F("fNEventSel","NEvent Selected for analysis",1,-0.5,0.5);
330   fHistList->Add(fNEventSel);
331   fPtAll = new TH1F("fPtAll","PtAll",fgkNPtBins, fgkPtMin, fgkPtMax);
332   fHistList->Add(fPtAll);
333   fPtSel = new TH1F("fPtSel","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
334   fHistList->Add(fPtSel);
335  
336   fPtAllminPtTPCvsPtAll = new TH2F("fPtAllminPtTPCvsPtAll","PtAllminPtTPCvsPtAll",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
337   fPtAllminPtTPCvsPtAll->SetXTitle("p_{t}^{Global}");
338   fPtAllminPtTPCvsPtAll->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
339   fHistList->Add(fPtAllminPtTPCvsPtAll);
340   
341   fPtAllminPtTPCvsPtAllNPointTPC = new TH3F("fPtAllminPtTPCvsPtAllNPointTPC","PtAllminPtTPCvsPtAllNPointTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,160,0.5,160.5);
342   fPtAllminPtTPCvsPtAllNPointTPC->SetXTitle("p_{t}^{Global}");
343   fPtAllminPtTPCvsPtAllNPointTPC->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
344   fPtAllminPtTPCvsPtAllNPointTPC->SetZTitle("N_{point,TPC}");
345   fHistList->Add(fPtAllminPtTPCvsPtAllNPointTPC);
346
347   fPtAllminPtTPCvsPtAllNPointTPCS = new TH3F("fPtAllminPtTPCvsPtAllNPointTPCS","PtAllminPtTPCvsPtAllNPointTPCS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,100,0.,1.);
348   fPtAllminPtTPCvsPtAllNPointTPCS->SetXTitle("p_{t}^{Global}");
349   fPtAllminPtTPCvsPtAllNPointTPCS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
350   fPtAllminPtTPCvsPtAllNPointTPCS->SetZTitle("N_{point,TPC}^{Shared}/N_{point,TPC}");
351   fHistList->Add(fPtAllminPtTPCvsPtAllNPointTPCS);
352
353   fPtAllminPtTPCvsPtAllDCAR = new TH3F("fPtAllminPtTPCvsPtAllDCAR","PtAllminPtTPCvsPtAllDCAR",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,80,-0.2,0.2);
354   fPtAllminPtTPCvsPtAllDCAR->SetXTitle("p_{t}^{Global}");
355   fPtAllminPtTPCvsPtAllDCAR->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
356   fPtAllminPtTPCvsPtAllDCAR->SetZTitle("DCA_{R}");
357   fHistList->Add(fPtAllminPtTPCvsPtAllDCAR);
358
359   fPtAllminPtTPCvsPtAllDCAZ = new TH3F("fPtAllminPtTPCvsPtAllDCAZ","PtAllminPtTPCvsPtAllDCAZ",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,80,-2.,2.);
360   fPtAllminPtTPCvsPtAllDCAZ->SetXTitle("p_{t}^{Global}");
361   fPtAllminPtTPCvsPtAllDCAZ->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
362   fPtAllminPtTPCvsPtAllDCAZ->SetZTitle("DCA_{Z}");
363   fHistList->Add(fPtAllminPtTPCvsPtAllDCAZ);
364
365   fPtAllminPtTPCvsPtAllPhi = new TH3F("fPtAllminPtTPCvsPtAllPhi","PtAllminPtTPCvsPtAllPhi",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkNPhiBins,kMinPhi,kMaxPhi);
366   fPtAllminPtTPCvsPtAllPhi->SetXTitle("p_{t}^{Global}");
367   fPtAllminPtTPCvsPtAllPhi->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
368   fPtAllminPtTPCvsPtAllPhi->SetZTitle("#phi");
369   fHistList->Add(fPtAllminPtTPCvsPtAllPhi);
370
371   fPtAllminPtTPCvsPtAllNPointITS = new TH3F("fPtAllminPtTPCvsPtAllNPointITS","PtAllminPtTPCvsPtAllNPointITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,9,-0.5,8.5);
372   fPtAllminPtTPCvsPtAllNPointITS->SetXTitle("p_{t}^{Global}");
373   fPtAllminPtTPCvsPtAllNPointITS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
374   fPtAllminPtTPCvsPtAllNPointITS->SetZTitle("N_{point,ITS}");
375   fHistList->Add(fPtAllminPtTPCvsPtAllNPointITS);
376   
377   fPtAllminPtTPCvsPtAllNSigmaToVertex = new TH3F("fPtAllminPtTPCvsPtAllNSigmaToVertex","PtAllminPtTPCvsPtAllNSigmaToVertex",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,0.,8.);
378   fPtAllminPtTPCvsPtAllNSigmaToVertex->SetXTitle("p_{t}^{Global}");
379   fPtAllminPtTPCvsPtAllNSigmaToVertex->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
380   fPtAllminPtTPCvsPtAllNSigmaToVertex->SetZTitle("N#sigma to vertex");
381   fHistList->Add(fPtAllminPtTPCvsPtAllNSigmaToVertex);
382
383   fPtAllminPtTPCvsPtAllChi2C = new TH3F("fPtAllminPtTPCvsPtAllChi2C","PtAllminPtTPCvsPtAllChi2C",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,0.,10.);
384   fPtAllminPtTPCvsPtAllChi2C->SetXTitle("p_{t}^{Global}");
385   fPtAllminPtTPCvsPtAllChi2C->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
386   fPtAllminPtTPCvsPtAllChi2C->SetZTitle("Constrained #chi^{2}");
387   fHistList->Add(fPtAllminPtTPCvsPtAllChi2C);
388
389   fPtAllminPtTPCvsPtAllRel1PtUncertainty = new TH3F("fPtAllminPtTPCvsPtAllRel1PtUncertainty","PtAllminPtTPCvsPtAllRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,30,0.,0.3);
390   fPtAllminPtTPCvsPtAllRel1PtUncertainty->SetXTitle("p_{t}^{Global}");
391   fPtAllminPtTPCvsPtAllRel1PtUncertainty->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
392   fPtAllminPtTPCvsPtAllRel1PtUncertainty->SetZTitle("Rel1PtUncertainty");
393   fHistList->Add(fPtAllminPtTPCvsPtAllRel1PtUncertainty);
394
395   fPtAllminPtTPCvsPtAllChi2PerNClusTPC = new TH3F("fPtAllminPtTPCvsPtAllChi2PerNClusTPC","PtAllminPtTPCvsPtAllChi2PerNClusTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
396   fPtAllminPtTPCvsPtAllChi2PerNClusTPC->SetXTitle("p_{t}^{Global}");
397   fPtAllminPtTPCvsPtAllChi2PerNClusTPC->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
398   fPtAllminPtTPCvsPtAllChi2PerNClusTPC->SetZTitle("#chi^{2}/(2*NClusTPC-5)");
399   fHistList->Add(fPtAllminPtTPCvsPtAllChi2PerNClusTPC);
400
401   fPtAllminPtTPCvsPtAllChi2PerNClusITS = new TH3F("fPtAllminPtTPCvsPtAllChi2PerNClusITS","PtAllminPtTPCvsPtAllChi2PerNClusITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
402   fPtAllminPtTPCvsPtAllChi2PerNClusITS->SetXTitle("p_{t}^{Global}");
403   fPtAllminPtTPCvsPtAllChi2PerNClusITS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
404   fPtAllminPtTPCvsPtAllChi2PerNClusITS->SetZTitle("#chi^{2}/(2*NClusITS-5)");
405   fHistList->Add(fPtAllminPtTPCvsPtAllChi2PerNClusITS);
406
407   fEtaPhiOutliers = new TH2F("fEtaPhiOutliers","PtAllminPtTPCvsPtAll",20, -1.,1.,fgkNPhiBins,kMinPhi,kMaxPhi);
408   fEtaPhiOutliers->SetXTitle("#eta");
409   fEtaPhiOutliers->SetYTitle("#phi");
410   fHistList->Add(fEtaPhiOutliers);
411
412   //Global vs ITSouter-TPCinner
413   fPtSelITSouter = new TH1F("fPtSelITSouter","PtSelITSouter",fgkNPtBins,fgkPtMin,fgkPtMax);
414   fHistList->Add(fPtSelITSouter);
415   
416   fPtITSouterminPtTPCvsPtAll = new TH2F("fPtITSouterminPtTPCvsPtAll","PtAllminPtTPCvsPtAll",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
417   fPtITSouterminPtTPCvsPtAll->SetXTitle("p_{t}^{Global}");
418   fPtITSouterminPtTPCvsPtAll->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
419   fHistList->Add(fPtITSouterminPtTPCvsPtAll);
420   
421   fPtITSouterminPtTPCvsPtAllNPointTPC = new TH3F("fPtITSouterminPtTPCvsPtAllNPointTPC","PtAllminPtTPCvsPtAllNPointTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,160,0.5,160.5);
422   fPtITSouterminPtTPCvsPtAllNPointTPC->SetXTitle("p_{t}^{Global}");
423   fPtITSouterminPtTPCvsPtAllNPointTPC->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
424   fPtITSouterminPtTPCvsPtAllNPointTPC->SetZTitle("N_{point,TPC}");
425   fHistList->Add(fPtITSouterminPtTPCvsPtAllNPointTPC);
426
427   fPtITSouterminPtTPCvsPtAllNPointTPCS = new TH3F("fPtITSouterminPtTPCvsPtAllNPointTPCS","PtAllminPtTPCvsPtAllNPointTPCS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,100,0.,1.);
428   fPtITSouterminPtTPCvsPtAllNPointTPCS->SetXTitle("p_{t}^{Global}");
429   fPtITSouterminPtTPCvsPtAllNPointTPCS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
430   fPtITSouterminPtTPCvsPtAllNPointTPCS->SetZTitle("N_{point,TPC}^{Shared}/N_{point,TPC}");
431   fHistList->Add(fPtITSouterminPtTPCvsPtAllNPointTPCS);
432
433   fPtITSouterminPtTPCvsPtAllDCAR = new TH3F("fPtITSouterminPtTPCvsPtAllDCAR","PtAllminPtTPCvsPtAllDCAR",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,80,-0.2,0.2);
434   fPtITSouterminPtTPCvsPtAllDCAR->SetXTitle("p_{t}^{Global}");
435   fPtITSouterminPtTPCvsPtAllDCAR->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
436   fPtITSouterminPtTPCvsPtAllDCAR->SetZTitle("DCA_{R}");
437   fHistList->Add(fPtITSouterminPtTPCvsPtAllDCAR);
438
439   fPtITSouterminPtTPCvsPtAllDCAZ = new TH3F("fPtITSouterminPtTPCvsPtAllDCAZ","PtAllminPtTPCvsPtAllDCAZ",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,80,-2.,2.);
440   fPtITSouterminPtTPCvsPtAllDCAZ->SetXTitle("p_{t}^{Global}");
441   fPtITSouterminPtTPCvsPtAllDCAZ->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
442   fPtITSouterminPtTPCvsPtAllDCAZ->SetZTitle("DCA_{Z}");
443   fHistList->Add(fPtITSouterminPtTPCvsPtAllDCAZ);
444
445   fPtITSouterminPtTPCvsPtAllPhi = new TH3F("fPtITSouterminPtTPCvsPtAllPhi","PtAllminPtTPCvsPtAllPhi",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkNPhiBins,kMinPhi,kMaxPhi);
446   fPtITSouterminPtTPCvsPtAllPhi->SetXTitle("p_{t}^{Global}");
447   fPtITSouterminPtTPCvsPtAllPhi->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
448   fPtITSouterminPtTPCvsPtAllPhi->SetZTitle("#phi");
449   fHistList->Add(fPtITSouterminPtTPCvsPtAllPhi);
450
451   fPtITSouterminPtTPCvsPtAllNPointITS = new TH3F("fPtITSouterminPtTPCvsPtAllNPointITS","PtAllminPtTPCvsPtAllNPointITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,9,-0.5,8.5);
452   fPtITSouterminPtTPCvsPtAllNPointITS->SetXTitle("p_{t}^{Global}");
453   fPtITSouterminPtTPCvsPtAllNPointITS->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
454   fPtITSouterminPtTPCvsPtAllNPointITS->SetZTitle("N_{point,ITS}");
455   fHistList->Add(fPtITSouterminPtTPCvsPtAllNPointITS);
456   
457   fPtITSouterminPtTPCvsPtAllNSigmaToVertex = new TH3F("fPtITSouterminPtTPCvsPtAllNSigmaToVertex","PtAllminPtTPCvsPtAllNSigmaToVertex",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,0.,8.);
458   fPtITSouterminPtTPCvsPtAllNSigmaToVertex->SetXTitle("p_{t}^{Global}");
459   fPtITSouterminPtTPCvsPtAllNSigmaToVertex->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
460   fPtITSouterminPtTPCvsPtAllNSigmaToVertex->SetZTitle("N#sigma to vertex");
461   fHistList->Add(fPtITSouterminPtTPCvsPtAllNSigmaToVertex);
462
463   fPtITSouterminPtTPCvsPtAllChi2C = new TH3F("fPtITSouterminPtTPCvsPtAllChi2C","PtAllminPtTPCvsPtAllChi2C",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,0.,10.);
464   fPtITSouterminPtTPCvsPtAllChi2C->SetXTitle("p_{t}^{Global}");
465   fPtITSouterminPtTPCvsPtAllChi2C->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
466   fPtITSouterminPtTPCvsPtAllChi2C->SetZTitle("Constrained #chi^{2}");
467   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2C);
468
469   fPtITSouterminPtTPCvsPtAllRel1PtUncertainty = new TH3F("fPtITSouterminPtTPCvsPtAllRel1PtUncertainty","PtAllminPtTPCvsPtAllRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,30,0.,0.3);
470   fPtITSouterminPtTPCvsPtAllRel1PtUncertainty->SetXTitle("p_{t}^{Global}");
471   fPtITSouterminPtTPCvsPtAllRel1PtUncertainty->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
472   fPtITSouterminPtTPCvsPtAllRel1PtUncertainty->SetZTitle("Rel1PtUncertainty");
473   fHistList->Add(fPtITSouterminPtTPCvsPtAllRel1PtUncertainty);
474
475   fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC","PtAllminPtTPCvsPtAllChi2PerNClusTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,35,0.,3.5);
476   fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC->SetXTitle("p_{t}^{Global}");
477   fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
478   fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC->SetZTitle("#chi^{2}/(2*NClusTPC-5)");
479   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC);
480
481   fPtITSouterminPtTPCvsPtAllChi2PerNClusITS = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITS","PtAllminPtTPCvsPtAllChi2PerNClusITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,35,0.,3.5);
482   fPtITSouterminPtTPCvsPtAllChi2PerNClusITS->SetXTitle("p_{t}^{Global}");
483   fPtITSouterminPtTPCvsPtAllChi2PerNClusITS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
484   fPtITSouterminPtTPCvsPtAllChi2PerNClusITS->SetZTitle("#chi^{2}/(2*NClusITS-5)");
485   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITS);
486
487   //As function of ITS layers
488   fPtITSouterminPtTPCvsPtAllITSLayer0 = new TH2F("fPtITSouterminPtTPCvsPtAllITSLayer0","PtAllminPtTPCvsPtAllITSLayer0",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
489   fPtITSouterminPtTPCvsPtAllITSLayer0->SetXTitle("p_{t}^{Global}");
490   fPtITSouterminPtTPCvsPtAllITSLayer0->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
491   fHistList->Add(fPtITSouterminPtTPCvsPtAllITSLayer0);
492
493   fPtITSouterminPtTPCvsPtAllITSLayer1 = new TH2F("fPtITSouterminPtTPCvsPtAllITSLayer1","PtAllminPtTPCvsPtAllITSLayer1",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
494   fPtITSouterminPtTPCvsPtAllITSLayer1->SetXTitle("p_{t}^{Global}");
495   fPtITSouterminPtTPCvsPtAllITSLayer1->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
496   fHistList->Add(fPtITSouterminPtTPCvsPtAllITSLayer1);
497
498   fPtITSouterminPtTPCvsPtAllITSLayer2 = new TH2F("fPtITSouterminPtTPCvsPtAllITSLayer2","PtAllminPtTPCvsPtAllITSLayer2",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
499   fPtITSouterminPtTPCvsPtAllITSLayer2->SetXTitle("p_{t}^{Global}");
500   fPtITSouterminPtTPCvsPtAllITSLayer2->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
501   fHistList->Add(fPtITSouterminPtTPCvsPtAllITSLayer2);
502
503   fPtITSouterminPtTPCvsPtAllITSLayer3 = new TH2F("fPtITSouterminPtTPCvsPtAllITSLayer3","PtAllminPtTPCvsPtAllITSLayer3",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
504   fPtITSouterminPtTPCvsPtAllITSLayer3->SetXTitle("p_{t}^{Global}");
505   fPtITSouterminPtTPCvsPtAllITSLayer3->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
506   fHistList->Add(fPtITSouterminPtTPCvsPtAllITSLayer3);
507
508   fPtITSouterminPtTPCvsPtAllITSLayer4 = new TH2F("fPtITSouterminPtTPCvsPtAllITSLayer4","PtAllminPtTPCvsPtAllITSLayer4",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
509   fPtITSouterminPtTPCvsPtAllITSLayer4->SetXTitle("p_{t}^{Global}");
510   fPtITSouterminPtTPCvsPtAllITSLayer4->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
511   fHistList->Add(fPtITSouterminPtTPCvsPtAllITSLayer4);
512
513   fPtITSouterminPtTPCvsPtAllITSLayer5 = new TH2F("fPtITSouterminPtTPCvsPtAllITSLayer5","PtAllminPtTPCvsPtAllITSLayer5",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
514   fPtITSouterminPtTPCvsPtAllITSLayer5->SetXTitle("p_{t}^{Global}");
515   fPtITSouterminPtTPCvsPtAllITSLayer5->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
516   fHistList->Add(fPtITSouterminPtTPCvsPtAllITSLayer5);
517
518   fPtITSouterminPtTPCvsPtAllNoSPD = new TH2F("fPtITSouterminPtTPCvsPtAllNoSPD","PtAllminPtTPCvsPtAllNoSPD",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
519   fPtITSouterminPtTPCvsPtAllNoSPD->SetXTitle("p_{t}^{Global}");
520   fPtITSouterminPtTPCvsPtAllNoSPD->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
521   fHistList->Add(fPtITSouterminPtTPCvsPtAllNoSPD);
522
523   fPtITSouterminPtTPCvsPtAllNoSDD = new TH2F("fPtITSouterminPtTPCvsPtAllNoSDD","PtAllminPtTPCvsPtAllNoSDD",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
524   fPtITSouterminPtTPCvsPtAllNoSDD->SetXTitle("p_{t}^{Global}");
525   fPtITSouterminPtTPCvsPtAllNoSDD->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
526   fHistList->Add(fPtITSouterminPtTPCvsPtAllNoSDD);
527
528   fPtITSouterminPtTPCvsPtAllNoSSD = new TH2F("fPtITSouterminPtTPCvsPtAllNoSSD","PtAllminPtTPCvsPtAllNoSSD",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
529   fPtITSouterminPtTPCvsPtAllNoSSD->SetXTitle("p_{t}^{Global}");
530   fPtITSouterminPtTPCvsPtAllNoSSD->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
531   fHistList->Add(fPtITSouterminPtTPCvsPtAllNoSSD);
532
533   //
534   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0 = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0","PtAllminPtTPCvsPtAllChi2PerNClusITSLayer0",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
535   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0->SetXTitle("p_{t}^{Global}");
536   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
537   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0->SetZTitle("#chi^{2}/NPointITS");
538   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0);
539
540   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1 = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1","PtAllminPtTPCvsPtAllChi2PerNClusITSLayer1",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
541   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1->SetXTitle("p_{t}^{Global}");
542   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
543   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1->SetZTitle("#chi^{2}/NPointITS");
544   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1);
545
546   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2 = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2","PtAllminPtTPCvsPtAllChi2PerNClusITSLayer2",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
547   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2->SetXTitle("p_{t}^{Global}");
548   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
549   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2->SetZTitle("#chi^{2}/NPointITS");
550   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2);
551
552   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3 = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3","PtAllminPtTPCvsPtAllChi2PerNClusITSLayer3",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
553   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3->SetXTitle("p_{t}^{Global}");
554   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
555   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3->SetZTitle("#chi^{2}/NPointITS");
556   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3);
557
558   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4 = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4","PtAllminPtTPCvsPtAllChi2PerNClusITSLayer4",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
559   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4->SetXTitle("p_{t}^{Global}");
560   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
561   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4->SetZTitle("#chi^{2}/NPointITS");
562   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4);
563
564   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5 = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5","PtAllminPtTPCvsPtAllChi2PerNClusITSLayer5",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
565   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5->SetXTitle("p_{t}^{Global}");
566   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
567   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5->SetZTitle("#chi^{2}/NPointITS");
568   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5);
569
570   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD","PtAllminPtTPCvsPtAllChi2PerNClusITSNoSPD",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
571   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD->SetXTitle("p_{t}^{Global}");
572   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
573   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD->SetZTitle("#chi^{2}/(2*NPointITS-5)");
574   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD);
575
576   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD","PtAllminPtTPCvsPtAllChi2PerNClusITSNoSDD",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
577   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD->SetXTitle("p_{t}^{Global}");
578   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
579   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD->SetZTitle("#chi^{2}/(2*NPointITS-5)");
580   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD);
581
582   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD = new TH3F("fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD","PtAllminPtTPCvsPtAllChi2PerNClusITSNoSSD",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkChi2PerClusBins,fgkChi2PerClusMin,fgkChi2PerClusMax);
583   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD->SetXTitle("p_{t}^{Global}");
584   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD->SetYTitle("(1/p_{t}^{ITSouter}-1/p_{t}^{TPCinner})/(1/p_{t}^{ITSouter})");
585   fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD->SetZTitle("#chi^{2}/(2*NPointITS-5)");
586   fHistList->Add(fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD);
587
588
589   //ITSrefit
590   fPtSelITS = new TH1F("fPtSelITSrefit","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
591   fHistListITS->Add(fPtSelITS);
592   
593   fPtITSminPtTPCvsPtITS = new TH2F("fPtITSminPtTPCvsPtITS","PtITSminPtTPCvsPtITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
594   fPtITSminPtTPCvsPtITS->SetXTitle("p_{t}^{Global}");
595   fPtITSminPtTPCvsPtITS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
596   fHistListITS->Add(fPtITSminPtTPCvsPtITS);
597   
598   fPtITSminPtTPCvsPtITSNPointTPC = new TH3F("fPtITSminPtTPCvsPtITSNPointTPC","PtITSminPtTPCvsPtITSNPointTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,160,0.5,160.5);
599   fPtITSminPtTPCvsPtITSNPointTPC->SetXTitle("p_{t}^{Global}");
600   fPtITSminPtTPCvsPtITSNPointTPC->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
601   fPtITSminPtTPCvsPtITSNPointTPC->SetZTitle("N_{point,TPC}");
602   fHistListITS->Add(fPtITSminPtTPCvsPtITSNPointTPC);
603
604   fPtITSminPtTPCvsPtITSNPointTPCS = new TH3F("fPtITSminPtTPCvsPtITSNPointTPCS","PtITSminPtTPCvsPtITSNPointTPCS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,100,0.,1.);
605   fPtITSminPtTPCvsPtITSNPointTPCS->SetXTitle("p_{t}^{Global}");
606   fPtITSminPtTPCvsPtITSNPointTPCS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
607   fPtITSminPtTPCvsPtITSNPointTPCS->SetZTitle("N_{point,TPC}^{Shared}/N_{point,TPC}");
608   fHistListITS->Add(fPtITSminPtTPCvsPtITSNPointTPCS);    
609
610   fPtITSminPtTPCvsPtITSDCAR = new TH3F("fPtITSminPtTPCvsPtITSDCAR","PtITSminPtTPCvsPtITSDCAR",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,80,-0.2,0.2);
611   fPtITSminPtTPCvsPtITSDCAR->SetXTitle("p_{t}^{Global}");
612   fPtITSminPtTPCvsPtITSDCAR->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
613   fPtITSminPtTPCvsPtITSDCAR->SetZTitle("DCA_{R}");
614   fHistListITS->Add(fPtITSminPtTPCvsPtITSDCAR);
615   
616   fPtITSminPtTPCvsPtITSDCAZ = new TH3F("fPtITSminPtTPCvsPtITSDCAZ","PtITSminPtTPCvsPtITSDCAZ",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,80,-2.,2.);
617   fPtITSminPtTPCvsPtITSDCAZ->SetXTitle("p_{t}^{Global}");
618   fPtITSminPtTPCvsPtITSDCAZ->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
619   fPtITSminPtTPCvsPtITSDCAZ->SetZTitle("DCA_{Z}");
620   fHistListITS->Add(fPtITSminPtTPCvsPtITSDCAZ);
621   
622   fPtITSminPtTPCvsPtITSPhi = new TH3F("fPtITSminPtTPCvsPtITSPhi","PtITSminPtTPCvsPtITSPhi",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkNPhiBins,kMinPhi,kMaxPhi);
623   fPtITSminPtTPCvsPtITSPhi->SetXTitle("p_{t}^{Global}");
624   fPtITSminPtTPCvsPtITSPhi->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
625   fPtITSminPtTPCvsPtITSPhi->SetZTitle("#phi");
626   fHistListITS->Add(fPtITSminPtTPCvsPtITSPhi);
627   
628   fPtITSminPtTPCvsPtITSNPointITS = new TH3F("fPtITSminPtTPCvsPtITSNPointITS","PtITSminPtTPCvsPtITSNPointITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,9,-0.5,8.5);
629   fPtITSminPtTPCvsPtITSNPointITS->SetXTitle("p_{t}^{Global}");
630   fPtITSminPtTPCvsPtITSNPointITS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
631   fPtITSminPtTPCvsPtITSNPointITS->SetZTitle("N_{point,ITS}");
632   fHistListITS->Add(fPtITSminPtTPCvsPtITSNPointITS); 
633   
634   fPtITSminPtTPCvsPtITSNSigmaToVertex = new TH3F("fPtITSminPtTPCvsPtITSNSigmaToVertex","PtITSminPtTPCvsPtITSNSigmaToVertex",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,0.,8.);
635   fPtITSminPtTPCvsPtITSNSigmaToVertex->SetXTitle("p_{t}^{Global}");
636   fPtITSminPtTPCvsPtITSNSigmaToVertex->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
637   fPtITSminPtTPCvsPtITSNSigmaToVertex->SetZTitle("N#sigma to vertex");
638   fHistListITS->Add(fPtITSminPtTPCvsPtITSNSigmaToVertex);
639
640   fPtITSminPtTPCvsPtITSChi2C = new TH3F("fPtITSminPtTPCvsPtITSChi2C","PtITSminPtTPCvsPtITSChi2C",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,0.,10.);
641   fPtITSminPtTPCvsPtITSChi2C->SetXTitle("p_{t}^{Global}");
642   fPtITSminPtTPCvsPtITSChi2C->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
643   fPtITSminPtTPCvsPtITSChi2C->SetZTitle("Constrained #chi^{2}");
644   fHistListITS->Add(fPtITSminPtTPCvsPtITSChi2C);
645
646   fPtITSminPtTPCvsPtITSRel1PtUncertainty = new TH3F("fPtITSminPtTPCvsPtITSRel1PtUncertainty","PtITSminPtTPCvsPtITSRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,30,0.,0.3);
647   fPtITSminPtTPCvsPtITSRel1PtUncertainty->SetXTitle("p_{t}^{Global}");
648   fPtITSminPtTPCvsPtITSRel1PtUncertainty->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
649   fPtITSminPtTPCvsPtITSRel1PtUncertainty->SetZTitle("Rel1PtUncertainty");
650   fHistListITS->Add(fPtITSminPtTPCvsPtITSRel1PtUncertainty);
651
652   fPtITSminPtTPCvsPtITSChi2PerNClusTPC = new TH3F("fPtITSminPtTPCvsPtITSChi2PerNClusTPC","PtITSminPtTPCvsPtITSChi2PerNClusTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,35,0.,3.5);
653   fPtITSminPtTPCvsPtITSChi2PerNClusTPC->SetXTitle("p_{t}^{Global}");
654   fPtITSminPtTPCvsPtITSChi2PerNClusTPC->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
655   fPtITSminPtTPCvsPtITSChi2PerNClusTPC->SetZTitle("#chi^{2}/(2*NClusTPC-5)");
656   fHistListITS->Add(fPtITSminPtTPCvsPtITSChi2PerNClusTPC);
657
658   fPtITSminPtTPCvsPtITSChi2PerNClusITS = new TH3F("fPtITSminPtTPCvsPtITSChi2PerNClusITS","PtITSminPtTPCvsPtITSChi2PerNClusITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,35,0.,3.5);
659   fPtITSminPtTPCvsPtITSChi2PerNClusITS->SetXTitle("p_{t}^{Global}");
660   fPtITSminPtTPCvsPtITSChi2PerNClusITS->SetYTitle("(1/p_{t}^{Global}-1/p_{t}^{TPC})/(1/p_{t}^{Global})");
661   fPtITSminPtTPCvsPtITSChi2PerNClusITS->SetZTitle("#chi^{2}/(2*NClusITS-5)");
662   fHistListITS->Add(fPtITSminPtTPCvsPtITSChi2PerNClusITS);
663
664   fPtRel1PtUncertaintyChi2PerClusTPC = new TH3F("fPtRel1PtUncertaintyChi2PerClusTPC","PtITSminPtTPCvsPtITSRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,30,0.,0.3,35,0.,3.5);
665   fPtRel1PtUncertaintyChi2PerClusTPC->SetXTitle("p_{t}^{global}");
666   fPtRel1PtUncertaintyChi2PerClusTPC->SetYTitle("Rel1PtUncertainty");
667   fPtRel1PtUncertaintyChi2PerClusTPC->SetZTitle("#chi^{2}/(2*N_{clusters}^{TPC}-5)");
668   fHistListITS->Add(fPtRel1PtUncertaintyChi2PerClusTPC);
669
670   fPtNPointTPCSChi2PerClusTPC = new TH3F("fPtNPointTPCSChi2PerClusTPC","PtITSminPtTPCvsPtITSNPointTPCS",fgkNPtBins, fgkPtMin,fgkPtMax,100,0.,1.,35,0.,3.5);
671   fPtNPointTPCSChi2PerClusTPC->SetXTitle("p_{t}^{global}");
672   fPtNPointTPCSChi2PerClusTPC->SetYTitle("N_{Point,TPC}^{Shared}/N_{Point,TPC}");
673   fPtNPointTPCSChi2PerClusTPC->SetZTitle("#chi^{2}/(2*N_{clusters}^{TPC}-5)");
674   fHistListITS->Add(fPtNPointTPCSChi2PerClusTPC);
675
676   fPtNPointTPCSRel1PtUncertainty = new TH3F("fPtNPointTPCSRel1PtUncertainty","PtITSminPtTPCvsPtITSNPointTPCS",fgkNPtBins, fgkPtMin,fgkPtMax,100,0.,1.,30,0.,0.3);
677   fPtNPointTPCSRel1PtUncertainty->SetXTitle("p_{t}^{global}");
678   fPtNPointTPCSRel1PtUncertainty->SetYTitle("N_{Point,TPC}^{Shared}/N_{Point,TPC}");
679   fPtNPointTPCSRel1PtUncertainty->SetZTitle("Rel1PtUncertainty");
680   fHistListITS->Add(fPtNPointTPCSRel1PtUncertainty);
681
682   fPtRel1PtUncertaintyChi2PerClusTPCSharedSel = new TH3F("fPtRel1PtUncertaintyChi2PerClusTPCSharedSel","PtITSminPtTPCvsPtITSRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,30,0.,0.3,35,0.,3.5);
683   fPtRel1PtUncertaintyChi2PerClusTPCSharedSel->SetXTitle("p_{t}^{global}");
684   fPtRel1PtUncertaintyChi2PerClusTPCSharedSel->SetYTitle("Rel1PtUncertainty");
685   fPtRel1PtUncertaintyChi2PerClusTPCSharedSel->SetZTitle("#chi^{2}/(2*N_{clusters}^{TPC}-5)");
686   fHistListITS->Add(fPtRel1PtUncertaintyChi2PerClusTPCSharedSel);
687
688   fPtAllTPC = new TH1F("fPtAllTPC","PtAll",fgkNPtBins, fgkPtMin, fgkPtMax);
689   fHistListTPC->Add(fPtAllTPC);
690   fPtSelTPC = new TH1F("fPtSelTPC","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
691   fHistListTPC->Add(fPtSelTPC);
692   fPtSelTPCITS = new TH1F("fPtSelTPCITS","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
693   fHistListTPC->Add(fPtSelTPCITS);
694
695   //****************************************************************************************************************//
696   //                                              Cosmic Candidates                                                 //
697   //****************************************************************************************************************//
698   fPtCosmicCandidates = new TH1F("fPtCosmicCandidates","fPtCosmicCandidates",fgkNPtBins, fgkPtMin, fgkPtMax);
699   fHistListCosmics->Add(fPtCosmicCandidates);  
700   fDeltaPtCosmicCandidates = new TH1F("fDeltaPtCosmicCandidates","fDeltaPtCosmicCandidates",fgkNPtBins, -50., 50.);
701   fHistListCosmics->Add(fDeltaPtCosmicCandidates);  
702   fDeltaPhi = new TH1F("fDeltaPhi","fDeltaPhi",fgkNPhiBins*2,-1*kMaxPhi,kMaxPhi);
703   fHistListCosmics->Add(fDeltaPhi);  
704   fDeltaEta = new TH1F("fDeltaEta","fDeltaEta",20, -2.,2.);
705   fHistListCosmics->Add(fDeltaEta);  
706
707   TH1::AddDirectory(oldStatus);   
708
709 }
710 //________________________________________________________________________
711 void AliPWG4HighPtQATPConly::Exec(Option_t *) {  
712   // Main loop
713   // Called for each event
714   AliDebug(2,Form(">> AliPWG4HighPtQATPConly::Exec \n"));  
715
716   // All events without selection
717   fNEventAll->Fill(0.);
718
719   if (!fESD) {
720     AliDebug(2,Form("ERROR: fESD not available"));
721     // Post output data
722      PostData(0, fHistList);
723      PostData(1, fHistListTPC);
724      PostData(2, fHistListITS);
725      PostData(3, fHistListCosmics);
726     return;
727   }
728
729   fESD->SetESDfriend(fESDfriend); //Attach the friend to the ESD
730   if (!fESDfriend) {
731     AliDebug(2,Form("ERROR: fESDfriend not available"));
732     // Post output data
733      PostData(0, fHistList);
734      PostData(1, fHistListTPC);
735      PostData(2, fHistListITS);
736      PostData(3, fHistListCosmics);
737     return;
738   }
739
740   Bool_t isSelected = ((AliInputEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler()))->IsEventSelected();
741   if(!isSelected) { //Select collison candidates
742     AliDebug(2,Form(" Trigger Selection: event REJECTED ... "));
743     // Post output data
744      PostData(0, fHistList);
745      PostData(1, fHistListTPC);
746      PostData(2, fHistListITS);
747      PostData(3, fHistListCosmics);
748     return;
749   }
750
751   //  AliMCEventHandler *eventHandler = dynamic_cast<AliMCEventHandler*> (AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
752 //   //  AliMCEventHandler* eventHandler = (AliMCEventHandler*)(AliAnalysisManager::GetAnalysisManager()->GetMCtruthEventHandler());
753   
754   AliStack* stack = 0x0;
755   AliMCEvent* mcEvent = 0x0;
756   
757   if(fMC) {
758     mcEvent = fMC;
759     if (!mcEvent) {
760       AliDebug(2,Form("ERROR: Could not retrieve MC event"));
761       PostData(0, fHistList);
762       PostData(1, fHistListTPC);
763       PostData(2, fHistListITS);
764       PostData(3, fHistListCosmics);
765       return;
766     }
767     
768     AliDebug(2,Form("MC particles: %d", mcEvent->GetNumberOfTracks()));
769     
770     stack = mcEvent->Stack();                //Particles Stack
771     
772     AliDebug(2,Form("MC particles stack: %d", stack->GetNtrack()));
773   }
774   
775
776   const AliESDVertex *vtx = fESD->GetPrimaryVertexTracks();
777   // Need vertex cut
778   if (vtx->GetNContributors() < 2) {
779     // Post output data
780     PostData(0, fHistList);
781     PostData(1, fHistListTPC);
782     PostData(2, fHistListITS);
783     PostData(3, fHistListCosmics);
784     return;
785   }
786
787   AliDebug(2,Form("Vertex title %s, status %d, nCont %d\n",vtx->GetTitle(), vtx->GetStatus(), vtx->GetNContributors()));
788   double primVtx[3];
789   vtx->GetXYZ(primVtx);
790   //  printf("primVtx: %g  %g  %g \n",primVtx[0],primVtx[1],primVtx[2]);
791   if(TMath::Sqrt(primVtx[0]*primVtx[0] + primVtx[1]*primVtx[1])>1. || TMath::Abs(primVtx[2]>10.)){
792     // Post output data
793     PostData(0, fHistList);
794     PostData(1, fHistListTPC);
795     PostData(2, fHistListITS);
796     PostData(3, fHistListCosmics);
797     return;
798   }
799   if(!fESD->GetNumberOfTracks() || fESD->GetNumberOfTracks()<2){ 
800     // Post output data
801     PostData(0, fHistList);
802     PostData(1, fHistListTPC);
803     PostData(2, fHistListITS);
804     PostData(3, fHistListCosmics);
805     return;
806   }
807   Int_t nTracks = fESD->GetNumberOfTracks();
808   AliDebug(2,Form("nTracks %d\n", nTracks));
809
810   if(!fTrackCuts) return;
811
812   // Selected events for analysis
813   fNEventSel->Fill(0.);
814
815   for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
816     
817     AliESDtrack *track = fESD->GetTrack(iTrack);
818     AliExternalTrackParam *trackTPC = (AliExternalTrackParam *)track->GetTPCInnerParam();
819     if(!track || !trackTPC) continue;
820     const AliESDfriendTrack* constfriendtrack = track->GetFriendTrack();
821     if (!constfriendtrack) { continue;}
822     AliESDfriendTrack friendtrack(*constfriendtrack);
823  
824     Float_t pt = track->Pt();
825     Float_t ptTPC = trackTPC->Pt();
826     Float_t phi = track->Phi();
827     Float_t dca2D, dcaZ;
828     track->GetImpactParameters(dca2D,dcaZ);
829     // Float_t dca2DTPC, dcaZTPC;
830     //track->GetImpactParametersTPC(dca2DTPC,dcaZTPC);
831     UChar_t itsMap = track->GetITSClusterMap();
832     Int_t nPointITS = 0;
833     for (Int_t i=0; i < 6; i++) {
834       if (itsMap & (1 << i))
835         nPointITS ++;
836     }
837     double mom[3];
838     track->GetPxPyPz(mom);
839     double momTPC[3];
840     trackTPC->GetPxPyPz(momTPC);
841     Float_t nSigmaToVertex = fTrackCuts->GetSigmaToVertex(track);// Calculates the number of sigma to the vertex for a track.
842     Float_t chi2C = track->GetConstrainedChi2();
843     Float_t relUncertainty1Pt = TMath::Sqrt(TMath::Abs(track->GetSigma1Pt2()))*pt;
844     Float_t chi2PerClusterTPC = -1.;
845     Float_t nClustersTPC = track->GetTPCNcls();
846     if(nClustersTPC>0.) chi2PerClusterTPC = track->GetTPCchi2()/(2.*nClustersTPC-5.);
847     Float_t chi2PerNPointITS = -1.;
848     if(nPointITS>3) chi2PerNPointITS = track->GetITSchi2()/(2.*(float)nPointITS-5.);
849
850     fPtAll->Fill(pt);
851     fPtAllTPC->Fill(ptTPC);
852
853
854     if (fTrackCuts->AcceptTrack(track)) {
855
856       Bool_t cosmic = kFALSE;    
857       if(pt>6.) { cosmic = IsCosmic(track,iTrack,6.); }
858       //    if(cosmic) continue;
859
860       fPtSel->Fill(pt);
861       fPtSelTPC->Fill(ptTPC);
862       if(ptTPC==0. || pt==0.) continue;
863       fPtAllminPtTPCvsPtAll->Fill(pt,(1./pt-1./ptTPC)/(1./pt) );
864       fPtAllminPtTPCvsPtAllNPointTPC->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nClustersTPC);
865       if(nClustersTPC>0.) fPtAllminPtTPCvsPtAllNPointTPCS->Fill(pt,(1./pt-1./ptTPC)/(1./pt),track->GetTPCnclsS()/nClustersTPC);
866       fPtAllminPtTPCvsPtAllDCAR->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dca2D);
867       fPtAllminPtTPCvsPtAllDCAZ->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dcaZ);
868       fPtAllminPtTPCvsPtAllPhi->Fill(pt,(1./pt-1./ptTPC)/(1./pt),phi);
869       fPtAllminPtTPCvsPtAllNPointITS->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nPointITS);
870       fPtAllminPtTPCvsPtAllNSigmaToVertex->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nSigmaToVertex);
871       fPtAllminPtTPCvsPtAllChi2C->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2C);
872       fPtAllminPtTPCvsPtAllRel1PtUncertainty->Fill(pt,(1./pt-1./ptTPC)/(1./pt),relUncertainty1Pt);
873       fPtAllminPtTPCvsPtAllChi2PerNClusTPC->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2PerClusterTPC);
874       if(nPointITS>3) fPtAllminPtTPCvsPtAllChi2PerNClusITS->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2PerNPointITS);
875       if(TMath::Abs((1./pt-1./ptTPC)/(1./pt))>0.8) fEtaPhiOutliers->Fill(track->Eta(),phi);
876       if (friendtrack.GetITSOut()) {
877         AliExternalTrackParam trackITSouter(*(friendtrack.GetITSOut())); 
878         Float_t ptITSouter = trackITSouter.Pt();
879         if(ptITSouter==0.) continue;
880         fPtSelITSouter->Fill(ptITSouter);
881         fPtITSouterminPtTPCvsPtAll->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter) );
882         fPtITSouterminPtTPCvsPtAllNPointTPC->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),nClustersTPC);
883         if(nClustersTPC>0.) fPtITSouterminPtTPCvsPtAllNPointTPCS->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),track->GetTPCnclsS()/nClustersTPC);
884         fPtITSouterminPtTPCvsPtAllDCAR->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),dca2D);
885         fPtITSouterminPtTPCvsPtAllDCAZ->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),dcaZ);
886         fPtITSouterminPtTPCvsPtAllPhi->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),phi);
887         fPtITSouterminPtTPCvsPtAllNPointITS->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),nPointITS);
888         fPtITSouterminPtTPCvsPtAllNSigmaToVertex->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),nSigmaToVertex);
889         fPtITSouterminPtTPCvsPtAllChi2C->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2C);
890         fPtITSouterminPtTPCvsPtAllRel1PtUncertainty->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),relUncertainty1Pt);
891         fPtITSouterminPtTPCvsPtAllChi2PerNClusTPC->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerClusterTPC);
892         if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITS->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
893         if(track->HasPointOnITSLayer(0)) {
894           fPtITSouterminPtTPCvsPtAllITSLayer0->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
895           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer0->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
896         }
897         if(!track->HasPointOnITSLayer(0) && track->HasPointOnITSLayer(1)) {
898           fPtITSouterminPtTPCvsPtAllITSLayer1->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
899           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer1->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
900         }
901         if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1) && track->HasPointOnITSLayer(2)) {
902           fPtITSouterminPtTPCvsPtAllITSLayer2->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
903           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer2->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
904         }
905         if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1) && !track->HasPointOnITSLayer(2) && track->HasPointOnITSLayer(3)) {
906           fPtITSouterminPtTPCvsPtAllITSLayer3->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
907           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer3->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
908         }
909         if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1) && !track->HasPointOnITSLayer(2) && !track->HasPointOnITSLayer(3) && track->HasPointOnITSLayer(4)) {
910           fPtITSouterminPtTPCvsPtAllITSLayer4->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
911           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer4->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
912         }
913         if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1) && !track->HasPointOnITSLayer(2) && !track->HasPointOnITSLayer(3) && !track->HasPointOnITSLayer(4) && track->HasPointOnITSLayer(5)) {
914           fPtITSouterminPtTPCvsPtAllITSLayer5->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
915           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSLayer5->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
916         }
917
918         if(!track->HasPointOnITSLayer(0) && !track->HasPointOnITSLayer(1)) {
919           fPtITSouterminPtTPCvsPtAllNoSPD->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
920           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSPD->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
921         }
922         if(!track->HasPointOnITSLayer(2) && !track->HasPointOnITSLayer(3)) {
923           fPtITSouterminPtTPCvsPtAllNoSDD->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
924           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSDD->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
925         }
926         if(!track->HasPointOnITSLayer(4) && !track->HasPointOnITSLayer(5)) {
927           fPtITSouterminPtTPCvsPtAllNoSSD->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter));
928           if(nPointITS>3) fPtITSouterminPtTPCvsPtAllChi2PerNClusITSNoSSD->Fill(pt,(1./ptITSouter-1./ptTPC)/(1./ptITSouter),chi2PerNPointITS);
929         }
930       }
931
932     }//fTrackCuts selection
933     
934     
935     //ITSrefit selection
936     if (fTrackCutsITS->AcceptTrack(track)) {
937       
938       fPtSelITS->Fill(pt);
939       fPtSelTPCITS->Fill(ptTPC);
940       fPtITSminPtTPCvsPtITS->Fill(pt,(1./pt-1./ptTPC)/(1./pt) );
941       fPtITSminPtTPCvsPtITSNPointTPC->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nClustersTPC);
942       if(nClustersTPC>0.) fPtITSminPtTPCvsPtITSNPointTPCS->Fill(pt,(1./pt-1./ptTPC)/(1./pt),track->GetTPCnclsS()/nClustersTPC);
943       fPtITSminPtTPCvsPtITSDCAR->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dca2D);
944       fPtITSminPtTPCvsPtITSDCAZ->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dcaZ);
945       fPtITSminPtTPCvsPtITSPhi->Fill(pt,(pt-ptTPC)/(pt),phi);
946       fPtITSminPtTPCvsPtITSNPointITS->Fill(pt,(pt-ptTPC)/(pt),nPointITS);
947       fPtITSminPtTPCvsPtITSNSigmaToVertex->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nSigmaToVertex);
948       fPtITSminPtTPCvsPtITSChi2C->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2C);
949       fPtITSminPtTPCvsPtITSRel1PtUncertainty->Fill(pt,(1./pt-1./ptTPC)/(1./pt),relUncertainty1Pt);
950       fPtITSminPtTPCvsPtITSChi2PerNClusTPC->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2PerClusterTPC);
951       if(nPointITS>3) fPtITSminPtTPCvsPtITSChi2PerNClusITS->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2PerNPointITS);
952
953       fPtRel1PtUncertaintyChi2PerClusTPC->Fill(pt,relUncertainty1Pt,chi2PerClusterTPC);
954       fPtNPointTPCSChi2PerClusTPC->Fill(pt,track->GetTPCnclsS()/nClustersTPC,chi2PerClusterTPC);
955       fPtNPointTPCSRel1PtUncertainty->Fill(pt,track->GetTPCnclsS()/nClustersTPC,relUncertainty1Pt);
956       if(track->GetTPCnclsS()/nClustersTPC>0.05) fPtRel1PtUncertaintyChi2PerClusTPCSharedSel->Fill(pt,relUncertainty1Pt,chi2PerClusterTPC);
957     }//fTrackCutsITS loop
958       
959   }//ESD track loop
960    
961   // Post output data
962   PostData(0, fHistList);
963   PostData(1, fHistListTPC);
964   PostData(2, fHistListITS);
965   PostData(3, fHistListCosmics);
966
967 }
968 //________________________________________________________________________
969 Bool_t AliPWG4HighPtQATPConly::IsCosmic(const AliESDtrack *track1 , Int_t trackNumber, Double_t ptMin)
970 {
971   Bool_t candidate = kFALSE;
972   Bool_t candidate2 = kFALSE;
973   if(!track1) return candidate;
974
975   Int_t nTracks = fESD->GetNumberOfTracks();
976  
977   for (Int_t iTrack = trackNumber+1; iTrack < nTracks; iTrack++) {
978     candidate2 = kFALSE;
979     AliESDtrack *track2 = fESD->GetTrack(iTrack);
980     if(!track2) continue;
981     if(!(fTrackCuts->AcceptTrack(track2))) continue;
982     if(track2->Pt()<ptMin) continue;
983     
984     //Check if same charge
985     if( (track1->GetSign()*track2->GetSign()) < 0. ) continue;
986     
987     //Check if back-to-back
988     Double_t mom1[3],mom2[3];
989     track1->GetPxPyPz(mom1);
990     track2->GetPxPyPz(mom2);
991     Double_t cosTheta = (mom1[0]*mom2[0]+mom1[1]*mom2[1]+mom1[2]*mom2[2])/( TMath::Sqrt(mom1[0]*mom1[0]+mom1[1]*mom1[1]+mom1[2]*mom1[2])*TMath::Sqrt(mom2[0]*mom2[0]+mom2[1]*mom2[1]+mom2[2]*mom2[2]) );
992     Double_t theta = TMath::ACos(cosTheta);
993     if(TMath::Abs(TMath::Pi()-theta)<fMaxCosmicAngle) { candidate = kTRUE; candidate2 = kTRUE;}
994
995     if(candidate2) {
996       fDeltaPtCosmicCandidates->Fill(track1->Pt()-track2->Pt());
997       fDeltaPhi->Fill(track1->Phi()-track2->Phi());
998       fDeltaEta->Fill(track1->Eta()-track2->Eta());
999     }
1000
1001   }
1002
1003   if(candidate) {
1004     fPtCosmicCandidates->Fill(track1->Pt());
1005   }
1006
1007    return candidate;
1008 }
1009
1010 //________________________________________________________________________
1011 void AliPWG4HighPtQATPConly::Terminate(Option_t *)
1012 {
1013
1014 }
1015
1016 #endif