StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
plotLYZ.C
1 //
3 // $Id: plotLYZ.C,v 1.4 2007/02/06 19:00:53 posk Exp $
4 //
5 // Plot histograms from flow.LeeYang.Zeros.root
6 //
7 // by Markus Oldenberg and Art Poskanzer
8 //
10 #include "TFile.h"
11 #include "TCanvas.h"
12 #include "TH1F.h"
13 #include "TLine.h"
14 #include "TMath.h"
15 #include <iostream.h>
16 #include <iomanip.h>
17 
18 void plotLYZ(TString fileName = "flow.hist.root", Bool_t firstPass = kFALSE, char* ext = "") {
19 
20  Bool_t mevSim, FTPC = kFALSE;
21  //Bool_t mevSim = kTRUE;
22  Bool_t FTPC = kTRUE;
23 
24  gROOT->SetStyle("Pub"); // set style
25  gStyle->SetFrameLineWidth(2);
26  gStyle->SetLineWidth(2);
27  gStyle->SetLineColor(kBlue);
28  gStyle->SetHistLineWidth(1);
29  gStyle->SetFuncWidth(2);
30 // gStyle->SetMarkerSize(2);
31 // gStyle->SetMarkerStyle(kFullCircle);
32 // gStyle->SetMarkerColor(kRed);
33  gROOT->ForceStyle();
34 
35  TFile *file = new TFile(fileName, "READ");
36 
37  Int_t maxSel = 2;
38  Int_t maxHar = 4; // 4
39  Int_t maxHarPlot = 2;
40  float ptMax = mevSim ? 2. : 6.;
41  float v1 = 4.0; // 4. for dir9 mevSim
42 
43  Float_t max, min, r0;
44 
45  // Get Ntheta
46  TString histNameNtheta("FlowLYZ_r0theta_Sel1_Har2");
47  TH1D* histNtheta = new TH1D;
48  histNtheta = (TH1D*)file->Get(histNameNtheta);
49  if (!histNtheta) {
50  cout << "### Can't find file " << histNameNtheta << endl;
51  return;
52  }
53  int maxTheta = histNtheta->GetNbinsX();
54  cout << "maxTheta = " << maxTheta << endl << endl;
55  if (maxTheta > 5) maxTheta = 5;
56  delete histNtheta;
57 
58  // Make the canvases
59  TCanvas *canMult = new TCanvas("MeanMultiplicity", "MeanMultiplicity");
60  TPaveLabel* run = new TPaveLabel(0.1,0.01,0.2,0.03,fileName.Data());
61  run->Draw();
62  TDatime now;
63  TPaveLabel* date = new TPaveLabel(0.7,0.01,0.9,0.03,now.AsString());
64  date->Draw();
65  TPad* graphPad = new TPad("Graphs","Graphs",0.01,0.05,0.97,0.95);
66  graphPad->Draw();
67 
68  TCanvas *canGtheta = new TCanvas("Gtheta", "Gtheta", 700, 900);
69  canGtheta->Divide(maxHarPlot*maxSel, maxTheta);
70 
71  TCanvas *canGthetaZoom = new TCanvas("Gtheta_Zoom", "Gtheta_Zoom", 700, 900);
72  canGthetaZoom->Divide(maxHarPlot*maxSel, maxTheta);
73 
74  TCanvas *can_r0theta = new TCanvas("r0", "r0");
75  can_r0theta->Divide(maxSel, maxHarPlot);
76 
77  TCanvas *can_vEta = new TCanvas("v(eta)", "v(eta)", 700, 900);
78  can_vEta->Divide(maxSel, maxHar);
79 
80  TCanvas *can_vPt = new TCanvas("v(pt)", "v(pt)", 700, 900);
81  can_vPt->Divide(maxSel, maxHar);
82 
83  TCanvas *can_vr0 = new TCanvas("vr0", "vr0");
84  can_vr0->Divide(maxSel);
85 
86  TCanvas *can_v = new TCanvas("v", "v");
87  can_v->Divide(maxSel);
88 
89  TCanvas *can_centX = new TCanvas("centX", "centX");
90  can_centX->Divide(maxSel, maxHarPlot);
91 
92  TCanvas *can_centY = new TCanvas("centY", "centY");
93  can_centY->Divide(maxSel, maxHarPlot);
94 
95  TCanvas *can_centQ = new TCanvas("centQ", "centQ");
96  can_centQ->Divide(maxSel, maxHarPlot);
97 
98  // Make the histograms
99  TH1D *histMult = new TH1D;
100  TH1F **hist_v = new TH1F*[maxSel];
101  TH1F **hist_vr0 = new TH1F*[maxSel];
102  TH1F **hist_v = new TH1F*[maxSel];
103  TH1F **hist_vEta = new TH1F*[maxSel*maxHar];
104  TH1F **hist_vPt = new TH1F*[maxSel*maxHar];
105  TH1D **hist_r0th = new TH1D*[maxSel*maxHarPlot];
106  TH1F **histG = new TH1F*[maxSel*maxHarPlot*maxTheta];
107  TH1F **histGZoom = new TH1F*[maxSel*maxHarPlot*maxTheta];
108  TH1D **hist_centX = new TH1D*[maxSel*maxHarPlot];
109  TH1D **hist_centY = new TH1D*[maxSel*maxHarPlot];
110  TH1D **hist_centQ = new TH1D*[maxSel*maxHarPlot];
111 
112  TLine **r0Line = new TLine*[maxSel*maxHarPlot*maxTheta];
113  TLine **r0LineZoom = new TLine*[maxSel*maxHarPlot*maxTheta];
114  TLine *ptZeroLine = new TLine(0., 0., ptMax, 0.);
115  TLine *pt5Line = new TLine(0., 5., 2., 5.);
116  TLine *ptV1Line = new TLine(0., v1, 2., v1);
117  TLine *etaZeroLine = new TLine(-1.5, 0., 1.5, 0.);
118  TLine *etaZeroLineFTPC = new TLine(-4.5, 0., 4.5, 0.);
119  TLine *GZeroLine = new TLine(0., 0., 0.35, 0.);
120  TLine *eta5Line = new TLine(-4.5, 5., 4.5, 5.);
121  TLine *etaV1LinePos = new TLine(0., v1, 4.5, v1);
122  TLine *etaV1LineNeg = new TLine(-4.5, -v1, 0., -v1);
123  TLine *vLine = new TLine(0.5, 0., maxHar+0.5, 0.);
124  TLine *recentZeroLine = new TLine(0.5, 0., 3.5, 0.);
125  TLine *recentedZeroLine = new TLine(0.5, 0., 2.5, 0.);
126 
127  // Multiplicity
128  TString histName("FlowLYZ_Mult");
129  cout << histName << endl;
130  graphPad->cd();
131  histMult = (TH1D*)file->Get(histName);
132  if (!mevSim) { histMult->Fit("gaus"); }
133  histMult->Draw();
134  Double_t entries = histMult->GetEntries();
135  TString* entriesChar = new TString("entries= ");
136  *entriesChar += (int)entries;
137  TLatex l;
138  l.SetNDC();
139  l.SetTextSize(0.05);
140  l.DrawLatex(0.65,0.8,entriesChar->Data());
141 
142  float _v, vErr;
143  for (Int_t sel = 0; sel < maxSel; sel++) {
144 
145  if (!firstPass) {
146  TString histName("FlowLYZ_v_Sel");
147  histName += sel+1;
148  cout << histName << endl;
149  hist_v[sel] = (TH1F*)file->Get(histName);
150  can_v->cd(sel+1);
151  hist_v[sel]->SetMinimum(0.);
152  hist_v[sel]->Draw();
153  vLine->Draw();
154  for (int j=1; j<=maxHar; j++) {
155  _v = hist_v[sel]->GetBinContent(j);
156  vErr = hist_v[sel]->GetBinError(j);
157  cout << setprecision(3) << "Sel = " << sel+1 << ": v" << j << " from pt = (" << _v <<
158  " +/- " << vErr << ") %" << endl;
159  }
160  }
161 
162  TString histName("FlowLYZ_vr0_Sel");
163  histName += sel+1;
164  cout << histName << endl;
165  hist_vr0[sel] = (TH1F*)file->Get(histName);
166  can_vr0->cd(sel+1);
167  hist_vr0[sel]->SetMinimum(0.);
168  hist_vr0[sel]->Draw();
169  for (int j=1; j<=maxHar; j++) {
170  _v = hist_vr0[sel]->GetBinContent(j);
171  vErr = hist_vr0[sel]->GetBinError(j);
172  cout << setprecision(3) << "Sel= " << sel+1 << ": v" << j << " from r0 = (" << _v <<
173  " +/- " << vErr << ") %" << endl;
174  }
175 
176  for (Int_t har = 0; har < maxHarPlot; har++) {
177  int n = sel + har;
178 
179  TString histName("FlowLYZ_r0theta_Sel");
180  histName += sel+1;
181  histName += "_Har";
182  histName += har+1;
183  cout << histName << endl;
184  hist_r0th[n] = (TH1D*)file->Get(histName);
185  can_r0theta->cd(1+sel+har*maxSel);
186  //hist_r0th[n]->Fit("pol0");
187  hist_r0th[n]->SetMinimum(0.);
188  hist_r0th[n]->Draw();
189  for (int th=1; th<=maxTheta; th++) {
190  _v = hist_r0th[n]->GetBinContent(th);
191  vErr = hist_r0th[n]->GetBinError(th);
192  if (TMath::IsNaN(vErr)) {
193  vErr = 0.;
194  hist_r0th[n]->SetBinError(th, 0.);
195  }
196  cout << setprecision(3) << "Sel=" << sel+1 << ", Har=" << har+1 <<": r0" << th << " = "
197  << _v << " +/- " << vErr << endl;
198  }
199 
200  if (!firstPass) {
201  TString histName("FlowCentX_Sel");
202  histName += sel+1;
203  histName += "_Har";
204  histName += har+1;
205  hist_centX[n] = (TH1D*)file->Get(histName);
206  if (hist_centX[n]) {
207  cout << histName << endl;
208  can_centX->cd(1+sel+har*maxSel);
209  hist_centX[n]->Draw();
210  recentZeroLine->Draw();
211  }
212 
213  TString histName("FlowCentY_Sel");
214  histName += sel+1;
215  histName += "_Har";
216  histName += har+1;
217  hist_centY[n] = (TH1D*)file->Get(histName);
218  if (hist_centY[n]) {
219  cout << histName << endl;
220  can_centY->cd(1+sel+har*maxSel);
221  hist_centY[n]->Draw();
222  recentZeroLine->Draw();
223  }
224  }
225 
226  TString histName("FlowQCent_Sel");
227  histName += sel+1;
228  histName += "_Har";
229  histName += har+1;
230  hist_centQ[n] = (TH1D*)file->Get(histName);
231  if (hist_centQ[n]) {
232  cout << histName << endl;
233  can_centQ->cd(1+sel+har*maxSel);
234  hist_centQ[n]->Draw();
235  recentedZeroLine->Draw();
236  }
237 
238  min = 0.00001;
239  float expan = 1.2;
240  for (Int_t theta = 0; theta < maxTheta; theta++) {
241  TString histName("FlowLYZ_Gtheta");
242  histName += theta;
243  histName += "_Sel";
244  histName += sel+1;
245  histName += "_Har";
246  histName += har+1;
247  cout << histName << endl;
248  histG[theta] = (TH1F*)file->Get(histName);
249  canGtheta->cd(1+har+sel*maxHarPlot+theta*maxSel*maxHarPlot);
250  TVirtualPad::Pad()->SetLogy();
251  histG[theta]->SetMinimum(min);
252  histG[theta]->DrawCopy("PH");
253  r0 = hist_r0th[n]->GetBinContent(theta+1);
254  r0Line[theta] = new TLine(r0, 0., r0, 1.);
255  r0Line[theta]->Draw();
256 
257  min = 0.000001;
258  max = 1.;
259  histGZoom[theta] = (TH1F*)file->Get(histName);
260  canGthetaZoom->cd(1+har+sel*maxHarPlot+theta*maxSel*maxHarPlot);
261  TVirtualPad::Pad()->SetLogy();
262  histGZoom[theta]->SetMaximum(max);
263  histGZoom[theta]->SetMinimum(min);
264  histGZoom[theta]->SetAxisRange(r0/expan, r0*expan, "X");
265  histGZoom[theta]->Draw("PH");
266  r0LineZoom[theta] = new TLine(r0, 0., r0, max);
267  r0LineZoom[theta]->Draw();
268 
269  }
270  }
271 
272  if (!firstPass) {
273  for (Int_t har = 0; har < maxHar; har++) {
274  int n = sel + har;
275 
276  TString histName("FlowLYZ_vEta_Sel");
277  histName += sel+1;
278  histName += "_Har";
279  histName += har+1;
280  cout << histName << endl;
281  hist_vEta[n] = (TH1D*)file->Get(histName);
282  can_vEta->cd(sel+1+har*maxSel);
283  hist_vEta[n]->SetMaximum(10.);
284  hist_vEta[n]->SetMinimum(-10.);
285  hist_vEta[n]->Draw("E");
286  if (FTPC) { etaZeroLineFTPC->Draw(); }
287  else { etaZeroLine->Draw(); }
288  if (mevSim) {
289  if (har==1) {
290  eta5Line->Draw();
291  } else if (har==0) {
292  etaV1LinePos->Draw();
293  etaV1LineNeg->Draw();
294  }
295  }
296 
297  TString histName("FlowLYZ_vPt_Sel");
298  histName += sel+1;
299  histName += "_Har";
300  histName += har+1;
301  cout << histName << endl;
302  hist_vPt[n] = (TH1D*)file->Get(histName);
303  can_vPt->cd(sel+1+har*maxSel);
304  if (mevSim) {
305  hist_vPt[n]->SetMaximum(10.);
306  hist_vPt[n]->SetMinimum(-10.);
307  } else {
308  hist_vPt[n]->SetMaximum(20.);
309  hist_vPt[n]->SetMinimum(-20.);
310  }
311  hist_vPt[n]->Draw("E");
312  ptZeroLine->Draw();
313  if (mevSim) {
314  if (har==1) {
315  pt5Line->Draw();
316  } else if (har==0) {
317  ptV1Line->Draw();
318  }
319  }
320  }
321  }
322  }
323 
324  if (strstr(ext,"ps")) {
325  canGtheta->SaveAs("FlowLYZ_Gtheta.ps");
326  canGthetaZoom->SaveAs("FlowLYZ_GthetaZoom.ps");
327  can_r0theta->SaveAs("FlowLYZ_r0.ps");
328  canMult->SaveAs("FlowLYZ_Mult.ps");
329  can_vEta->SaveAs("FlowLYZ_vEta.ps");
330  can_vPt->SaveAs("FlowLYZ_vPt.ps");
331  can_vr0->SaveAs("FlowLYZ_vro.ps");
332  can_v->SaveAs("FlowLYZ_v.ps");
333  can_centX->SaveAs("FlowCentX.ps");
334  can_centY->SaveAs("FlowCentY.ps");
335  can_centQ->SaveAs("FlowCentQ.ps");
336  } else if (strstr(ext,"gif")) {
337  canGtheta->SaveAs("FlowLYZ_Gtheta.gif");
338  canGthetaZoom->SaveAs("FlowLYZ_GthetaZoom.gif");
339  can_r0theta->SaveAs("FlowLYZ_r0.gif");
340  canMult->SaveAs("FlowLYZ_Mult.gif");
341  can_vEta->SaveAs("FlowLYZ_vEta.gif");
342  can_vPt->SaveAs("FlowLYZ_vPt.gif");
343  can_vr0->SaveAs("FlowLYZ_vro.gif");
344  can_v->SaveAs("FlowLYZ_v.gif");
345  }
346 
347  return;
348 }
349 
351 //
352 // $Log: plotLYZ.C,v $
353 // Revision 1.4 2007/02/06 19:00:53 posk
354 // In Lee Yang Zeros method, introduced recentering of Q vector.
355 // Reactivated eta symmetry cut.
356 //
357 // Revision 1.2 2006/03/22 22:02:14 posk
358 // Updates to macros.
359 //
360 //