AliPhysics  27da6ee (27da6ee)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliFlowCommonHist.cxx
Go to the documentation of this file.
1 /**************************************************************************
2  * Copyright(c) 1998-1999, ALICE Experiment at CERN, All rights reserved. *
3  * *
4  * Author: The ALICE Off-line Project. *
5  * Contributors are mentioned in the code where appropriate. *
6  * *
7  * Permission to use, copy, modify and distribute this software and its *
8  * documentation strictly for non-commercial purposes is hereby granted *
9  * without fee, provided that the above copyright notice appears in all *
10  * copies and that both the copyright notice and this permission notice *
11  * appear in the supporting documentation. The authors make no claims *
12  * about the suitability of this software for any purpose. It is *
13  * provided "as is" without express or implied warranty. *
14  **************************************************************************/
15 
16 #include "Riostream.h" //needed as include
17 #include "AliFlowCommonConstants.h" //needed as include
18 #include "AliFlowCommonHist.h"
19 #include "AliFlowEventSimple.h"
20 #include "AliFlowTrackSimple.h"
21 
22 #include "TString.h"
23 #include "TProfile.h"
24 #include "TMath.h" //needed as include
25 #include "TList.h"
26 #include "TH2F.h"
27 #include "AliFlowVector.h"
28 #include "TBrowser.h"
29 
30 class TH1F;
31 class TH1D;
32 
33 // AliFlowCommonHist:
34 //
35 // Description: Class to organise common histograms for Flow Analysis
36 //
37 // authors: N. van der Kolk (kolk@nikhef.nl), A. Bilandzic (anteb@nikhef.nl), RS
38 
39 
40 using std::cout;
41 using std::endl;
43 
44 //-----------------------------------------------------------------------
45 
47  TNamed(),
48  fBookOnlyBasic(kFALSE),
49  fHistMultRP(NULL),
50  fHistMultPOI(NULL),
51  fHistMultPOIvsRP(NULL),
52  fHistPtRP(NULL),
53  fHistPtPOI(NULL),
54  fHistPtSub0(NULL),
55  fHistPtSub1(NULL),
56  fHistPhiRP(NULL),
57  fHistPhiPOI(NULL),
58  fHistPhiSub0(NULL),
59  fHistPhiSub1(NULL),
60  fHistEtaRP(NULL),
61  fHistEtaPOI(NULL),
62  fHistEtaSub0(NULL),
63  fHistEtaSub1(NULL),
64  fHistPhiEtaRP(NULL),
65  fHistPhiEtaPOI(NULL),
66  fHistProMeanPtperBin(NULL),
67  fHistWeightvsPhi(NULL),
68  fHistQ(NULL),
69  fHistAngleQ(NULL),
70  fHistAngleQSub0(NULL),
71  fHistAngleQSub1(NULL),
72  fHarmonic(NULL),
73  fRefMultVsNoOfRPs(NULL),
74  fHistRefMult(NULL),
75  fHistMassPOI(NULL),
76  fHistList(NULL)
77 {
78 
79  //default constructor
80 
81 }
82 
84  TNamed(),
85  fBookOnlyBasic(a.fBookOnlyBasic),
86  fHistMultRP(new TH1F(*a.fHistMultRP)),
87  fHistMultPOI(new TH1F(*a.fHistMultPOI)),
88  fHistMultPOIvsRP(new TH2F(*a.fHistMultPOIvsRP)),
89  fHistPtRP(new TH1F(*a.fHistPtRP)),
90  fHistPtPOI(new TH1F(*a.fHistPtPOI)),
91  fHistPtSub0(new TH1F(*a.fHistPtSub0)),
92  fHistPtSub1(new TH1F(*a.fHistPtSub1)),
93  fHistPhiRP(new TH1F(*a.fHistPhiRP)),
94  fHistPhiPOI(new TH1F(*a.fHistPhiPOI)),
95  fHistPhiSub0(new TH1F(*a.fHistPhiSub0)),
96  fHistPhiSub1(new TH1F(*a.fHistPhiSub1)),
97  fHistEtaRP(new TH1F(*a.fHistEtaRP)),
98  fHistEtaPOI(new TH1F(*a.fHistEtaPOI)),
99  fHistEtaSub0(new TH1F(*a.fHistEtaSub0)),
100  fHistEtaSub1(new TH1F(*a.fHistEtaSub1)),
101  fHistPhiEtaRP(new TH2F(*a.fHistPhiEtaRP)),
102  fHistPhiEtaPOI(new TH2F(*a.fHistPhiEtaPOI)),
103  fHistProMeanPtperBin(new TProfile(*a.fHistProMeanPtperBin)),
104  fHistWeightvsPhi(new TH2F(*a.fHistWeightvsPhi)),
105  fHistQ(new TH1F(*a.fHistQ)),
106  fHistAngleQ(new TH1F(*a.fHistAngleQ)),
107  fHistAngleQSub0(new TH1F(*a.fHistAngleQSub0)),
108  fHistAngleQSub1(new TH1F(*a.fHistAngleQSub1)),
109  fHarmonic(new TProfile(*a.fHarmonic)),
110  fRefMultVsNoOfRPs(new TProfile(*a.fRefMultVsNoOfRPs)),
111  fHistRefMult(new TH1F(*a.fHistRefMult)),
112  fHistMassPOI(new TH2F(*a.fHistMassPOI)),
113  fHistList(NULL)
114 {
115  // copy constructor
116 
117  fHistList = new TList();
118  fHistList-> Add(fHistMultRP);
119  fHistList-> Add(fHistMultPOI);
121  fHistList-> Add(fHistPtRP);
122  fHistList-> Add(fHistPtPOI);
123  if(!fBookOnlyBasic){fHistList-> Add(fHistPtSub0);}
124  if(!fBookOnlyBasic){fHistList-> Add(fHistPtSub1);}
125  fHistList-> Add(fHistPhiRP);
126  fHistList-> Add(fHistPhiPOI);
128  if(!fBookOnlyBasic){fHistList-> Add(fHistPhiSub1);}
129  fHistList-> Add(fHistEtaRP);
130  fHistList-> Add(fHistEtaPOI);
137  fHistList-> Add(fHarmonic);
139  fHistList-> Add(fHistRefMult);
140  fHistList-> Add(fHistMassPOI);
141  if(!fBookOnlyBasic){fHistList-> Add(fHistQ);}
142  if(!fBookOnlyBasic){fHistList-> Add(fHistAngleQ);}
145  // TListIter next = TListIter(a.fHistList);
146 
147 }
148 
149 
150 //-----------------------------------------------------------------------
151 
152  AliFlowCommonHist::AliFlowCommonHist(const char *anInput, const char *title, Bool_t bookOnlyBasic):
153  TNamed(anInput,title),
154  fBookOnlyBasic(bookOnlyBasic),
155  fHistMultRP(NULL),
156  fHistMultPOI(NULL),
157  fHistMultPOIvsRP(NULL),
158  fHistPtRP(NULL),
159  fHistPtPOI(NULL),
160  fHistPtSub0(NULL),
161  fHistPtSub1(NULL),
162  fHistPhiRP(NULL),
163  fHistPhiPOI(NULL),
164  fHistPhiSub0(NULL),
165  fHistPhiSub1(NULL),
166  fHistEtaRP(NULL),
167  fHistEtaPOI(NULL),
168  fHistEtaSub0(NULL),
169  fHistEtaSub1(NULL),
170  fHistPhiEtaRP(NULL),
171  fHistPhiEtaPOI(NULL),
172  fHistProMeanPtperBin(NULL),
173  fHistWeightvsPhi(NULL),
174  fHistQ(NULL),
175  fHistAngleQ(NULL),
176  fHistAngleQSub0(NULL),
177  fHistAngleQSub1(NULL),
178  fHarmonic(NULL),
179  fRefMultVsNoOfRPs(NULL),
180  fHistRefMult(NULL),
181  fHistMassPOI(NULL),
182  fHistList(NULL)
183 {
184 
185  //constructor creating histograms
186  Int_t iNbinsMult = AliFlowCommonConstants::GetMaster()->GetNbinsMult();
187  Int_t iNbinsPt = AliFlowCommonConstants::GetMaster()->GetNbinsPt();
188  Int_t iNbinsPhi = AliFlowCommonConstants::GetMaster()->GetNbinsPhi();
189  Int_t iNbinsEta = AliFlowCommonConstants::GetMaster()->GetNbinsEta();
190  Int_t iNbinsQ = AliFlowCommonConstants::GetMaster()->GetNbinsQ();
191  Int_t iNbinsMass = AliFlowCommonConstants::GetMaster()->GetNbinsMass();
192  TString sName;
193 
194  Double_t dMultMin = AliFlowCommonConstants::GetMaster()->GetMultMin();
195  Double_t dMultMax = AliFlowCommonConstants::GetMaster()->GetMultMax();
196  Double_t dPtMin = AliFlowCommonConstants::GetMaster()->GetPtMin();
197  Double_t dPtMax = AliFlowCommonConstants::GetMaster()->GetPtMax();
198  Double_t dPhiMin = AliFlowCommonConstants::GetMaster()->GetPhiMin();
199  Double_t dPhiMax = AliFlowCommonConstants::GetMaster()->GetPhiMax();
200  Double_t dEtaMin = AliFlowCommonConstants::GetMaster()->GetEtaMin();
201  Double_t dEtaMax = AliFlowCommonConstants::GetMaster()->GetEtaMax();
202  Double_t dQMin = AliFlowCommonConstants::GetMaster()->GetQMin();
203  Double_t dQMax = AliFlowCommonConstants::GetMaster()->GetQMax();
204  Double_t dHistWeightvsPhiMin = AliFlowCommonConstants::GetMaster()->GetHistWeightvsPhiMin();
205  Double_t dHistWeightvsPhiMax = AliFlowCommonConstants::GetMaster()->GetHistWeightvsPhiMax();
206  Double_t dMassMin = AliFlowCommonConstants::GetMaster()->GetMassMin();
207  Double_t dMassMax = AliFlowCommonConstants::GetMaster()->GetMassMax();
208 
209 
210  // cout<<"The settings for the common histograms are as follows:"<<endl;
211  // cout<<"Multiplicity: "<<iNbinsMult<<" bins between "<<dMultMin<<" and "<<dMultMax<<endl;
212  // cout<<"Pt: "<<iNbinsPt<<" bins between "<<dPtMin<<" and "<<dPtMax<<endl;
213  // cout<<"Phi: "<<iNbinsPhi<<" bins between "<<dPhiMin<<" and "<<dPhiMax<<endl;
214  // cout<<"Eta: "<<iNbinsEta<<" bins between "<<dEtaMin<<" and "<<dEtaMax<<endl;
215  // cout<<"Q: "<<iNbinsQ<<" bins between "<<dQMin<<" and "<<dQMax<<endl;
216  // cout<<"Mass: "<<iNbinsMass<<" bins between "<<dMassMin<<" and "<<dMassMax<<endl;
217 
218  //Multiplicity
219  sName = "Control_Flow_MultRP_";
220  sName +=anInput;
221  fHistMultRP = new TH1F(sName.Data(), sName.Data(),iNbinsMult, dMultMin, dMultMax);
222  fHistMultRP ->SetXTitle("Multiplicity for RP selection");
223  fHistMultRP ->SetYTitle("Counts");
224 
225  sName = "Control_Flow_MultPOI_";
226  sName +=anInput;
227  fHistMultPOI = new TH1F(sName.Data(), sName.Data(),iNbinsMult, dMultMin, dMultMax);
228  fHistMultPOI ->SetXTitle("Multiplicity for POI selection");
229  fHistMultPOI ->SetYTitle("Counts");
230 
231  if(!fBookOnlyBasic){
232  sName = "Control_Flow_MultPOIvsRP_";
233  sName +=anInput;
234  fHistMultPOIvsRP = new TH2F(sName.Data(), sName.Data(),iNbinsMult, dMultMin, dMultMax,100, dMultMin, dMultMax);
235  fHistMultPOIvsRP->SetXTitle("Multiplicity for RP selection");
236  fHistMultPOIvsRP->SetYTitle("Multiplicity for POI selection");
237  }
238 
239  //Pt
240  sName = "Control_Flow_PtRP_";
241  sName +=anInput;
242  fHistPtRP = new TH1F(sName.Data(), sName.Data(),iNbinsPt, dPtMin, dPtMax);
243  fHistPtRP ->SetXTitle("P_{t} (GeV/c) for RP selection");
244  fHistPtRP ->SetYTitle("Counts");
245 
246  sName = "Control_Flow_PtPOI_";
247  sName +=anInput;
248  fHistPtPOI = new TH1F(sName.Data(), sName.Data(),iNbinsPt, dPtMin, dPtMax);
249  //binning has to be the same as for fHistProVPt! use to get Nprime!
250  fHistPtPOI ->SetXTitle("P_{t} (GeV/c) for POI selection");
251  fHistPtPOI ->SetYTitle("Counts");
252 
253  if(!fBookOnlyBasic){
254  sName = "Control_Flow_PtSub0_";
255  sName +=anInput;
256  fHistPtSub0 = new TH1F(sName.Data(), sName.Data(),iNbinsPt, dPtMin, dPtMax);
257  fHistPtSub0 ->SetXTitle("P_{t} (GeV/c) for Subevent 0 selection");
258  fHistPtSub0 ->SetYTitle("Counts");
259  }
260 
261  if(!fBookOnlyBasic){
262  sName = "Control_Flow_PtSub1_";
263  sName +=anInput;
264  fHistPtSub1 = new TH1F(sName.Data(), sName.Data(),iNbinsPt, dPtMin, dPtMax);
265  fHistPtSub1 ->SetXTitle("P_{t} (GeV/c) for Subevent 1 selection");
266  fHistPtSub1 ->SetYTitle("Counts");
267  }
268 
269  //Phi
270  sName = "Control_Flow_PhiRP_";
271  sName +=anInput;
272  fHistPhiRP = new TH1F(sName.Data(), sName.Data(),iNbinsPhi, dPhiMin, dPhiMax);
273  fHistPhiRP ->SetXTitle("#phi for RP selection");
274  fHistPhiRP ->SetYTitle("Counts");
275 
276  sName = "Control_Flow_PhiPOI_";
277  sName +=anInput;
278  fHistPhiPOI = new TH1F(sName.Data(), sName.Data(),iNbinsPhi, dPhiMin, dPhiMax);
279  fHistPhiPOI ->SetXTitle("#phi for POI selection");
280  fHistPhiPOI ->SetYTitle("Counts");
281 
282  if(!fBookOnlyBasic){
283  sName = "Control_Flow_PhiSub0_";
284  sName +=anInput;
285  fHistPhiSub0 = new TH1F(sName.Data(), sName.Data(),iNbinsPhi, dPhiMin, dPhiMax);
286  fHistPhiSub0 ->SetXTitle("#phi for Subevent 0 selection");
287  fHistPhiSub0 ->SetYTitle("Counts");
288  }
289 
290  if(!fBookOnlyBasic){
291  sName = "Control_Flow_PhiSub1_";
292  sName +=anInput;
293  fHistPhiSub1 = new TH1F(sName.Data(), sName.Data(),iNbinsPhi, dPhiMin, dPhiMax);
294  fHistPhiSub1 ->SetXTitle("#phi for Subevent 1 selection");
295  fHistPhiSub1 ->SetYTitle("Counts");
296  }
297 
298  //Eta
299  sName = "Control_Flow_EtaRP_";
300  sName +=anInput;
301  fHistEtaRP = new TH1F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax);
302  fHistEtaRP ->SetXTitle("#eta for RP selection");
303  fHistEtaRP ->SetYTitle("Counts");
304 
305  sName = "Control_Flow_EtaPOI_";
306  sName +=anInput;
307  fHistEtaPOI = new TH1F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax);
308  fHistEtaPOI ->SetXTitle("#eta for POI selection");
309  fHistEtaPOI ->SetYTitle("Counts");
310 
311  if(!fBookOnlyBasic){
312  sName = "Control_Flow_EtaSub0_";
313  sName +=anInput;
314  fHistEtaSub0 = new TH1F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax);
315  fHistEtaSub0 ->SetXTitle("#eta for Subevent 0 selection");
316  fHistEtaSub0 ->SetYTitle("Counts");
317  }
318 
319  if(!fBookOnlyBasic){
320  sName = "Control_Flow_EtaSub1_";
321  sName +=anInput;
322  fHistEtaSub1 = new TH1F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax);
323  fHistEtaSub1 ->SetXTitle("#eta for Subevent 1 selection");
324  fHistEtaSub1 ->SetYTitle("Counts");
325  }
326 
327  if(!fBookOnlyBasic){
328  //Phi vs Eta
329  sName = "Control_Flow_PhiEtaRP_";
330  sName +=anInput;
331  fHistPhiEtaRP = new TH2F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax, iNbinsPhi, dPhiMin, dPhiMax);
332  fHistPhiEtaRP ->SetXTitle("#eta");
333  fHistPhiEtaRP ->SetYTitle("#phi");
334  }
335 
336  if(!fBookOnlyBasic){
337  sName = "Control_Flow_PhiEtaPOI_";
338  sName +=anInput;
339  fHistPhiEtaPOI = new TH2F(sName.Data(), sName.Data(),iNbinsEta, dEtaMin, dEtaMax, iNbinsPhi, dPhiMin, dPhiMax);
340  fHistPhiEtaPOI ->SetXTitle("#eta");
341  fHistPhiEtaPOI ->SetYTitle("#phi");
342  }
343 
344  //Mean Pt per pt bin
345  sName = "Control_FlowPro_MeanPtperBin_";
346  sName +=anInput;
347  fHistProMeanPtperBin = new TProfile(sName.Data(), sName.Data(),iNbinsPt,dPtMin,dPtMax);
348  fHistProMeanPtperBin ->SetXTitle("P_{t} (GeV/c) ");
349  fHistProMeanPtperBin ->SetYTitle("<P_{t}> (GeV/c) ");
350 
351 
352  if(!fBookOnlyBasic){
353  //Particle weight
354  sName = "Control_Flow_WeightvsPhi_";
355  sName +=anInput;
356  fHistWeightvsPhi = new TH2F(sName.Data(), sName.Data(), iNbinsPhi, dPhiMin, dPhiMax, 300, dHistWeightvsPhiMin, dHistWeightvsPhiMax);
357  fHistWeightvsPhi ->SetXTitle("#phi");
358  fHistWeightvsPhi ->SetYTitle("weight");
359  }
360 
361  if(!fBookOnlyBasic){
362  //Q vector
363  sName = "Control_Flow_Q_";
364  sName +=anInput;
365  fHistQ = new TH1F(sName.Data(), sName.Data(),iNbinsQ, dQMin, dQMax);
366  fHistQ ->SetXTitle("Q_{vector}/Mult");
367  fHistQ ->SetYTitle("Counts");
368  }
369 
370  if(!fBookOnlyBasic){
371  //Angle of Q vector
372  sName = "Control_Flow_AngleQ_";
373  sName +=anInput;
374  fHistAngleQ = new TH1F(sName.Data(), sName.Data(),72, 0., TMath::Pi());
375  fHistAngleQ ->SetXTitle("Angle of Q_{vector}");
376  fHistAngleQ ->SetYTitle("Counts");
377  }
378 
379  if(!fBookOnlyBasic){
380  sName = "Control_Flow_AngleQSub0_";
381  sName +=anInput;
382  fHistAngleQSub0 = new TH1F(sName.Data(), sName.Data(),72, 0., TMath::Pi());
383  fHistAngleQSub0 ->SetXTitle("Angle of Q_{vector} for Subevent 0");
384  fHistAngleQSub0 ->SetYTitle("Counts");
385  }
386 
387  if(!fBookOnlyBasic){
388  sName = "Control_Flow_AngleQSub1_";
389  sName +=anInput;
390  fHistAngleQSub1 = new TH1F(sName.Data(), sName.Data(),72, 0., TMath::Pi());
391  fHistAngleQSub1 ->SetXTitle("Angle of Q_{vector} for Subevent 1");
392  fHistAngleQSub1 ->SetYTitle("Counts");
393  }
394 
395  //harmonic
396  sName = "Control_Flow_Harmonic_";
397  sName +=anInput;
398  fHarmonic = new TProfile(sName.Data(),sName.Data(),1,0,1);
399  fHarmonic ->SetYTitle("harmonic");
400 
401  //<reference multiplicity> versus # of RPs
402  sName = "Reference_Multiplicity_Vs_Number_Of_RPs_";
403  sName +=anInput;
404  fRefMultVsNoOfRPs = new TProfile(sName.Data(),sName.Data(),iNbinsMult, dMultMin, dMultMax);
405  fRefMultVsNoOfRPs->SetYTitle("<reference multiplicity>");
406  fRefMultVsNoOfRPs->SetXTitle("# of RPs");
407 
408  //reference multiplicity
409  sName = "Control_Flow_Ref_Mult_";
410  sName +=anInput;
411  fHistRefMult = new TH1F(sName.Data(), sName.Data(),iNbinsMult, dMultMin, dMultMax);
412  fHistRefMult->SetXTitle("Reference multiplicity");
413  fHistRefMult->SetYTitle("Counts");
414 
415  //mass for POI selection
416  sName = "Control_Flow_Mass_POI";
417  sName +=anInput;
418  fHistMassPOI = new TH2F(sName.Data(), sName.Data(), iNbinsMass, dMassMin, dMassMax,
419  iNbinsPt, dPtMin, dPtMax);
420  Double_t mbWidth = 0.;
421  if(iNbinsMass != 0)
422  {
423  mbWidth = (dMassMax-dMassMin)/iNbinsMass*1000.;
424  }
425  fHistMassPOI->SetXTitle("Mass (GeV/c^{2})");
426  fHistMassPOI->SetYTitle("P_{t} (GeV/c)");
427  fHistMassPOI->SetZTitle( Form("Counts/(%.2f MeV/c^{2})",mbWidth));
428 
429  //list of histograms if added here also add in copy constructor
430  fHistList = new TList();
431  fHistList-> Add(fHistMultRP);
432  fHistList-> Add(fHistMultPOI);
434  fHistList-> Add(fHistPtRP);
435  fHistList-> Add(fHistPtPOI);
436  if(!fBookOnlyBasic){fHistList-> Add(fHistPtSub0);}
437  if(!fBookOnlyBasic){fHistList-> Add(fHistPtSub1);}
438  fHistList-> Add(fHistPhiRP);
439  fHistList-> Add(fHistPhiPOI);
442  fHistList-> Add(fHistEtaRP);
443  fHistList-> Add(fHistEtaPOI);
444  if(!fBookOnlyBasic){fHistList-> Add(fHistEtaSub0);}
446  if(!fBookOnlyBasic){fHistList-> Add(fHistPhiEtaRP);}
450  fHistList-> Add(fHarmonic);
452  fHistList-> Add(fHistRefMult);
453  fHistList-> Add(fHistMassPOI);
454  if(!fBookOnlyBasic){fHistList-> Add(fHistQ);}
455  if(!fBookOnlyBasic){fHistList-> Add(fHistAngleQ);}
458 
459 }
460 
461 
462 //-----------------------------------------------------------------------
463 
465 {
466  //deletes histograms
467  delete fHistMultRP;
468  delete fHistMultPOI;
469  delete fHistMultPOIvsRP;
470  delete fHistPtRP;
471  delete fHistPtPOI;
472  delete fHistPtSub0;
473  delete fHistPtSub1;
474  delete fHistPhiRP;
475  delete fHistPhiPOI;
476  delete fHistPhiSub0;
477  delete fHistPhiSub1;
478  delete fHistEtaRP;
479  delete fHistEtaPOI;
480  delete fHistEtaSub0;
481  delete fHistEtaSub1;
482  delete fHistPhiEtaRP;
483  delete fHistPhiEtaPOI;
484  delete fHistProMeanPtperBin;
485  delete fHistWeightvsPhi;
486  delete fHistQ;
487  delete fHistAngleQ;
488  delete fHistAngleQSub0;
489  delete fHistAngleQSub1;
490  delete fHarmonic;
491  delete fRefMultVsNoOfRPs;
492  delete fHistRefMult;
493  delete fHistMassPOI;
494  delete fHistList;
495 }
496 
497 
498 //-----------------------------------------------------------------------
499 
500 Bool_t AliFlowCommonHist::FillControlHistograms(AliFlowEventSimple* anEvent,TList *weightsList, Bool_t usePhiWeights, Bool_t usePtWeights, Bool_t useEtaWeights)
501 {
502  //Fills the control histograms
503  if (!anEvent){
504  cout<<"##### FillControlHistograms: FlowEvent pointer null"<<endl;
505  return kFALSE;
506  }
507 
508  //track datamembers
509  Double_t dPt = 0.;
510  Double_t dPhi = 0.;
511  Double_t dEta = 0.;
512  Double_t dWeight = 1.;
513 
514  //weights used for corrections
515  Double_t dWPhi = 1.;
516  Double_t dWPt = 1.;
517  Double_t dWEta = 1.;
518 
519  TH1F *phiWeights = NULL;
520  TH1F *phiWeightsSub0 = NULL;
521  TH1F *phiWeightsSub1 = NULL;
522  TH1D *ptWeights = NULL;
523  TH1D *etaWeights = NULL;
524 
525  Int_t nBinsPhi = 0;
526  Int_t nBinsPhiSub0 = 0;
527  Int_t nBinsPhiSub1 = 0;
528  Double_t dBinWidthPt = 0.;
529  Double_t dPtMin = 0.;
530  Double_t dBinWidthEta = 0.;
531  Double_t dEtaMin = 0.;
532 
533  if(weightsList)
534  {
535  if(usePhiWeights)
536  {
537  phiWeights = dynamic_cast<TH1F *>(weightsList->FindObject("phi_weights"));
538  if(phiWeights) nBinsPhi = phiWeights->GetNbinsX();
539  phiWeightsSub0 = dynamic_cast<TH1F *>(weightsList->FindObject("phi_weights_sub0"));
540  if(phiWeightsSub0) nBinsPhiSub0 = phiWeightsSub0->GetNbinsX();
541  phiWeightsSub1 = dynamic_cast<TH1F *>(weightsList->FindObject("phi_weights_sub1"));
542  if(phiWeightsSub1) nBinsPhiSub1 = phiWeightsSub1->GetNbinsX();
543  }
544  if(usePtWeights)
545  {
546  ptWeights = dynamic_cast<TH1D *>(weightsList->FindObject("pt_weights"));
547  if(ptWeights)
548  {
549  dBinWidthPt = ptWeights->GetBinWidth(1); // assuming that all bins have the same width
550  dPtMin = (ptWeights->GetXaxis())->GetXmin();
551  }
552  }
553  if(useEtaWeights)
554  {
555  etaWeights = dynamic_cast<TH1D *>(weightsList->FindObject("eta_weights"));
556  if(etaWeights)
557  {
558  dBinWidthEta = etaWeights->GetBinWidth(1); // assuming that all bins have the same width
559  dEtaMin = (etaWeights->GetXaxis())->GetXmin();
560  }
561  }
562  } // end of if(weightsList)
563 
564 
565 
566  //fill the histograms
567  AliFlowVector vQ = anEvent->GetQ(2, weightsList, usePhiWeights, usePtWeights, useEtaWeights);
568  //weight by the Multiplicity
569  Double_t dQX = 0.;
570  Double_t dQY = 0.;
571  if (vQ.GetMult()!=0) {
572  dQX = vQ.X()/vQ.GetMult();
573  dQY = vQ.Y()/vQ.GetMult();
574  }
575  vQ.Set(dQX,dQY);
576  if(!fBookOnlyBasic){fHistQ->Fill(vQ.Mod());}
577  if(!fBookOnlyBasic){fHistAngleQ->Fill(vQ.Phi()/2);}
578 
579  AliFlowVector* vQSub = new AliFlowVector[2];
580  anEvent->Get2Qsub(vQSub, 2, weightsList, usePhiWeights, usePtWeights, useEtaWeights);
581  AliFlowVector vQa = vQSub[0];
582  AliFlowVector vQb = vQSub[1];
583  if(!fBookOnlyBasic){fHistAngleQSub0->Fill(vQa.Phi()/2);}
584  if(!fBookOnlyBasic){fHistAngleQSub1->Fill(vQb.Phi()/2);}
585 
586  Double_t dMultRP = 0.;
587  Double_t dMultPOI = 0.;
588 
589  Int_t iNumberOfTracks = anEvent->NumberOfTracks();
590  AliFlowTrackSimple* pTrack = NULL;
591 
592  for (Int_t i=0;i<iNumberOfTracks;i++) {
593  pTrack = anEvent->GetTrack(i);
594  if (pTrack ) {
595  dWeight = pTrack->Weight();
596  dPt = pTrack->Pt();
597  dPhi = pTrack->Phi();
598  if (dPhi<0.) dPhi+=2*TMath::Pi();
599  dEta = pTrack->Eta();
600 
601  //weights are only used for the RP selection
602  if (pTrack->InRPSelection()){
603  // determine Phi weight:
604  if(phiWeights && nBinsPhi) {
605  dWPhi = phiWeights->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhi/TMath::TwoPi())));
606  }
607  // determine v'(pt) weight:
608  if(ptWeights && dBinWidthPt) {
609  dWPt=ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-dPtMin)/dBinWidthPt)));
610  }
611  // determine v'(eta) weight:
612  if(etaWeights && dBinWidthEta) {
613  dWEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-dEtaMin)/dBinWidthEta)));
614  }
615 
616  //the total weight is the product
617  Double_t dW = dWeight*dWPhi*dWPt*dWEta;
618 
619  //pt
620  fHistPtRP->Fill(dPt,dW);
621  //phi
622  fHistPhiRP->Fill(dPhi,dW);
623  //eta
624  fHistEtaRP->Fill(dEta,dW);
625  //eta vs phi
626  if(!fBookOnlyBasic){fHistPhiEtaRP->Fill(dEta,dPhi,dW);}
627  //weight vs phi
628  if(!fBookOnlyBasic){fHistWeightvsPhi->Fill(dPhi,dW);}
629  //count
630  dMultRP += dW;
631  }
632  if (pTrack->InRPSelection() && pTrack->InSubevent(0)) {
633  // determine Phi weight:
634  if(phiWeightsSub0 && nBinsPhiSub0){
635  dWPhi = phiWeightsSub0->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhiSub0/TMath::TwoPi())));
636  }
637  // determine v'(pt) weight:
638  if(ptWeights && dBinWidthPt) {
639  dWPt=ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-dPtMin)/dBinWidthPt)));
640  }
641  // determine v'(eta) weight:
642  if(etaWeights && dBinWidthEta) {
643  dWEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-dEtaMin)/dBinWidthEta)));
644  }
645 
646  //the total weight is the product
647  Double_t dW = dWeight*dWPhi*dWPt*dWEta;
648 
649  //pt
650  if(!fBookOnlyBasic){fHistPtSub0 ->Fill(dPt,dW);}
651  //phi
652  if(!fBookOnlyBasic){fHistPhiSub0 ->Fill(dPhi,dW);}
653  //eta
654  if(!fBookOnlyBasic){fHistEtaSub0 ->Fill(dEta,dW);}
655  }
656  if (pTrack->InRPSelection() && pTrack->InSubevent(1)) {
657  // determine Phi weight:
658  if(phiWeightsSub1 && nBinsPhiSub1){
659  dWPhi = phiWeightsSub1->GetBinContent(1+(Int_t)(TMath::Floor(dPhi*nBinsPhiSub1/TMath::TwoPi())));
660  }
661  // determine v'(pt) weight:
662  if(ptWeights && dBinWidthPt) {
663  dWPt=ptWeights->GetBinContent(1+(Int_t)(TMath::Floor((dPt-dPtMin)/dBinWidthPt)));
664  }
665  // determine v'(eta) weight:
666  if(etaWeights && dBinWidthEta) {
667  dWEta=etaWeights->GetBinContent(1+(Int_t)(TMath::Floor((dEta-dEtaMin)/dBinWidthEta)));
668  }
669 
670  //the total weight is the product
671  Double_t dW = dWeight*dWPhi*dWPt*dWEta;
672 
673  //pt
674  if(!fBookOnlyBasic){fHistPtSub1 -> Fill(dPt,dW);}
675  //phi
676  if(!fBookOnlyBasic){fHistPhiSub1 -> Fill(dPhi,dW);}
677  //eta
678  if(!fBookOnlyBasic){fHistEtaSub1 -> Fill(dEta,dW);}
679  }
680  if (pTrack->InPOISelection()){
681 
682  Double_t dW = dWeight; //no pt, phi or eta weights
683 
684  //pt
685  fHistPtPOI ->Fill(dPt,dW);
686  //phi
687  fHistPhiPOI ->Fill(dPhi,dW);
688  //eta
689  fHistEtaPOI ->Fill(dEta,dW);
690  //eta vs phi
691  if(!fBookOnlyBasic){fHistPhiEtaPOI->Fill(dEta,dPhi,dW);}
692  //mean pt
693  fHistProMeanPtperBin ->Fill(dPt,dPt,dW);
694  //mass
695  fHistMassPOI->Fill(pTrack->Mass(),dPt,dW);
696  //count
697  dMultPOI += dW;
698  }
699  } //track
700  } //loop over tracks
701 
702  fHistMultRP->Fill(dMultRP);
703  fHistMultPOI->Fill(dMultPOI);
704  if(!fBookOnlyBasic){fHistMultPOIvsRP->Fill(dMultRP,dMultPOI);}
705 
706  //<reference multiplicity> versus # of RPs:
707  fRefMultVsNoOfRPs->Fill(dMultRP+0.5,anEvent->GetReferenceMultiplicity(),1.);
708 
709  //reference multiplicity:
710  fHistRefMult->Fill(anEvent->GetReferenceMultiplicity());
711 
712  return kTRUE;
713 }
714 
715 //-----------------------------------------------------------------------
716 
718 {
719  //get entries in bin aBin from fHistPtRP
720  Double_t dEntries = fHistPtRP->GetBinContent(aBin);
721 
722  return dEntries;
723 
724 }
725 
726 //-----------------------------------------------------------------------
727 
729 {
730  //get entries in bin aBin from fHistPtPOI
731  Double_t dEntries = fHistPtPOI->GetBinContent(aBin);
732 
733  return dEntries;
734 
735 }
736 
737 //-----------------------------------------------------------------------
738 
740 {
741  //get entries in bin aBin from fHistPtRP
742  Double_t dEntries = fHistEtaRP->GetBinContent(aBin);
743 
744  return dEntries;
745 
746 }
747 
748 //-----------------------------------------------------------------------
749 
751 {
752  //get entries in bin aBin from fHistPtPOI
753  Double_t dEntries = fHistEtaPOI->GetBinContent(aBin);
754 
755  return dEntries;
756 
757 }
758 
759 //-----------------------------------------------------------------------
760 
761 Double_t AliFlowCommonHist::GetMeanPt(Int_t aBin)
762 {
763  //Get entry from bin aBin from fHistProMeanPtperBin
764  Double_t dMeanPt = fHistProMeanPtperBin->GetBinContent(aBin);
765 
766  return dMeanPt;
767 
768 }
769 
770 
771 //-----------------------------------------------------------------------
772  Double_t AliFlowCommonHist::Merge(TCollection *aList)
773 {
774  //merge fuction
775  //cout<<"entering merge function"<<endl;
776  if (!aList) return 0;
777  if (aList->IsEmpty()) return 0; //no merging is needed
778 
779  Int_t iCount = 0;
780  TIter next(aList); // list is supposed to contain only objects of the same type as this
781  AliFlowCommonHist *toMerge=NULL;
782  while ((toMerge=(AliFlowCommonHist*)next())) {
783  TList tempList;
784  tempList.Add(toMerge->GetHistList());
785  fHistList->Merge(&tempList);
786  iCount++;
787  }
788 
789  //cout<<"Merged"<<endl;
790  return (double)iCount;
791 
792 }
793 
794 void AliFlowCommonHist::Print(Option_t *option) const
795 {
796  // -*-*-*-*-*Print some global quantities for this histogram collection class *-*-*-*-*-*-*-*
797  // ===============================================
798  // printf( "TH1.Print Name = %s, Entries= %d, Total sum= %g\n",GetName(),Int_t(fEntries),GetSumOfWeights());
799  printf( "Class.Print Name = %s, Histogram list:\n",GetName());
800 
801  if (fHistList) {
802  fHistList->Print(option);
803  }
804  else
805  {
806  printf( "Empty histogram list \n");
807  }
808 }
809 
810 //-----------------------------------------------------------------------
811  void AliFlowCommonHist::Browse(TBrowser *b)
812 {
813 
814  if (!b) return;
815  if (fHistList) b->Add(fHistList,"AliFlowCommonHistList");
816 }
817 
818 
819 
820 
Double_t GetHistWeightvsPhiMin() const
virtual AliFlowVector GetQ(Int_t n=2, TList *weightsList=NULL, Bool_t usePhiWeights=kFALSE, Bool_t usePtWeights=kFALSE, Bool_t useEtaWeights=kFALSE)
Double_t GetEntriesInEtaBinPOI(Int_t iBin)
const char * title
Definition: MakeQAPdf.C:26
Double_t GetMeanPt(Int_t iBin)
AliFlowTrackSimple * GetTrack(Int_t i)
Bool_t InSubevent(Int_t i) const
Double_t GetEntriesInPtBinPOI(Int_t iBin)
Double_t GetMult() const
Definition: AliFlowVector.h:46
virtual Double_t Merge(TCollection *aList)
ClassImp(AliFlowCommonHist) AliFlowCommonHist
Double_t Mass() const
void Browse(TBrowser *b)
Double_t GetHistWeightvsPhiMax() const
Bool_t InRPSelection() const
Bool_t FillControlHistograms(AliFlowEventSimple *anEvent, TList *weightsList=NULL, Bool_t usePhiWeights=kFALSE, Bool_t usePtWeights=kFALSE, Bool_t useEtaWeights=kFALSE)
Double_t Phi() const
void Print(Option_t *option="") const
Double_t GetEntriesInPtBinRP(Int_t iBin)
static AliFlowCommonConstants * GetMaster()
Double_t GetEntriesInEtaBinRP(Int_t iBin)
virtual void Get2Qsub(AliFlowVector *Qarray, Int_t n=2, TList *weightsList=NULL, Bool_t usePhiWeights=kFALSE, Bool_t usePtWeights=kFALSE, Bool_t useEtaWeights=kFALSE)
Int_t GetReferenceMultiplicity() const
Double_t Pt() const
TProfile * fRefMultVsNoOfRPs
Double_t Weight() const
TProfile * fHistProMeanPtperBin
Bool_t InPOISelection(Int_t poiType=1) const
Double_t Eta() const
Int_t NumberOfTracks() const