]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FLOW/AliFlowLYZHist1.cxx
bug in StepManager fixed, thanks to Andreas
[u/mrichter/AliRoot.git] / PWG2 / FLOW / AliFlowLYZHist1.cxx
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
34 ClassImp(AliFlowLYZHist1)
35
36 //-----------------------------------------------------------------------
37
38   AliFlowLYZHist1::AliFlowLYZHist1():
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,const char *aTitle):
51     TNamed(anInput,aTitle),
52     fHistGtheta(0),
53     fHistProReGtheta(0),
54     fHistProImGtheta(0),
55     fHistList(NULL)
56 {
57
58   //constructor creating histograms 
59   Int_t iNbins = AliFlowLYZConstants::kNbins;
60   Double_t dMin = AliFlowLYZConstants::fgMin;
61   Double_t dMax = AliFlowLYZConstants::fgMax;
62   TString title, name;
63  
64   
65   //fHistGtheta
66   name = "First_Flow_Gtheta";
67   name +=theta;
68   name +="_LYZ";
69   title = "First_Flow_Gtheta";
70   title +=theta;
71   title +="_LYZ";
72   fHistGtheta = new TH1D(name.Data(),title.Data(),iNbins,dMin,dMax);  
73   fHistGtheta->SetXTitle("r");
74   fHistGtheta->SetYTitle("|G^{#theta}(ir)|^{2}");
75   
76   //fHistProReGtheta
77   name = "First_FlowPro_ReGtheta";
78   name +=theta;
79   name +="_LYZ";
80   title = "First_FlowPro_ReGtheta";
81   title +=theta;
82   title +="_LYZ";
83   fHistProReGtheta = new TProfile(name.Data(),title.Data(),iNbins,dMin,dMax);
84   fHistProReGtheta->SetXTitle("r");
85   fHistProReGtheta->SetYTitle("Re G^{#theta}(ir)");
86   
87   //fHistProImGtheta
88   name = "First_FlowPro_ImGtheta";
89   name +=theta;
90   name +="_LYZ";
91   title = "First_FlowPro_ImGtheta";
92   title +=theta;
93   title +="_LYZ";
94   fHistProImGtheta = new TProfile(name.Data(),title.Data(),iNbins,dMin,dMax);
95   fHistProImGtheta->SetXTitle("r");
96   fHistProImGtheta->SetYTitle("Im G^{#theta}(ir)");
97
98   //list of histograms
99   fHistList = new TList();
100   fHistList-> Add(fHistGtheta);
101   fHistList-> Add(fHistProReGtheta);
102   fHistList-> Add(fHistProImGtheta);
103       
104 }
105   
106 //----------------------------------------------------------------------- 
107
108 AliFlowLYZHist1::~AliFlowLYZHist1()
109 {
110   //deletes histograms
111   delete fHistGtheta;
112   delete fHistProReGtheta;
113   delete fHistProImGtheta;
114   delete fHistList;
115 }
116
117 //----------------------------------------------------------------------- 
118
119 void AliFlowLYZHist1::Fill(Double_t f, TComplex c)
120 {
121   //fill the histograms
122
123   fHistProReGtheta->Fill(f, c.Re());
124   fHistProImGtheta->Fill(f, c.Im());
125 }
126
127 //----------------------------------------------------------------------- 
128
129 TH1D* AliFlowLYZHist1::FillGtheta()
130 {
131   //fill the fHistGtheta histograms
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 dBin = fHistGtheta->GetBinCenter(bin); 
137           Double_t dRe = fHistProReGtheta->GetBinContent(bin);
138           Double_t dIm = fHistProImGtheta->GetBinContent(bin);
139           TComplex cGtheta(dRe,dIm);
140           //fill fHistGtheta with the modulus squared of cGtheta
141           fHistGtheta->Fill(dBin,cGtheta.Rho2());
142         }
143
144   return fHistGtheta;
145 }
146
147 //----------------------------------------------------------------------- 
148
149 Double_t AliFlowLYZHist1::GetR0()
150 {
151   //find the first minimum of the square of the modulus of Gtheta 
152
153   Int_t iNbins = fHistGtheta->GetNbinsX();
154   Double_t dR0 = 0.; 
155
156   for (Int_t b=2;b<iNbins;b++)
157     {
158       Double_t dG0 = fHistGtheta->GetBinContent(b);
159       Double_t dGnext = fHistGtheta->GetBinContent(b+1);
160       Double_t dGnextnext = fHistGtheta->GetBinContent(b+2);
161       
162       if (dGnext > dG0 && dGnextnext > dG0)
163         {
164           Double_t dGlast = fHistGtheta->GetBinContent(b-1);
165           Double_t dXlast = fHistGtheta->GetBinCenter(b-1);
166           Double_t dX0 = fHistGtheta->GetBinCenter(b);
167           Double_t dXnext = fHistGtheta->GetBinCenter(b+1);
168
169           dR0 = dX0 - ((dX0-dXlast)*(dX0-dXlast)*(dG0-dGnext) - (dX0-dXnext)*(dX0-dXnext)*(dG0-dGlast))/
170             (2.*((dX0-dXlast)*(dG0-dGnext) - (dX0-dXnext)*(dG0-dGlast))); //interpolated minimum
171           
172           break; //stop loop if minimum is found
173         } //if
174
175     }//b
176
177       
178   return dR0;
179 }
180    
181 //----------------------------------------------------------------------- 
182 Double_t AliFlowLYZHist1::GetBinCenter(Int_t i)
183 {
184   //gets bincenter of histogram
185   Double_t dR = fHistGtheta->GetBinCenter(i);
186   return dR;
187 }
188
189 //----------------------------------------------------------------------- 
190
191 Int_t AliFlowLYZHist1::GetNBins()
192 {
193   //gets iNbins
194   Int_t iNbins = fHistGtheta->GetNbinsX();
195   return iNbins;
196 }
197
198 //----------------------------------------------------------------------- 
199  Double_t AliFlowLYZHist1::Merge(TCollection *aList)
200 {
201   //merge fuction
202   if (!aList) return 0;
203   if (aList->IsEmpty()) return 0; //no merging is needed
204
205   Int_t iCount = 0;
206   TIter next(aList); // list is supposed to contain only objects of the same type as this
207   AliFlowLYZHist1 *toMerge;
208   // make a temporary list
209   TList *pTemp = new TList();
210   while ((toMerge=(AliFlowLYZHist1*)next())) {
211     pTemp->Add(toMerge->GetHistList()); 
212     iCount++;
213   }
214   // Now call merge for fHistList providing temp list
215   fHistList->Merge(pTemp);
216   // Cleanup
217   delete pTemp;
218     
219   return (double)iCount;
220     
221 }
222
223 //----------------------------------------------------------------------- 
224 void AliFlowLYZHist1::Print(Option_t *option) const
225 {
226   //   -*-*-*-*-*Print some global quantities for this histogram collection class *-*-*-*-*-*-*-*
227   //             ===============================================
228   //   printf( "TH1.Print Name  = %s, Entries= %d, Total sum= %g\n",GetName(),Int_t(fEntries),GetSumOfWeights());
229   printf( "Class.Print Name = %s, Histogram list:\n",GetName());
230
231   if (fHistList) {  
232     fHistList->Print(option);
233   }
234   else
235     {
236       printf( "Empty histogram list \n");
237     }
238 }
239
240 //----------------------------------------------------------------------- 
241  void AliFlowLYZHist1::Browse(TBrowser *b)
242 {
243
244   if (!b) return;
245   if (fHistList) b->Add(fHistList,"AliFlowLYZHist1List");
246 }
247
248
249