]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG1/AliComparisonRes.cxx
New output added from analyze function (Jacek)
[u/mrichter/AliRoot.git] / PWG1 / AliComparisonRes.cxx
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.
9 //
10 // Author: J.Otwinowski 04/02/2008 
11 //------------------------------------------------------------------------------
12
13 /*
14  
15   // after running comparison task, read the file, and get component
16   gROOT->LoadMacro("$ALICE_ROOT/PWG1/Macros/LoadMyLibs.C");
17   LoadMyLibs();
18
19   TFile f("Output.root");
20   AliComparisonRes * compObj = (AliComparisonRes*)f.Get("AliComparisonRes");
21
22   // analyse comparison data
23   compObj->Analyse();
24
25   // the output histograms/graphs will be stored in the folder "folderRes" 
26   compObj->GetAnalysisFolder()->ls("*");
27
28   // user can save whole comparison object (or only folder with anlysed histograms) 
29   // in the seperate output file (e.g.)
30   TFile fout("Analysed_Res.root","recreate");
31   compObj->Write(); // compObj->GetAnalysisFolder()->Write();
32   fout.Close();
33
34 */
35
36 #include <iostream>
37
38 #include "TFile.h"
39 #include "TCint.h"
40 #include "TH3F.h"
41 #include "TH2F.h"
42 #include "TF1.h"
43 #include "TProfile.h"
44 #include "TProfile2D.h"
45 #include "TGraph2D.h"
46 #include "TCanvas.h"
47 #include "TGraph.h"
48
49 #include "AliESDEvent.h"   
50 #include "AliESD.h"
51 #include "AliESDfriend.h"
52 #include "AliESDfriendTrack.h"
53 #include "AliESDVertex.h"
54 #include "AliRecInfoCuts.h" 
55 #include "AliMCInfoCuts.h" 
56 #include "AliLog.h" 
57 #include "AliTracker.h" 
58
59 #include "AliMathBase.h"
60 #include "AliTreeDraw.h" 
61
62 #include "AliMCInfo.h" 
63 #include "AliESDRecInfo.h" 
64 #include "AliComparisonRes.h" 
65
66 using namespace std;
67
68 ClassImp(AliComparisonRes)
69
70 //_____________________________________________________________________________
71 AliComparisonRes::AliComparisonRes():
72   AliComparisonObject("AliComparisonRes"),
73
74   // Resolution 
75   fPtResolLPT(0),        // pt resolution - low pt
76   fPtResolHPT(0),        // pt resolution - high pt 
77   fPtPullLPT(0),         // pt resolution - low pt
78   fPtPullHPT(0),         // pt resolution - high pt 
79   fPhiResolTan(0),       //-> angular resolution 
80   fTanResolTan(0),       //-> angular resolution
81   fPhiPullTan(0),        //-> angular resolution
82   fTanPullTan(0),        //-> angular resolution
83
84   //
85   // Resolution constrained param
86   //
87   fCPhiResolTan(0),  // angular resolution -  constrained
88   fCTanResolTan(0),  // angular resolution -  constrained
89   fCPtResolTan(0),   // pt resolution      -  constrained
90   fCPhiPullTan(0),   // angular resolution -  constrained
91   fCTanPullTan(0),   // angular resolution -  constrained
92   fCPtPullTan(0),    // pt resolution      -  constrained
93
94   //
95   // Parametrisation histograms
96   //
97
98   f1Pt2ResolS1PtTPC(0),
99   f1Pt2ResolS1PtTPCITS(0),
100   fYResolS1PtTPC(0),
101   fYResolS1PtTPCITS(0),
102   fZResolS1PtTPC(0),
103   fZResolS1PtTPCITS(0),
104   fPhiResolS1PtTPC(0),
105   fPhiResolS1PtTPCITS(0),
106   fThetaResolS1PtTPC(0),
107   fThetaResolS1PtTPCITS(0),
108
109   // constrained
110   fC1Pt2ResolS1PtTPC(0),
111   fC1Pt2ResolS1PtTPCITS(0),
112   fCYResolS1PtTPC(0),
113   fCYResolS1PtTPCITS(0),
114   fCZResolS1PtTPC(0),
115   fCZResolS1PtTPCITS(0),
116   fCPhiResolS1PtTPC(0),
117   fCPhiResolS1PtTPCITS(0),
118   fCThetaResolS1PtTPC(0),
119   fCThetaResolS1PtTPCITS(0),
120
121   // vertex
122   fVertex(0),
123  
124   // Cuts 
125   fCutsRC(0),  
126   fCutsMC(0),  
127
128   // histogram folder 
129   fAnalysisFolder(0)
130 {
131   Init();
132   
133   // vertex (0,0,0)
134   fVertex = new AliESDVertex();
135   fVertex->SetXv(0.0);
136   fVertex->SetYv(0.0);
137   fVertex->SetZv(0.0);
138 }
139
140 //_____________________________________________________________________________
141 AliComparisonRes::~AliComparisonRes(){
142   
143   // Resolution histograms
144   if(fPtResolLPT) delete  fPtResolLPT; fPtResolLPT=0;     
145   if(fPtResolHPT) delete  fPtResolHPT; fPtResolHPT=0;    
146   if(fPtPullLPT)  delete  fPtPullLPT;  fPtPullLPT=0;    
147   if(fPtPullHPT)  delete  fPtPullHPT;  fPtPullHPT=0;   
148
149   if(fPhiResolTan)  delete  fPhiResolTan; fPhiResolTan=0;   
150   if(fTanResolTan)  delete  fTanResolTan; fTanResolTan=0;   
151   if(fPhiPullTan)  delete  fPhiPullTan; fPhiPullTan=0;   
152   if(fTanPullTan)  delete  fTanPullTan; fTanPullTan=0;   
153
154   // Resolution histograms (constrained param)
155   if(fCPhiResolTan) delete fCPhiResolTan; fCPhiResolTan=0;
156   if(fCTanResolTan) delete fCTanResolTan; fCTanResolTan=0;
157   if(fCPtResolTan)  delete fCPtResolTan;  fCPtResolTan=0; 
158   if(fCPhiPullTan)  delete fCPhiPullTan;  fCPhiPullTan=0;
159   if(fCTanPullTan)  delete fCTanPullTan;  fCTanPullTan=0;
160   if(fCPtPullTan)   delete fCPtPullTan;   fCPtPullTan=0;
161
162   // Parametrisation histograms
163   // 
164   if(f1Pt2ResolS1PtTPC) delete f1Pt2ResolS1PtTPC; f1Pt2ResolS1PtTPC=0;
165   if(f1Pt2ResolS1PtTPCITS) delete f1Pt2ResolS1PtTPCITS; f1Pt2ResolS1PtTPCITS=0;
166   if(fYResolS1PtTPC) delete fYResolS1PtTPC; fYResolS1PtTPC=0;
167   if(fYResolS1PtTPCITS) delete fYResolS1PtTPCITS; fYResolS1PtTPCITS=0;
168   if(fZResolS1PtTPC) delete fZResolS1PtTPC; fZResolS1PtTPC=0;
169   if(fZResolS1PtTPCITS) delete fZResolS1PtTPCITS; fZResolS1PtTPCITS=0;
170   if(fPhiResolS1PtTPC) delete fPhiResolS1PtTPC; fPhiResolS1PtTPC=0;
171   if(fPhiResolS1PtTPCITS) delete fPhiResolS1PtTPCITS; fPhiResolS1PtTPCITS=0;
172   if(fThetaResolS1PtTPC) delete fThetaResolS1PtTPC; fThetaResolS1PtTPC=0;
173   if(fThetaResolS1PtTPCITS) delete fThetaResolS1PtTPCITS; fThetaResolS1PtTPCITS=0;
174
175   // constrained
176   if(fC1Pt2ResolS1PtTPC) delete fC1Pt2ResolS1PtTPC; fC1Pt2ResolS1PtTPC=0;
177   if(fC1Pt2ResolS1PtTPCITS) delete fC1Pt2ResolS1PtTPCITS; fC1Pt2ResolS1PtTPCITS=0;
178   if(fCYResolS1PtTPC) delete fCYResolS1PtTPC; fCYResolS1PtTPC=0;
179   if(fCYResolS1PtTPCITS) delete fCYResolS1PtTPCITS; fCYResolS1PtTPCITS=0;
180   if(fCZResolS1PtTPC) delete fCZResolS1PtTPC; fCZResolS1PtTPC=0;
181   if(fCZResolS1PtTPCITS) delete fCZResolS1PtTPCITS; fCZResolS1PtTPCITS=0;
182   if(fCPhiResolS1PtTPC) delete fCPhiResolS1PtTPC; fCPhiResolS1PtTPC=0;
183   if(fCPhiResolS1PtTPCITS) delete fCPhiResolS1PtTPCITS; fCPhiResolS1PtTPCITS=0;
184   if(fCThetaResolS1PtTPC) delete fCThetaResolS1PtTPC; fCThetaResolS1PtTPC=0;
185   if(fCThetaResolS1PtTPCITS) delete fCThetaResolS1PtTPCITS; fCThetaResolS1PtTPCITS=0;
186
187   if(fVertex) delete fVertex; fVertex=0;
188
189   if(fAnalysisFolder) delete fAnalysisFolder; fAnalysisFolder=0;
190 }
191
192 //_____________________________________________________________________________
193 void AliComparisonRes::Init(){
194
195   // Init histograms
196   fCPhiResolTan = new TH2F("CPhiResolTan","CPhiResolTan",50, -2,2,200,-0.025,0.025);   
197   fCPhiResolTan->SetXTitle("tan(#theta)");
198   fCPhiResolTan->SetYTitle("#Delta#phi");
199
200   fCTanResolTan = new TH2F("CTanResolTan","CTanResolTan",50, -2,2,200,-0.025,0.025);
201   fCTanResolTan->SetXTitle("tan(#theta)");
202   fCTanResolTan->SetYTitle("#Delta#theta");
203
204   fCPtResolTan=new TH2F("CPtResol","CPtResol",50, -2,2,200,-0.2,0.2);    
205   fCPtResolTan->SetXTitle("Tan(#theta)");
206   fCPtResolTan->SetYTitle("#Deltap_{t}/p_{t}");
207
208   fCPhiPullTan = new TH2F("CPhiPullTan","CPhiPullTan",50, -2,2,200,-5,5);   
209   fCPhiPullTan->SetXTitle("Tan(#theta)");
210   fCPhiPullTan->SetYTitle("#Delta#phi/#Sigma");
211
212   fCTanPullTan = new TH2F("CTanPullTan","CTanPullTan",50, -2,2,200,-5,5);
213   fCTanPullTan->SetXTitle("Tan(#theta)");
214   fCTanPullTan->SetYTitle("#Delta#theta/#Sigma");
215
216   fCPtPullTan=new TH2F("CPtPull","CPtPull",50, -2,2,200,-5,5);    
217   fCPtPullTan->SetXTitle("Tan(#theta)");
218   fCPtPullTan->SetYTitle("(1/mcp_{t}-1/p_{t})/#Sigma");
219
220   fPtResolLPT = new TH2F("Pt_resol_lpt","pt resol",10, 0.1,3,200,-0.2,0.2);
221   fPtResolLPT->SetXTitle("p_{t}");
222   fPtResolLPT->SetYTitle("#Deltap_{t}/p_{t}");
223
224   fPtResolHPT = new TH2F("Pt_resol_hpt","pt resol",10, 2,100,200,-0.3,0.3);  
225   fPtResolHPT->SetXTitle("p_{t}");
226   fPtResolHPT->SetYTitle("#Deltap_{t}/p_{t}");
227
228   fPtPullLPT = new TH2F("Pt_pull_lpt","pt pull",10, 0.1,3,200,-6,6);
229   fPtPullLPT->SetXTitle("p_{t}");
230   fPtPullLPT->SetYTitle("#Deltap_{t}/#Sigma");
231
232   fPtPullHPT = new TH2F("Pt_pull_hpt","pt pull",10,2,100,200,-6,6);  
233   fPtPullHPT->SetXTitle("p_{t}");
234   fPtPullHPT->SetYTitle("#Deltap_{t}/#Sigma");
235   
236   fPhiResolTan = new TH2F("PhiResolTan","PhiResolTan",50, -2,2,200,-0.025,0.025);   
237   fPhiResolTan->SetXTitle("tan(#theta)");
238   fPhiResolTan->SetYTitle("#Delta#phi");
239
240   fTanResolTan = new TH2F("TanResolTan","TanResolTan",50, -2,2,200,-0.025,0.025);
241   fTanResolTan->SetXTitle("tan(#theta)");
242   fTanResolTan->SetYTitle("#Delta#theta");
243
244   fPhiPullTan = new TH2F("PhiPullTan","PhiPullTan",50, -2,2,200,-5,5);   
245   fPhiPullTan->SetXTitle("Tan(#theta)");
246   fPhiPullTan->SetYTitle("#Delta#phi/#Sigma");
247
248   fTanPullTan = new TH2F("TanPullTan","TanPullTan",50, -2,2,200,-5,5);
249   fTanPullTan->SetXTitle("Tan(#theta)");
250   fTanPullTan->SetYTitle("#Delta#theta/#Sigma");
251
252   //
253   // Parametrisation histograms
254   // 
255
256   f1Pt2ResolS1PtTPC = new TH2F("f1Pt2ResolS1PtTPC","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs sqrt(1/pt))",100,0,3,200,-0.010,0.010);  
257   f1Pt2ResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
258   f1Pt2ResolS1PtTPC->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
259
260   f1Pt2ResolS1PtTPCITS = new TH2F("f1Pt2ResolS1PtTPCITS","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs sqrt(1/pt))",100,0,3,200,-0.010,0.010);  
261   f1Pt2ResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
262   f1Pt2ResolS1PtTPCITS->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
263
264   fYResolS1PtTPC = new TH2F("fYResolS1PtTPC","fYResolS1PtTPC",100, 0,3,200,-1.0,1.0);   
265   fYResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
266   fYResolS1PtTPC->SetYTitle("#DeltaY");
267
268   fYResolS1PtTPCITS = new TH2F("fYResolS1PtTPCITS","fYResolS1PtTPCITS",100, 0,3,200,-0.05,0.05);   
269   fYResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
270   fYResolS1PtTPCITS->SetYTitle("#DeltaY");
271
272   fZResolS1PtTPC = new TH2F("fZResolS1PtTPC","fZResolS1PtTPC",100, 0,3,200,-1.0,1.0);   
273   fZResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
274   fZResolS1PtTPC->SetYTitle("#DeltaZ");
275
276   fZResolS1PtTPCITS = new TH2F("fZResolS1PtTPCITS","fZResolS1PtTPCITS",100, 0,3,200,-0.05,0.05);   
277   fZResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
278   fZResolS1PtTPCITS->SetYTitle("#DeltaZ");
279
280   fPhiResolS1PtTPC = new TH2F("fPhiResolS1PtTPC","fPhiResolS1PtTPC",100, 0,3,200,-0.025,0.025);   
281   fPhiResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
282   fPhiResolS1PtTPC->SetYTitle("#Delta#phi");
283
284   fPhiResolS1PtTPCITS = new TH2F("fPhiResolS1PtTPCITS","fPhiResolS1PtTPCITS",100, 0,3,200,-0.01,0.01);   
285   fPhiResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
286   fPhiResolS1PtTPCITS->SetYTitle("#Delta#phi");
287
288   fThetaResolS1PtTPC = new TH2F("fThetaResolS1PtTPC","fThetaResolS1PtTPC",100, 0,3,200,-0.025,0.025);   
289   fThetaResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
290   fThetaResolS1PtTPC->SetYTitle("#Delta#theta");
291
292   fThetaResolS1PtTPCITS = new TH2F("fThetaResolS1PtTPCITS","fThetaResolS1PtTPCITS",100, 0,3,200,-0.01,0.01);   
293   fThetaResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
294   fThetaResolS1PtTPCITS->SetYTitle("#Delta#theta");
295   
296   // constrained
297   fC1Pt2ResolS1PtTPC = new TH2F("fC1Pt2ResolS1PtTPC","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs 1/pt)",100,0,3,200,-0.010,0.010);  
298   fC1Pt2ResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
299   fC1Pt2ResolS1PtTPC->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
300
301   fC1Pt2ResolS1PtTPCITS = new TH2F("fC1Pt2ResolS1PtTPCITS","(1/mcpt-1/pt)/(1+1/mcpt)^2 vs 1/pt)",100,0,3,200,-0.010,0.010);  
302   fC1Pt2ResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
303   fC1Pt2ResolS1PtTPCITS->SetYTitle("(1/mcp_{t}-1/p_{t})/(1+1/mcp_{t})^2)");
304
305   fCYResolS1PtTPC = new TH2F("fCYResolS1PtTPC","fCYResolS1PtTPC",100, 0,3,200,-1.0,1.0);   
306   fCYResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
307   fCYResolS1PtTPC->SetYTitle("#DeltaY");
308
309   fCYResolS1PtTPCITS = new TH2F("fCYResolS1PtTPCITS","fCYResolS1PtTPCITS",100, 0,3,200,-0.01,0.01);   
310   fCYResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
311   fCYResolS1PtTPCITS->SetYTitle("#DeltaY");
312
313   fCZResolS1PtTPC = new TH2F("fCZResolS1PtTPC","fCZResolS1PtTPC",100, 0,3,200,-1.0,1.0);   
314   fCZResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
315   fCZResolS1PtTPC->SetYTitle("#DeltaZ");
316
317   fCZResolS1PtTPCITS = new TH2F("fCZResolS1PtTPCITS","fCZResolS1PtTPCITS",100, 0,3,200,-0.025,0.025);   
318   fCZResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
319   fCZResolS1PtTPCITS->SetYTitle("#DeltaZ");
320
321   fCPhiResolS1PtTPC = new TH2F("fCPhiResolS1PtTPC","fCPhiResolS1PtTPC",100, 0,3,200,-0.025,0.025);   
322   fCPhiResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
323   fCPhiResolS1PtTPC->SetYTitle("#Delta#phi");
324
325   fCPhiResolS1PtTPCITS = new TH2F("fCPhiResolS1PtTPCITS","fCPhiResolS1PtTPCITS",100, 0,3,200,-0.003,0.003);   
326   fCPhiResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
327   fCPhiResolS1PtTPCITS->SetYTitle("#Delta#phi");
328
329   fCThetaResolS1PtTPC = new TH2F("fCThetaResolS1PtTPC","fCThetaResolS1PtTPC",100, 0,3,200,-0.025,0.025);   
330   fCThetaResolS1PtTPC->SetXTitle("#sqrt{1/mcp_{t}}");
331   fCThetaResolS1PtTPC->SetYTitle("#Delta#theta");
332
333   fCThetaResolS1PtTPCITS = new TH2F("fCThetaResolS1PtTPCITS","fCThetaResolS1PtTPCITS",100, 0,3,200,-0.005,0.005);   
334   fCThetaResolS1PtTPCITS->SetXTitle("#sqrt{1/mcp_{t}}");
335   fCThetaResolS1PtTPCITS->SetYTitle("#Delta#theta");
336
337   // Init cuts 
338   if(!fCutsMC) 
339     AliDebug(AliLog::kError, "ERROR: Cannot find AliMCInfoCuts object");
340   if(!fCutsRC) 
341     AliDebug(AliLog::kError, "ERROR: Cannot find AliRecInfoCuts object");
342
343   // init folder
344   fAnalysisFolder = CreateFolder("folderRes","Analysis Resolution Folder");
345 }
346
347 //_____________________________________________________________________________
348 void AliComparisonRes::Process(AliMCInfo* infoMC, AliESDRecInfo *infoRC)
349 {
350   // Fill resolution comparison information 
351   AliExternalTrackParam *track = 0;
352   Double_t field      = AliTracker::GetBz(); // nominal Bz field [kG]
353   Double_t kMaxD      = 123456.0; // max distance
354
355   Int_t clusterITS[200];
356   Double_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
357
358   Float_t deltaPt, pullPt, deltaPhi,pullPhi, deltaTan, pullTan, delta1Pt2, deltaY1Pt, deltaZ1Pt, deltaPhi1Pt, deltaTheta1Pt; 
359   Float_t deltaPtTPC, pullPtTPC, deltaPhiTPC, deltaTanTPC, delta1Pt2TPC, deltaY1PtTPC, deltaZ1PtTPC, deltaPhi1PtTPC, deltaTheta1PtTPC; 
360
361   Float_t mcpt = infoMC->GetParticle().Pt();
362   Float_t s1mcpt = TMath::Sqrt(1./infoMC->GetParticle().Pt());
363   Float_t tantheta = TMath::Tan(infoMC->GetParticle().Theta()-TMath::Pi()*0.5);
364
365   // distance to Prim. vertex 
366   const Double_t* dv = infoMC->GetVDist(); 
367   Bool_t isPrim = TMath::Sqrt(dv[0]*dv[0] + dv[1]*dv[1])<fCutsMC->GetMaxR() && TMath::Abs(dv[2])<fCutsMC->GetMaxVz();
368
369   // Check selection cuts
370   if (fCutsMC->IsPdgParticle(TMath::Abs(infoMC->GetParticle().GetPdgCode())) == kFALSE) return; 
371   if (!isPrim) return;
372   if (infoRC->GetStatus(1)!=3) return; // TPC refit
373   if (!infoRC->GetESDtrack()) return;  
374   if (infoRC->GetESDtrack()->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return;
375
376   // calculate and set prim. vertex
377   fVertex->SetXv( infoMC->GetParticle().Vx() - dv[0] );
378   fVertex->SetYv( infoMC->GetParticle().Vy() - dv[1] );
379   fVertex->SetZv( infoMC->GetParticle().Vz() - dv[2] );
380
381   deltaPt= (mcpt-infoRC->GetESDtrack()->Pt())/mcpt;  
382   pullPt= (1/mcpt-infoRC->GetESDtrack()->OneOverPt())/TMath::Sqrt(infoRC->GetESDtrack()->GetSigma1Pt2());  
383   deltaPhi = TMath::ATan2(infoRC->GetESDtrack()->Py(),infoRC->GetESDtrack()->Px())-
384                      TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
385   pullPhi = deltaPhi/TMath::Sqrt(infoRC->GetESDtrack()->GetSigmaSnp2()); 
386
387   deltaTan = TMath::ATan2(infoRC->GetESDtrack()->Pz(),infoRC->GetESDtrack()->Pt())-
388                      TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
389   pullTan = deltaPhi/TMath::Sqrt(infoRC->GetESDtrack()->GetSigmaSnp2());
390
391   delta1Pt2 = (1/mcpt-infoRC->GetESDtrack()->OneOverPt())/TMath::Power(1+1/mcpt,2);       
392   deltaY1Pt = (infoMC->GetParticle().Vy()-infoRC->GetESDtrack()->GetY()) / (0.2+1/mcpt);
393   deltaZ1Pt = (infoMC->GetParticle().Vz()-infoRC->GetESDtrack()->GetZ()) / (0.2+1/mcpt);
394   deltaPhi1Pt = deltaPhi   / (0.1+1/mcpt);
395   deltaTheta1Pt = deltaTan / (0.1+1/mcpt);
396
397   // calculate track parameters at vertex
398   const AliExternalTrackParam *innerTPC =  0;
399   if ((innerTPC = infoRC->GetESDtrack()->GetTPCInnerParam()) != 0)
400   {
401     if ((track = new AliExternalTrackParam(*infoRC->GetESDtrack()->GetTPCInnerParam())) != 0 )
402     {
403       Bool_t bDCAStatus = track->PropagateToDCA(fVertex,field,kMaxD,dca,cov);
404
405       // Fill parametrisation histograms (only TPC track)
406       if(bDCAStatus) 
407           {
408                         deltaPtTPC= (mcpt-innerTPC->Pt())/mcpt;  
409                         pullPtTPC= (1/mcpt-innerTPC->OneOverPt())/TMath::Sqrt(innerTPC->GetSigma1Pt2());  
410                         deltaPhiTPC = TMath::ATan2(innerTPC->Py(),innerTPC->Px())-
411                                                                 TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
412
413                         deltaTanTPC = TMath::ATan2(innerTPC->Pz(),innerTPC->Pt())-
414                                                                 TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
415
416                         delta1Pt2TPC = (1/mcpt-innerTPC->OneOverPt())/TMath::Power(1+1/mcpt,2);       
417                         deltaY1PtTPC= (infoMC->GetParticle().Vy()-innerTPC->GetY()) / (0.2+1/mcpt);
418                         deltaZ1PtTPC = (infoMC->GetParticle().Vz()-innerTPC->GetZ()) / (0.2+1/mcpt);
419                         deltaPhi1PtTPC = deltaPhiTPC   / (0.1+1/mcpt);
420                         deltaTheta1PtTPC = deltaTanTPC / (0.1+1/mcpt);
421
422                         f1Pt2ResolS1PtTPC->Fill(s1mcpt,delta1Pt2TPC);
423                         fYResolS1PtTPC->Fill(s1mcpt,deltaY1PtTPC);
424                         fZResolS1PtTPC->Fill(s1mcpt,deltaZ1PtTPC);
425                         fPhiResolS1PtTPC->Fill(s1mcpt,deltaPhi1PtTPC);
426                         fThetaResolS1PtTPC->Fill(s1mcpt,deltaTheta1PtTPC);
427           }
428           delete track;
429     }
430   }
431
432   // TPC and ITS (nb. of clusters >2) in the system
433   if(infoRC->GetESDtrack()->GetITSclusters(clusterITS)>2) 
434   {
435       f1Pt2ResolS1PtTPCITS->Fill(s1mcpt,delta1Pt2);
436       fYResolS1PtTPCITS->Fill(s1mcpt,deltaY1Pt);
437       fZResolS1PtTPCITS->Fill(s1mcpt,deltaZ1Pt);
438       fPhiResolS1PtTPCITS->Fill(s1mcpt,deltaPhi1Pt);
439       fThetaResolS1PtTPCITS->Fill(s1mcpt,deltaTheta1Pt);
440   }
441
442   // Fill histograms
443   fPtResolLPT->Fill(mcpt,deltaPt);
444   fPtResolHPT->Fill(mcpt,deltaPt);
445   fPtPullLPT->Fill(mcpt,pullPt);
446   fPtPullHPT->Fill(mcpt,pullPt);  
447
448   fPhiResolTan->Fill(tantheta,deltaPhi);
449   fPhiPullTan->Fill(tantheta,pullPhi);
450   fTanResolTan->Fill(tantheta,deltaTan);
451   fTanPullTan->Fill(tantheta,pullTan);
452
453 }
454
455 //_____________________________________________________________________________
456 void AliComparisonRes::ProcessConstrained(AliMCInfo* infoMC, AliESDRecInfo *infoRC)
457 {
458   // Fill resolution comparison information (constarained parameters) 
459   //
460   AliExternalTrackParam *track = 0;
461   Double_t field      = AliTracker::GetBz(); // nominal Bz field [kG]
462   Double_t kMaxD      = 123456.0; // max distance
463
464   Int_t clusterITS[200];
465   Double_t dca[2], cov[3]; // dca_xy, dca_z, sigma_xy, sigma_xy_z, sigma_z
466  
467   Float_t deltaPt, pullPt, deltaPhi, pullPhi, deltaTan, pullTan, delta1Pt2, deltaY1Pt, deltaZ1Pt, deltaPhi1Pt, deltaTheta1Pt; 
468   Float_t deltaPtTPC, pullPtTPC, deltaPhiTPC, deltaTanTPC, delta1Pt2TPC, deltaY1PtTPC, deltaZ1PtTPC, deltaPhi1PtTPC, deltaTheta1PtTPC; 
469
470   Float_t mcpt = infoMC->GetParticle().Pt();
471   Float_t s1mcpt = TMath::Sqrt(1./infoMC->GetParticle().Pt());
472   Float_t tantheta = TMath::Tan(infoMC->GetParticle().Theta()-TMath::Pi()*0.5);
473
474   // distance to Prim. vertex 
475   const Double_t* dv = infoMC->GetVDist(); 
476   Bool_t isPrim = TMath::Sqrt(dv[0]*dv[0] + dv[1]*dv[1])<fCutsMC->GetMaxR() && TMath::Abs(dv[2])<fCutsMC->GetMaxVz();
477   
478   // Check selection cuts
479   if (fCutsMC->IsPdgParticle(TMath::Abs(infoMC->GetParticle().GetPdgCode())) == kFALSE) return; 
480   if (!isPrim) return;
481   if (infoRC->GetStatus(1)!=3) return;
482   if (!infoRC->GetESDtrack()) return;  
483   if (infoRC->GetESDtrack()->GetTPCNcls()<fCutsRC->GetMinNClustersTPC()) return;
484   if (!infoRC->GetESDtrack()->GetConstrainedParam()) return;
485
486 // calculate and set prim. vertex
487   fVertex->SetXv( infoMC->GetParticle().Vx() - dv[0] );
488   fVertex->SetYv( infoMC->GetParticle().Vy() - dv[1] );
489   fVertex->SetZv( infoMC->GetParticle().Vz() - dv[2] );
490
491   // constrained parameters resolution
492   const AliExternalTrackParam * cparam = infoRC->GetESDtrack()->GetConstrainedParam();
493   deltaPt= (mcpt-cparam->Pt())/mcpt;  
494   pullPt= (1/mcpt-cparam->OneOverPt())/TMath::Sqrt(cparam->GetSigma1Pt2());          
495   deltaPhi = TMath::ATan2(cparam->Py(),cparam->Px())-
496                      TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
497   pullPhi = deltaPhi/TMath::Sqrt(cparam->GetSigmaSnp2()); 
498   deltaTan = TMath::ATan2(cparam->Pz(),cparam->Pt())-TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
499   pullTan = deltaPhi/TMath::Sqrt(cparam->GetSigmaSnp2()); 
500
501
502   delta1Pt2 = (1/mcpt-cparam->OneOverPt())/TMath::Power(1+1/mcpt,2);       
503
504   deltaY1Pt = (infoMC->GetParticle().Vy()-cparam->GetY()) / (0.2+1/mcpt);
505   deltaZ1Pt = (infoMC->GetParticle().Vz()-cparam->GetZ()) / (0.2+1/mcpt);
506   deltaPhi1Pt = deltaPhi   / (0.1+1/mcpt);
507   deltaTheta1Pt = deltaTan / (0.1+1/mcpt);
508
509   // calculate track parameters at vertex
510   const AliExternalTrackParam *innerTPC =  0;
511   if ((innerTPC = infoRC->GetESDtrack()->GetTPCInnerParam()) != 0)
512   {
513     if ((track = new AliExternalTrackParam(*infoRC->GetESDtrack()->GetTPCInnerParam())) != 0 )
514     {
515       Bool_t bDCAStatus = track->PropagateToDCA(fVertex,field,kMaxD,dca,cov);
516
517       // Fill parametrisation histograms (only TPC track)
518       if(bDCAStatus) 
519           {
520                   deltaPtTPC= (mcpt-innerTPC->Pt())/mcpt;  
521                   pullPtTPC= (1/mcpt-innerTPC->OneOverPt())/TMath::Sqrt(innerTPC->GetSigma1Pt2());  
522                   deltaPhiTPC = TMath::ATan2(innerTPC->Py(),innerTPC->Px())-
523                                                                 TMath::ATan2(infoMC->GetParticle().Py(),infoMC->GetParticle().Px());
524
525                   deltaTanTPC = TMath::ATan2(innerTPC->Pz(),innerTPC->Pt())-
526                                                                 TMath::ATan2(infoMC->GetParticle().Pz(),infoMC->GetParticle().Pt());
527
528                   delta1Pt2TPC = (1/mcpt-innerTPC->OneOverPt())/TMath::Power(1+1/mcpt,2);       
529                   deltaY1PtTPC= (infoMC->GetParticle().Vy()-innerTPC->GetY()) / (0.2+1/mcpt);
530                   deltaZ1PtTPC = (infoMC->GetParticle().Vz()-innerTPC->GetZ()) / (0.2+1/mcpt);
531                   deltaPhi1PtTPC = deltaPhiTPC   / (0.1+1/mcpt);
532                   deltaTheta1PtTPC = deltaTanTPC / (0.1+1/mcpt);
533
534           fC1Pt2ResolS1PtTPC->Fill(s1mcpt,delta1Pt2TPC);
535           fCYResolS1PtTPC->Fill(s1mcpt,deltaY1PtTPC);
536           fCZResolS1PtTPC->Fill(s1mcpt,deltaZ1PtTPC);
537           fCPhiResolS1PtTPC->Fill(s1mcpt,deltaPhi1PtTPC);
538           fCThetaResolS1PtTPC->Fill(s1mcpt,deltaTheta1PtTPC);
539           }
540           delete track;
541     }
542   }
543
544  // TPC and ITS (nb. of clusters >2) in the system
545   if(infoRC->GetESDtrack()->GetITSclusters(clusterITS)>2) 
546   {
547       fC1Pt2ResolS1PtTPCITS->Fill(s1mcpt,delta1Pt2);
548       fCYResolS1PtTPCITS->Fill(s1mcpt,deltaY1Pt);
549       fCZResolS1PtTPCITS->Fill(s1mcpt,deltaZ1Pt);
550       fCPhiResolS1PtTPCITS->Fill(s1mcpt,deltaPhi1Pt);
551       fCThetaResolS1PtTPCITS->Fill(s1mcpt,deltaTheta1Pt);
552   }
553
554   // Fill histograms
555   fCPtResolTan->Fill(tantheta,deltaPt);
556   fCPtPullTan->Fill(tantheta,pullPt);
557   fCPhiResolTan->Fill(tantheta,deltaPhi);
558   fCPhiPullTan->Fill(tantheta,pullPhi);
559   fCTanResolTan->Fill(tantheta,deltaTan);
560   fCTanPullTan->Fill(tantheta,pullTan);
561 }
562  
563 //_____________________________________________________________________________
564 void AliComparisonRes::Exec(AliMCInfo* infoMC, AliESDRecInfo *infoRC){
565   
566   // Process comparison information 
567   Process(infoMC,infoRC);
568   ProcessConstrained(infoMC,infoRC);
569 }
570
571 //_____________________________________________________________________________
572 TH1F* AliComparisonRes::MakeResol(TH2F * his, Int_t integ, Bool_t type){
573   // Create resolution histograms
574  
575   TH1F *hisr, *hism;
576   if (!gPad) new TCanvas;
577   hisr = AliTreeDraw::CreateResHistoI(his,&hism,integ);
578   if (type) return hism;
579   else 
580     return hisr;
581 }
582
583 //_____________________________________________________________________________
584 void AliComparisonRes::Analyse(){
585   // Analyse comparison information and store output histograms
586   // in the folder "folderRes"
587   //
588  
589   TH1::AddDirectory(kFALSE);
590
591   AliComparisonRes * comp=this;
592   TH1F *hiss=0;
593   TObjArray *aFolderObj = new TObjArray;
594
595   // write results in the folder 
596
597   TCanvas * c = new TCanvas("Phi resol Tan","Phi resol Tan");
598   c->cd();
599   //
600   hiss = comp->MakeResol(comp->fCPtResolTan,1,0);
601   hiss->SetXTitle("Tan(#theta)");
602   hiss->SetYTitle("#sigmap_{t}/p_{t}");
603   hiss->Draw(); 
604   hiss->SetName("CptResolTan");
605   
606   aFolderObj->Add(hiss);
607
608   //
609   hiss = comp->MakeResol(comp->fPtResolLPT,1,0);
610   hiss->SetXTitle("p_{t} (GeV)");
611   hiss->SetYTitle("#sigmap_{t}/p_{t}");
612   hiss->Draw(); 
613   hiss->SetName("PtResolPt");
614   
615   aFolderObj->Add(hiss);
616
617   //
618   hiss = comp->MakeResol(comp->fPtResolLPT,1,1);
619   hiss->SetXTitle("p_{t} (GeV)");
620   hiss->SetYTitle("mean #Deltap_{t}/p_{t}");
621   hiss->Draw(); 
622   hiss->SetName("PtMeanPt");
623   
624   aFolderObj->Add(hiss);
625
626   //
627   hiss = comp->MakeResol(comp->fPhiResolTan,1,0);
628   hiss->SetXTitle("Tan(#theta)");
629   hiss->SetYTitle("#sigma#phi (rad)");
630   hiss->Draw();
631   hiss->SetName("PhiResolTan");
632   
633   aFolderObj->Add(hiss);
634   //
635   hiss = comp->MakeResol(comp->fTanResolTan,1,0);
636   hiss->SetXTitle("Tan(#theta)");
637   hiss->SetYTitle("#sigma#theta (rad)");
638   hiss->Draw();
639   hiss->SetName("ThetaResolTan");
640   
641   aFolderObj->Add(hiss);
642
643   //
644   hiss = comp->MakeResol(comp->fPhiResolTan,1,1);
645   hiss->SetXTitle("Tan(#theta)");
646   hiss->SetYTitle("mean #Delta#phi (rad)");
647   hiss->Draw();
648   hiss->SetName("PhiMeanTan");
649   
650   aFolderObj->Add(hiss);
651   //
652   hiss = comp->MakeResol(comp->fTanResolTan,1,1);
653   hiss->SetXTitle("Tan(#theta)");
654   hiss->SetYTitle("mean #Delta#theta (rad)");
655   hiss->Draw();
656   hiss->SetName("ThetaMeanTan");
657   
658   aFolderObj->Add(hiss);
659
660   //
661   hiss = comp->MakeResol(comp->fCPhiResolTan,1,0);
662   hiss->SetXTitle("Tan(#theta)");
663   hiss->SetYTitle("#sigma#phi (rad)");
664   hiss->Draw();
665   hiss->SetName("CPhiResolTan");
666   
667   aFolderObj->Add(hiss);
668   //
669   hiss = comp->MakeResol(comp->fCTanResolTan,1,0);
670   hiss->SetXTitle("Tan(#theta)");
671   hiss->SetYTitle("#sigma#theta (rad)");
672   hiss->Draw();
673   hiss->SetName("CThetaResolTan");
674   
675   aFolderObj->Add(hiss);
676   //
677   hiss = comp->MakeResol(comp->fCPtPullTan,1,0);
678   hiss->SetXTitle("Tan(#theta)");
679   hiss->SetYTitle("1/mcp_{t}-1/p_{t}/#Sigma(1/p_{t})");
680   hiss->Draw();
681   hiss->SetName("CptPullTan");
682   
683   aFolderObj->Add(hiss);
684   //
685   hiss = comp->MakeResol(comp->fC1Pt2ResolS1PtTPC,1,0);
686   hiss->SetXTitle("#sqrt(1/mcp_{t})");
687   hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
688   hiss->Draw();
689   hiss->SetName("C1Pt2ResolS1PtTPC");
690   
691   aFolderObj->Add(hiss);
692
693   hiss = comp->MakeResol(comp->fC1Pt2ResolS1PtTPCITS,1,0);
694   hiss->SetXTitle("#sqrt(1/mcp_{t})");
695   hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
696   hiss->Draw();
697   hiss->SetName("C1Pt2ResolS1PtTPCITS");
698   
699   aFolderObj->Add(hiss);
700   //
701   hiss = comp->MakeResol(comp->fCYResolS1PtTPC,1,0);
702   hiss->SetXTitle("#sqrt(1/mcp_{t})");
703   hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
704   hiss->Draw();
705   hiss->SetName("CYResolS1PtTPC");
706   
707   aFolderObj->Add(hiss);
708
709   hiss = comp->MakeResol(comp->fCYResolS1PtTPCITS,1,0);
710   hiss->SetXTitle("#sqrt(1/mcp_{t})");
711   hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
712   hiss->Draw();
713   hiss->SetName("CYResolS1PtTPCITS");
714   
715   aFolderObj->Add(hiss);
716   //
717   hiss = comp->MakeResol(comp->fCZResolS1PtTPC,1,0);
718   hiss->SetXTitle("#sqrt(1/mcp_{t})");
719   hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
720   hiss->Draw();
721   hiss->SetName("CZResolS1PtTPC");
722   
723   aFolderObj->Add(hiss);
724
725   hiss = comp->MakeResol(comp->fCZResolS1PtTPCITS,1,0);
726   hiss->SetXTitle("#sqrt(1/mcp_{t})");
727   hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
728   hiss->Draw();
729   hiss->SetName("CZResolS1PtTPCITS");
730   
731   aFolderObj->Add(hiss);
732   //
733   hiss = comp->MakeResol(comp->fCPhiResolS1PtTPC,1,0);
734   hiss->SetXTitle("#sqrt(1/mcp_{t})");
735   hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
736   hiss->Draw();
737   hiss->SetName("CPhiResolS1PtTPC");
738   
739   aFolderObj->Add(hiss);
740
741   hiss = comp->MakeResol(comp->fCPhiResolS1PtTPCITS,1,0);
742   hiss->SetXTitle("#sqrt(1/mcp_{t})");
743   hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
744   hiss->Draw();
745   hiss->SetName("CPhiResolS1PtTPCITS");
746   
747   aFolderObj->Add(hiss);
748   //
749   hiss = comp->MakeResol(comp->fCThetaResolS1PtTPC,1,0);
750   hiss->SetXTitle("#sqrt(1/mcp_{t})");
751   hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
752   hiss->Draw();
753   hiss->SetName("CThetaResolS1PtTPC");
754   
755   aFolderObj->Add(hiss);
756
757   hiss = comp->MakeResol(comp->fCThetaResolS1PtTPCITS,1,0);
758   hiss->SetXTitle("#sqrt(1/mcp_{t})");
759   hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
760   hiss->Draw();
761   hiss->SetName("CThetaResolS1PtTPCITS");
762   
763   aFolderObj->Add(hiss);
764
765   //
766   hiss = comp->MakeResol(comp->f1Pt2ResolS1PtTPC,1,0);
767   hiss->SetXTitle("#sqrt(1/mcp_{t})");
768   hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
769   hiss->Draw();
770   hiss->SetName("OnePt2ResolS1PtTPC");
771   
772   aFolderObj->Add(hiss);
773
774   hiss = comp->MakeResol(comp->f1Pt2ResolS1PtTPCITS,1,0);
775   hiss->SetXTitle("#sqrt(1/mcp_{t})");
776   hiss->SetYTitle("1/mcp_{t}-1/p_{t}/(1+1/p_{t})^2");
777   hiss->Draw();
778   hiss->SetName("OnePt2ResolS1PtTPCITS");
779   
780   aFolderObj->Add(hiss);
781   //
782   hiss = comp->MakeResol(comp->fYResolS1PtTPC,1,0);
783   hiss->SetXTitle("#sqrt(1/mcp_{t})");
784   hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
785   hiss->Draw();
786   hiss->SetName("YResolS1PtTPC");
787   
788   aFolderObj->Add(hiss);
789
790   hiss = comp->MakeResol(comp->fYResolS1PtTPCITS,1,0);
791   hiss->SetXTitle("#sqrt(1/mcp_{t})");
792   hiss->SetYTitle("(mcy-y)/(0.2+1/mcp_{t})");
793   hiss->Draw();
794   hiss->SetName("YResolS1PtTPCITS");
795   
796   aFolderObj->Add(hiss);
797   //
798   hiss = comp->MakeResol(comp->fZResolS1PtTPC,1,0);
799   hiss->SetXTitle("#sqrt(1/mcp_{t})");
800   hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
801   hiss->Draw();
802   hiss->SetName("ZResolS1PtTPC");
803   
804   aFolderObj->Add(hiss);
805
806   hiss = comp->MakeResol(comp->fZResolS1PtTPCITS,1,0);
807   hiss->SetXTitle("#sqrt(1/mcp_{t})");
808   hiss->SetYTitle("(mcz-z)/(0.2+1/mcp_{t})");
809   hiss->Draw();
810   hiss->SetName("ZResolS1PtTPCITS");
811   
812   aFolderObj->Add(hiss);
813   //
814   hiss = comp->MakeResol(comp->fPhiResolS1PtTPC,1,0);
815   hiss->SetXTitle("#sqrt(1/mcp_{t})");
816   hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
817   hiss->Draw();
818   hiss->SetName("PhiResolS1PtTPC");
819   
820   aFolderObj->Add(hiss);
821
822   hiss = comp->MakeResol(comp->fPhiResolS1PtTPCITS,1,0);
823   hiss->SetXTitle("#sqrt(1/mcp_{t})");
824   hiss->SetYTitle("(mc#phi-#phi)/(0.1+1/mcp_{t})");
825   hiss->Draw();
826   hiss->SetName("PhiResolS1PtTPCITS");
827   
828   aFolderObj->Add(hiss);
829   //
830   hiss = comp->MakeResol(comp->fThetaResolS1PtTPC,1,0);
831   hiss->SetXTitle("#sqrt(1/mcp_{t})");
832   hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
833   hiss->Draw();
834   hiss->SetName("ThetaResolS1PtTPC");
835   
836   aFolderObj->Add(hiss);
837
838   hiss = comp->MakeResol(comp->fThetaResolS1PtTPCITS,1,0);
839   hiss->SetXTitle("#sqrt(1/mcp_{t})");
840   hiss->SetYTitle("(mc#theta-#theta)/(0.1+1/mcp_{t})");
841   hiss->Draw();
842   hiss->SetName("ThetaResolS1PtTPCITS");
843   
844   aFolderObj->Add(hiss);
845
846   // export objects to analysis folder
847   fAnalysisFolder = ExportToFolder(aFolderObj);
848
849   // delete only TObjArray
850   if(aFolderObj) delete aFolderObj;
851 }
852
853 //_____________________________________________________________________________
854 TFolder* AliComparisonRes::ExportToFolder(TObjArray * array) 
855 {
856   // recreate folder avery time and export objects to new one
857   //
858   AliComparisonRes * comp=this;
859   TFolder *folder = comp->GetAnalysisFolder();
860
861   TString name, title;
862   TFolder *newFolder = 0;
863   Int_t i = 0;
864   Int_t size = array->GetSize();
865
866   if(folder) { 
867      // get name and title from old folder
868      name = folder->GetName();  
869      title = folder->GetTitle();  
870
871          // delete old one
872      delete folder;
873
874          // create new one
875      newFolder = CreateFolder(name.Data(),title.Data());
876      newFolder->SetOwner();
877
878          // add objects to folder
879      while(i < size) {
880            newFolder->Add(array->At(i));
881            i++;
882          }
883   }
884
885 return newFolder;
886 }
887
888 //_____________________________________________________________________________
889 Long64_t AliComparisonRes::Merge(TCollection* list) 
890 {
891   // Merge list of objects (needed by PROOF)
892
893   if (!list)
894   return 0;
895
896   if (list->IsEmpty())
897   return 1;
898
899   TIterator* iter = list->MakeIterator();
900   TObject* obj = 0;
901
902   // collection of generated histograms
903   Int_t count=0;
904   while((obj = iter->Next()) != 0) 
905   {
906   AliComparisonRes* entry = dynamic_cast<AliComparisonRes*>(obj);
907   if (entry == 0) continue; 
908
909   fPtResolLPT->Add(entry->fPtResolLPT);
910   fPtResolHPT->Add(entry->fPtResolHPT);
911   fPtPullLPT->Add(entry->fPtPullLPT);
912   fPtPullHPT->Add(entry->fPtPullHPT);
913   fPhiResolTan->Add(entry->fPhiResolTan);
914   fTanResolTan->Add(entry->fTanResolTan);
915   fPhiPullTan->Add(entry->fPhiPullTan);
916   fTanPullTan->Add(entry->fTanPullTan);
917
918   // Histograms for 1/pt parameterisation
919   f1Pt2ResolS1PtTPC->Add(entry->f1Pt2ResolS1PtTPC);
920   fYResolS1PtTPC->Add(entry->fYResolS1PtTPC);
921   fZResolS1PtTPC->Add(entry->fZResolS1PtTPC);
922   fPhiResolS1PtTPC->Add(entry->fPhiResolS1PtTPC);
923   fThetaResolS1PtTPC->Add(entry->fThetaResolS1PtTPC);
924
925   f1Pt2ResolS1PtTPCITS->Add(entry->f1Pt2ResolS1PtTPCITS);
926   fYResolS1PtTPCITS->Add(entry->fYResolS1PtTPCITS);
927   fZResolS1PtTPCITS->Add(entry->fZResolS1PtTPCITS);
928   fPhiResolS1PtTPCITS->Add(entry->fPhiResolS1PtTPCITS);
929   fThetaResolS1PtTPCITS->Add(entry->fThetaResolS1PtTPCITS);
930
931   // Resolution histograms (constrained param)
932   fCPhiResolTan->Add(entry->fCPhiResolTan);
933   fCTanResolTan->Add(entry->fCTanResolTan);
934   fCPtResolTan->Add(entry->fCPtResolTan);
935   fCPhiPullTan->Add(entry->fCPhiPullTan);
936   fCTanPullTan->Add(entry->fCTanPullTan);
937   fCPtPullTan->Add(entry->fCPtPullTan);
938
939   //  Histograms for 1/pt parameterisation (constrained)
940   fC1Pt2ResolS1PtTPC->Add(entry->fC1Pt2ResolS1PtTPC);
941   fCYResolS1PtTPC->Add(entry->fCYResolS1PtTPC);
942   fCZResolS1PtTPC->Add(entry->fCZResolS1PtTPC);
943   fCPhiResolS1PtTPC->Add(entry->fCPhiResolS1PtTPC);
944   fCThetaResolS1PtTPC->Add(entry->fCThetaResolS1PtTPC);
945
946   fC1Pt2ResolS1PtTPCITS->Add(entry->fC1Pt2ResolS1PtTPCITS);
947   fCYResolS1PtTPCITS->Add(entry->fCYResolS1PtTPCITS);
948   fCZResolS1PtTPCITS->Add(entry->fCZResolS1PtTPCITS);
949   fCPhiResolS1PtTPCITS->Add(entry->fCPhiResolS1PtTPCITS);
950   fCThetaResolS1PtTPCITS->Add(entry->fCThetaResolS1PtTPCITS);
951
952   count++;
953   }
954
955 return count;
956 }
957
958 //_____________________________________________________________________________
959 TFolder* AliComparisonRes::CreateFolder(TString name,TString title) { 
960 // create folder for analysed histograms
961 //
962 TFolder *folder = 0;
963   folder = new TFolder(name.Data(),title.Data());
964
965   return folder;
966 }