1 //------------------------------------------------------------------------------
2 // Implementation of AliComparisonRes class. It keeps information from
3 // comparison of reconstructed and MC particle tracks. In addtion,
4 // it keeps selection cuts used during comparison. The comparison
5 // information is stored in the ROOT histograms. Analysis of these
6 // histograms can be done by using Analyse() class function. The result of
7 // the analysis (histograms/graphs) are stored in the folder which is
8 // a data member of AliComparisonRes.
10 // Author: J.Otwinowski 04/02/2008
11 //------------------------------------------------------------------------------
15 // after running comparison task, read the file, and get component
16 gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
19 TFile f("Output.root");
20 //AliComparisonRes * compObj = (AliComparisonRes*)f.Get("AliComparisonRes");
21 AliComparisonRes * compObj = (AliComparisonRes*)cOutput->FindObject("AliComparisonRes");
23 // analyse comparison data
26 // the output histograms/graphs will be stored in the folder "folderRes"
27 compObj->GetAnalysisFolder()->ls("*");
29 // user can save whole comparison object (or only folder with anlysed histograms)
30 // in the seperate output file (e.g.)
31 TFile fout("Analysed_Res.root","recreate");
32 compObj->Write(); // compObj->GetAnalysisFolder()->Write();
45 #include "TProfile2D.h"
50 #include "AliESDEvent.h"
52 #include "AliESDfriend.h"
53 #include "AliESDfriendTrack.h"
54 #include "AliESDVertex.h"
55 #include "AliRecInfoCuts.h"
56 #include "AliMCInfoCuts.h"
58 #include "AliTracker.h"
60 #include "AliMathBase.h"
61 #include "AliTreeDraw.h"
63 #include "AliMCInfo.h"
64 #include "AliESDRecInfo.h"
65 #include "AliComparisonRes.h"
69 ClassImp(AliComparisonRes)
71 //_____________________________________________________________________________
72 AliComparisonRes::AliComparisonRes():
73 AliComparisonObject("AliComparisonRes"),
76 fMCVertex(0), //-> MC primary vertex
77 fRecVertex(0), //-> Reconstructed primary vertex
96 fPhiResolTanTPCITS(0),
97 fTanResolTanTPCITS(0),
102 // Resolution constrained param
104 fCPhiResolTan(0), // angular resolution - constrained
105 fCTanResolTan(0), // angular resolution - constrained
106 fCPtResolTan(0), // pt resolution - constrained
107 fCPhiPullTan(0), // angular resolution - constrained
108 fCTanPullTan(0), // angular resolution - constrained
109 fCPtPullTan(0), // pt resolution - constrained
112 // Parametrisation histograms
115 f1Pt2ResolS1PtTPC(0),
116 f1Pt2ResolS1PtTPCITS(0),
118 fYResolS1PtTPCITS(0),
120 fZResolS1PtTPCITS(0),
122 fPhiResolS1PtTPCITS(0),
123 fThetaResolS1PtTPC(0),
124 fThetaResolS1PtTPCITS(0),
127 fC1Pt2ResolS1PtTPC(0),
128 fC1Pt2ResolS1PtTPCITS(0),
130 fCYResolS1PtTPCITS(0),
132 fCZResolS1PtTPCITS(0),
133 fCPhiResolS1PtTPC(0),
134 fCPhiResolS1PtTPCITS(0),
135 fCThetaResolS1PtTPC(0),
136 fCThetaResolS1PtTPCITS(0),
148 //_____________________________________________________________________________
149 AliComparisonRes::~AliComparisonRes(){
152 if(fMCVertex) delete fMCVertex; fMCVertex=0;
153 if(fRecVertex) delete fRecVertex; fRecVertex=0;
155 // Global observables
157 if(fPhiTanPtTPC) delete fPhiTanPtTPC; fPhiTanPtTPC=0;
158 if(fPhiTanPtTPCITS) delete fPhiTanPtTPCITS; fPhiTanPtTPCITS=0;
160 // Resolution histograms
161 if(fPtResolTPC) delete fPtResolTPC; fPtResolTPC=0;
162 if(fPtPullTPC) delete fPtPullTPC; fPtPullTPC=0;
163 if(fPhiResolTanTPC) delete fPhiResolTanTPC; fPhiResolTanTPC=0;
164 if(fTanResolTanTPC) delete fTanResolTanTPC; fTanResolTanTPC=0;
165 if(fPhiPullTanTPC) delete fPhiPullTanTPC; fPhiPullTanTPC=0;
166 if(fTanPullTanTPC) delete fTanPullTanTPC; fTanPullTanTPC=0;
168 if(fPtResolTPCITS) delete fPtResolTPCITS; fPtResolTPCITS=0;
169 if(fPtPullTPCITS) delete fPtPullTPCITS; fPtPullTPCITS=0;
170 if(fPhiResolTanTPCITS) delete fPhiResolTanTPCITS; fPhiResolTanTPCITS=0;
171 if(fTanResolTanTPCITS) delete fTanResolTanTPCITS; fTanResolTanTPCITS=0;
172 if(fPhiPullTanTPCITS) delete fPhiPullTanTPCITS; fPhiPullTanTPCITS=0;
173 if(fTanPullTanTPCITS) delete fTanPullTanTPCITS; fTanPullTanTPCITS=0;
175 // Resolution histograms (constrained param)
176 if(fCPhiResolTan) delete fCPhiResolTan; fCPhiResolTan=0;
177 if(fCTanResolTan) delete fCTanResolTan; fCTanResolTan=0;
178 if(fCPtResolTan) delete fCPtResolTan; fCPtResolTan=0;
179 if(fCPhiPullTan) delete fCPhiPullTan; fCPhiPullTan=0;
180 if(fCTanPullTan) delete fCTanPullTan; fCTanPullTan=0;
181 if(fCPtPullTan) delete fCPtPullTan; fCPtPullTan=0;
183 // Parametrisation histograms
185 if(f1Pt2ResolS1PtTPC) delete f1Pt2ResolS1PtTPC; f1Pt2ResolS1PtTPC=0;
186 if(f1Pt2ResolS1PtTPCITS) delete f1Pt2ResolS1PtTPCITS; f1Pt2ResolS1PtTPCITS=0;
187 if(fYResolS1PtTPC) delete fYResolS1PtTPC; fYResolS1PtTPC=0;
188 if(fYResolS1PtTPCITS) delete fYResolS1PtTPCITS; fYResolS1PtTPCITS=0;
189 if(fZResolS1PtTPC) delete fZResolS1PtTPC; fZResolS1PtTPC=0;
190 if(fZResolS1PtTPCITS) delete fZResolS1PtTPCITS; fZResolS1PtTPCITS=0;
191 if(fPhiResolS1PtTPC) delete fPhiResolS1PtTPC; fPhiResolS1PtTPC=0;
192 if(fPhiResolS1PtTPCITS) delete fPhiResolS1PtTPCITS; fPhiResolS1PtTPCITS=0;
193 if(fThetaResolS1PtTPC) delete fThetaResolS1PtTPC; fThetaResolS1PtTPC=0;
194 if(fThetaResolS1PtTPCITS) delete fThetaResolS1PtTPCITS; fThetaResolS1PtTPCITS=0;
197 if(fC1Pt2ResolS1PtTPC) delete fC1Pt2ResolS1PtTPC; fC1Pt2ResolS1PtTPC=0;
198 if(fC1Pt2ResolS1PtTPCITS) delete fC1Pt2ResolS1PtTPCITS; fC1Pt2ResolS1PtTPCITS=0;
199 if(fCYResolS1PtTPC) delete fCYResolS1PtTPC; fCYResolS1PtTPC=0;
200 if(fCYResolS1PtTPCITS) delete fCYResolS1PtTPCITS; fCYResolS1PtTPCITS=0;
201 if(fCZResolS1PtTPC) delete fCZResolS1PtTPC; fCZResolS1PtTPC=0;
202 if(fCZResolS1PtTPCITS) delete fCZResolS1PtTPCITS; fCZResolS1PtTPCITS=0;
203 if(fCPhiResolS1PtTPC) delete fCPhiResolS1PtTPC; fCPhiResolS1PtTPC=0;
204 if(fCPhiResolS1PtTPCITS) delete fCPhiResolS1PtTPCITS; fCPhiResolS1PtTPCITS=0;
205 if(fCThetaResolS1PtTPC) delete fCThetaResolS1PtTPC; fCThetaResolS1PtTPC=0;
206 if(fCThetaResolS1PtTPCITS) delete fCThetaResolS1PtTPCITS; fCThetaResolS1PtTPCITS=0;
208 if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
211 //_____________________________________________________________________________
212 void AliComparisonRes::Init(){
214 // Control histograms
215 fMCVertex = new TH3F("MCVertex","MC vertex Xv:Yv:Zv",100,-0.05,0.05,100,-0.05,0.05,100,-10,10);
216 fMCVertex->SetXTitle("Xv [cm]");
217 fMCVertex->SetYTitle("Yv [cm]");
218 fMCVertex->SetZTitle("Zv [cm]");
220 fRecVertex = new TH3F("RecVertex","Rec prim. vertex Xv:Yv:Zv",100,-0.01,0.01,100,-0.01,0.01,100,-0.5,0.5);
221 fRecVertex->SetXTitle("Xv [cm]");
222 fRecVertex->SetYTitle("Yv [cm]");
223 fRecVertex->SetZTitle("Zv [cm]");
226 fPhiTanPtTPC = new TH3F("PhiTanPtTPC","phi vs tan#theta vs pt - TPC only",200,-3.5,3.5,50,-2,2,100,0.1,3.0);
227 fPhiTanPtTPC->SetXTitle("#phi [rad]");
228 fPhiTanPtTPC->SetYTitle("tan#theta");
229 fPhiTanPtTPC->SetZTitle("p_{t} [GeV]");
231 fPhiTanPtTPCITS = new TH3F("PhiTanPtTPCITS","phi vs tan#theta vs pt - TPC+ITS",200,-3.5,3.5,50,-2,2,100,0.1,3.0);
232 fPhiTanPtTPCITS->SetXTitle("#phi [rad]");
233 fPhiTanPtTPCITS->SetYTitle("tan#theta");
234 fPhiTanPtTPCITS->SetZTitle("p_{t} [GeV]");
237 fPtResolTPC = new TH2F("PtResolTPC","pt resol",10, 0.1,3,200,-0.2,0.2);
238 fPtResolTPC->SetXTitle("p_{t} [GeV]");
239 fPtResolTPC->SetYTitle("#Deltap_{t}/p_{t}");
241 fPtPullTPC = new TH2F("fPtPullTPC","pt pull",10, 0.1,3,200,-6,6);
242 fPtPullTPC->SetXTitle("p_{t} [GeV]");
243 fPtPullTPC->SetYTitle("#Deltap_{t}/#Sigma");
245 fPhiResolTanTPC = new TH2F("PhiResolTanTPC","PhiResolTanTPC",50, -2,2,200,-0.025,0.025);
246 fPhiResolTanTPC->SetXTitle("tan(#theta)");
247 fPhiResolTanTPC->SetYTitle("#Delta#phi [rad]");
249 fTanResolTanTPC = new TH2F("TanResolTanTPC","TanResolTanTPC",50, -2,2,200,-0.025,0.025);
250 fTanResolTanTPC->SetXTitle("tan(#theta)");
251 fTanResolTanTPC->SetYTitle("#Delta#theta [rad]");
253 fPhiPullTanTPC = new TH2F("PhiPullTanTPC","PhiPullTanTPC",50, -2,2,200,-5,5);
254 fPhiPullTanTPC->SetXTitle("Tan(#theta)");
255 fPhiPullTanTPC->SetYTitle("#Delta#phi/#Sigma");
257 fTanPullTanTPC = new TH2F("TanPullTanTPC","TanPullTanTPC",50, -2,2,200,-5,5);
258 fTanPullTanTPC->SetXTitle("Tan(#theta)");
259 fTanPullTanTPC->SetYTitle("#Delta#theta/#Sigma");
261 fPtResolTPCITS = new TH2F("PtResolTPCITS","pt resol",10, 0.1,3,200,-0.2,0.2);
262 fPtResolTPCITS->SetXTitle("p_{t} [GeV]");
263 fPtResolTPCITS->SetYTitle("#Deltap_{t}/p_{t}");
265 fPtPullTPCITS = new TH2F("fPtPullTPCITS","pt pull",10, 0.1,3,200,-6,6);
266 fPtPullTPCITS->SetXTitle("p_{t} [GeV]");
267 fPtPullTPCITS->SetYTitle("#Deltap_{t}/#Sigma");
269 fPhiResolTanTPCITS = new TH2F("PhiResolTanTPCITS","PhiResolTanTPCITS",50, -2,2,200,-0.025,0.025);
270 fPhiResolTanTPCITS->SetXTitle("tan(#theta)");
271 fPhiResolTanTPCITS->SetYTitle("#Delta#phi [rad]");
273 fTanResolTanTPCITS = new TH2F("TanResolTanTPCITS","TanResolTanTPCITS",50, -2,2,200,-0.025,0.025);
274 fTanResolTanTPCITS->SetXTitle("tan(#theta)");
275 fTanResolTanTPCITS->SetYTitle("#Delta#theta [rad]");
277 fPhiPullTanTPCITS = new TH2F("PhiPullTanTPCITS","PhiPullTanTPCITS",50, -2,2,200,-5,5);
278 fPhiPullTanTPCITS->SetXTitle("Tan(#theta)");
279 fPhiPullTanTPCITS->SetYTitle("#Delta#phi/#Sigma");
281 fTanPullTanTPCITS = new TH2F("TanPullTanTPCITS","TanPullTanTPCITS",50, -2,2,200,-5,5);
282 fTanPullTanTPCITS->SetXTitle("Tan(#theta)");
283 fTanPullTanTPCITS->SetYTitle("#Delta#theta/#Sigma");
285 // Resolution constrained param
286 fCPhiResolTan = new TH2F("CPhiResolTan","CPhiResolTan",50, -2,2,200,-0.025,0.025);
287 fCPhiResolTan->SetXTitle("tan(#theta)");
288 fCPhiResolTan->SetYTitle("#Delta#phi [rad]");
290 fCTanResolTan = new TH2F("CTanResolTan","CTanResolTan",50, -2,2,200,-0.025,0.025);
291 fCTanResolTan->SetXTitle("tan(#theta)");
292 fCTanResolTan->SetYTitle("#Delta#theta [rad]");
294 fCPtResolTan=new TH2F("CPtResol","CPtResol",50, -2,2,200,-0.2,0.2);
295 fCPtResolTan->SetXTitle("Tan(#theta)");
296 fCPtResolTan->SetYTitle("#Deltap_{t}/p_{t}");
298 fCPhiPullTan = new TH2F("CPhiPullTan","CPhiPullTan",50, -2,2,200,-5,5);
299 fCPhiPullTan->SetXTitle("Tan(#theta)");
300 fCPhiPullTan->SetYTitle("#Delta#phi/#Sigma");
302 fCTanPullTan = new TH2F("CTanPullTan","CTanPullTan",50, -2,2,200,-5,5);
303 fCTanPullTan->SetXTitle("Tan(#theta)");
304 fCTanPullTan->SetYTitle("#Delta#theta/#Sigma");
306 fCPtPullTan=new TH2F("CPtPull","CPtPull",50, -2,2,200,-5,5);
307 fCPtPullTan->SetXTitle("Tan(#theta)");
308 fCPtPullTan->SetYTitle("(1/mcp_{t}-1/p_{t})/#Sigma");
311 // Parametrisation histograms
314 f1Pt2ResolS1PtTPC = new TH2F("f1Pt2ResolS1PtTPC","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs sqrt(1/pt))",100,0,3,200,-0.010,0.010);
315 f1Pt2ResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
316 f1Pt2ResolS1PtTPC->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
318 f1Pt2ResolS1PtTPCITS = new TH2F("f1Pt2ResolS1PtTPCITS","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs sqrt(1/pt))",100,0,3,200,-0.010,0.010);
319 f1Pt2ResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
320 f1Pt2ResolS1PtTPCITS->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
322 fYResolS1PtTPC = new TH2F("fYResolS1PtTPC","fYResolS1PtTPC",100, 0,3,200,-1.0,1.0);
323 fYResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
324 fYResolS1PtTPC->SetYTitle("#DeltaY");
326 fYResolS1PtTPCITS = new TH2F("fYResolS1PtTPCITS","fYResolS1PtTPCITS",100, 0,3,200,-0.05,0.05);
327 fYResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
328 fYResolS1PtTPCITS->SetYTitle("#DeltaY");
330 fZResolS1PtTPC = new TH2F("fZResolS1PtTPC","fZResolS1PtTPC",100, 0,3,200,-1.0,1.0);
331 fZResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
332 fZResolS1PtTPC->SetYTitle("#DeltaZ");
334 fZResolS1PtTPCITS = new TH2F("fZResolS1PtTPCITS","fZResolS1PtTPCITS",100, 0,3,200,-0.05,0.05);
335 fZResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
336 fZResolS1PtTPCITS->SetYTitle("#DeltaZ");
338 fPhiResolS1PtTPC = new TH2F("fPhiResolS1PtTPC","fPhiResolS1PtTPC",100, 0,3,200,-0.025,0.025);
339 fPhiResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
340 fPhiResolS1PtTPC->SetYTitle("#Delta#phi");
342 fPhiResolS1PtTPCITS = new TH2F("fPhiResolS1PtTPCITS","fPhiResolS1PtTPCITS",100, 0,3,200,-0.01,0.01);
343 fPhiResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
344 fPhiResolS1PtTPCITS->SetYTitle("#Delta#phi");
346 fThetaResolS1PtTPC = new TH2F("fThetaResolS1PtTPC","fThetaResolS1PtTPC",100, 0,3,200,-0.025,0.025);
347 fThetaResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
348 fThetaResolS1PtTPC->SetYTitle("#Delta#theta");
350 fThetaResolS1PtTPCITS = new TH2F("fThetaResolS1PtTPCITS","fThetaResolS1PtTPCITS",100, 0,3,200,-0.01,0.01);
351 fThetaResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
352 fThetaResolS1PtTPCITS->SetYTitle("#Delta#theta");
355 fC1Pt2ResolS1PtTPC = new TH2F("fC1Pt2ResolS1PtTPC","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs 1/pt)",100,0,3,200,-0.010,0.010);
356 fC1Pt2ResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
357 fC1Pt2ResolS1PtTPC->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
359 fC1Pt2ResolS1PtTPCITS = new TH2F("fC1Pt2ResolS1PtTPCITS","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs 1/pt)",100,0,3,200,-0.010,0.010);
360 fC1Pt2ResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
361 fC1Pt2ResolS1PtTPCITS->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
363 fCYResolS1PtTPC = new TH2F("fCYResolS1PtTPC","fCYResolS1PtTPC",100, 0,3,200,-1.0,1.0);
364 fCYResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
365 fCYResolS1PtTPC->SetYTitle("#DeltaY");
367 fCYResolS1PtTPCITS = new TH2F("fCYResolS1PtTPCITS","fCYResolS1PtTPCITS",100, 0,3,200,-0.01,0.01);
368 fCYResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
369 fCYResolS1PtTPCITS->SetYTitle("#DeltaY");
371 fCZResolS1PtTPC = new TH2F("fCZResolS1PtTPC","fCZResolS1PtTPC",100, 0,3,200,-1.0,1.0);
372 fCZResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
373 fCZResolS1PtTPC->SetYTitle("#DeltaZ");
375 fCZResolS1PtTPCITS = new TH2F("fCZResolS1PtTPCITS","fCZResolS1PtTPCITS",100, 0,3,200,-0.025,0.025);
376 fCZResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
377 fCZResolS1PtTPCITS->SetYTitle("#DeltaZ");
379 fCPhiResolS1PtTPC = new TH2F("fCPhiResolS1PtTPC","fCPhiResolS1PtTPC",100, 0,3,200,-0.025,0.025);
380 fCPhiResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
381 fCPhiResolS1PtTPC->SetYTitle("#Delta#phi");
383 fCPhiResolS1PtTPCITS = new TH2F("fCPhiResolS1PtTPCITS","fCPhiResolS1PtTPCITS",100, 0,3,200,-0.003,0.003);
384 fCPhiResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
385 fCPhiResolS1PtTPCITS->SetYTitle("#Delta#phi");
387 fCThetaResolS1PtTPC = new TH2F("fCThetaResolS1PtTPC","fCThetaResolS1PtTPC",100, 0,3,200,-0.025,0.025);
388 fCThetaResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
389 fCThetaResolS1PtTPC->SetYTitle("#Delta#theta");
391 fCThetaResolS1PtTPCITS = new TH2F("fCThetaResolS1PtTPCITS","fCThetaResolS1PtTPCITS",100, 0,3,200,-0.005,0.005);
392 fCThetaResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
393 fCThetaResolS1PtTPCITS->SetYTitle("#Delta#theta");
397 AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
399 AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
402 fAnalysisFolder = CreateFolder("folderRes","Analysis Resolution Folder");
405 //_____________________________________________________________________________
406 void AliComparisonRes::Process(AliMCInfo* infoMC, AliESDRecInfo *infoRC)
408 // Fill resolution comparison information
409 AliExternalTrackParam *track = 0;
410 Double_t field = AliTracker::GetBz(); // nominal Bz field [kG]
411 Double_t kMaxD = 123456.0; // max distance
412 AliESDVertex vertexMC; // MC primary vertex
414 Int_t clusterITS[200];
415 Double_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
417 Float_t deltaPt, pullPt, deltaPhi,pullPhi, deltaTan, pullTan, delta1Pt2, deltaY1Pt, deltaZ1Pt, deltaPhi1Pt, deltaTheta1Pt;
418 Float_t deltaPtTPC, pullPtTPC, deltaPhiTPC, pullPhiTPC, deltaTanTPC, pullTanTPC, delta1Pt2TPC, deltaY1PtTPC, deltaZ1PtTPC, deltaPhi1PtTPC, deltaTheta1PtTPC;
420 Float_t mcpt = infoMC->GetParticle().Pt();
421 Float_t s1mcpt = TMath::Sqrt(1./infoMC->GetParticle().Pt());
422 Float_t tantheta = TMath::Tan(infoMC->GetParticle().Theta()-TMath::Pi()*0.5);
423 Float_t eta = infoMC->GetParticle().Eta();
425 // distance to Prim. vertex
426 const Double_t* dv = infoMC->GetVDist();
427 Bool_t isPrim = TMath::Sqrt(dv[0]*dv[0] + dv[1]*dv[1])<fCutsMC->GetMaxR() && TMath::Abs(dv[2])<fCutsMC->GetMaxVz();
429 // Check selection cuts
430 if (fCutsMC->IsPdgParticle(TMath::Abs(infoMC->GetParticle().GetPdgCode())) == kFALSE) return;
432 if (infoRC->GetStatus(1)!=3) return; // TPC refit
433 if (!infoRC->GetESDtrack()) return;
434 if (infoRC->GetESDtrack()->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return;
436 // calculate and set prim. vertex
437 vertexMC.SetXv( infoMC->GetParticle().Vx() - dv[0] );
438 vertexMC.SetYv( infoMC->GetParticle().Vy() - dv[1] );
439 vertexMC.SetZv( infoMC->GetParticle().Vz() - dv[2] );
441 // Fill MC vertex histo
442 fMCVertex->Fill(vertexMC.GetXv(),vertexMC.GetYv(),vertexMC.GetZv());
445 // Fill Rec vertex histo
446 //fRecVertex->Fill(infoRC->GetESDvertex()->GetXv(),infoRC->GetESDvertex()->GetYv(),infoRC->GetESDvertex()->GetZv());
447 //printf("Rec vertex xv %f, yv %f, zv %f \n",infoRC->GetESDvertex()->GetXv(),infoRC->GetESDvertex()->GetYv(),infoRC->GetESDvertex()->GetZv());
449 deltaPt= (mcpt-infoRC->GetESDtrack()->Pt())/mcpt;
450 pullPt= (1/mcpt-infoRC->GetESDtrack()->OneOverPt())/TMath::Sqrt(infoRC->GetESDtrack()->GetSigma1Pt2());
451 deltaPhi = TMath::ATan2(infoRC->GetESDtrack()->Py(),infoRC->GetESDtrack()->Px())-
452 TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
453 pullPhi = deltaPhi/TMath::Sqrt(infoRC->GetESDtrack()->GetSigmaSnp2());
455 deltaTan = TMath::ATan2(infoRC->GetESDtrack()->Pz(),infoRC->GetESDtrack()->Pt())-
456 TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
457 pullTan = deltaTan/TMath::Sqrt(infoRC->GetESDtrack()->GetSigmaTgl2());
459 delta1Pt2 = (1/mcpt-infoRC->GetESDtrack()->OneOverPt())/TMath::Power(1+1/mcpt,2);
460 deltaY1Pt = (infoMC->GetParticle().Vy()-infoRC->GetESDtrack()->GetY()) / (0.2+1/mcpt);
461 deltaZ1Pt = (infoMC->GetParticle().Vz()-infoRC->GetESDtrack()->GetZ()) / (0.2+1/mcpt);
462 deltaPhi1Pt = deltaPhi / (0.1+1/mcpt);
463 deltaTheta1Pt = deltaTan / (0.1+1/mcpt);
465 // calculate track parameters at vertex
466 const AliExternalTrackParam *innerTPC = 0;
467 if ((innerTPC = infoRC->GetESDtrack()->GetTPCInnerParam()) != 0)
469 if ((track = new AliExternalTrackParam(*infoRC->GetESDtrack()->GetTPCInnerParam())) != 0 )
471 Bool_t bDCAStatus = track->PropagateToDCA(&vertexMC,field,kMaxD,dca,cov);
473 // Fill parametrisation histograms (only TPC track)
477 deltaPtTPC= (mcpt-innerTPC->Pt())/mcpt;
478 pullPtTPC= (1/mcpt-innerTPC->OneOverPt())/TMath::Sqrt(innerTPC->GetSigma1Pt2());
479 deltaPhiTPC = TMath::ATan2(innerTPC->Py(),innerTPC->Px())-
480 TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
481 pullPhiTPC = deltaPhiTPC/TMath::Sqrt(innerTPC->GetSigmaSnp2());
483 deltaTanTPC = TMath::ATan2(innerTPC->Pz(),innerTPC->Pt())-
484 TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
485 pullTanTPC = deltaTanTPC/TMath::Sqrt(infoRC->GetESDtrack()->GetSigmaTgl2());
487 delta1Pt2TPC = (1/mcpt-innerTPC->OneOverPt())/TMath::Power(1+1/mcpt,2);
488 deltaY1PtTPC= (infoMC->GetParticle().Vy()-innerTPC->GetY()) / (0.2+1/mcpt);
489 deltaZ1PtTPC = (infoMC->GetParticle().Vz()-innerTPC->GetZ()) / (0.2+1/mcpt);
490 deltaPhi1PtTPC = deltaPhiTPC / (0.1+1/mcpt);
491 deltaTheta1PtTPC = deltaTanTPC / (0.1+1/mcpt);
493 fPhiTanPtTPC->Fill(TMath::ATan2(innerTPC->Py(),innerTPC->Px()),
494 TMath::ATan2(innerTPC->Pz(),innerTPC->Pt()),
497 f1Pt2ResolS1PtTPC->Fill(s1mcpt,delta1Pt2TPC);
498 fYResolS1PtTPC->Fill(s1mcpt,deltaY1PtTPC);
499 fZResolS1PtTPC->Fill(s1mcpt,deltaZ1PtTPC);
500 fPhiResolS1PtTPC->Fill(s1mcpt,deltaPhi1PtTPC);
501 fThetaResolS1PtTPC->Fill(s1mcpt,deltaTheta1PtTPC);
503 fPtResolTPC->Fill(mcpt,deltaPtTPC);
504 fPtPullTPC->Fill(mcpt,pullPtTPC);
505 fPhiResolTanTPC->Fill(tantheta,deltaPhiTPC);
506 fPhiPullTanTPC->Fill(tantheta,pullPhiTPC);
507 fTanResolTanTPC->Fill(tantheta,deltaTanTPC);
508 fTanPullTanTPC->Fill(tantheta,pullTanTPC);
514 // TPC and ITS clusters in the system
515 if(infoRC->GetESDtrack()->GetITSclusters(clusterITS)>fCutsRC->GetMinNClustersITS())
517 fPhiTanPtTPCITS->Fill(TMath::ATan2(infoRC->GetESDtrack()->Py(),infoRC->GetESDtrack()->Px()),
518 TMath::ATan2(infoRC->GetESDtrack()->Pz(),infoRC->GetESDtrack()->Pt()),
519 infoRC->GetESDtrack()->Pt());
521 f1Pt2ResolS1PtTPCITS->Fill(s1mcpt,delta1Pt2);
522 fYResolS1PtTPCITS->Fill(s1mcpt,deltaY1Pt);
523 fZResolS1PtTPCITS->Fill(s1mcpt,deltaZ1Pt);
524 fPhiResolS1PtTPCITS->Fill(s1mcpt,deltaPhi1Pt);
525 fThetaResolS1PtTPCITS->Fill(s1mcpt,deltaTheta1Pt);
527 fPtResolTPCITS->Fill(mcpt,deltaPt);
528 fPtPullTPCITS->Fill(mcpt,pullPt);
530 fPhiResolTanTPCITS->Fill(tantheta,deltaPhi);
531 fPhiPullTanTPCITS->Fill(tantheta,pullPhi);
532 fTanResolTanTPCITS->Fill(tantheta,deltaTan);
533 fTanPullTanTPCITS->Fill(tantheta,pullTan);
537 //_____________________________________________________________________________
538 void AliComparisonRes::ProcessConstrained(AliMCInfo* infoMC, AliESDRecInfo *infoRC)
540 // Fill resolution comparison information (constarained parameters)
542 AliExternalTrackParam *track = 0;
543 Double_t field = AliTracker::GetBz(); // nominal Bz field [kG]
544 Double_t kMaxD = 123456.0; // max distance
545 AliESDVertex vertexMC; // MC primary vertex
547 Int_t clusterITS[200];
548 Double_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
550 Float_t deltaPt, pullPt, deltaPhi, pullPhi, deltaTan, pullTan, delta1Pt2, deltaY1Pt, deltaZ1Pt, deltaPhi1Pt, deltaTheta1Pt;
551 Float_t deltaPtTPC, pullPtTPC, deltaPhiTPC, deltaTanTPC, delta1Pt2TPC, deltaY1PtTPC, deltaZ1PtTPC, deltaPhi1PtTPC, deltaTheta1PtTPC;
553 Float_t mcpt = infoMC->GetParticle().Pt();
554 Float_t s1mcpt = TMath::Sqrt(1./infoMC->GetParticle().Pt());
555 Float_t tantheta = TMath::Tan(infoMC->GetParticle().Theta()-TMath::Pi()*0.5);
557 // distance to Prim. vertex
558 const Double_t* dv = infoMC->GetVDist();
559 Bool_t isPrim = TMath::Sqrt(dv[0]*dv[0] + dv[1]*dv[1])<fCutsMC->GetMaxR() && TMath::Abs(dv[2])<fCutsMC->GetMaxVz();
561 // Check selection cuts
563 if (fCutsMC->IsPdgParticle(TMath::Abs(infoMC->GetParticle().GetPdgCode())) == kFALSE) return;
565 if (infoRC->GetStatus(1)!=3) return;
566 if (!infoRC->GetESDtrack()) return;
567 if (infoRC->GetESDtrack()->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return;
568 if (!infoRC->GetESDtrack()->GetConstrainedParam()) return;
570 // calculate and set prim. vertex
571 vertexMC.SetXv( infoMC->GetParticle().Vx() - dv[0] );
572 vertexMC.SetYv( infoMC->GetParticle().Vy() - dv[1] );
573 vertexMC.SetZv( infoMC->GetParticle().Vz() - dv[2] );
575 // constrained parameters resolution
576 const AliExternalTrackParam * cparam = infoRC->GetESDtrack()->GetConstrainedParam();
577 deltaPt= (mcpt-cparam->Pt())/mcpt;
578 pullPt= (1/mcpt-cparam->OneOverPt())/TMath::Sqrt(cparam->GetSigma1Pt2());
579 deltaPhi = TMath::ATan2(cparam->Py(),cparam->Px())-
580 TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
581 pullPhi = deltaPhi/TMath::Sqrt(cparam->GetSigmaSnp2());
582 deltaTan = TMath::ATan2(cparam->Pz(),cparam->Pt())-TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
583 pullTan = deltaPhi/TMath::Sqrt(cparam->GetSigmaSnp2());
586 delta1Pt2 = (1/mcpt-cparam->OneOverPt())/TMath::Power(1+1/mcpt,2);
588 deltaY1Pt = (infoMC->GetParticle().Vy()-cparam->GetY()) / (0.2+1/mcpt);
589 deltaZ1Pt = (infoMC->GetParticle().Vz()-cparam->GetZ()) / (0.2+1/mcpt);
590 deltaPhi1Pt = deltaPhi / (0.1+1/mcpt);
591 deltaTheta1Pt = deltaTan / (0.1+1/mcpt);
593 // calculate track parameters at vertex
594 const AliExternalTrackParam *innerTPC = 0;
595 if ((innerTPC = infoRC->GetESDtrack()->GetTPCInnerParam()) != 0)
597 if ((track = new AliExternalTrackParam(*infoRC->GetESDtrack()->GetTPCInnerParam())) != 0 )
599 Bool_t bDCAStatus = track->PropagateToDCA(&vertexMC,field,kMaxD,dca,cov);
601 // Fill parametrisation histograms (only TPC track)
604 deltaPtTPC= (mcpt-innerTPC->Pt())/mcpt;
605 pullPtTPC= (1/mcpt-innerTPC->OneOverPt())/TMath::Sqrt(innerTPC->GetSigma1Pt2());
606 deltaPhiTPC = TMath::ATan2(innerTPC->Py(),innerTPC->Px())-
607 TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
609 deltaTanTPC = TMath::ATan2(innerTPC->Pz(),innerTPC->Pt())-
610 TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
612 delta1Pt2TPC = (1/mcpt-innerTPC->OneOverPt())/TMath::Power(1+1/mcpt,2);
613 deltaY1PtTPC= (infoMC->GetParticle().Vy()-innerTPC->GetY()) / (0.2+1/mcpt);
614 deltaZ1PtTPC = (infoMC->GetParticle().Vz()-innerTPC->GetZ()) / (0.2+1/mcpt);
615 deltaPhi1PtTPC = deltaPhiTPC / (0.1+1/mcpt);
616 deltaTheta1PtTPC = deltaTanTPC / (0.1+1/mcpt);
618 fC1Pt2ResolS1PtTPC->Fill(s1mcpt,delta1Pt2TPC);
619 fCYResolS1PtTPC->Fill(s1mcpt,deltaY1PtTPC);
620 fCZResolS1PtTPC->Fill(s1mcpt,deltaZ1PtTPC);
621 fCPhiResolS1PtTPC->Fill(s1mcpt,deltaPhi1PtTPC);
622 fCThetaResolS1PtTPC->Fill(s1mcpt,deltaTheta1PtTPC);
628 // TPC and ITS in the system
629 if(infoRC->GetESDtrack()->GetITSclusters(clusterITS)>fCutsRC->GetMinNClustersITS())
631 fC1Pt2ResolS1PtTPCITS->Fill(s1mcpt,delta1Pt2);
632 fCYResolS1PtTPCITS->Fill(s1mcpt,deltaY1Pt);
633 fCZResolS1PtTPCITS->Fill(s1mcpt,deltaZ1Pt);
634 fCPhiResolS1PtTPCITS->Fill(s1mcpt,deltaPhi1Pt);
635 fCThetaResolS1PtTPCITS->Fill(s1mcpt,deltaTheta1Pt);
638 fCPtResolTan->Fill(tantheta,deltaPt);
639 fCPtPullTan->Fill(tantheta,pullPt);
640 fCPhiResolTan->Fill(tantheta,deltaPhi);
641 fCPhiPullTan->Fill(tantheta,pullPhi);
642 fCTanResolTan->Fill(tantheta,deltaTan);
643 fCTanPullTan->Fill(tantheta,pullTan);
647 //_____________________________________________________________________________
648 void AliComparisonRes::Exec(AliMCInfo* infoMC, AliESDRecInfo *infoRC){
650 // Process comparison information
651 Process(infoMC,infoRC);
652 ProcessConstrained(infoMC,infoRC);
655 //_____________________________________________________________________________
656 TH1F* AliComparisonRes::MakeResol(TH2F * his, Int_t integ, Bool_t type){
657 // Create resolution histograms
660 if (!gPad) new TCanvas;
661 hisr = AliTreeDraw::CreateResHistoI(his,&hism,integ);
662 if (type) return hism;
667 //_____________________________________________________________________________
668 void AliComparisonRes::Analyse(){
669 // Analyse comparison information and store output histograms
670 // in the folder "folderRes"
673 TH1::AddDirectory(kFALSE);
675 AliComparisonRes * comp=this;
677 TObjArray *aFolderObj = new TObjArray;
679 // write results in the folder
681 TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
684 hiss = comp->MakeResol(comp->fPtResolTPC,1,0);
685 hiss->SetXTitle("p_{t} (GeV)");
686 hiss->SetYTitle("#sigmap_{t}/p_{t}");
688 hiss->SetName("PtResolPtTPC");
690 aFolderObj->Add(hiss);
693 hiss = comp->MakeResol(comp->fPtResolTPC,1,1);
694 hiss->SetXTitle("p_{t} (GeV)");
695 hiss->SetYTitle("mean #Deltap_{t}/p_{t}");
697 hiss->SetName("PtMeanPtTPC");
699 aFolderObj->Add(hiss);
702 hiss = comp->MakeResol(comp->fPhiResolTanTPC,1,0);
703 hiss->SetXTitle("Tan(#theta)");
704 hiss->SetYTitle("#sigma#phi (rad)");
706 hiss->SetName("PhiResolTanTPC");
708 aFolderObj->Add(hiss);
711 hiss = comp->MakeResol(comp->fPhiResolTanTPC,1,1);
712 hiss->SetXTitle("Tan(#theta)");
713 hiss->SetYTitle("mean #Delta#phi (rad)");
715 hiss->SetName("PhiMeanTanTPC");
718 aFolderObj->Add(hiss);
720 hiss = comp->MakeResol(comp->fTanResolTanTPC,1,0);
721 hiss->SetXTitle("Tan(#theta)");
722 hiss->SetYTitle("#sigma#theta (rad)");
724 hiss->SetName("ThetaResolTanTPC");
726 aFolderObj->Add(hiss);
729 hiss = comp->MakeResol(comp->fTanResolTanTPC,1,1);
730 hiss->SetXTitle("Tan(#theta)");
731 hiss->SetYTitle("mean #Delta#theta (rad)");
733 hiss->SetName("ThetaMeanTanTPC");
735 aFolderObj->Add(hiss);
739 hiss = comp->MakeResol(comp->fPtResolTPCITS,1,0);
740 hiss->SetXTitle("p_{t}");
741 hiss->SetYTitle("#sigmap_{t}/p_{t}");
743 hiss->SetName("PtResolPtTPCITS");
745 aFolderObj->Add(hiss);
748 hiss = comp->MakeResol(comp->fPtResolTPCITS,1,1);
749 hiss->SetXTitle("p_{t}");
750 hiss->SetYTitle("mean #Deltap_{t}/p_{t}");
752 hiss->SetName("PtMeanPtTPCITS");
754 aFolderObj->Add(hiss);
757 hiss = comp->MakeResol(comp->fPhiResolTanTPCITS,1,0);
758 hiss->SetXTitle("Tan(#theta)");
759 hiss->SetYTitle("#sigma#phi (rad)");
761 hiss->SetName("PhiResolTanTPCITS");
763 aFolderObj->Add(hiss);
766 hiss = comp->MakeResol(comp->fPhiResolTanTPCITS,1,1);
767 hiss->SetXTitle("Tan(#theta)");
768 hiss->SetYTitle("mean #Delta#phi (rad)");
770 hiss->SetName("PhiMeanTanTPCITS");
772 aFolderObj->Add(hiss);
775 hiss = comp->MakeResol(comp->fTanResolTanTPCITS,1,0);
776 hiss->SetXTitle("Tan(#theta)");
777 hiss->SetYTitle("#sigma#theta (rad)");
779 hiss->SetName("ThetaResolTanTPCITS");
781 aFolderObj->Add(hiss);
784 hiss = comp->MakeResol(comp->fTanResolTanTPCITS,1,1);
785 hiss->SetXTitle("Tan(#theta)");
786 hiss->SetYTitle("mean #Delta#theta (rad)");
788 hiss->SetName("ThetaMeanTanTPCITS");
790 aFolderObj->Add(hiss);
794 hiss = comp->MakeResol(comp->fCPtResolTan,1,0);
795 hiss->SetXTitle("Tan(#theta)");
796 hiss->SetYTitle("#sigmap_{t}/p_{t}");
798 hiss->SetName("CptResolTan");
800 aFolderObj->Add(hiss);
802 hiss = comp->MakeResol(comp->fCPhiResolTan,1,0);
803 hiss->SetXTitle("Tan(#theta)");
804 hiss->SetYTitle("#sigma#phi (rad)");
806 hiss->SetName("CPhiResolTan");
808 aFolderObj->Add(hiss);
810 hiss = comp->MakeResol(comp->fCTanResolTan,1,0);
811 hiss->SetXTitle("Tan(#theta)");
812 hiss->SetYTitle("#sigma#theta (rad)");
814 hiss->SetName("CThetaResolTan");
816 aFolderObj->Add(hiss);
818 hiss = comp->MakeResol(comp->fCPtPullTan,1,0);
819 hiss->SetXTitle("Tan(#theta)");
820 hiss->SetYTitle("1/mcp_{t}-1/p_{t}/#Sigma(1/p_{t})");
822 hiss->SetName("CptPullTan");
824 aFolderObj->Add(hiss);
826 hiss = comp->MakeResol(comp->fC1Pt2ResolS1PtTPC,1,0);
827 hiss->SetXTitle("#sqrt(1/mcp_{t})");
828 hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
830 hiss->SetName("C1Pt2ResolS1PtTPC");
832 aFolderObj->Add(hiss);
834 hiss = comp->MakeResol(comp->fC1Pt2ResolS1PtTPCITS,1,0);
835 hiss->SetXTitle("#sqrt(1/mcp_{t})");
836 hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
838 hiss->SetName("C1Pt2ResolS1PtTPCITS");
840 aFolderObj->Add(hiss);
842 hiss = comp->MakeResol(comp->fCYResolS1PtTPC,1,0);
843 hiss->SetXTitle("#sqrt(1/mcp_{t})");
844 hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
846 hiss->SetName("CYResolS1PtTPC");
848 aFolderObj->Add(hiss);
850 hiss = comp->MakeResol(comp->fCYResolS1PtTPCITS,1,0);
851 hiss->SetXTitle("#sqrt(1/mcp_{t})");
852 hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
854 hiss->SetName("CYResolS1PtTPCITS");
856 aFolderObj->Add(hiss);
858 hiss = comp->MakeResol(comp->fCZResolS1PtTPC,1,0);
859 hiss->SetXTitle("#sqrt(1/mcp_{t})");
860 hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
862 hiss->SetName("CZResolS1PtTPC");
864 aFolderObj->Add(hiss);
866 hiss = comp->MakeResol(comp->fCZResolS1PtTPCITS,1,0);
867 hiss->SetXTitle("#sqrt(1/mcp_{t})");
868 hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
870 hiss->SetName("CZResolS1PtTPCITS");
872 aFolderObj->Add(hiss);
874 hiss = comp->MakeResol(comp->fCPhiResolS1PtTPC,1,0);
875 hiss->SetXTitle("#sqrt(1/mcp_{t})");
876 hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
878 hiss->SetName("CPhiResolS1PtTPC");
880 aFolderObj->Add(hiss);
882 hiss = comp->MakeResol(comp->fCPhiResolS1PtTPCITS,1,0);
883 hiss->SetXTitle("#sqrt(1/mcp_{t})");
884 hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
886 hiss->SetName("CPhiResolS1PtTPCITS");
888 aFolderObj->Add(hiss);
890 hiss = comp->MakeResol(comp->fCThetaResolS1PtTPC,1,0);
891 hiss->SetXTitle("#sqrt(1/mcp_{t})");
892 hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
894 hiss->SetName("CThetaResolS1PtTPC");
896 aFolderObj->Add(hiss);
898 hiss = comp->MakeResol(comp->fCThetaResolS1PtTPCITS,1,0);
899 hiss->SetXTitle("#sqrt(1/mcp_{t})");
900 hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
902 hiss->SetName("CThetaResolS1PtTPCITS");
904 aFolderObj->Add(hiss);
907 hiss = comp->MakeResol(comp->f1Pt2ResolS1PtTPC,1,0);
908 hiss->SetXTitle("#sqrt(1/mcp_{t})");
909 hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
911 hiss->SetName("OnePt2ResolS1PtTPC");
913 aFolderObj->Add(hiss);
915 hiss = comp->MakeResol(comp->f1Pt2ResolS1PtTPCITS,1,0);
916 hiss->SetXTitle("#sqrt(1/mcp_{t})");
917 hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
919 hiss->SetName("OnePt2ResolS1PtTPCITS");
921 aFolderObj->Add(hiss);
923 hiss = comp->MakeResol(comp->fYResolS1PtTPC,1,0);
924 hiss->SetXTitle("#sqrt(1/mcp_{t})");
925 hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
927 hiss->SetName("YResolS1PtTPC");
929 aFolderObj->Add(hiss);
931 hiss = comp->MakeResol(comp->fYResolS1PtTPCITS,1,0);
932 hiss->SetXTitle("#sqrt(1/mcp_{t})");
933 hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
935 hiss->SetName("YResolS1PtTPCITS");
937 aFolderObj->Add(hiss);
939 hiss = comp->MakeResol(comp->fZResolS1PtTPC,1,0);
940 hiss->SetXTitle("#sqrt(1/mcp_{t})");
941 hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
943 hiss->SetName("ZResolS1PtTPC");
945 aFolderObj->Add(hiss);
947 hiss = comp->MakeResol(comp->fZResolS1PtTPCITS,1,0);
948 hiss->SetXTitle("#sqrt(1/mcp_{t})");
949 hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
951 hiss->SetName("ZResolS1PtTPCITS");
953 aFolderObj->Add(hiss);
955 hiss = comp->MakeResol(comp->fPhiResolS1PtTPC,1,0);
956 hiss->SetXTitle("#sqrt(1/mcp_{t})");
957 hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
959 hiss->SetName("PhiResolS1PtTPC");
961 aFolderObj->Add(hiss);
963 hiss = comp->MakeResol(comp->fPhiResolS1PtTPCITS,1,0);
964 hiss->SetXTitle("#sqrt(1/mcp_{t})");
965 hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
967 hiss->SetName("PhiResolS1PtTPCITS");
969 aFolderObj->Add(hiss);
971 hiss = comp->MakeResol(comp->fThetaResolS1PtTPC,1,0);
972 hiss->SetXTitle("#sqrt(1/mcp_{t})");
973 hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
975 hiss->SetName("ThetaResolS1PtTPC");
977 aFolderObj->Add(hiss);
979 hiss = comp->MakeResol(comp->fThetaResolS1PtTPCITS,1,0);
980 hiss->SetXTitle("#sqrt(1/mcp_{t})");
981 hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
983 hiss->SetName("ThetaResolS1PtTPCITS");
985 aFolderObj->Add(hiss);
987 // export objects to analysis folder
988 fAnalysisFolder = ExportToFolder(aFolderObj);
990 // delete only TObjArray
991 if(aFolderObj) delete aFolderObj;
994 //_____________________________________________________________________________
995 TFolder* AliComparisonRes::ExportToFolder(TObjArray * array)
997 // recreate folder avery time and export objects to new one
999 AliComparisonRes * comp=this;
1000 TFolder *folder = comp->GetAnalysisFolder();
1002 TString name, title;
1003 TFolder *newFolder = 0;
1005 Int_t size = array->GetSize();
1008 // get name and title from old folder
1009 name = folder->GetName();
1010 title = folder->GetTitle();
1016 newFolder = CreateFolder(name.Data(),title.Data());
1017 newFolder->SetOwner();
1019 // add objects to folder
1021 newFolder->Add(array->At(i));
1029 //_____________________________________________________________________________
1030 Long64_t AliComparisonRes::Merge(TCollection* list)
1032 // Merge list of objects (needed by PROOF)
1037 if (list->IsEmpty())
1040 TIterator* iter = list->MakeIterator();
1043 // collection of generated histograms
1045 while((obj = iter->Next()) != 0)
1047 AliComparisonRes* entry = dynamic_cast<AliComparisonRes*>(obj);
1048 if (entry == 0) continue;
1050 fMCVertex->Add(entry->fMCVertex);
1051 fRecVertex->Add(entry->fRecVertex);
1053 fPhiTanPtTPC->Add(entry->fPhiTanPtTPC);
1054 fPhiTanPtTPCITS->Add(entry->fPhiTanPtTPCITS);
1056 fPtResolTPC->Add(entry->fPtResolTPC);
1057 fPtPullTPC->Add(entry->fPtPullTPC);
1058 fPhiResolTanTPC->Add(entry->fPhiResolTanTPC);
1059 fTanResolTanTPC->Add(entry->fTanResolTanTPC);
1060 fPhiPullTanTPC->Add(entry->fPhiPullTanTPC);
1061 fTanPullTanTPC->Add(entry->fTanPullTanTPC);
1063 fPtResolTPCITS->Add(entry->fPtResolTPCITS);
1064 fPtPullTPCITS->Add(entry->fPtPullTPCITS);
1065 fPhiResolTanTPCITS->Add(entry->fPhiResolTanTPCITS);
1066 fTanResolTanTPCITS->Add(entry->fTanResolTanTPCITS);
1067 fPhiPullTanTPCITS->Add(entry->fPhiPullTanTPCITS);
1068 fTanPullTanTPCITS->Add(entry->fTanPullTanTPCITS);
1070 // Histograms for 1/pt parameterisation
1071 f1Pt2ResolS1PtTPC->Add(entry->f1Pt2ResolS1PtTPC);
1072 fYResolS1PtTPC->Add(entry->fYResolS1PtTPC);
1073 fZResolS1PtTPC->Add(entry->fZResolS1PtTPC);
1074 fPhiResolS1PtTPC->Add(entry->fPhiResolS1PtTPC);
1075 fThetaResolS1PtTPC->Add(entry->fThetaResolS1PtTPC);
1077 f1Pt2ResolS1PtTPCITS->Add(entry->f1Pt2ResolS1PtTPCITS);
1078 fYResolS1PtTPCITS->Add(entry->fYResolS1PtTPCITS);
1079 fZResolS1PtTPCITS->Add(entry->fZResolS1PtTPCITS);
1080 fPhiResolS1PtTPCITS->Add(entry->fPhiResolS1PtTPCITS);
1081 fThetaResolS1PtTPCITS->Add(entry->fThetaResolS1PtTPCITS);
1083 // Resolution histograms (constrained param)
1084 fCPhiResolTan->Add(entry->fCPhiResolTan);
1085 fCTanResolTan->Add(entry->fCTanResolTan);
1086 fCPtResolTan->Add(entry->fCPtResolTan);
1087 fCPhiPullTan->Add(entry->fCPhiPullTan);
1088 fCTanPullTan->Add(entry->fCTanPullTan);
1089 fCPtPullTan->Add(entry->fCPtPullTan);
1091 // Histograms for 1/pt parameterisation (constrained)
1092 fC1Pt2ResolS1PtTPC->Add(entry->fC1Pt2ResolS1PtTPC);
1093 fCYResolS1PtTPC->Add(entry->fCYResolS1PtTPC);
1094 fCZResolS1PtTPC->Add(entry->fCZResolS1PtTPC);
1095 fCPhiResolS1PtTPC->Add(entry->fCPhiResolS1PtTPC);
1096 fCThetaResolS1PtTPC->Add(entry->fCThetaResolS1PtTPC);
1098 fC1Pt2ResolS1PtTPCITS->Add(entry->fC1Pt2ResolS1PtTPCITS);
1099 fCYResolS1PtTPCITS->Add(entry->fCYResolS1PtTPCITS);
1100 fCZResolS1PtTPCITS->Add(entry->fCZResolS1PtTPCITS);
1101 fCPhiResolS1PtTPCITS->Add(entry->fCPhiResolS1PtTPCITS);
1102 fCThetaResolS1PtTPCITS->Add(entry->fCThetaResolS1PtTPCITS);
1110 //_____________________________________________________________________________
1111 TFolder* AliComparisonRes::CreateFolder(TString name,TString title) {
1112 // create folder for analysed histograms
1114 TFolder *folder = 0;
1115 folder = new TFolder(name.Data(),title.Data());