]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG4/JetTasks/AliHistogramsUE.cxx
Updates to run with deltas (L. Cunqueiro)
[u/mrichter/AliRoot.git] / PWG4 / JetTasks / AliHistogramsUE.cxx
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  *                                                                        *
4  * Author: A.Abrahantes, E.Lopez, S.Vallero                               *
5  * Version 1.0                                                            *
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 /* $Id:$ */
17 ////////////////////////////////////////////////
18 //--------------------------------------------- 
19 // Class  to handle histograms for UE analysis
20 //---------------------------------------------
21 ////////////////////////////////////////////////
22 #include <TROOT.h>
23 #include <TBranch.h>
24 #include <TCanvas.h>
25 #include <TChain.h>
26 #include <TFile.h>
27 #include <TH1F.h>
28 #include <TH1I.h>
29 #include <TH2F.h>
30 #include <TLegend.h>
31 #include <TList.h>
32 #include <TLorentzVector.h>
33 #include <TMath.h>
34 #include <TProfile.h>
35 #include <TRandom.h>
36 #include <TStyle.h>
37 #include <TSystem.h>
38 #include <TTree.h>
39 #include <TVector3.h>
40
41 #include "AliHistogramsUE.h"
42 #include "AliAnalysisTaskUE.h"
43 #include "AliAnalyseUE.h"
44 #include "AliCFContainer.h"
45 #include "AliCFManager.h"
46 #include "AliLog.h"
47
48 ClassImp( AliHistogramsUE)
49
50 //____________________________________________________________________
51 AliHistogramsUE:: AliHistogramsUE():
52 TObject(),
53 fBinsPtInHist(0),
54 fMinJetPtInHist(0.),
55 fMaxJetPtInHist(0.),
56 fTrackEtaCut(0.),
57 fJet1EtaCut(0.),
58 fListOfHistos(0x0),  
59 fhNJets(0x0),
60 fhEleadingPt(0x0),
61 fhMinRegPtDist(0x0),
62 fhRegionMultMin(0x0),
63 fhMinRegAvePt(0x0), 
64 fhMinRegSumPt(0x0),            
65 fhMinRegMaxPtPart(0x0),
66 fhMinRegSumPtvsMult(0x0),
67 fhdNdEtaPhiDist(0x0),        
68 fhdNdEtaPhiDistMC(0x0),        
69 fhFullRegPartPtDistVsEt(0x0), 
70 fhFullRegPartPtDistVsEtMC(0x0), 
71 fhTransRegPartPtDistVsEt(0x0),
72 fhTransRegPartPtDistVsEtMC(0x0),
73 fhRegionSumPtMaxVsEt(0x0),
74 fhRegionMultMax(0x0),         
75 fhRegionMultMaxVsEt(0x0),     
76 fhRegionSumPtMinVsEt(0x0), //fhRegionMultMin(0x0),         
77 fhRegionMultMinVsEt(0x0),
78 fhRegionAveSumPtVsEt(0x0),
79 fhRegionDiffSumPtVsEt(0x0),
80 fhRegionAvePartPtMaxVsEt(0x0),
81 fhRegionAvePartPtMinVsEt(0x0),
82 fhRegionMaxPartPtMaxVsEt(0x0),
83 fhRegForwardMult(0x0),
84 fhRegForwardSumPtvsMult(0x0),
85 fhRegBackwardMult(0x0),
86 fhRegBackwardSumPtvsMult(0x0),
87 fhRegForwardPartPtDistVsEt(0x0),
88 fhRegForwardPartPtDistVsEtMC(0x0),
89 fhRegBackwardPartPtDistVsEt(0x0),
90 fhRegBackwardPartPtDistVsEtMC(0x0),
91 fhRegTransMult(0x0),
92 fhRegTransSumPtVsMult(0x0),
93 fhMinRegSumPtJetPtBin(0x0),
94 fhMaxRegSumPtJetPtBin(0x0),
95 fhVertexMult(0x0),
96 fh1Xsec(0x0),
97 fh1Trials(0x0),
98 fhDCAxy(0x0),
99 fhDCAxyPrimary(0x0)
100 //For Corrections
101 {
102   // Default constructor
103
104 }
105
106 //____________________________________________________________________
107 AliHistogramsUE:: AliHistogramsUE(TList *list):
108 TObject(),
109 fBinsPtInHist(0),
110 fMinJetPtInHist(0.),
111 fMaxJetPtInHist(0.),
112 fTrackEtaCut(0.),
113 fJet1EtaCut(0.),
114 fListOfHistos(0x0),  
115 fhNJets(0x0),
116 fhEleadingPt(0x0),
117 fhMinRegPtDist(0x0),
118 fhRegionMultMin(0x0),
119 fhMinRegAvePt(0x0), 
120 fhMinRegSumPt(0x0),            
121 fhMinRegMaxPtPart(0x0),
122 fhMinRegSumPtvsMult(0x0),
123 fhdNdEtaPhiDist(0x0),        
124 fhdNdEtaPhiDistMC(0x0),        
125 fhFullRegPartPtDistVsEt(0x0), 
126 fhFullRegPartPtDistVsEtMC(0x0), 
127 fhTransRegPartPtDistVsEt(0x0),
128 fhTransRegPartPtDistVsEtMC(0x0),
129 fhRegionSumPtMaxVsEt(0x0),
130 fhRegionMultMax(0x0),         
131 fhRegionMultMaxVsEt(0x0),     
132 fhRegionSumPtMinVsEt(0x0), //fhRegionMultMin(0x0),         
133 fhRegionMultMinVsEt(0x0),
134 fhRegionAveSumPtVsEt(0x0),
135 fhRegionDiffSumPtVsEt(0x0),
136 fhRegionAvePartPtMaxVsEt(0x0),
137 fhRegionAvePartPtMinVsEt(0x0),
138 fhRegionMaxPartPtMaxVsEt(0x0),
139 fhRegForwardMult(0x0),
140 fhRegForwardSumPtvsMult(0x0),
141 fhRegBackwardMult(0x0),
142 fhRegBackwardSumPtvsMult(0x0),
143 fhRegForwardPartPtDistVsEt(0x0),
144 fhRegForwardPartPtDistVsEtMC(0x0),
145 fhRegBackwardPartPtDistVsEt(0x0),
146 fhRegBackwardPartPtDistVsEtMC(0x0),
147 fhRegTransMult(0x0),
148 fhRegTransSumPtVsMult(0x0),
149 fhMinRegSumPtJetPtBin(0x0),
150 fhMaxRegSumPtJetPtBin(0x0),
151 fhVertexMult(0x0),
152 fh1Xsec(0x0),
153 fh1Trials(0x0),
154 fhDCAxy(0x0),
155 fhDCAxyPrimary(0x0)
156 //For Corrections
157 {
158   // Constructor, initialize members from list
159         fhNJets = (TH1F*)list->FindObject("hNJets"); 
160         fhEleadingPt = (TH1F*)list->FindObject("hEleadingPt");
161         fhMinRegPtDist = (TH1F*)list->FindObject("hMinRegPtDist");
162         fhRegionMultMin = (TH1F*)list->FindObject("hRegionMultMin");
163         fhMinRegAvePt = (TH1F*)list->FindObject("hMinRegAvePt");
164         fhMinRegSumPt = (TH1F*)list->FindObject("hMinRegSumPt");
165         fhMinRegMaxPtPart = (TH1F*)list->FindObject("hMinRegMaxPtPart");
166         fhMinRegSumPtvsMult = (TH1F*)list->FindObject("hMinRegSumPtvsMult");
167         fhdNdEtaPhiDist = (TH2F*)list->FindObject("hdNdEtaPhiDist");  
168         fhdNdEtaPhiDistMC = (TH2F*)list->FindObject("hdNdEtaPhiDistMC");  
169         fhFullRegPartPtDistVsEt = (TH2F*)list->FindObject("hFullRegPartPtDistVsEt");
170         fhFullRegPartPtDistVsEtMC = (TH2F*)list->FindObject("hFullRegPartPtDistVsEtMC");
171         fhTransRegPartPtDistVsEt = (TH2F*)list->FindObject("hTransRegPartPtDistVsEt");
172         fhTransRegPartPtDistVsEtMC = (TH2F*)list->FindObject("hTransRegPartPtDistVsEtMC");
173         fhRegionSumPtMaxVsEt = (TH1F*)list->FindObject("hRegionSumPtMaxVsEt"); 
174         fhRegionMultMax = (TH1I*)list->FindObject("hRegionMultMax");
175         fhRegionMultMaxVsEt = (TH1F*)list->FindObject("hRegionMultMaxVsEt");
176         fhRegionSumPtMinVsEt = (TH1F*)list->FindObject("hRegionSumPtMinVsEt");
177         fhRegionMultMinVsEt = (TH1F*)list->FindObject("hRegionMultMinVsEt");
178         fhRegionAveSumPtVsEt = (TH1F*)list->FindObject("hRegionAveSumPtVsEt");
179         fhRegionDiffSumPtVsEt = (TH1F*)list->FindObject("hRegionDiffSumPtVsEt");
180         fhRegionAvePartPtMaxVsEt = (TH1F*)list->FindObject("hRegionAvePartPtMaxVsEt");
181         fhRegionAvePartPtMinVsEt = (TH1F*)list->FindObject("hRegionAvePartPtMinVsEt");
182         fhRegionMaxPartPtMaxVsEt = (TH1F*)list->FindObject("hRegionMaxPartPtMaxVsEt");
183         fhRegForwardMult = (TH2F*)list->FindObject("hRegForwardMult"); 
184         fhRegForwardSumPtvsMult = (TH2F*)list->FindObject("hRegForwardSumPtvsMult");
185         fhRegBackwardMult = (TH2F*)list->FindObject("hRegBackwardMult");
186         fhRegBackwardSumPtvsMult = (TH2F*)list->FindObject("hRegBackwardSumPtvsMult");
187         fhRegForwardPartPtDistVsEt = (TH2F*)list->FindObject("hRegForwardPartPtDistVsEt");
188         fhRegForwardPartPtDistVsEtMC = (TH2F*)list->FindObject("hRegForwardPartPtDistVsEtMC");
189         fhRegBackwardPartPtDistVsEt = (TH2F*)list->FindObject("hRegBackwardPartPtDistVsEt");
190         fhRegBackwardPartPtDistVsEtMC = (TH2F*)list->FindObject("hRegBackwardPartPtDistVsEtMC");
191         fhRegTransMult = (TH2F*)list->FindObject("hRegTransMult");
192         fhRegTransSumPtVsMult = (TH2F*)list->FindObject("hRegTransSumPtVsMult");
193         fhMinRegSumPtJetPtBin = (TH2F*)list->FindObject("hMinRegSumPtJetPtBin");
194         fhMaxRegSumPtJetPtBin = (TH2F*)list->FindObject("hMaxRegSumPtJetPtBin");
195         fhVertexMult = (TH1F*)list->FindObject("hVertexMult");
196         fh1Xsec = (TProfile*)list->FindObject("h1Xsec");
197         fh1Trials = (TH1F*)list->FindObject("h1Trials"); 
198         fhDCAxy = (TH2F*)list->FindObject("hDCAxy");
199         fhDCAxyPrimary = (TH2F*)list->FindObject("hDCAxyPrimary");
200 }
201 //____________________________________________________________________
202 AliHistogramsUE:: AliHistogramsUE(const AliHistogramsUE & original):
203 TObject(),
204 fBinsPtInHist(original.fBinsPtInHist),
205 fMinJetPtInHist(original.fMinJetPtInHist),
206 fMaxJetPtInHist(original.fMaxJetPtInHist),
207 fTrackEtaCut(original.fTrackEtaCut),
208 fJet1EtaCut(original.fJet1EtaCut),
209 fListOfHistos(original.fListOfHistos),  
210 fhNJets(original.fhNJets),
211 fhEleadingPt(original.fhEleadingPt),
212 fhMinRegPtDist(original.fhMinRegPtDist),
213 fhRegionMultMin(original.fhRegionMultMin),
214 fhMinRegAvePt(original.fhMinRegAvePt), 
215 fhMinRegSumPt(original.fhMinRegSumPt),            
216 fhMinRegMaxPtPart(original.fhMinRegMaxPtPart),
217 fhMinRegSumPtvsMult(original.fhMinRegSumPtvsMult),
218 fhdNdEtaPhiDist(original.fhdNdEtaPhiDist),        
219 fhdNdEtaPhiDistMC(original.fhdNdEtaPhiDistMC),        
220 fhFullRegPartPtDistVsEt(original.fhFullRegPartPtDistVsEt), 
221 fhFullRegPartPtDistVsEtMC(original.fhFullRegPartPtDistVsEtMC), 
222 fhTransRegPartPtDistVsEt(original.fhTransRegPartPtDistVsEt),
223 fhTransRegPartPtDistVsEtMC(original.fhTransRegPartPtDistVsEtMC),
224 fhRegionSumPtMaxVsEt(original.fhRegionSumPtMaxVsEt),
225 fhRegionMultMax(original.fhRegionMultMax),         
226 fhRegionMultMaxVsEt(original.fhRegionMultMaxVsEt),     
227 fhRegionSumPtMinVsEt(original.fhRegionSumPtMinVsEt),         
228 fhRegionMultMinVsEt(original.fhRegionMultMinVsEt),
229 fhRegionAveSumPtVsEt(original.fhRegionAveSumPtVsEt),
230 fhRegionDiffSumPtVsEt(original.fhRegionDiffSumPtVsEt),
231 fhRegionAvePartPtMaxVsEt(original.fhRegionAvePartPtMaxVsEt),
232 fhRegionAvePartPtMinVsEt(original.fhRegionAvePartPtMinVsEt),
233 fhRegionMaxPartPtMaxVsEt(original.fhRegionMaxPartPtMaxVsEt),
234 fhRegForwardMult(original.fhRegForwardMult),
235 fhRegForwardSumPtvsMult(original.fhRegForwardSumPtvsMult),
236 fhRegBackwardMult(original.fhRegBackwardMult),
237 fhRegBackwardSumPtvsMult(original.fhRegBackwardSumPtvsMult),
238 fhRegForwardPartPtDistVsEt(original.fhRegForwardPartPtDistVsEt),
239 fhRegForwardPartPtDistVsEtMC(original.fhRegForwardPartPtDistVsEtMC),
240 fhRegBackwardPartPtDistVsEt(original.fhRegBackwardPartPtDistVsEt),
241 fhRegBackwardPartPtDistVsEtMC(original.fhRegBackwardPartPtDistVsEtMC),
242 fhRegTransMult(original.fhRegTransMult),
243 fhRegTransSumPtVsMult(original.fhRegTransSumPtVsMult),
244 fhMinRegSumPtJetPtBin(original.fhMinRegSumPtJetPtBin),
245 fhMaxRegSumPtJetPtBin(original.fhMaxRegSumPtJetPtBin),
246 fhVertexMult(original.fhVertexMult),
247 fh1Xsec(original.fh1Xsec),
248 fh1Trials(original.fh1Trials),
249 fhDCAxy(original.fhDCAxy),
250 fhDCAxyPrimary(original.fhDCAxyPrimary)
251 //For Corrections
252 {
253
254   // Copy constructor
255
256 }
257
258
259 //______________________________________________________________
260 AliHistogramsUE & AliHistogramsUE::operator = (const AliHistogramsUE & /*source*/)
261 {
262
263   // assignment operator
264   return *this;
265
266 }
267
268 //______________________________________________________________
269 TObjArray* AliHistogramsUE::CreateCanvas(const Int_t ncanv){
270         
271         // Create canvas for plotting 
272         printf("Creating %d canvas ... \n",ncanv);
273         TObjArray *arr=new TObjArray;
274         TString name;
275         for(Int_t i=0;i<ncanv;i++ ){
276         name=Form("Canvas %d",i);
277         TCanvas *c=new TCanvas(name,name);
278         c->SetFillColor(0);
279         gStyle->SetOptStat(0);
280         gStyle->SetOptTitle(0);
281         arr->Add(c);
282         }
283
284 return arr;
285
286 }
287
288 //____________________________________________________________________
289 void  AliHistogramsUE::CreateHistograms(TList *list,Int_t bins, Double_t min, Double_t max, Double_t etacut)
290 {
291
292   // Create all histograms necessary for UE analysis 
293   fBinsPtInHist = bins;
294   fMinJetPtInHist = min;
295   fMaxJetPtInHist = max;
296   fTrackEtaCut= etacut;
297
298   //Number of reconstructed clusters  
299   fhNJets = new TH1F("hNJets", "Number of clusters",  20, 0, 20);
300   fhNJets->SetXTitle("Number of reconstructed clusters");
301   fhNJets->SetYTitle("#");
302   fhNJets->Sumw2();
303   list->Add( fhNJets );                 // At(0) 
304   
305   //pT distribution of leading clusters
306   fhEleadingPt  = new TH1F("hEleadingPt",   "Leading cluster p_{T}",  bins, min,   max);
307   fhEleadingPt->SetXTitle("p_{T} of  cluster (GeV/c)");
308   fhEleadingPt->SetYTitle("1/N_{ev} dN/dp_{T} (|#eta|<0.5)");
309   fhEleadingPt->Sumw2();
310   list->Add( fhEleadingPt );            // At(1)
311   
312   //Track pT distribution in MIN zone
313   fhMinRegPtDist = new TH1F("hMinRegPtDist",   "p_{T} distribution in MIN zone",  50,0.,20.);
314   fhMinRegPtDist->SetXTitle("Track p_{T} (GeV/c)");
315   fhMinRegPtDist->SetYTitle("dN/dp_{T}");
316   fhMinRegPtDist->Sumw2();
317   list->Add( fhMinRegPtDist );          // At(2)
318   
319   //Multiplicity in MIN zone
320   fhRegionMultMin = new TH1F("hRegionMultMin",      "N_{ch}^{90, min}",  21, -0.5,   20.5);
321   fhRegionMultMin->SetXTitle("N_{ch tracks}");
322   fhRegionMultMin->Sumw2();
323   list->Add( fhRegionMultMin );         // At(3)            
324   
325   //Average pT in MIN region
326   fhMinRegAvePt = new TH1F("hMinRegAvePt", "#LTp_{T}#GT",  50, 0.,   20.);
327   fhMinRegAvePt->SetXTitle("p_{T} (GeV/c)");
328   fhMinRegAvePt->Sumw2();
329   list->Add( fhMinRegAvePt );           // At(4)
330   
331   //Sum pT in MIN region
332   fhMinRegSumPt = new TH1F("hMinRegSumPt", "#Sigma p_{T} ",  50, 0.,   20.);
333   fhMinRegSumPt->SetYTitle("Ed^{3}N_{tracks}/dp^{3} (c^{3}/GeV^{2})");  
334   fhMinRegSumPt->SetXTitle("#Sigma p_{T} (GeV/c)");
335   fhMinRegSumPt->Sumw2();
336   list->Add( fhMinRegSumPt );           // At(5)
337   
338   //Track with maximum pT in MIN region
339   fhMinRegMaxPtPart = new TH1F("hMinRegMaxPtPart", "max(p_{T})|_{event} ",  50, 0.,   20.);
340   fhMinRegMaxPtPart->SetYTitle("Ed^{3}N_{tracks}/dp^{3} (c^{3}/GeV^{2})");  
341   fhMinRegMaxPtPart->SetXTitle("p_{T} (GeV/c)");
342   fhMinRegMaxPtPart->Sumw2();
343   list->Add( fhMinRegMaxPtPart );       // At(6)
344   
345   //Sum pT vs. multiplicity in MIN region
346   fhMinRegSumPtvsMult = new TH1F("hMinRegSumPtvsMult", "#Sigma p_{T} vs. Multiplicity ",  21, -0.5,   20.5);
347   fhMinRegSumPtvsMult->SetYTitle("#Sigma p_{T} (GeV/c)");  
348   fhMinRegSumPtvsMult->SetXTitle("N_{charge}");
349   fhMinRegSumPtvsMult->Sumw2();
350   list->Add( fhMinRegSumPtvsMult );     // At(7);
351   
352   //Phi correlation track-cluster vs. leading cluster pT 
353   fhdNdEtaPhiDist  = new TH2F("hdNdEtaPhiDist", Form("Charge particle density |#eta|<%3.1f vs #Delta#phi", fTrackEtaCut),62, 0.,   2.*TMath::Pi(), bins, min, max);
354   fhdNdEtaPhiDist->SetXTitle("#Delta#phi");
355   fhdNdEtaPhiDist->SetYTitle("Leading cluster p_{T}");
356   fhdNdEtaPhiDist->Sumw2();
357   list->Add( fhdNdEtaPhiDist );        // At(8)
358   
359   //Can be used to get track pT distribution for different cluster pT bins (full region)
360   fhFullRegPartPtDistVsEt = new TH2F("hFullRegPartPtDistVsEt", Form( "dN/dp_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut),100,0.,50., bins, min, max);
361   fhFullRegPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
362   fhFullRegPartPtDistVsEt->SetXTitle("p_{T}");
363   fhFullRegPartPtDistVsEt->Sumw2();
364   list->Add( fhFullRegPartPtDistVsEt );  // At(9) 
365   
366   //Can be used to get part pT distribution for different cluster pT bins (transverse region)
367   fhTransRegPartPtDistVsEt = new TH2F("hTransRegPartPtDistVsEt", Form( "dN/dp_{T} in tranvese regions |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut),100,0.,50., bins, min,   max);
368   fhTransRegPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
369   fhTransRegPartPtDistVsEt->SetXTitle("p_{T}");
370   fhTransRegPartPtDistVsEt->Sumw2();
371   list->Add( fhTransRegPartPtDistVsEt );  // At(10)
372   
373   //Sum pT in MAX region vs. leading-cluster pT
374   fhRegionSumPtMaxVsEt = new TH1F("hRegionSumPtMaxVsEt",  "P_{T}^{90, max} vs Leading cluster p_{T}",  bins, min,   max);
375   fhRegionSumPtMaxVsEt->SetXTitle("p_{T} (GeV/c)");
376   fhRegionSumPtMaxVsEt->Sumw2();
377   list->Add( fhRegionSumPtMaxVsEt );      // At(11)
378   
379
380   //Sum pT in MIN region vs. leading-cluster pT
381   fhRegionSumPtMinVsEt = new TH1F("hRegionSumPtMinVsEt",   "P_{T}^{90, min} vs Leading cluster p_{T}",  bins, min,   max);
382   fhRegionSumPtMinVsEt->SetXTitle("p_{T} (GeV/c)");
383   fhRegionSumPtMinVsEt->Sumw2();
384   list->Add( fhRegionSumPtMinVsEt );      // At(12)
385  
386   //Multiplicity in MAX region
387   fhRegionMultMax = new TH1I("hRegionMultMax",      "N_{ch}^{90, max}",  21, -0.5,   20.5);
388   fhRegionMultMax->SetXTitle("N_{ch tracks}");
389   fhRegionMultMax->Sumw2();
390   list->Add( fhRegionMultMax );           // At(13)
391   
392   //Multiplicity in MAX region vs. leading-cluster pT
393   fhRegionMultMaxVsEt = new TH1F("hRegionMultMaxVsEt",  "N_{ch}^{90, max} vs Leading cluster p_{T}",  bins, min,   max);
394   fhRegionMultMaxVsEt->SetXTitle("p_{T} (GeV/c)");
395   fhRegionMultMaxVsEt->Sumw2();
396   list->Add( fhRegionMultMaxVsEt );       // At(14)
397   
398   //Multiplicity in MIN region vs. leading-cluster pT
399   fhRegionMultMinVsEt = new TH1F("hRegionMultMinVsEt",  "N_{ch}^{90, min} vs Leading cluster p_{T}",  bins, min,   max);
400   fhRegionMultMinVsEt->SetXTitle("p_{T} (GeV/c)");
401   fhRegionMultMinVsEt->Sumw2();
402   list->Add( fhRegionMultMinVsEt );      // At(15)
403  
404   //Average sum pT in TRANSVERSE(MIN+MAX) region vs. leading-cluster pT 
405   fhRegionAveSumPtVsEt = new TH1F("hRegionAveSumPtVsEt", "(P_{T}^{90, max} + P_{T}^{90, min})/2 vs Leading cluster p_{T}",  bins, min,   max);
406   fhRegionAveSumPtVsEt->SetXTitle("p_{T} (GeV/c)");
407   fhRegionAveSumPtVsEt->Sumw2();
408   list->Add( fhRegionAveSumPtVsEt );     // At(16)
409   
410   //Difference sum pT (MAX-MIN) vs.  leading-cluster pT
411   fhRegionDiffSumPtVsEt= new TH1F("hRegionDiffSumPtVsEt", "(P_{T}^{90, max} - P_{T}^{90, min}) vs Leading cluster p_{T}",  bins, min,   max);
412   fhRegionDiffSumPtVsEt->SetXTitle("P_{T} (GeV/c)");
413   fhRegionDiffSumPtVsEt->Sumw2();
414   list->Add( fhRegionDiffSumPtVsEt );    // At(17)
415   
416   //Average track pT in MAX region vs. leading-cluster pT
417   fhRegionAvePartPtMaxVsEt = new TH1F("hRegionAvePartPtMaxVsEt", "#LTp_{T}#GT^{90, max} vs Leading cluster p_{T}",  bins, min,   max);
418   fhRegionAvePartPtMaxVsEt->SetXTitle("p_{T} (GeV/c)");
419   fhRegionAvePartPtMaxVsEt->Sumw2();
420   list->Add( fhRegionAvePartPtMaxVsEt );  // At(18)
421   
422   //Average track pT in MIN region vs. leading-cluster pT
423   fhRegionAvePartPtMinVsEt = new TH1F("hRegionAvePartPtMinVsEt", "#LTp_{T}#GT^{90, min} vs Leading cluster p_{T}",  bins, min,   max);
424   fhRegionAvePartPtMinVsEt->SetXTitle("p_{T} (GeV/c)");
425   fhRegionAvePartPtMinVsEt->Sumw2();
426   list->Add( fhRegionAvePartPtMinVsEt );   // At(19)
427   
428   //Maximum track pT in MAX region vs. leading-cluster pT
429   fhRegionMaxPartPtMaxVsEt = new TH1F("hRegionMaxPartPtMaxVsEt", "max(p_{T})^{90} vs Leading cluster p_{T}",  bins, min,   max);
430   fhRegionMaxPartPtMaxVsEt->SetXTitle("p_{T} (GeV/c)");
431   fhRegionMaxPartPtMaxVsEt->Sumw2();
432   list->Add( fhRegionMaxPartPtMaxVsEt );    // At(20)
433   
434   //Multiplicity in FORWARD region
435   fhRegForwardMult = new TH2F("hRegForwardMult", "N_{ch}^{forward}", bins, min, max, 21, -0.5,   20.5);
436   fhRegForwardMult->SetXTitle("N_{ch tracks}");
437   fhRegForwardMult->Sumw2();
438   list->Add( fhRegForwardMult );           // At(25)
439   
440   //Sum pT in FORWARD region vs. multiplicity
441   fhRegForwardSumPtvsMult = new TH2F("hRegForwardSumPtvsMult", "Forward #Sigma p_{T} vs. Multiplicity ", bins, min, max, 21, -0.5,   20.5);
442   fhRegForwardSumPtvsMult->SetYTitle("#Sigma p_{T} (GeV/c)");  
443   fhRegForwardSumPtvsMult->SetXTitle("N_{charge}");
444   fhRegForwardSumPtvsMult->Sumw2();
445   list->Add( fhRegForwardSumPtvsMult );     // At(26);
446   
447   //Multiplicity in BACKWARD region
448   fhRegBackwardMult = new TH2F("hRegBackwardMult", "N_{ch}^{backward}", bins, min, max, 21, -0.5,   20.5);
449   fhRegBackwardMult->SetXTitle("N_{ch tracks}");
450   fhRegBackwardMult->Sumw2();
451   list->Add( fhRegBackwardMult );           // At(27)
452  
453   //Sum pT in BACKWARD region vs. multiplicity
454   fhRegBackwardSumPtvsMult = new TH2F("hRegBackwardSumPtvsMult", "Backward #Sigma p_{T} vs. Multiplicity ", bins, min, max, 21, -0.5,   20.5);
455   fhRegBackwardSumPtvsMult->SetYTitle("#Sigma p_{T} (GeV/c)");  
456   fhRegBackwardSumPtvsMult->SetXTitle("N_{charge}");
457   fhRegBackwardSumPtvsMult->Sumw2();
458   list->Add( fhRegBackwardSumPtvsMult );     // At(28);
459   
460   //Track pT distribution in FORWARD region vs. leading-cluster pT 
461   fhRegForwardPartPtDistVsEt = new TH2F("hRegForwardPartPtDistVsEt", Form( "dN/dP_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut), 100,0.,50., bins, min, max);
462   fhRegForwardPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
463   fhRegForwardPartPtDistVsEt->SetXTitle("p_{T} (GeV/c)");
464   fhRegForwardPartPtDistVsEt->Sumw2();
465   list->Add( fhRegForwardPartPtDistVsEt );  // At(29) 
466   
467   //Track pT distribution in BACKWARD region vs. leading-cluster pT 
468   fhRegBackwardPartPtDistVsEt = new TH2F("hRegBackwardPartPtDistVsEt", Form( "dN/dP_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut), 100,0.,50., bins, min, max);
469   fhRegBackwardPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
470   fhRegBackwardPartPtDistVsEt->SetXTitle("p_{T}");
471   fhRegBackwardPartPtDistVsEt->Sumw2();
472   list->Add( fhRegBackwardPartPtDistVsEt );  // At(30) 
473   
474   //Multiplicity in TRANSVERSE (MIN+MAX) region
475   fhRegTransMult  = new TH2F("hRegTransMult", "N_{ch}^{transv}", bins, min, max, 21, -0.5,   20.5);
476   fhRegTransMult->SetXTitle("N_{ch tracks}");
477   fhRegTransMult->Sumw2();
478   list->Add( fhRegTransMult );           // At(31)
479   
480   //Sum pT in TRANSVERSE (MIN+MAX) region vs. multiplicity
481   fhRegTransSumPtVsMult = new TH2F("hRegTransSumPtVsMult", "Transverse #Sigma p_{T} vs. Multiplicity ",bins, min, max, 21, -0.5,   20.5);
482   fhRegTransSumPtVsMult->SetYTitle("#Sigma p_{T} (GeV/c)");  
483   fhRegTransSumPtVsMult->SetXTitle("N_{charge}");
484   fhRegTransSumPtVsMult->Sumw2();
485   list->Add( fhRegTransSumPtVsMult );     // At(32);
486   
487   //Sum pT in MIN region per cluster pT bin
488   fhMinRegSumPtJetPtBin = new TH2F("hMinRegSumPtJetPtBin", "Transverse Min Reg #Sigma p_{T} per cluster pT bin",bins, min, max, 50, 0.,   20.);
489   fhMinRegSumPtJetPtBin->SetXTitle("Leading cluster p_{T}");
490   fhMinRegSumPtJetPtBin->Sumw2();
491   list->Add( fhMinRegSumPtJetPtBin );           // At(33)
492   
493   //Sum pT in MAX region per cluster pT bin
494   fhMaxRegSumPtJetPtBin = new TH2F("hMaxRegSumPtJetPtBin",      "Transverse Max Reg #Sigma p_{T} per cluster pT bin", bins, min, max, 50, 0.,   20.);
495   fhMaxRegSumPtJetPtBin->SetXTitle("Leading cluster p_{T}");
496   fhMaxRegSumPtJetPtBin->Sumw2();
497   list->Add( fhMaxRegSumPtJetPtBin );           // At(34)
498   
499   //Multiplicity in main vertex
500   fhVertexMult = new TH1F("hVertexMult",      "Multiplicity in Main Vertex", 81, -0.5 , 80.5);
501   fhVertexMult->SetXTitle("Main Vertex Multiplicity");
502   fhVertexMult->Sumw2();
503   list->Add( fhVertexMult ); //At(35)
504   
505   fh1Xsec = new TProfile("h1Xsec","xsec from pyxsec.root",1,0,1); 
506   fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
507   fh1Xsec->Sumw2();
508   list->Add( fh1Xsec );            //At(36)
509   
510   fh1Trials = new TH1F("h1Trials","trials from pyxsec.root",1,0,1);
511   fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
512   fh1Trials->Sumw2();
513   list->Add( fh1Trials ); //At(37)
514   
515   
516   fListOfHistos=list; 
517 }
518
519 //____________________________________________________________________
520 void AliHistogramsUE::CreateCorrectionsContainer(AliCFManager* cfman,Int_t bins, Double_t min, Double_t max, Double_t etacut, Double_t jetetacut){
521   
522   fBinsPtInHist = bins;
523   fMinJetPtInHist = min;
524   fMaxJetPtInHist = max;
525   fTrackEtaCut= etacut;
526   fJet1EtaCut = jetetacut;
527   //Define some constant
528   const Double_t minpT=fMinJetPtInHist; 
529   const Double_t maxpT=fMaxJetPtInHist;
530
531   const Double_t mineta=(-1.*fJet1EtaCut);
532   const Double_t maxeta=fJet1EtaCut;
533         
534   const Double_t minprocess=-0.5; 
535   const Double_t maxprocess=9.5;
536
537   const Double_t mindeltaeta=-5.; 
538   const Double_t maxdeltaeta=5.;
539
540   const Double_t mindeltaphi=0.;
541   const Double_t maxdeltaphi=7.;
542   
543   const Double_t  minradius=0.; 
544   const Double_t  maxradius=10.;
545
546   //Define sensitive variables
547   UInt_t ipT       = 0;     //leading track pT
548   UInt_t ieta      = 1;     //leading track eta 
549   UInt_t iprocess  = 2;     //process type (ND,DD,SD)  
550   UInt_t ipTMC     = 3;     //MC leading track pT  
551   UInt_t ietaMC    = 4;     //MC leading track eta 
552   UInt_t ideltaeta = 5;     //leading track eta  reco-MC
553   UInt_t ideltaphi = 6;     //leading track phi  reco-MC 
554   UInt_t iradius   = 7;     //leading track radius reco-MC 
555   
556   //Set-up grid
557   UInt_t nstep = 6;
558   const Int_t nvar = 8;
559   const Int_t nbinspT = fBinsPtInHist;
560   const Int_t nbinseta = 10;
561   const Int_t nbinsprocess = 10;
562   const Int_t nbinsdeltaeta = 20;
563   const Int_t nbinsdeltaphi = 20; 
564   const Int_t nbinsradius = 20; 
565  
566
567   Int_t iBin[nvar];
568   iBin[0] = nbinspT;
569   iBin[1] = nbinseta;
570   iBin[2] = nbinsprocess;
571   iBin[3] = nbinspT;
572   iBin[4] = nbinseta;
573   iBin[5] = nbinsdeltaeta;
574   iBin[6] = nbinsdeltaphi;
575   iBin[7] = nbinsradius;
576   
577   //lower bounds
578   Double_t *binLimpT=new Double_t[nbinspT+1];
579   for (Int_t i=0; i<=nbinspT; i++) binLimpT[i]=(Double_t)minpT + (maxpT-minpT)/nbinspT*(Double_t)i ;
580   
581   Double_t *binLimeta=new Double_t[nbinseta+1];
582   for (Int_t i=0; i<=nbinseta; i++) binLimeta[i]=(Double_t)mineta + (maxeta-mineta)/nbinseta*(Double_t)i ;
583   
584   Double_t *binLimprocess=new Double_t[nbinsprocess+1];
585   for (Int_t i=0; i<=nbinsprocess; i++) binLimprocess[i]=(Double_t)minprocess + (maxprocess-minprocess)/nbinsprocess*(Double_t)i ;
586   
587   Double_t *binLimdeltaeta=new Double_t[nbinsdeltaeta+1];
588   for (Int_t i=0; i<=nbinsdeltaeta; i++) binLimdeltaeta[i]=(Double_t)mindeltaeta + (maxdeltaeta-mindeltaeta)/nbinsdeltaeta*(Double_t)i ;
589   
590   Double_t *binLimdeltaphi=new Double_t[nbinsdeltaphi+1];
591   for (Int_t i=0; i<=nbinsdeltaphi; i++) binLimdeltaphi[i]=(Double_t)mindeltaphi + (maxdeltaphi-mindeltaphi)/nbinsdeltaphi*(Double_t)i ;
592   
593   Double_t *binLimradius=new Double_t[nbinsradius+1];
594   for (Int_t i=0; i<=nbinsradius; i++) binLimradius[i]=(Double_t)minradius + (maxradius-minradius)/nbinsradius*(Double_t)i ;
595   
596
597   //Container
598   AliCFContainer * container = new AliCFContainer("container1", "EventSelection",nstep,nvar,iBin);
599   container->SetBinLimits(ipT,binLimpT);
600   container->SetBinLimits(ieta,binLimeta);
601   container->SetBinLimits(iprocess,binLimprocess);
602   container->SetBinLimits(ipTMC,binLimpT);
603   container->SetBinLimits(ietaMC,binLimeta);
604   container->SetBinLimits(ideltaeta,binLimdeltaeta);
605   container->SetBinLimits(ideltaphi,binLimdeltaphi);
606   container->SetBinLimits(iradius,binLimradius);
607
608   container->SetVarTitle(ipT,"Leading track p_{T} (reco.)");
609   container->SetVarTitle(ieta,"Leading track #eta (reco.)");
610   container->SetVarTitle(iprocess,"Process");
611   container->SetVarTitle(ipTMC,"Leading track p_{T} (true)");
612   container->SetVarTitle(ietaMC,"Leading track #eta (true)");
613   container->SetVarTitle(ideltaeta,"Leading track #Delta #eta (reco.-true)");
614   container->SetVarTitle(ideltaphi,"Leading track #Delta #phi (reco.-true)");
615   container->SetVarTitle(iradius,"Leading track R (reco.-true)");
616
617   //set steps
618   container->SetStepTitle(0,"Triggered");
619   container->SetStepTitle(1,"Pass physics selection");
620   container->SetStepTitle(2,"Pass primary vertex cuts");
621   container->SetStepTitle(3,"Required analysis topology ");
622   container->SetStepTitle(4,"Leading track p_{T} > 1 GeV/c");
623   container->SetStepTitle(5,"Leading track correctly identified");
624   
625   cfman->SetEventContainer(container);
626   
627 }
628
629
630
631 //____________________________________________________________________
632 void  AliHistogramsUE::CreateHistogramsCorrections(TList *list,Int_t bins, Double_t min, Double_t max, Double_t etacut)
633 {
634
635   // Create all histograms necessary for UE corrections
636   fBinsPtInHist = bins;
637   fMinJetPtInHist = min;
638   fMaxJetPtInHist = max;
639   fTrackEtaCut= etacut;
640
641   //Number of reconstructed clusters  
642   fhNJets = new TH1F("hNJets", "Number of clusters",  20, 0, 20);
643   fhNJets->SetXTitle("Number of reconstructed clusters");
644   fhNJets->SetYTitle("#");
645   fhNJets->Sumw2();
646   list->Add( fhNJets );                  
647   
648   //Cross-section from MC
649   fh1Xsec = new TProfile("h1Xsec","xsec from pyxsec.root",1,0,1); 
650   fh1Xsec->GetXaxis()->SetBinLabel(1,"<#sigma>");
651   fh1Xsec->Sumw2();
652   list->Add( fh1Xsec );   
653
654   //Number of trials from MC
655   fh1Trials = new TH1F("h1Trials","trials from pyxsec.root",1,0,1);
656   fh1Trials->GetXaxis()->SetBinLabel(1,"#sum{ntrials}");
657   fh1Trials->Sumw2();
658   list->Add( fh1Trials ); 
659
660   //FOR TRACK EFFICIENCY 
661   //Phi correlation track-cluster vs. leading cluster pT 
662   fhdNdEtaPhiDist  = new TH2F("hdNdEtaPhiDist", Form("Charge particle density |#eta|<%3.1f vs #Delta#phi", fTrackEtaCut),62, 0.,   2.*TMath::Pi(), bins, min, max);
663   fhdNdEtaPhiDist->SetXTitle("#Delta#phi");
664   fhdNdEtaPhiDist->SetYTitle("Leading cluster p_{T}");
665   fhdNdEtaPhiDist->Sumw2();
666   list->Add( fhdNdEtaPhiDist );        
667   //idem for MC true 
668   fhdNdEtaPhiDistMC  = new TH2F("hdNdEtaPhiDistMC", Form("Charge particle density |#eta|<%3.1f vs #Delta#phi", fTrackEtaCut),62, 0.,   2.*TMath::Pi(), bins, min, max);
669   fhdNdEtaPhiDistMC->SetXTitle("#Delta#phi");
670   fhdNdEtaPhiDistMC->SetYTitle("Leading cluster p_{T}");
671   fhdNdEtaPhiDistMC->Sumw2();
672   list->Add( fhdNdEtaPhiDistMC );      
673   
674
675   //Can be used to get track pT distribution for different cluster pT bins (full region)
676   fhFullRegPartPtDistVsEt = new TH2F("hFullRegPartPtDistVsEt", Form( "dN/dp_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut),100,0.,50., bins, min, max);
677   fhFullRegPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
678   fhFullRegPartPtDistVsEt->SetXTitle("p_{T}");
679   fhFullRegPartPtDistVsEt->Sumw2();
680   list->Add( fhFullRegPartPtDistVsEt );   
681   //idem for MC true
682   fhFullRegPartPtDistVsEtMC = new TH2F("hFullRegPartPtDistVsEtMC", Form( "dN/dp_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut),100,0.,50., bins, min, max);
683   fhFullRegPartPtDistVsEtMC->SetYTitle("Leading cluster p_{T}");
684   fhFullRegPartPtDistVsEtMC->SetXTitle("p_{T}");
685   fhFullRegPartPtDistVsEtMC->Sumw2();
686   list->Add( fhFullRegPartPtDistVsEtMC );   
687
688
689   //Can be used to get part pT distribution for different cluster pT bins (transverse region)
690   fhTransRegPartPtDistVsEt = new TH2F("hTransRegPartPtDistVsEt", Form( "dN/dp_{T} in tranvese regions |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut),100,0.,50., bins, min,   max);
691   fhTransRegPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
692   fhTransRegPartPtDistVsEt->SetXTitle("p_{T}");
693   fhTransRegPartPtDistVsEt->Sumw2();
694   list->Add( fhTransRegPartPtDistVsEt );  
695   //idem for MC true
696   fhTransRegPartPtDistVsEtMC = new TH2F("hTransRegPartPtDistVsEtMC", Form( "dN/dp_{T} in tranvese regions |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut),100,0.,50., bins, min,   max);
697   fhTransRegPartPtDistVsEtMC->SetYTitle("Leading cluster p_{T}");
698   fhTransRegPartPtDistVsEtMC->SetXTitle("p_{T}");
699   fhTransRegPartPtDistVsEtMC->Sumw2();
700   list->Add( fhTransRegPartPtDistVsEtMC ); 
701
702
703   //Track pT distribution in FORWARD region vs. leading-cluster pT 
704   fhRegForwardPartPtDistVsEt = new TH2F("hRegForwardPartPtDistVsEt", Form( "dN/dP_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut), 100,0.,50., bins, min, max);
705   fhRegForwardPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
706   fhRegForwardPartPtDistVsEt->SetXTitle("p_{T} (GeV/c)");
707   fhRegForwardPartPtDistVsEt->Sumw2();
708   list->Add( fhRegForwardPartPtDistVsEt );   
709   //idem for MC true 
710   fhRegForwardPartPtDistVsEtMC = new TH2F("hRegForwardPartPtDistVsEtMC", Form( "dN/dP_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut), 100,0.,50., bins, min, max);
711   fhRegForwardPartPtDistVsEtMC->SetYTitle("Leading cluster p_{T}");
712   fhRegForwardPartPtDistVsEtMC->SetXTitle("p_{T} (GeV/c)");
713   fhRegForwardPartPtDistVsEtMC->Sumw2();
714   list->Add( fhRegForwardPartPtDistVsEtMC );   
715   
716   //Track pT distribution in BACKWARD region vs. leading-cluster pT 
717   fhRegBackwardPartPtDistVsEt = new TH2F("hRegBackwardPartPtDistVsEt", Form( "dN/dP_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut), 100,0.,50., bins, min, max);
718   fhRegBackwardPartPtDistVsEt->SetYTitle("Leading cluster p_{T}");
719   fhRegBackwardPartPtDistVsEt->SetXTitle("p_{T}");
720   fhRegBackwardPartPtDistVsEt->Sumw2();
721   list->Add( fhRegBackwardPartPtDistVsEt );   
722   //idem for MC true 
723   fhRegBackwardPartPtDistVsEtMC = new TH2F("hRegBackwardPartPtDistVsEtMC", Form( "dN/dP_{T} |#eta|<%3.1f vs Leading cluster p_{T}", fTrackEtaCut), 100,0.,50., bins, min, max);
724   fhRegBackwardPartPtDistVsEtMC->SetYTitle("Leading cluster p_{T}");
725   fhRegBackwardPartPtDistVsEtMC->SetXTitle("p_{T}");
726   fhRegBackwardPartPtDistVsEtMC->Sumw2();
727   list->Add( fhRegBackwardPartPtDistVsEtMC );  
728
729   //FOR DCA DISTRIBUTION 
730   fhDCAxy = new TH2F("hDCAxy","Transverse DCA vs. leading track p_{T} ",50,-5.,5.,bins,min,max); 
731   fhDCAxy->SetYTitle("Leading cluster p_{T}");
732   fhDCAxy->SetXTitle("d_{0}");
733   fhDCAxy->Sumw2();
734   list->Add( fhDCAxy );
735   
736   //Fill only if reconstructed points back to a true primary  
737   fhDCAxyPrimary = new TH2F("hDCAxyPrimary","Transverse DCA vs. leading track p_{T} (primaries)",50,-5.,5.,bins,min,max);
738   fhDCAxyPrimary->SetYTitle("Leading cluster p_{T}");
739   fhDCAxyPrimary->SetXTitle("d_{0}");
740   fhDCAxyPrimary->Sumw2();
741   list->Add( fhDCAxyPrimary );  
742  
743   fListOfHistos=list; 
744 }
745
746 //____________________________________________________________________
747 void AliHistogramsUE::DrawUE(Int_t debug){
748
749     // To draw histograms at the end of task running 
750     // Normalize histos to region area TODO: 
751     // Normalization done at Analysis, taking into account 
752     // area variations on per-event basis (cone case)
753    
754     //HIGH WARNING!!!!!: DO NOT SCALE ANY OF THE ORIGINAL HISTOGRAMS
755     //MAKE A COPY, DRAW IT, And later sacale that copy. CAF Issue!!!!!
756     
757     Int_t binsPtInHist = fhEleadingPt->GetNbinsX();
758     Double_t minJetPtInHist = fhEleadingPt->GetXaxis()->GetBinLowEdge(1);
759     Double_t maxJetPtInHist = fhEleadingPt->GetXaxis()->GetBinUpEdge(binsPtInHist);
760      
761     //Sum pT
762     TCanvas* c1 = new TCanvas("c1",Form("sumPt dist (%s)", GetTitle()),60,60,1100,700);
763     c1->Divide(2,2);
764     c1->cd(1);
765     TH1F *h1r = new TH1F("hRegionEtvsSumPtMax" , "", binsPtInHist,  minJetPtInHist, maxJetPtInHist);
766     //TH1F *h1r = new TH1F();
767     h1r->Divide(fhRegionSumPtMaxVsEt,fhEleadingPt,1,1);
768     //h1r->Scale( areafactor );
769     h1r->SetMarkerStyle(20);
770     h1r->SetXTitle("P_{T} of Leading Jet (GeV/c)");
771     h1r->SetYTitle("P_{T}^{90, max}");
772     h1r->DrawCopy("p");
773     
774     c1->cd(2);
775     TH1F *h2r = new TH1F("hRegionEtvsSumPtMin" , "", binsPtInHist,  minJetPtInHist, maxJetPtInHist);
776     h2r->Divide(fhRegionSumPtMinVsEt,fhEleadingPt,1,1);
777     //h2r->Scale( areafactor );
778     h2r->SetMarkerStyle(20);
779     h2r->SetXTitle("P_{T} of Leading Jet (GeV/c)");
780     h2r->SetYTitle("P_{T}^{90, min}");
781     h2r->DrawCopy("p");
782     
783     c1->cd(3);
784     TH1F *h4r = new TH1F("hRegionEtvsDiffPt" , "", binsPtInHist,  minJetPtInHist, maxJetPtInHist);
785     //TH1F *h41r = new TH1F("hRegForwvsDiffPt" , "", fbinsPtInHist,  fMinJetPtInHist, fMaxJetPtInHist);
786     //TH1F *h42r = new TH1F("hRegBackvsDiffPt" , "", fbinsPtInHist,  fMinJetPtInHist, fMaxJetPtInHist);
787     //h41r->Divide(fhRegForwardSumPtVsEt,fhEleadingPt,1,1);
788     //h42r->Divide(fhRegBackwardSumPtVsEt,fhEleadingPt,1,1);
789     h4r->Divide(fhRegionAveSumPtVsEt,fhEleadingPt,1,1);
790     //h4r->Scale(2.); // make average
791     //h4r->Scale( areafactor );
792     h4r->SetYTitle("#DeltaP_{T}^{90}");
793     h4r->SetXTitle("P_{T} of Leading Jet (GeV/c)");
794     h4r->SetMarkerStyle(20);
795     h4r->DrawCopy("p");
796
797     c1->cd(4);
798     TH1F *h5r = new TH1F("hRegionMultMaxVsEtleading",   "",  binsPtInHist, minJetPtInHist,   maxJetPtInHist);
799     TH1F *h6r = new TH1F("hRegionMultMinVsEtleading",   "",  binsPtInHist, minJetPtInHist,   maxJetPtInHist);
800     h5r->Divide(fhRegionMultMaxVsEt,fhEleadingPt,1,1);
801     h6r->Divide(fhRegionMultMinVsEt,fhEleadingPt,1,1);
802     //h5r->Scale( areafactor );
803     h5r->SetYTitle("N_{Tracks}^{90}");
804     h5r->SetXTitle("P_{T} of Leading Jet (GeV/c)");
805     h5r->SetMarkerStyle(20);
806     h5r->DrawCopy("p");
807     h6r->SetMarkerStyle(21);
808     h6r->SetMarkerColor(2);
809     h6r->SetYTitle("N_{Tracks}^{90}");
810     h6r->SetXTitle("P_{T} of Leading Jet (GeV/c)");
811     h6r->DrawCopy("p same");
812     c1->Update();
813     
814     //Get Normalization
815     //Double_t xsec = fh1Xsec->GetBinContent(1);
816     Double_t xsec = fh1Xsec->GetBinContent(1);
817     Double_t ntrials = fh1Trials->GetBinContent(1);
818     Double_t normFactor = xsec/ntrials;
819     if(debug > 1)Printf("xSec %f nTrials %f Norm %f \n",xsec,ntrials,normFactor);
820     
821     
822     //Jet pT distribution
823     TCanvas* c2 = new TCanvas("c2","Jet Pt dist",160,160,1200,800);
824     TH1 * copy = 0;
825     c2->Divide(2,2);
826     c2->cd(1);
827     fhEleadingPt->SetMarkerStyle(20);
828     fhEleadingPt->SetMarkerColor(2);
829     //if( normFactor > 0.) fhEleadingPt->Scale(normFactor);
830     //fhEleadingPt->Draw("p");
831     copy = fhEleadingPt->DrawCopy("p");
832     if( normFactor > 0.) copy->Scale(normFactor);
833     gPad->SetLogy();
834     
835     c2->cd(2);
836     Int_t xbin1 = fhdNdEtaPhiDist->GetYaxis()->FindFixBin(minJetPtInHist);
837     Int_t xbin2 = fhdNdEtaPhiDist->GetYaxis()->FindFixBin(maxJetPtInHist);
838     TH1D* dNdEtaPhiDistAllJets = fhdNdEtaPhiDist->ProjectionX("dNdEtaPhiDistAllJets",xbin1,xbin2);
839     dNdEtaPhiDistAllJets->SetMarkerStyle(20);
840     dNdEtaPhiDistAllJets->SetMarkerColor(2);
841     dNdEtaPhiDistAllJets->DrawCopy("p");
842     gPad->SetLogy();
843     
844     c2->cd(3);      
845     fhNJets->DrawCopy();
846     //c2->cd(4);      
847     //fhValidRegion->DrawCopy("p");
848     //fhTransRegPartPtDist  = (TH1F*)fListOfHistos->At(2);
849     //fhRegionMultMin          = (TH1F*)fListOfHistos->At(3);
850     //fhMinRegAvePt     = (TH1F*)fListOfHistos->At(4);
851     //fhMinRegSumPt     = (TH1F*)fListOfHistos->At(5);   
852     //fhMinRegMaxPtPart     = (TH1F*)fListOfHistos->At(6);
853     //fhMinRegSumPtvsMult   = (TH1F*)fListOfHistos->At(7);
854     c2->Update(); 
855     
856     //pT distributions
857     TCanvas* c3 = new TCanvas("c3"," pT dist",160,160,1200,800);
858     c3->Divide(2,2);
859     c3->cd(1);
860      //fhTransRegPartPtDist->SetMarkerStyle(20);
861      //fhTransRegPartPtDist->SetMarkerColor(2); 
862      //fhTransRegPartPtDist->Scale(areafactor/fhTransRegPartPtDist->GetEntries());
863      //fhTransRegPartPtDist->DrawCopy("p");
864      //gPad->SetLogy();
865      
866
867     c3->cd(2); 
868     fhMinRegSumPt->SetMarkerStyle(20);
869     fhMinRegSumPt->SetMarkerColor(2);  
870     //fhMinRegSumPt->Scale(areafactor);
871     fhMinRegSumPt->DrawCopy("p");
872     gPad->SetLogy();
873     
874     c3->cd(3);
875     fhMinRegAvePt->SetMarkerStyle(20);
876     fhMinRegAvePt->SetMarkerColor(2);  
877     //fhMinRegAvePt->Scale(areafactor);
878     fhMinRegAvePt->DrawCopy("p");
879     gPad->SetLogy();
880     
881     c3->cd(4);
882     TH1F *h7r = new TH1F("hRegionMultMinVsMult",   "",  21, -0.5,   20.5);
883     h7r->Divide(fhMinRegSumPtvsMult,fhRegionMultMin,1,1);
884     h7r->SetMarkerStyle(20);
885     h7r->SetMarkerColor(2);   
886     h7r->DrawCopy("p");
887     c3->Update();
888     
889     
890
891     //Save canvas
892     c1->SaveAs("c1.pdf");
893     AliInfo("Canvas 1 saved");
894     c2->SaveAs("c2.pdf");
895     AliInfo("Canvas 2 saved");
896     c3->SaveAs("c3.pdf");
897     AliInfo("Canvas 3 saved");
898   
899 }
900
901 //____________________________________________________________________
902 void AliHistogramsUE::FillHistogram(const char* name, Double_t fillX){
903   
904   // Fill 1D histogram with double
905   ((TH1F*)fListOfHistos->FindObject(name))->Fill(fillX);
906
907 }
908
909 //____________________________________________________________________
910 void AliHistogramsUE::FillHistogram(const char* name, Int_t fillX){
911   
912   // Fill 1D histogram with integer
913   ((TH1F*)fListOfHistos->FindObject(name))->Fill(fillX);
914
915 }
916
917 //____________________________________________________________________
918 void AliHistogramsUE::FillHistogram(const char* name, Double_t fillX, Double_t fillY){
919   
920  // Case of TH1F with weight or TH2F w/o weight
921  TObject *obj = fListOfHistos->FindObject(name);
922  if (obj->InheritsFrom("TH1F")){
923         ((TH1F*)fListOfHistos->FindObject(name))->Fill(fillX, fillY);
924  } else {
925         ((TH2F*)fListOfHistos->FindObject(name))->Fill(fillX, fillY);
926         }
927
928
929 }
930
931 //____________________________________________________________________
932 void AliHistogramsUE::FillHistogram(const char* name, Double_t fillX, Double_t fillY, Double_t weight){
933
934   // Fill 2D histogram with double and weight   
935   ((TH2F*)fListOfHistos->FindObject(name))->Fill(fillX, fillY, weight);
936
937 }
938
939 //____________________________________________________________________
940 void AliHistogramsUE::FillHistogram(const char* name, Double_t fillX, Int_t fillY, Double_t weight){
941   
942   // Fill 2D histogram with integer and weight   
943   ((TH2F*)fListOfHistos->FindObject(name))->Fill(fillX, fillY, weight);
944
945 }
946
947 //____________________________________________________________________
948 TObjArray*  AliHistogramsUE::GetHistosForPlotting(TString data, TString branches){
949
950         // Instance filled histos for plotting purpose
951         printf("Creating histograms ... \n");
952                 
953         printf("Reading file: %s\n",data.Data());
954         
955         // Read input files ----------------------------------------- 
956         TFile *fdata = new TFile(data.Data());
957         TDirectoryFile *ddata[20];
958         TList *ldata[20];
959
960         TObjArray  *arrb=branches.Tokenize(";");
961         TIter next(arrb);
962         TObject *o=0;
963         Int_t br=0;
964         while ( (o=next()) ){
965                 ddata[br] = (TDirectoryFile*)fdata->Get(Form("PWG4_UE%s",o->GetName()));
966                 if(!ddata[br]) printf("ERROR: No histo dir found! \n");
967                 ldata[br] = (TList*)ddata[br]->Get(Form("histosUE%s",o->GetName()));
968                 printf("Reading branch: %s\n",o->GetName());
969                 if(!ldata[br]) printf("ERROR: No histo list found! \n");
970                 br++;
971         }
972         
973         TObjArray *arr=new TObjArray;
974
975         TH1F *hjets[20]     = {0,};       // accepted leading jets
976         TH1F *hnjets[20]    = {0,};       // number of accepted jets
977         TH2F *hetaphi[20]   = {0,};       // delta-phi particle-jet correlation
978         TH2F *hptfull[20]   = {0,};       // particle pT all regions vs. jet pT 
979         TH2F *hpttransv[20] = {0,};       // particle pT transv. regions vs. jet pT 
980         TH1F *hmax[20]      = {0,};       // sum pT in MAX region
981         TH1F *hmin[20]      = {0,};       // sum pT in MIN region
982         TH1F *hmultmax[20]  = {0,};       // multiplicity in MAX region
983         TH1F *hmultmin[20]  = {0,};       // multiplicity in MIN region
984
985         for (Int_t i =0; i<br; i++){
986
987                 //Number of jets --------------------------------------------
988                 hnjets[i] = (TH1F*) ldata[i]->FindObject("fhNJets"); 
989                 hnjets[i]->GetXaxis()->SetTitle("Number of jets");
990                 hnjets[i]->GetYaxis()->SetTitle("1/n_{ev} dN/dN_{jets}");
991                 hnjets[i]->SetMarkerStyle(20); 
992                 hnjets[i]->SetMarkerColor(1+i);
993
994
995                 //Leading jets ----------------------------------------------
996                 hjets[i] = (TH1F*) ldata[i]->FindObject("hEleadingPt"); 
997                 hjets[i]->GetXaxis()->SetTitle("P_{T} (GeV/c)");
998                 hjets[i]->GetYaxis()->SetTitle("1/n_{ev} dN/dp_{T} (|#eta<0.5|)");
999                 hjets[i]->SetMarkerStyle(20); 
1000                 hjets[i]->SetMarkerColor(1+i);
1001                 hjets[i]->SetMinimum(0.1);
1002                 hjets[i]->SetMaximum(1000.);
1003
1004
1005                 //Transverse Region MAX -------------------------------------
1006                 hmax[i] = (TH1F*) ldata[i]->FindObject("hRegionSumPtMaxVsEt");
1007                 if (!hmax[i])AliInfo("Histo not found!!!");
1008                 hmax[i]->GetXaxis()->SetTitle("Leading jet P_{T} (GeV/c)");
1009                 hmax[i]->GetYaxis()->SetTitle("P_{T}^{90,max} (GeV/c)");
1010                 hmax[i]->SetMarkerStyle(20);
1011                 hmax[i]->SetMarkerColor(1+i);
1012                 hmax[i]->Divide(hjets[i]); // normalize for jet spectrum
1013                 hmax[i]->SetMaximum(5.);
1014
1015
1016                 //Transverse Region MIN -------------------------------------
1017                 hmin[i] = (TH1F*) ldata[i]->FindObject("hRegionSumPtMinVsEt");
1018                 hmin[i]->GetXaxis()->SetTitle("Leading jet P_{T} (GeV/c)");
1019                 hmin[i]->GetYaxis()->SetTitle("P_{T}^{90,min} (GeV/c)");
1020                 hmin[i]->SetMarkerStyle(20);
1021                 hmin[i]->SetMarkerColor(1+i);
1022                 hmin[i]->SetMaximum(3.);
1023                 hmin[i]->Divide(hjets[i]); // normalize for jet spectrum
1024
1025
1026                 //Multiplicity MAX ------------------------------------------
1027                 hmultmax[i] = (TH1F*) ldata[i]->FindObject("hRegionMultMaxVsEt");
1028                 hmultmax[i]->GetXaxis()->SetTitle("Leading Jet P_{T} (GeV/c)");
1029                 hmultmax[i]->GetYaxis()->SetTitle("N_{ch}^{90,max}");
1030                 hmultmax[i]->SetMarkerStyle(20);
1031                 hmultmax[i]->SetMarkerColor(1+i);
1032                 hmultmax[i]->SetMaximum(10.);
1033                 hmultmax[i]->Divide(hjets[i]); // normalize for jet spectrum
1034
1035
1036                 //Multiplicity MIN ------------------------------------------
1037                 hmultmin[i] = (TH1F*) ldata[i]->FindObject("hRegionMultMinVsEt");
1038                 hmultmin[i]->GetXaxis()->SetTitle("Leading Jet P_{T} (GeV/c)");
1039                 hmultmin[i]->GetYaxis()->SetTitle("N_{ch}^{90,min}");
1040                 hmultmin[i]->SetMarkerStyle(20);
1041                 hmultmin[i]->SetMarkerColor(1+i);
1042                 hmultmin[i]->SetMaximum(3.);
1043                 hmultmin[i]->Divide(hjets[i]); // normalize for jet spectrum
1044             
1045
1046                 // Phi-correlation with leading jet --------------------------
1047                 hetaphi[i] = (TH2F*) ldata[i]->FindObject("hdNdEtaPhiDist");
1048                 hetaphi[i]->GetXaxis()->SetTitle("#Delta #phi (w.r.t. leading jet)");
1049                 hetaphi[i]->SetMarkerStyle(20);
1050             
1051
1052                 // pT distribution in full region vs. jet pT --------------------------
1053                 hptfull[i] = (TH2F*) ldata[i]->FindObject("hFullRegPartPtDistVsEt");
1054                 hptfull[i]->GetYaxis()->SetTitle("Leading Jet P_{T} (GeV/c)");
1055                 hptfull[i]->GetXaxis()->SetTitle("Track P_{T} (GeV/c)");
1056             
1057
1058                 // pT distribution in transv region vs. jet pT --------------------------
1059                 hpttransv[i] = (TH2F*) ldata[i]->FindObject("hTransRegPartPtDistVsEt");
1060                 hpttransv[i]->GetYaxis()->SetTitle("Leading Jet P_{T} (GeV/c)");
1061                 hpttransv[i]->GetXaxis()->SetTitle("Track P_{T} (GeV/c)");
1062             
1063
1064                 // Return Histos
1065                 arr->Add(hnjets[i]);      //at 0 number of jets
1066                 arr->Add(hjets[i]);     //at 1 leading jets
1067                 arr->Add(hmax[i]);      //at 2 sum pT MAX
1068                 arr->Add(hmin[i]);      //at 3 sumpT MIN
1069                 arr->Add(hmultmax[i]);   //at 4 multiplicity MAX 
1070                 arr->Add(hmultmin[i]);   //at 5 multiplicity MIN
1071                 arr->Add(hetaphi[i]);     //at 6 phi correlation
1072                 arr->Add(hptfull[i]);     //at 7 pT distr in full region
1073                 arr->Add(hpttransv[i]);   //at 8 pT distr in transv region
1074
1075         }
1076         return arr;
1077 }
1078
1079 //_______________________________________________________________________________________
1080 void AliHistogramsUE::SetStyle(){
1081
1082   // Set plotting style
1083   gPad->SetFrameFillColor(0);
1084   gPad->SetFillColor(0);
1085   gPad->SetBorderSize(2);
1086   gPad->SetGridy();
1087   gPad->SetFrameBorderMode(0);
1088   //gStyle->SetOptStat(0);
1089   gStyle->SetOptTitle(0);
1090
1091 }
1092
1093
1094 //____________________________________________________________________
1095 TList* AliHistogramsUE::GetHistograms(){
1096
1097   // Return list of relevant histograms 
1098   return fListOfHistos;
1099
1100 }
1101
1102 //____________________________________________________________________
1103 void AliHistogramsUE::PlotBranchesUE(TString file, TString branches, Double_t minJetProjection){
1104
1105   // Function to be called by external macro to plot analysis from different jet branches 
1106   
1107   //Count the branches 
1108   TObjArray  *arr=branches.Tokenize(";");
1109   TIter next(arr);
1110   TObject *o=0;
1111   Int_t br=0;
1112   while ( (o=next()) ){
1113     if(o)br++;
1114   }
1115   
1116
1117   //Canvas
1118   TObjArray *arrC=CreateCanvas(9);
1119  
1120   //Histograms 
1121   TObjArray *arrH=GetHistosForPlotting(file.Data(),branches.Data());
1122   TH1F *hnjets[20] = {0,};
1123   TH1F *hjets[20]  = {0,};
1124   TH1F *hmax[20]   = {0,};
1125   TH1F *hmin[20]   = {0,};
1126   TH1F *hmultmax[20]  = {0,};
1127   TH1F *hmultmin[20]  = {0,};
1128   TH2F *hetaphi[20]   = {0,};
1129   TH2F *hptfull[20]   = {0,};
1130   TH2F *hpttransv[20] = {0,};
1131
1132   for (Int_t i= 0; i<br; i++){
1133         hnjets[i]=((TH1F*)arrH->At(9*i));
1134         hjets[i]=((TH1F*)arrH->At(1+(9*i)));
1135         hmax[i]=((TH1F*)arrH->At(2+(9*i)));
1136         hmin[i]=((TH1F*)arrH->At(3+(9*i)));
1137         hmultmax[i]=((TH1F*)arrH->At(4+(9*i)));
1138         hmultmin[i]=((TH1F*)arrH->At(5+(9*i)));
1139         hetaphi[i]=((TH2F*)arrH->At(6+(9*i))); 
1140         hptfull[i]=((TH2F*)arrH->At(7+(9*i))); 
1141         hpttransv[i]=((TH2F*)arrH->At(8+(9*i))); 
1142         }
1143
1144   //Define jet-pT range in projections
1145   Int_t binstartprojection = hetaphi[0]->GetYaxis()->FindFixBin(minJetProjection); //be careful...  
1146   Int_t binstopprojection = hetaphi[0]->GetYaxis()->GetNbins();
1147
1148   Double_t normphi[20]; 
1149
1150   for (Int_t i= 0; i<br; i++){
1151         normphi[i]=hjets[i]->Integral(binstartprojection,binstopprojection);    
1152         hnjets[i]->Scale(1./(hnjets[i]->GetBinWidth(1)*hnjets[i]->GetEntries()));
1153         hjets[i]->Scale(1./(hjets[i]->GetBinWidth(1)*hnjets[i]->GetEntries()));
1154         }
1155
1156   //LEGENDS
1157   //Legend jets
1158   TLegend *leg=new TLegend(0.5,0.6,0.89,0.89);  // for jet spectrum
1159   leg->SetFillColor(0);
1160   leg->SetHeader("Jet Finders:");
1161   //Legend density
1162   TLegend *legd = new TLegend(0.1077586,0.6016949,0.4971264,0.8919492,NULL,"brNDC");
1163   legd->SetFillColor(0);
1164   legd->SetHeader("Jet Finders:");
1165   //Legend pT distributions
1166   TLegend *legp = new TLegend(0.1364943,0.1292373,0.5258621,0.4194915,NULL,"brNDC");
1167   legp->SetFillColor(0);
1168   legp->SetHeader("Jet Finders:");
1169
1170   arr=branches.Tokenize(";");
1171   TIter next1(arr);
1172   o=0;
1173   Int_t brleg=0;
1174   while ( (o=next1()) ){
1175         leg->AddEntry(hjets[brleg],Form("UE%s",o->GetName()),"p");
1176         legd->AddEntry(hjets[brleg],Form("UE%s",o->GetName()),"p");
1177         legp->AddEntry(hjets[brleg],Form("UE%s",o->GetName()),"p");
1178         brleg++;
1179         }
1180
1181   //1) NUMBER OF CLUSTERS 
1182   TCanvas *c0=((TCanvas*)arrC->At(0)); 
1183   c0->cd();
1184   SetStyle();
1185   gPad->SetLogy();
1186   hnjets[0]->SetMaximum(1.4);
1187   hnjets[0]->Draw("E1");
1188   for (Int_t i= 1; i<br; i++){
1189         hnjets[i]->Draw("same");
1190         }
1191   leg->Draw("same");
1192
1193   //2) LEADING CLUSTERS pT 
1194   TCanvas *c1=((TCanvas*)arrC->At(1)); 
1195   c1->cd();
1196   SetStyle();
1197   gPad->SetLogy();
1198   hjets[0]->Draw("E1");
1199   hjets[0]->SetMaximum(1.4);
1200   for (Int_t i= 1; i<br; i++){
1201         hjets[i]->Draw("same");
1202         }
1203   leg->Draw("same");
1204
1205   //3) SUM-pT IN MAX REGION
1206   TCanvas *c2=((TCanvas*)arrC->At(2)); 
1207   c2->cd();
1208   SetStyle();
1209   hmax[0]->Draw("E1");
1210   for (Int_t i= 1; i<br; i++){
1211         hmax[i]->Draw("same");
1212         }
1213   leg->Draw("same");
1214
1215   //4) SUM-pT IN MIN REGION
1216   TCanvas *c3=((TCanvas*)arrC->At(3)); 
1217   c3->cd();
1218   SetStyle();
1219   hmin[0]->GetYaxis()->SetRangeUser(0.,2.5);
1220   hmin[0]->Draw("E1");
1221   for (Int_t i= 1; i<br; i++){
1222         hmin[i]->Draw("same");
1223         }
1224   leg->Draw("same");
1225
1226   //5) MULTIPLICITY IN MAX REGION
1227   TCanvas *c4=((TCanvas*)arrC->At(4));
1228   c4->cd();
1229   SetStyle();
1230   hmultmax[0]->GetYaxis()->SetRangeUser(0.,5.);
1231   hmultmax[0]->Draw("E1");
1232   for (Int_t i= 1; i<br; i++){
1233         hmultmax[i]->Draw("same");
1234         }
1235   leg->Draw("same");
1236
1237   //6) MULTIPLICITY IN MIN REGION
1238   TCanvas *c5=((TCanvas*)arrC->At(5));
1239   c5->cd();
1240   SetStyle();
1241   hmultmin[0]->GetYaxis()->SetRangeUser(0.,2.5);
1242   hmultmin[0]->Draw("E1");
1243   for (Int_t i= 1; i<br; i++){
1244         hmultmin[i]->Draw("same");
1245         }
1246   leg->Draw("same");
1247
1248   //7) JET-TRACK CORRELATION
1249   TCanvas *c6=((TCanvas*)arrC->At(6));
1250   c6->cd();
1251   SetStyle();
1252   gPad->SetLogy();
1253   TH1D *tmpetaphi[20];
1254   for (Int_t i= 0; i<br; i++){
1255         tmpetaphi[i]=new TH1D();
1256         tmpetaphi[i]=hetaphi[i]->ProjectionX(Form("data%d",i),binstartprojection);
1257         tmpetaphi[i]->GetYaxis()->SetTitle("1/n_{lj} dN/d#Delta #phi");
1258         tmpetaphi[i]->Scale(1./(hetaphi[i]->GetBinWidth(1)*normphi[i]));
1259         tmpetaphi[i]->SetMarkerColor(i+1);
1260         tmpetaphi[i]->GetXaxis()->SetLimits(-3.*TMath::Pi()/2.,TMath::Pi()/2.);
1261         tmpetaphi[i]->GetYaxis()->SetRangeUser(0.5,20.);
1262         if (i==0) tmpetaphi[i]->Draw("E1");
1263         else tmpetaphi[i]->Draw("same");
1264         // evaluate mean multiplicity in transverse regions
1265         Double_t err1=0.;
1266         Double_t err2=0.;
1267         Double_t err3=0.;
1268         Double_t mean=(tmpetaphi[i]->IntegralAndError(1,5,err1)+tmpetaphi[i]->IntegralAndError(27,36,err2)+tmpetaphi[i]->IntegralAndError(58,62,err3))/(20.);
1269         err1=TMath::Sqrt(err1*err1+err2*err2+err3*err3)/20.;
1270         Printf("Branch: %d  MeanTransvMult: %f err: %f",i,mean,err1);
1271         }
1272   legd->Draw("same");
1273
1274   //8) TRACK-pT DISTRIBUTION IN FULL REGION
1275   TCanvas *c7=((TCanvas*)arrC->At(7));
1276   c7->cd();
1277   SetStyle();
1278   gPad->SetLogy();
1279   gPad->SetLogx();
1280   gPad->SetTitle("Full region (projection X)");
1281   TH1D *tmpfull[20];
1282   for (Int_t i= 0; i<br; i++){
1283         tmpfull[i]=new TH1D();
1284         tmpfull[i]=hptfull[i]->ProjectionX(Form("full%d",i),binstartprojection);
1285         tmpfull[i]->GetYaxis()->SetTitle("1/n_{lj} dN/dp_{T}");
1286         tmpfull[i]->Scale(1./(tmpfull[i]->GetBinWidth(1)*normphi[i]));
1287         tmpfull[i]->SetMarkerStyle(20);
1288         tmpfull[i]->SetMarkerColor(i+1);
1289         tmpfull[i]->SetMaximum(100.);
1290         if (i==0)tmpfull[i]->Draw("E1");
1291         else tmpfull[i]->Draw("same"); 
1292         }
1293   legp->Draw("same");
1294
1295   //9) TRACK-pT DISTRIBUTION IN TRANSVERSE (MIN+MAX) REGION
1296   TCanvas *c8=((TCanvas*)arrC->At(8));
1297   c8->cd();
1298   SetStyle();
1299   gPad->SetLogy();
1300   gPad->SetLogx();
1301   gPad->SetTitle("Transverse regions (projection X)");
1302   TH1D *tmptransv[20];
1303   for (Int_t i= 0; i<br; i++){
1304         tmptransv[i]=new TH1D();
1305         tmptransv[i]=hpttransv[i]->ProjectionX(Form("transv%d",i),binstartprojection);
1306         tmptransv[i]->GetYaxis()->SetTitle("1/n_{lj} dN/dp_{T}");
1307         tmptransv[i]->Scale(1./(tmptransv[i]->GetBinWidth(1)*normphi[i]));
1308         tmptransv[i]->SetMarkerStyle(20);
1309         tmptransv[i]->SetMarkerColor(i+1);
1310         tmptransv[i]->SetMaximum(10.);
1311         if (i==0)tmptransv[i]->Draw("E1");
1312         else tmptransv[i]->Draw("same"); 
1313         }
1314   legp->Draw("same");
1315 }