AliPhysics  71e3bc7 (71e3bc7)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
AliFlowLYZHist1.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 
17 //#include "Riostream.h" //in case one wants to make print statements
18 #include "TProfile.h"
19 #include "TString.h"
20 #include "TComplex.h"
21 #include "TList.h"
22 #include "AliFlowLYZConstants.h"
23 #include "AliFlowLYZHist1.h"
24 #include "TBrowser.h"
25 
26 class TH1D;
27 
28 // Class to organize the histograms in the first run
29 // in the Lee Yang Zeros Flow analysis.
30 // Also contains methods to find the first minimum R0
31 // of the generating function.
32 // author: N. van der Kolk (kolk@nikhef.nl)
33 
35 
36 //-----------------------------------------------------------------------
37 
39  TNamed(),
40  fHistGtheta(0),
41  fHistProReGtheta(0),
42  fHistProImGtheta(0),
43  fHistList(NULL)
44 {
45  //default constructor
46 }
47 
48 //-----------------------------------------------------------------------
49 
50  AliFlowLYZHist1::AliFlowLYZHist1(Int_t theta, const char *anInput, Bool_t useSum):
51  TNamed(anInput,anInput),
52  fHistGtheta(0),
53  fHistProReGtheta(0),
54  fHistProImGtheta(0),
55  fHistList(NULL)
56 {
57 
58  //constructor creating histograms
62  Double_t dMin = 0.;
63 
64  TString name, addlast;
65 
66  if (useSum) { addlast = "LYZSUM"; }
67  else { addlast = "LYZPROD"; }
68 
69  //fHistGtheta
70  name = "First_Flow_Gtheta";
71  name +=theta;
72  name +=addlast;
73  if (useSum) { fHistGtheta = new TH1D(name.Data(),name.Data(),iNbins,dMin,dMaxSUM); }
74  else { fHistGtheta = new TH1D(name.Data(),name.Data(),iNbins,dMin,dMaxPROD); }
75  fHistGtheta->SetXTitle("r");
76  fHistGtheta->SetYTitle("|G^{#theta}(ir)|^{2}");
77 
78  //fHistProReGtheta
79  name = "First_FlowPro_ReGtheta";
80  name +=theta;
81  name +=addlast;
82  if (useSum) { fHistProReGtheta = new TProfile(name.Data(),name.Data(),iNbins,dMin,dMaxSUM); }
83  else { fHistProReGtheta = new TProfile(name.Data(),name.Data(),iNbins,dMin,dMaxPROD); }
84  fHistProReGtheta->SetXTitle("r");
85  fHistProReGtheta->SetYTitle("Re G^{#theta}(ir)");
86 
87  //fHistProImGtheta
88  name = "First_FlowPro_ImGtheta";
89  name +=theta;
90  name +=addlast;
91  if (useSum) { fHistProImGtheta = new TProfile(name.Data(),name.Data(),iNbins,dMin,dMaxSUM); }
92  else { fHistProImGtheta = new TProfile(name.Data(),name.Data(),iNbins,dMin,dMaxPROD); }
93  fHistProImGtheta->SetXTitle("r");
94  fHistProImGtheta->SetYTitle("Im G^{#theta}(ir)");
95 
96  //list of histograms
97  fHistList = new TList();
98  fHistList-> Add(fHistGtheta);
101 
102 }
103 
104 //-----------------------------------------------------------------------
105 
107 {
108  //deletes histograms
109  delete fHistGtheta;
110  delete fHistProReGtheta;
111  delete fHistProImGtheta;
112  delete fHistList;
113 }
114 
115 //-----------------------------------------------------------------------
116 
118 {
119  //fill the histograms
120 
121  fHistProReGtheta->Fill(f, c.Re());
122  fHistProImGtheta->Fill(f, c.Im());
123 }
124 
125 //-----------------------------------------------------------------------
126 
128 {
129  //method called in Finish() of AliFlowAnalysisWithLeeYangZeroes
130  //fills the fHistGtheta histograms
131 
132  Int_t iNbins = fHistGtheta->GetNbinsX();
133  for (Int_t bin=1;bin<=iNbins;bin++)
134  {
135  //get bincentre of bins in histogram
136  Double_t dRe = fHistProReGtheta->GetBinContent(bin);
137  Double_t dIm = fHistProImGtheta->GetBinContent(bin);
138  TComplex cGtheta(dRe,dIm);
139  //fill fHistGtheta with the modulus squared of cGtheta
140  //to avoid errors when using a merged outputfile use SetBinContent() and not Fill()
141  fHistGtheta->SetBinContent(bin,cGtheta.Rho2());
142  fHistGtheta->SetBinError(bin,0.0);
143  }
144 
145  return fHistGtheta;
146 }
147 
148 //-----------------------------------------------------------------------
149 
151 {
152  //find the first minimum of the square of the modulus of Gtheta
153 
154  Int_t iNbins = fHistGtheta->GetNbinsX();
155  Double_t dR0 = 0.;
156 
157  for (Int_t b=2;b<iNbins;b++)
158  {
159  Double_t dG0 = fHistGtheta->GetBinContent(b);
160  Double_t dGnext = fHistGtheta->GetBinContent(b+1);
161  Double_t dGnextnext = fHistGtheta->GetBinContent(b+2);
162 
163  if (dGnext > dG0 && dGnextnext > dG0)
164  {
165  Double_t dGlast = fHistGtheta->GetBinContent(b-1);
166  Double_t dXlast = fHistGtheta->GetBinCenter(b-1);
167  Double_t dX0 = fHistGtheta->GetBinCenter(b);
168  Double_t dXnext = fHistGtheta->GetBinCenter(b+1);
169 
170  dR0 = dX0 - ((dX0-dXlast)*(dX0-dXlast)*(dG0-dGnext) - (dX0-dXnext)*(dX0-dXnext)*(dG0-dGlast))/
171  (2.*((dX0-dXlast)*(dG0-dGnext) - (dX0-dXnext)*(dG0-dGlast))); //parabolic interpolated minimum
172 
173  break; //stop loop if minimum is found
174  } //if
175 
176  }//b
177 
178 
179  return dR0;
180 }
181 
182 //-----------------------------------------------------------------------
184 {
185  //gets bincenter of histogram
186  Double_t dR = fHistGtheta->GetBinCenter(i);
187  return dR;
188 }
189 
190 //-----------------------------------------------------------------------
191 
193 {
194  //gets iNbins
195  Int_t iNbins = fHistGtheta->GetNbinsX();
196  return iNbins;
197 }
198 
199 //-----------------------------------------------------------------------
201 {
202  //merge fuction
203  if (!aList) return 0;
204  if (aList->IsEmpty()) return 0; //no merging is needed
205 
206  Int_t iCount = 0;
207  TIter next(aList); // list is supposed to contain only objects of the same type as this
208  AliFlowLYZHist1 *toMerge;
209  // make a temporary list
210  TList *pTemp = new TList();
211  while ((toMerge=(AliFlowLYZHist1*)next())) {
212  pTemp->Add(toMerge->GetHistList());
213  iCount++;
214  }
215  // Now call merge for fHistList providing temp list
216  fHistList->Merge(pTemp);
217  // Cleanup
218  delete pTemp;
219 
220  return (double)iCount;
221 
222 }
223 
224 //-----------------------------------------------------------------------
225 void AliFlowLYZHist1::Print(Option_t *option) const
226 {
227  // -*-*-*-*-*Print some global quantities for this histogram collection class *-*-*-*-*-*-*-*
228  // ===============================================
229  // printf( "TH1.Print Name = %s, Entries= %d, Total sum= %g\n",GetName(),Int_t(fEntries),GetSumOfWeights());
230  printf( "Class.Print Name = %s, Histogram list:\n",GetName());
231 
232  if (fHistList) {
233  fHistList->Print(option);
234  }
235  else
236  {
237  printf( "Empty histogram list \n");
238  }
239 }
240 
241 //-----------------------------------------------------------------------
242  void AliFlowLYZHist1::Browse(TBrowser *b)
243 {
244 
245  if (!b) return;
246  if (fHistList) b->Add(fHistList,"AliFlowLYZHist1List");
247 }
248 
249 
250 
void Fill(Double_t f, TComplex c)
double Double_t
Definition: External.C:58
virtual ~AliFlowLYZHist1()
TProfile * fHistProImGtheta
void Print(Option_t *option="") const
void Browse(TBrowser *b)
Double_t GetBinCenter(Int_t i)
TCanvas * c
Definition: TestFitELoss.C:172
static AliFlowLYZConstants * GetMaster()
int Int_t
Definition: External.C:63
TList * GetHistList()
Double_t GetMaxPROD() const
ClassImp(AliFlowLYZHist1) AliFlowLYZHist1
Definition: External.C:212
virtual Double_t Merge(TCollection *aList)
TProfile * fHistProReGtheta
const char Option_t
Definition: External.C:48
bool Bool_t
Definition: External.C:53
Double_t GetMaxSUM() const