]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliPWG4HighPtQATPConly.cxx
changes to run also on AOD MC
[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 #include "AliPWG4HighPtQATPConly.h"
26
27 #include "TH1.h"
28 #include "TH2.h"
29 #include "TH3.h"
30 #include "TList.h"
31 #include "TChain.h"
32 #include "TH3F.h"
33 #include "AliAnalysisManager.h"
34 #include "AliESDInputHandler.h"
35 #include "AliESDtrack.h"
36 #include "AliESDtrackCuts.h"
37 #include "AliExternalTrackParam.h"
38
39 using namespace std; //required for resolving the 'cout' symbol
40
41 ClassImp(AliPWG4HighPtQATPConly)
42
43 AliPWG4HighPtQATPConly::AliPWG4HighPtQATPConly(): AliAnalysisTask("AliPWG4HighPtQATPConly", ""), 
44   fESD(0), 
45   fTrackCuts(0), 
46   fTrackCutsITS(0),
47   fNEvent(0), // just to avoid warnings, inititialized in InitPointers too
48   fPtAll(0),  //
49   fPtSel(0),  //
50   fPtAllminPtTPCvsPtAll(0),
51   fPtAllminPtTPCvsPtAllNPointTPC(0),
52   fPtAllminPtTPCvsPtAllDCAR(0),
53   fPtAllminPtTPCvsPtAllDCAZ(0),
54   fPtAllminPtTPCvsPtAllPhi(0),
55   fPtAllminPtTPCvsPtAllNPointITS(0),
56   fPtAllminPtTPCvsPtAllNSigmaToVertex(0),
57   fPtAllminPtTPCvsPtAllChi2C(0),
58   fPtAllminPtTPCvsPtAllRel1PtUncertainty(0),
59   fHistList(0),
60   fPtAllTPC(0),
61   fPtSelTPC(0),
62   fPtSelTPCITS(0),
63   fHistListTPC(0),
64   fPtSelITS(0),
65   fPtITSminPtTPCvsPtITS(0),
66   fPtITSminPtTPCvsPtITSNPointTPC(0),
67   fPtITSminPtTPCvsPtITSDCAR(0),
68   fPtITSminPtTPCvsPtITSDCAZ(0),
69   fPtITSminPtTPCvsPtITSPhi(0),
70   fPtITSminPtTPCvsPtITSNPointITS(0),
71   fPtITSminPtTPCvsPtITSNSigmaToVertex(0),
72   fPtITSminPtTPCvsPtITSChi2C(0),
73   fPtITSminPtTPCvsPtITSRel1PtUncertainty(0),
74   fHistListITS(0)
75 {
76   InitHistPointers();
77 }
78 //________________________________________________________________________
79 AliPWG4HighPtQATPConly::AliPWG4HighPtQATPConly(const char *name): 
80   AliAnalysisTask(name, ""), 
81   fESD(0),
82   fTrackCuts(),
83   fTrackCutsITS(),
84   fNEvent(0), // just to avoid warnings, inititialized in InitPointers too
85   fPtAll(0),  //
86   fPtSel(0),  // 
87   fPtAllminPtTPCvsPtAll(0),
88   fPtAllminPtTPCvsPtAllNPointTPC(0),
89   fPtAllminPtTPCvsPtAllDCAR(0),
90   fPtAllminPtTPCvsPtAllDCAZ(0),
91   fPtAllminPtTPCvsPtAllPhi(0),
92   fPtAllminPtTPCvsPtAllNPointITS(0),
93   fPtAllminPtTPCvsPtAllNSigmaToVertex(0),
94   fPtAllminPtTPCvsPtAllChi2C(0),
95   fPtAllminPtTPCvsPtAllRel1PtUncertainty(0),
96   fHistList(0),
97   fPtAllTPC(0),
98   fPtSelTPC(0),
99   fPtSelTPCITS(0),
100   fHistListTPC(0),
101   fPtSelITS(0),
102   fPtITSminPtTPCvsPtITS(0),
103   fPtITSminPtTPCvsPtITSNPointTPC(0),
104   fPtITSminPtTPCvsPtITSDCAR(0),
105   fPtITSminPtTPCvsPtITSDCAZ(0),
106   fPtITSminPtTPCvsPtITSPhi(0),
107   fPtITSminPtTPCvsPtITSNPointITS(0),
108   fPtITSminPtTPCvsPtITSNSigmaToVertex(0),
109   fPtITSminPtTPCvsPtITSChi2C(0),
110   fPtITSminPtTPCvsPtITSRel1PtUncertainty(0),
111   fHistListITS(0)
112 {
113   //
114   // Constructor. Initialization of Inputs and Outputs
115   //
116   Info("AliPWG4HighPtQATPConly","Calling Constructor");
117   // Input slot #0 works with a TChain ESD
118   DefineInput(0, TChain::Class());
119   // Output slot #0 writes into a TList
120   DefineOutput(0, TList::Class());
121   // Output slot #1 writes into a TList
122   DefineOutput(1, TList::Class());
123   // Output slot #2 writes into a TList
124   DefineOutput(2, TList::Class());
125   InitHistPointers();
126   //  TH1::AddDirectory(kFALSE);
127   //  TH2::AddDirectory(kFALSE);
128   //  TH3::AddDirectory(kFALSE);
129   }
130
131 //________________________________________________________________________
132 void AliPWG4HighPtQATPConly::ConnectInputData(Option_t *) 
133 {
134   // Connect ESD here
135   // Called once
136   printf(">> AliPWG4HighPtQATPConly::ConnectInputData \n");
137   TTree* tree = dynamic_cast<TTree*> (GetInputData(0));
138   if (!tree) {
139     Printf("ERROR: Could not read chain from input slot 0");
140   } else {
141     
142     AliESDInputHandler *esdH = dynamic_cast<AliESDInputHandler*> (AliAnalysisManager::GetAnalysisManager()->GetInputEventHandler());
143     
144     if (!esdH) {
145       Printf("ERROR: Could not get ESDInputHandler");
146     } else
147       fESD = esdH->GetEvent();
148   }
149 }
150
151 //________________________________________________________________________
152 void AliPWG4HighPtQATPConly::InitHistPointers() {
153   //Initialize histograms
154   fNEvent = 0;
155   fPtAll = 0;
156   fPtSel = 0;
157   //TPC only histos compared to global tracking
158   fPtAllTPC = 0;
159   fPtSelTPC = 0;
160   fPtAllminPtTPCvsPtAll = 0;
161   fPtAllminPtTPCvsPtAllNPointTPC = 0;
162   fPtAllminPtTPCvsPtAllDCAR = 0;
163   fPtAllminPtTPCvsPtAllDCAZ = 0;
164   fPtAllminPtTPCvsPtAllPhi = 0;
165   fPtAllminPtTPCvsPtAllNPointITS = 0;
166   fPtAllminPtTPCvsPtAllNSigmaToVertex = 0;
167   fPtAllminPtTPCvsPtAllChi2C = 0;
168   fPtAllminPtTPCvsPtAllRel1PtUncertainty = 0;
169   //ITSrefit histos compared to TPConly tracks
170   fPtSelITS = 0;
171   fPtITSminPtTPCvsPtITS = 0;
172   fPtITSminPtTPCvsPtITSNPointTPC = 0;
173   fPtITSminPtTPCvsPtITSDCAR = 0;
174   fPtITSminPtTPCvsPtITSDCAZ = 0;
175   fPtITSminPtTPCvsPtITSPhi = 0;
176   fPtITSminPtTPCvsPtITSNPointITS = 0;
177   fPtITSminPtTPCvsPtITSNSigmaToVertex = 0;
178   fPtITSminPtTPCvsPtITSChi2C = 0;
179   fPtITSminPtTPCvsPtITSRel1PtUncertainty = 0;
180 }
181 //________________________________________________________________________
182 void AliPWG4HighPtQATPConly::CreateOutputObjects() {
183   //Create output objects
184   printf(">> AliPWG4HighPtQATPConly::CreateOutputObjects \n");
185   OpenFile(0);
186   fHistList = new TList();
187   OpenFile(1);
188   fHistListTPC = new TList();
189   OpenFile(2);
190   fHistListITS = new TList();
191
192   Int_t fgkNPhiBins=18;
193   Float_t kMinPhi = 0.;
194   Float_t kMaxPhi = 2.*TMath::Pi();
195   
196   Int_t fgkNPtBins=98;
197   Float_t fgkPtMin=2.;
198   Float_t fgkPtMax=100.;
199   Int_t fgkResPtBins=40;
200
201   fNEvent = new TH1F("fNEvent","NEvent",1,-0.5,0.5);
202   fHistList->Add(fNEvent);
203   fPtAll = new TH1F("fPtAll","PtAll",fgkNPtBins, fgkPtMin, fgkPtMax);
204   fHistList->Add(fPtAll);
205   fPtSel = new TH1F("fPtSel","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
206   fHistList->Add(fPtSel);
207   
208   fPtAllminPtTPCvsPtAll = new TH2F("fPtAllminPtTPCvsPtAll","PtAllminPtTPCvsPtAll",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
209   fPtAllminPtTPCvsPtAll->SetXTitle("p_{t}^{All}");
210   fPtAllminPtTPCvsPtAll->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
211   fHistList->Add(fPtAllminPtTPCvsPtAll);
212   
213   fPtAllminPtTPCvsPtAllNPointTPC = new TH3F("fPtAllminPtTPCvsPtAllNPointTPC","PtAllminPtTPCvsPtAllNPointTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,160,0.5,160.5);
214   fPtAllminPtTPCvsPtAllNPointTPC->SetXTitle("p_{t}^{All}");
215   fPtAllminPtTPCvsPtAllNPointTPC->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
216   fPtAllminPtTPCvsPtAllNPointTPC->SetZTitle("N_{point,TPC}");
217   fHistList->Add(fPtAllminPtTPCvsPtAllNPointTPC);
218
219   fPtAllminPtTPCvsPtAllDCAR = new TH3F("fPtAllminPtTPCvsPtAllDCAR","PtAllminPtTPCvsPtAllDCAR",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,-1.,1.);
220   fPtAllminPtTPCvsPtAllDCAR->SetXTitle("p_{t}^{All}");
221   fPtAllminPtTPCvsPtAllDCAR->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
222   fPtAllminPtTPCvsPtAllDCAR->SetZTitle("DCA_{R}");
223   fHistList->Add(fPtAllminPtTPCvsPtAllDCAR);
224
225   fPtAllminPtTPCvsPtAllDCAZ = new TH3F("fPtAllminPtTPCvsPtAllDCAZ","PtAllminPtTPCvsPtAllDCAZ",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,-2.,2.);
226   fPtAllminPtTPCvsPtAllDCAZ->SetXTitle("p_{t}^{All}");
227   fPtAllminPtTPCvsPtAllDCAZ->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
228   fPtAllminPtTPCvsPtAllDCAZ->SetZTitle("DCA_{Z}");
229   fHistList->Add(fPtAllminPtTPCvsPtAllDCAZ);
230
231   fPtAllminPtTPCvsPtAllPhi = new TH3F("fPtAllminPtTPCvsPtAllPhi","PtAllminPtTPCvsPtAllPhi",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkNPhiBins,kMinPhi,kMaxPhi);
232   fPtAllminPtTPCvsPtAllPhi->SetXTitle("p_{t}^{All}");
233   fPtAllminPtTPCvsPtAllPhi->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
234   fPtAllminPtTPCvsPtAllPhi->SetZTitle("#phi");
235   fHistList->Add(fPtAllminPtTPCvsPtAllPhi);
236
237   fPtAllminPtTPCvsPtAllNPointITS = new TH3F("fPtAllminPtTPCvsPtAllNPointITS","PtAllminPtTPCvsPtAllNPointITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,9,-0.5,8.5);
238   fPtAllminPtTPCvsPtAllNPointITS->SetXTitle("p_{t}^{All}");
239   fPtAllminPtTPCvsPtAllNPointITS->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
240   fPtAllminPtTPCvsPtAllNPointITS->SetZTitle("N_{point,ITS}}");
241   fHistList->Add(fPtAllminPtTPCvsPtAllNPointITS);
242   
243   fPtAllminPtTPCvsPtAllNSigmaToVertex = new TH3F("fPtAllminPtTPCvsPtAllNSigmaToVertex","PtAllminPtTPCvsPtAllNSigmaToVertex",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,0.,8.);
244   fPtAllminPtTPCvsPtAllNSigmaToVertex->SetXTitle("p_{t}^{All}");
245   fPtAllminPtTPCvsPtAllNSigmaToVertex->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
246   fPtAllminPtTPCvsPtAllNSigmaToVertex->SetZTitle("N#sigma to vertex");
247   fHistList->Add(fPtAllminPtTPCvsPtAllNSigmaToVertex);
248
249   fPtAllminPtTPCvsPtAllChi2C = new TH3F("fPtAllminPtTPCvsPtAllChi2C","PtAllminPtTPCvsPtAllChi2C",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,0.,10.);
250   fPtAllminPtTPCvsPtAllChi2C->SetXTitle("p_{t}^{All}");
251   fPtAllminPtTPCvsPtAllChi2C->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
252   fPtAllminPtTPCvsPtAllChi2C->SetZTitle("Constrained #chi^{2}");
253   fHistList->Add(fPtAllminPtTPCvsPtAllChi2C);
254
255   fPtAllminPtTPCvsPtAllRel1PtUncertainty = new TH3F("fPtAllminPtTPCvsPtAllRel1PtUncertainty","PtAllminPtTPCvsPtAllRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,30,0.,0.3);
256   fPtAllminPtTPCvsPtAllRel1PtUncertainty->SetXTitle("p_{t}^{All}");
257   fPtAllminPtTPCvsPtAllRel1PtUncertainty->SetYTitle("(1/p_{t}^{All}-1/p_{t}^{TPC})/(1/p_{t}^{All})");
258   fPtAllminPtTPCvsPtAllRel1PtUncertainty->SetZTitle("Rel1PtUncertainty");
259   fHistList->Add(fPtAllminPtTPCvsPtAllRel1PtUncertainty);
260
261   //ITSrefit
262   fPtSelITS = new TH1F("fPtSelITSrefit","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
263   fHistListITS->Add(fPtSelITS);
264   
265   fPtITSminPtTPCvsPtITS = new TH2F("fPtITSminPtTPCvsPtITS","PtITSminPtTPCvsPtITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.);
266   fPtITSminPtTPCvsPtITS->SetXTitle("p_{t}^{ITS}");
267   fPtITSminPtTPCvsPtITS->SetYTitle("(1/p_{t}^{ITS}-1/p_{t}^{TPC})/(1/p_{t}^{ITS})");
268   fHistListITS->Add(fPtITSminPtTPCvsPtITS);
269   
270   fPtITSminPtTPCvsPtITSNPointTPC = new TH3F("fPtITSminPtTPCvsPtITSNPointTPC","PtITSminPtTPCvsPtITSNPointTPC",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,160,0.5,160.5);
271   fPtITSminPtTPCvsPtITSNPointTPC->SetXTitle("p_{t}^{ITSrefit}");
272   fPtITSminPtTPCvsPtITSNPointTPC->SetYTitle("(1/p_{t}^{ITSrefit}-1/p_{t}^{TPC})/(1/p_{t}^{ITSrefit})");
273   fPtITSminPtTPCvsPtITSNPointTPC->SetZTitle("N_{point,TPC}");
274   fHistListITS->Add(fPtITSminPtTPCvsPtITSNPointTPC);
275     
276   fPtITSminPtTPCvsPtITSDCAR = new TH3F("fPtITSminPtTPCvsPtITSDCAR","PtITSminPtTPCvsPtITSDCAR",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,-1.,1.);
277   fPtITSminPtTPCvsPtITSDCAR->SetXTitle("p_{t}^{ITSrefit}");
278   fPtITSminPtTPCvsPtITSDCAR->SetYTitle("(1/p_{t}^{ITSrefit}-1/p_{t}^{TPC})/(1/p_{t}^{ITSrefit})");
279   fPtITSminPtTPCvsPtITSDCAR->SetZTitle("DCA_{R}");
280   fHistListITS->Add(fPtITSminPtTPCvsPtITSDCAR);
281   
282   fPtITSminPtTPCvsPtITSDCAZ = new TH3F("fPtITSminPtTPCvsPtITSDCAZ","PtITSminPtTPCvsPtITSDCAZ",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,-2.,2.);
283   fPtITSminPtTPCvsPtITSDCAZ->SetXTitle("p_{t}^{ITSrefit}");
284   fPtITSminPtTPCvsPtITSDCAZ->SetYTitle("(1/p_{t}^{ITSrefit}-1/p_{t}^{TPC})/(1/p_{t}^{ITSrefit})");
285   fPtITSminPtTPCvsPtITSDCAZ->SetZTitle("DCA_{Z}");
286   fHistListITS->Add(fPtITSminPtTPCvsPtITSDCAZ);
287   
288   fPtITSminPtTPCvsPtITSPhi = new TH3F("fPtITSminPtTPCvsPtITSPhi","PtITSminPtTPCvsPtITSPhi",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,fgkNPhiBins,kMinPhi,kMaxPhi);
289   fPtITSminPtTPCvsPtITSPhi->SetXTitle("p_{t}^{ITSrefit}");
290   fPtITSminPtTPCvsPtITSPhi->SetYTitle("(1/p_{t}^{ITSrefit}-1/p_{t}^{TPC})/(1/p_{t}^{ITSrefit})");
291   fPtITSminPtTPCvsPtITSPhi->SetZTitle("#phi");
292   fHistListITS->Add(fPtITSminPtTPCvsPtITSPhi);
293   
294   fPtITSminPtTPCvsPtITSNPointITS = new TH3F("fPtITSminPtTPCvsPtITSNPointITS","PtITSminPtTPCvsPtITSNPointITS",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,9,-0.5,8.5);
295   fPtITSminPtTPCvsPtITSNPointITS->SetXTitle("p_{t}^{ITSrefit}");
296   fPtITSminPtTPCvsPtITSNPointITS->SetYTitle("(1/p_{t}^{ITSrefit}-1/p_{t}^{TPC})/(1/p_{t}^{ITSrefit})");
297   fPtITSminPtTPCvsPtITSNPointITS->SetZTitle("N_{point,ITS}}");
298   fHistListITS->Add(fPtITSminPtTPCvsPtITSNPointITS); 
299   
300   fPtITSminPtTPCvsPtITSNSigmaToVertex = new TH3F("fPtITSminPtTPCvsPtITSNSigmaToVertex","PtITSminPtTPCvsPtITSNSigmaToVertex",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,40,0.,8.);
301   fPtITSminPtTPCvsPtITSNSigmaToVertex->SetXTitle("p_{t}^{ITS}");
302   fPtITSminPtTPCvsPtITSNSigmaToVertex->SetYTitle("(1/p_{t}^{ITS}-1/p_{t}^{TPC})/(1/p_{t}^{ITS})");
303   fPtITSminPtTPCvsPtITSNSigmaToVertex->SetZTitle("N#sigma to vertex");
304   fHistListITS->Add(fPtITSminPtTPCvsPtITSNSigmaToVertex);
305
306   fPtITSminPtTPCvsPtITSChi2C = new TH3F("fPtITSminPtTPCvsPtITSChi2C","PtITSminPtTPCvsPtITSChi2C",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,20,0.,10.);
307   fPtITSminPtTPCvsPtITSChi2C->SetXTitle("p_{t}^{ITS}");
308   fPtITSminPtTPCvsPtITSChi2C->SetYTitle("(1/p_{t}^{ITS}-1/p_{t}^{TPC})/(1/p_{t}^{ITS})");
309   fPtITSminPtTPCvsPtITSChi2C->SetZTitle("Constrained #chi^{2}");
310   fHistListITS->Add(fPtITSminPtTPCvsPtITSChi2C);
311
312   fPtITSminPtTPCvsPtITSRel1PtUncertainty = new TH3F("fPtITSminPtTPCvsPtITSRel1PtUncertainty","PtITSminPtTPCvsPtITSRel1PtUncertainty",fgkNPtBins, fgkPtMin,fgkPtMax,fgkResPtBins,-1.,1.,30,0.,0.3);
313   fPtITSminPtTPCvsPtITSRel1PtUncertainty->SetXTitle("p_{t}^{ITS}");
314   fPtITSminPtTPCvsPtITSRel1PtUncertainty->SetYTitle("(1/p_{t}^{ITS}-1/p_{t}^{TPC})/(1/p_{t}^{ITS})");
315   fPtITSminPtTPCvsPtITSRel1PtUncertainty->SetZTitle("Rel1PtUncertainty");
316   fHistListITS->Add(fPtITSminPtTPCvsPtITSRel1PtUncertainty);
317
318   fPtAllTPC = new TH1F("fPtAllTPC","PtAll",fgkNPtBins, fgkPtMin, fgkPtMax);
319   fHistListTPC->Add(fPtAllTPC);
320   fPtSelTPC = new TH1F("fPtSelTPC","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
321   fHistListTPC->Add(fPtSelTPC);
322   fPtSelTPCITS = new TH1F("fPtSelTPCITS","PtSel",fgkNPtBins, fgkPtMin, fgkPtMax);
323   fHistListTPC->Add(fPtSelTPCITS);
324   
325 }
326 //________________________________________________________________________
327 void AliPWG4HighPtQATPConly::Exec(Option_t *) {  
328   // Main loop
329   // Called for each event
330   printf(">> AliPWG4HighPtQATPConly::Exec \n");
331   
332  if (!fESD) {
333     Printf("ERROR: fESD not available");
334     return;
335   }
336
337   const AliESDVertex *vtx = fESD->GetPrimaryVertex();
338
339   // Need vertex cut
340   if (vtx->GetNContributors() < 2)
341     return;
342
343   printf("Vertex title %s, status %d, nCont %d\n",vtx->GetTitle(), vtx->GetStatus(), vtx->GetNContributors());
344   // Need to keep track of evts without vertex
345   fNEvent->Fill(0.);
346
347   if(!fESD->GetNumberOfTracks() || fESD->GetNumberOfTracks()<2) return;
348   Int_t nTracks = fESD->GetNumberOfTracks();
349   printf("nTracks %d\n", nTracks);
350   
351   for (Int_t iTrack = 0; iTrack < nTracks; iTrack++) {
352     
353     AliESDtrack *track = fESD->GetTrack(iTrack);
354     AliExternalTrackParam *trackTPC = (AliExternalTrackParam *)track->GetTPCInnerParam();
355     AliESDtrack *trackTPConly = fTrackCuts->GetTPCOnlyTrack(fESD,iTrack); 
356     if(!track || !trackTPC || !trackTPConly) continue;
357
358     Float_t pt = track->Pt();
359     Float_t ptTPC = trackTPC->Pt();
360     Float_t phi = track->Phi();
361     Float_t dca2D, dcaZ;
362     track->GetImpactParameters(dca2D,dcaZ);
363     // Float_t dca2DTPC, dcaZTPC;
364     //track->GetImpactParametersTPC(dca2DTPC,dcaZTPC);
365     UChar_t itsMap = track->GetITSClusterMap();
366     Int_t nPointITS = 0;
367     for (Int_t i=0; i < 6; i++) {
368       if (itsMap & (1 << i))
369         nPointITS ++;
370     }
371     Float_t nSigmaToVertex = fTrackCuts->GetSigmaToVertex(track);// Calculates the number of sigma to the vertex for a track.
372     Float_t chi2C = track->GetConstrainedChi2();
373     // Float_t relUncertainty1Pt = TMath::Sqrt(extCov[14])*pt;
374     Float_t relUncertainty1Pt = TMath::Sqrt(track->GetSigma1Pt2())*pt;
375
376     fPtAll->Fill(pt);
377     fPtAllTPC->Fill(ptTPC);
378     
379     if (fTrackCuts->AcceptTrack(track)) {
380
381       fPtSel->Fill(pt);
382       
383       fPtSelTPC->Fill(ptTPC);
384       fPtAllminPtTPCvsPtAll->Fill(pt,(1./pt-1./ptTPC)/(1./pt) );
385       fPtAllminPtTPCvsPtAllNPointTPC->Fill(pt,(1./pt-1./ptTPC)/(1./pt),track->GetTPCNcls());
386       fPtAllminPtTPCvsPtAllDCAR->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dca2D);
387       fPtAllminPtTPCvsPtAllDCAZ->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dcaZ);
388       fPtAllminPtTPCvsPtAllPhi->Fill(pt,(1./pt-1./ptTPC)/(1./pt),phi);
389       fPtAllminPtTPCvsPtAllNPointITS->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nPointITS);
390       fPtAllminPtTPCvsPtAllNSigmaToVertex->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nSigmaToVertex);
391       fPtAllminPtTPCvsPtAllChi2C->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2C);
392       fPtAllminPtTPCvsPtAllRel1PtUncertainty->Fill(pt,(1./pt-1./ptTPC)/(1./pt),relUncertainty1Pt);
393     }//fTrackCuts selection
394     
395     
396     //ITSrefit selection
397     if (fTrackCutsITS->AcceptTrack(track)) {
398       
399       fPtSelITS->Fill(pt);
400       fPtSelTPCITS->Fill(ptTPC);
401       fPtITSminPtTPCvsPtITS->Fill(pt,(1./pt-1./ptTPC)/(1./pt) );
402       fPtITSminPtTPCvsPtITSNPointTPC->Fill(pt,(1./pt-1./ptTPC)/(1./pt),track->GetTPCNcls());
403       fPtITSminPtTPCvsPtITSDCAR->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dca2D);
404       fPtITSminPtTPCvsPtITSDCAZ->Fill(pt,(1./pt-1./ptTPC)/(1./pt),dcaZ);
405       fPtITSminPtTPCvsPtITSPhi->Fill(pt,(pt-ptTPC)/(pt),phi);
406       fPtITSminPtTPCvsPtITSNPointITS->Fill(pt,(pt-ptTPC)/(pt),nPointITS);
407       fPtITSminPtTPCvsPtITSNSigmaToVertex->Fill(pt,(1./pt-1./ptTPC)/(1./pt),nSigmaToVertex);
408       fPtITSminPtTPCvsPtITSChi2C->Fill(pt,(1./pt-1./ptTPC)/(1./pt),chi2C);
409       fPtITSminPtTPCvsPtITSRel1PtUncertainty->Fill(pt,(1./pt-1./ptTPC)/(1./pt),relUncertainty1Pt);
410     }//fTrackCutsITS loop
411       
412   }//ESD track loop
413    
414   // Post output data
415   PostData(0, fHistList);
416   PostData(1, fHistListTPC);
417   PostData(2, fHistListITS);
418
419 }
420 //________________________________________________________________________
421 void AliPWG4HighPtQATPConly::Terminate(Option_t *)
422 {
423   printf("->AliPWG4HighPtQATPConly::Terminate \n");
424
425   printf("-<AliPWG4HighPtQATPConly::Terminate \n");
426 }