AliPhysics  1909eaa (1909eaa)
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
SummarydNdetaDrawer.C
Go to the documentation of this file.
1 #include "SummaryDrawer.C"
2 #include <TPaveText.h>
3 #include <TMultiGraph.h>
4 
6 {
7 public:
8  enum {
9  kForward = 0x01,
10  kCentral = 0x02,
11  kSums = 0x04,
12  kResults = 0x08,
13  kMinBias = 0x10,
14  kMC = 0x80,
15  kNormal = 0x0F
16  };
18  : SummaryDrawer()
19  {}
20  const char* ColName(const char* prefix, bool results=false)
21  {
22  return Form("%sdNdeta%s", prefix, results ? "Results" : "Sums");
23  }
24  //____________________________________________________________________
25  void Run(const char* fname="forward_dndeta.root", UShort_t flags=kNormal)
26  {
27  // --- Open the file -----------------------------------------------
28  TString filename(fname);
29  TFile* file = TFile::Open(filename.Data(), "READ");
30  if (!file) {
31  Error("Run", "Failed to open \"%s\"", filename.Data());
32  return;
33  }
34  // Options
35  Bool_t forward = flags & kForward;
36  Bool_t central = flags & kCentral;
37  Bool_t sums = flags & kSums;
38  Bool_t results = flags & kResults;
39  Bool_t onlyMB = flags & kMinBias;
40  Bool_t mc = flags & kMC;
41  fPause = flags & kPause;
42 
43  // --- Force MB for pp ---------------------------------------------
44  UShort_t sys = 0;
45  TCollection* c = GetCollection(file, ColName("Forward", false));
46  GetParameter(c, "sys", sys);
47  if (sys == 1) {
48  // onlyMB = true;
49  // Info("Run", "Found sys==1 -> Forcing MB");
50  }
51 
52  // --- Test of MC --------------------------------------------------
53  TCollection* mcC = GetCollection(file, ColName("MCTruth"), false);
54  if (mcC) {
55  TCollection* mcAll = GetCollection(mcC, "all");
56  if (mcAll && GetObject(mcAll, "sum")) {
57  Info("Run", "Found MC truth output");
58  mc = true;
59  }
60  }
61  // --- Make our canvas ---------------------------------------------
62  TString pdfName(filename);
63  pdfName.ReplaceAll(".root", ".pdf");
64  CreateCanvas(pdfName, flags & kLandscape);
65 
66  // --- Make a Title page -------------------------------------------
67  DrawTitlePage(file, mc, onlyMB);
68 
69  // --- Do each sub-algorithm ---------------------------------------
70  THStack* rF = 0;
71  if (forward && sums) DrawSums(file, "Forward", onlyMB);
72  if (forward && results) rF = DrawRes(file, "Forward", onlyMB);
73 
74  THStack* rC = 0;
75  if (central && sums) DrawSums(file, "Central", onlyMB);
76  if (central && results) rC = DrawRes(file, "Central", onlyMB);
77 
78  THStack* rM = 0;
79  if (mc && sums) DrawSums(file, "MCTruth", onlyMB);
80  if (mc && results) rM = DrawRes(file, "MCTruth", onlyMB);
81 
82  if (rC && rF && results) DrawBoth(rC, rF, rM, file, onlyMB);
83 
84  CloseCanvas();
85  }
86 
87 protected:
88  //____________________________________________________________________
89  TAxis* GetCentAxis(const TCollection* parent, Bool_t verbose=false)
90  {
91  TObject* cO = GetObject(parent, "centAxis", verbose);
92  TAxis* cA = 0;
93  if (!cO) return 0;
94 
95  if (cO->IsA()->InheritsFrom(TAxis::Class()))
96  cA = static_cast<TAxis*>(cO);
97  else if (cO->IsA()->InheritsFrom(TH1::Class())) {
98  TH1* cH = static_cast<TH1*>(cO);
99  cA = cH->GetXaxis();
100  }
101  // if (cA) cA->Dump();
102  if (!cA || !cA->GetXbins() || !cA->GetXbins()->GetArray() ||
103  cA->GetXmin() > cA->GetXmax()) return 0;
104  return cA;
105  }
106  //____________________________________________________________________
108  const TString& base,
109  Float_t cLow,
110  Float_t cHigh,
111  TString& title)
112  {
113  TString folder;
114  title = TString::Format("%s %s: ", base.Data(), title.Data());
115  if (cLow < 0 || cHigh < 0 || cLow >= cHigh) {
116  folder = "all";
117  title.Append("All selected events");
118  }
119  else {
120  folder.Form("cent%03dd%02d_%03dd%02d",
121  Int_t(cLow), Int_t(cLow*100) % 100,
122  Int_t(cHigh), Int_t(cHigh*100) % 100);
123  title.Append(Form("%6.2f%% - %6.2f%%", cLow, cHigh));
124  }
125 
126  return GetCollection(sums, folder);
127 
128  }
129  //____________________________________________________________________
130  void DrawTitlePage(TFile* file, Bool_t mc, Bool_t onlyMB)
131  {
132  TCollection* c = GetCollection(file, ColName("Forward", true));
133 
134  fBody->cd();
135 
136  Double_t y = .9;
137  TLatex* ltx = new TLatex(.5, y,
138  "#frac{1}{#it{N}}#kern[.1]"
139  "{#frac{d#it{N_{ch}}}{d#it{#eta}}}");
140  ltx->SetTextSize(0.07);
141  ltx->SetTextFont(42);
142  ltx->SetTextAlign(22);
143  ltx->SetNDC();
144  ltx->Draw();
145  y -= .075;
146 
147  if (mc) {
148  ltx = new TLatex(.5, y, "Simulation input");
149  ltx->SetNDC();
150  ltx->SetTextAlign(23);
151  ltx->SetTextFont(42);
152  ltx->SetTextSize(.03);
153  ltx->Draw();
154  y -= .035;
155  }
156  if (onlyMB) {
157  ltx = new TLatex(.5, y, "No centrality");
158  ltx->SetNDC();
159  ltx->SetTextAlign(23);
160  ltx->SetTextFont(42);
161  ltx->SetTextSize(.03);
162  ltx->Draw();
163  y -= .035;
164  }
165 
166  DrawResTitle(c, y, onlyMB);
167 
168  PrintCanvas("1/N dN/d#eta");
169  }
170  //____________________________________________________________________
171  void DrawBoth(THStack* rC, THStack* rF, THStack* rM,
172  TFile* file, Bool_t onlyMB)
173  {
174  fBody->cd();
175  Double_t y1 = fLandscape ? 0 : .3;
176  Double_t x2 = fLandscape ? .7 : 1;
177  Double_t x1 = fLandscape ? x2 : 0;
178  Double_t y2 = fLandscape ? 1 : y1;
179  TPad* p1 = new TPad("p1", "p1", 0, y1, x2, 1, 0, 0);
180  TPad* p2 = new TPad("p2", "p2", x1, 0, 1, y2, 0, 0);
181 
182  fBody->cd();
183  p1->Draw();
184  p1->cd();
185 
186  TH1* h = 0;
187  if (rM) {
188  TIter nextM(rM->GetHists());
189  while ((h = static_cast<TH1*>(nextM()))) rC->Add(h);
190  }
191 
192  TIter next(rF->GetHists());
193  while ((h = static_cast<TH1*>(next()))) rC->Add(h);
194 
195 
196  rC->Draw("nostack");
197 
198  TCollection* fS = GetCollection(file, ColName("Forward", false));
199  UShort_t sys, sNN;
200  ULong_t trigger;
201  GetParameter(fS, "sNN", sNN);
202  GetParameter(fS, "sys", sys);
203  GetParameter(fS, "trigger", trigger);
204  TAxis* centAxis = GetCentAxis(fS);
205  UShort_t cLow = centAxis && !onlyMB ? centAxis->GetXmin() : 0;
206  UShort_t cHigh = centAxis && !onlyMB ? centAxis->GetXmax() : 100;
207 
208  // CompileScript("OtherData.C", "", "RefData", false);
209 
210  // If we have V0AND trigger, get NSD other data
211  TMultiGraph* other = 0;
212 #if 0
213  if (!centAxis) {
214  Int_t oT = (trigger == 0x2000) ? 0x4 : trigger;
215  TString oC = Form("RefData::GetData(%hu,%hu,%hu,%hu,%hu,0xF)",
216  sys, sNN, oT, cLow, cHigh);
217  other = reinterpret_cast<TMultiGraph*>(gROOT->ProcessLine(oC));
218  }
219  else {
220  other = new TMultiGraph("other", "");
221  Int_t nCent = centAxis->GetNbins();
222  for (Int_t i = 1; i <= nCent; i++) {
223  TString oC = Form("RefData::GetData(%hu,%hu,%hu,%hu,%hu,0xF)",
224  sys, sNN, 0, UShort_t(centAxis->GetBinLowEdge(i)),
225  UShort_t(centAxis->GetBinUpEdge(i)));
226  TMultiGraph* oM =
227  reinterpret_cast<TMultiGraph*>(gROOT->ProcessLine(oC));
228  if (oM) other->Add(oM);
229  }
230  }
231 #endif
232  if (other) {
233  // p1->Clear();
234  // other->Draw("ap");
235  // Double_t oMax = other->GetHistogram()->GetMaximum();
236  // Double_t rMax = rC->GetMaximum("nostack");
237  // other->SetMaximum(1.2*TMath::Max(oMax, rMax));
238  // rC->Draw("same nostack");
239  TObject* g = 0;
240  TIter nextG(other->GetListOfGraphs());
241  while ((g = nextG())) {
242  // Printf("Drawing %s/%s", g->GetName(), g->GetTitle());
243  g->DrawClone("same p");
244  }
245  }
246 
247  fBody->cd();
248  p2->Draw();
249  p2->cd();
250 
251 
252  TLegend* l = new TLegend(0.01, 0.1, 0.99, 0.99,
253  onlyMB || !centAxis ? "" : "Centralities");
254  l->SetNColumns(fLandscape ? 1 : 2);
255  l->SetFillStyle(0);
256  l->SetBorderSize(0);
257  CleanStack(rC, l, onlyMB ? 0 : centAxis);
258  TString seen;
259  if (other) {
260  TIter nextG(other->GetListOfGraphs());
261  TObject* g = 0;
262  while ((g = nextG())) {
263  if (seen.Index(g->GetTitle()) != kNPOS) continue;
264  seen.Append(Form("|%s", g->GetTitle()));
265  TLegendEntry* e = l->AddEntry("dummy", g->GetTitle(), "p");
266  TGraph* gg = static_cast<TGraph*>(g);
267  e->SetMarkerStyle(gg->GetMarkerStyle());
268  e->SetMarkerSize(gg->GetMarkerSize());
269  e->SetMarkerColor(kBlack);
270  }
271  }
272  l->Draw();
273 
274  PrintCanvas("Both");
275  }
276  //____________________________________________________________________
277  void DrawSums(TDirectory* top, const TString& base, bool onlyMB)
278  {
279  TCollection* c = GetCollection(top, ColName(base));
280  if (!c) return;
281 
282  TAxis* centAxis = (onlyMB ? 0 : GetCentAxis(c));
283  if (centAxis && centAxis->GetNbins() < 1) centAxis = 0;
284 
285  Int_t txtPad = 0;
286  Double_t xSave = fParVal->GetX();
287  Double_t size = 0.05;
288  fParVal->SetX(.45);
289  Double_t y = .8;
290 
291  if (!onlyMB && centAxis) {
292  size = 0.03;
293  fBody->Divide(1, 2);
294  txtPad = 1;
295 
296  fBody->cd(1);
297  for (Int_t i = 1; i <= centAxis->GetNbins(); i++) {
298  DrawParameter(y, (i == 1 ? "Centrality classes" : ""),
299  Form("%3d%% - %3d%%",
300  Int_t(centAxis->GetBinLowEdge(i)),
301  Int_t(centAxis->GetBinUpEdge(i))), size);
302  }
303 
304  TH1* cent = GetH1(c, "cent");
305  cent->GetXaxis()->SetRangeUser(0,100);
306  cent->SetFillColor(kRed+1);
307  cent->SetFillStyle(3002);
308  cent->SetXTitle("Centrality [%]");
309  cent->SetYTitle("Events");
310  cent->SetMaximum(1.3*cent->GetMaximum());
311  TH1* centAcc = GetH1(c, "centAcc");
312  centAcc->SetFillStyle(3002);
313 
314  TLatex* overUnder = new TLatex(0.15, .88,
315  Form("#splitline{<0: %d}{>100: %d}",
316  int(cent->GetBinContent(1)),
317  int(cent->GetBinContent(102))));
318  overUnder->SetTextColor(kRed+1);
319  overUnder->SetNDC();
320  overUnder->SetTextAlign(13);
321  overUnder->SetTextFont(42);
322  TLatex* overUnderAcc = new TLatex(0.3, .88,
323  Form("#splitline{<0: %d}{>100: %d}",
324  int(centAcc->GetBinContent(1)),
325  int(centAcc->GetBinContent(102))));
326  overUnderAcc->SetTextColor(kGreen+1);
327  overUnderAcc->SetNDC();
328  overUnderAcc->SetTextAlign(13);
329  overUnderAcc->SetTextFont(42);
330 
331  DrawInPad(fBody, 2, cent);
332  DrawInPad(fBody, 2, centAcc, "same", kLegend);
333  DrawInPad(fBody, 2, overUnder, "same");
334  DrawInPad(fBody, 2, overUnderAcc, "same");
335 
336  }
337  fBody->cd(txtPad);
338 
339  UShort_t sys, sNN, scheme;
340  ULong_t trigger;
341  GetParameter(c, "sNN", sNN);
342  GetParameter(c, "sys", sys);
343  GetParameter(c, "scheme", scheme);
344  GetParameter(c, "trigger", trigger);
345 
346  TString schemeString;
347  if (scheme == 0) schemeString = "1/N_{accepted}";
348  if (scheme & 0x1) schemeString.Append("1/#epsilon_{V}1/#epsilon_{T}");
349  if (scheme & 0x2) schemeString.Append("Shape ");
350  if (scheme & 0x4) schemeString.Append("A+C-E ");
351  if (scheme & 0x8) schemeString.Append("#epsilon_{T,MC} ");
352  if (scheme & 0x10) schemeString.Append("0-bin");
353 
354  TString trigString; TriggerString(trigger, trigString);
355  TString sysString; SysString(sys, sysString);
356  TString sNNString; SNNString(sNN, sNNString);
357 
358  DrawParameter(y, "Collision system", sysString, size);
359  DrawParameter(y, "#sqrt{s_{NN}}", sNNString, size);
360  DrawParameter(y, "Normalization scheme", schemeString,size);
361  DrawParameter(y, "Triggers", trigString, size);
362 
363  fParVal->SetX(xSave);
364 
365  PrintCanvas(Form("%s sums", base.Data()));
366 
367  Int_t cLow = centAxis ? centAxis->GetXmin() : 0;
368  Int_t cHigh = centAxis ? -centAxis->GetXmax() : -100;
369  DrawCentSum(c, base, cLow, cHigh);
370  if (onlyMB || !centAxis) return;
371 
372  for (Int_t i = 1; i <= centAxis->GetNbins(); i++)
373  DrawCentSum(c, base, centAxis->GetBinLowEdge(i),
374  centAxis->GetBinUpEdge(i));
375  }
376  //____________________________________________________________________
377  void DrawCentSum(const TCollection* sums, const TString& base,
378  Float_t cLow, Float_t cHigh)
379  {
380  // Info("DrawCentSum", "Drawing centrality sum [%d,%d] in %s (%s)",
381  // cLow, cHigh, sums->GetName(), base.Data());
382  TString title("sums");
383  TCollection* c = GetCentCollection(sums, base, cLow, cHigh, title);
384  if (!c) return;
385 
386  TH2* bin = GetH2(c, "sum");
387  TH2* bin0 = GetH2(c, "sum0");
388  TH1* type = GetH1(c, "events");
389  TH1* trig = GetH1(c, "triggers");
390  TH1* stat = GetH1(c, "status");
391  if (!bin0 || !bin || !trig || !type) return;
392 
393  type->SetFillStyle(3001);
394  type->SetFillColor(kGreen+1);
395 
396  fBody->Divide(2, 2);
397 
398  DrawInPad(fBody, 1, trig, "HIST TEXT");
399  DrawInPad(fBody, 2, type, "HIST TEXT");
400  DrawInPad(fBody, 3, bin, "colz");
401 
402  if (bin0->GetEntries() <= 0) {
403  DrawInPad(fBody, 4, stat, "HIST TEXT");
404  // fBody->cd(4);
405  // TLatex* l = new TLatex(0.5, 0.5, "No 0-bin events");
406  // l->SetNDC();
407  // l->SetTextAlign(22);
408  // l->Draw();
409  }
410  else
411  DrawInPad(fBody, 4, bin0, "colz");
412 
413  PrintCanvas(title);
414  }
415  //____________________________________________________________________
417  {
418  Double_t xSave = fParVal->GetX();
419  Double_t size = 0.05;
420  fParVal->SetX(.5);
421  // Double_t y = .9;
422  TAxis* centAxis = GetCentAxis(c);
423  if (!onlyMB && centAxis) {
424  size = 0.03;
425  for (Int_t i = 1; i <= centAxis->GetNbins(); i++) {
426  DrawParameter(y, (i == 1 ? "Centrality classes" : ""),
427  Form("%6.2f%% - %6.2f%%",
428  centAxis->GetBinLowEdge(i),
429  centAxis->GetBinUpEdge(i)), size);
430  }
431  }
432  TObject* oSNN = GetObject(c, "sNN");
433  TString tSNN; SNNString(oSNN->GetUniqueID(), tSNN);
434  TObject* oTrg = GetObject(c,"trigger");
435  DrawParameter(y, "Collision system", GetObject(c, "sys")->GetTitle(), size);
436  DrawParameter(y, "#sqrt{s_{NN}}",tSNN, size);
437  DrawParameter(y, "Trigger",(oTrg ? oTrg->GetTitle() : "?"), size);
438  TObject* oscheme = GetObject(c,"scheme");
439  TString scheme = oscheme ? oscheme->GetTitle() : "";
440  if (scheme.IsNull()) scheme = "1/N_{accepted}";
441  DrawParameter(y, "Normalization scheme", scheme, size);
442 
443  Double_t epsT = 0, epsT0 = 0;
444  GetParameter(c, "triggerEff", epsT);
445  GetParameter(c, "triggerEff0", epsT0);
446  DrawParameter(y, "#epsilon_{T}", Form("%5.3f", epsT), size);
447  DrawParameter(y, "#epsilon_{T,zero bin}", Form("%5.3f", epsT0), size);
448  Double_t deltaIP =0;
449  GetParameter(c, "deltaIP", deltaIP);
450  DrawParameter(y, "IP #delta_{xy}", Form("%5.3fmm", deltaIP), size);
451 
452  TObject* options = GetObject(c, "options");
453  TString opts(options->GetTitle());
454  TObjArray* tokens = opts.Tokenize(",");
455  TObjString* opt = 0;;
456  TIter oNext(tokens);
457  Bool_t first = true;
458  while ((opt = static_cast<TObjString*>(oNext()))) {
459  DrawParameter(y, (first ? "options" : ""),
460  opt->String().Strip(TString::kBoth), size);
461  first = false;
462  }
463 
464  fParVal->SetX(xSave);
465  }
466 
467  //____________________________________________________________________
468  THStack* DrawRes(TDirectory* top, const TString& base, Bool_t onlyMB)
469  {
470  // Info("DrawRes", "Drawing results for %s", base.Data());
471  TCollection* c = GetCollection(top, ColName(base, true));
472  if (!c) return 0;
473  TCollection* s = GetCollection(top, ColName(base, false));
474  // if (!s) return 0;
475 
476  fBody->cd();
477  Double_t y = .9;
478  DrawResTitle(c, y, onlyMB);
479  PrintCanvas(Form("%s results", base.Data()));
480 
481  Int_t sys = GetObject(c, "sys")->GetUniqueID();
482  TH1* emp = GetH1(c, "empirical");
483  TF1* dc = static_cast<TF1*>(GetObject(c,"deltaCorr"));
484  TF1* vw = static_cast<TF1*>(GetObject(s,"ipZw"));
485  TProfile* sc = static_cast<TProfile*>(GetObject(s,"sumVsC"));
486  if (vw) vw->SetRange(-4,6);
487  Int_t nPad = 0;
488  if (emp) nPad++;
489  if (dc) nPad++;
490  if (vw) nPad++;
491  if (sc) nPad++;
492  if (nPad > 0) {
493  fBody->Divide(nPad,1);
494  Int_t iPad = 1;
495  if (emp) DrawInPad(fBody, iPad++, emp, "", 0, "Empirical");
496  if (dc) DrawInPad(fBody, iPad++, dc, "", 0, "\\hbox{IP} \\delta_{xy}");
497  if (vw) {
498  DrawInPad(fBody, iPad++, vw, "", 0, "\\hbox{IP} \\delta_{z}");
499  Double_t y = .95;
500  DrawParameter(y, "#mu_{Z}", Form("%5.3f", vw->GetParameter(0)));
501  DrawParameter(y, "#sigma_{Z}", Form("%5.3f", vw->GetParameter(1)));
502  DrawParameter(y, "#mu_{Z,ref}", Form("%5.3f", vw->GetParameter(2)));
503  DrawParameter(y, "#sigma_{Z,ref}", Form("%5.3f", vw->GetParameter(3)));
504  }
505  if (sc) DrawInPad(fBody, iPad, sc, "", 0, "#LT#Sigma signal#GT");
506  PrintCanvas(Form("%s results - corrections", base.Data()));
507  }
508 
509  TAxis* centAxis = (onlyMB ? 0 : GetCentAxis(c));
510  if (centAxis && centAxis->GetNbins() < 1) centAxis = 0;
511 
512  THStack* dndeta_ = GetStack(c, "dndeta");
513  if (!dndeta_ || !dndeta_->GetHists() ||
514  dndeta_->GetHists()->GetEntries() < 0) return 0;
515 
516  TLegend* l = new TLegend(0.1, 0.1, 0.9, 0.9,
517  onlyMB || !centAxis? "" : "Centralities");
518  l->SetNColumns(fLandscape ? 1 : 2);
519  l->SetFillStyle(0);
520  l->SetBorderSize(0);
521  THStack* dndeta = CleanStack(dndeta_, l, centAxis);
522 
523  THStack* dndetaEmp = GetStack(c, "dndetaEmp");
524  if (!dndetaEmp || !dndetaEmp->GetHists() ||
525  dndetaEmp->GetHists()->GetEntries() < 0) dndetaEmp = 0;
526 
527  THStack* leftRight = GetStack(c, "leftRight");
528  if (!leftRight || !leftRight->GetHists() ||
529  leftRight->GetHists()->GetEntries() < 0) leftRight = 0;
530 
531  if (!onlyMB) {
532  Double_t y1 = fLandscape ? 0 : .3;
533  Double_t x2 = fLandscape ? .7 : 1;
534  Double_t x1 = fLandscape ? x2 : 0;
535  Double_t y2 = fLandscape ? 1 : y1;
536  TPad* p1 = new TPad("p1", "p1", 0, y1, x2, 1, 0, 0);
537  TPad* p2 = new TPad("p2", "p2", x1, 0, 1, y2, 0, 0);
538  fBody->cd();
539  p1->Draw();
540  p1->cd();
541  fBody->cd();
542  p2->Draw();
543  p2->cd();
544  p1->Divide(1,3,0,0);
545 
546  // fBody->Divide(1, 3, 0, 0);
547 
548  DrawInPad(p2, 0, l, "");
549  DrawInPad(p1, 1, dndeta, "nostack", 0,
550  "\\mathrm{d}N_{\\mathrm{ch}}/\\mathrm{d}\\eta|"
551  "_{\\mathrm{incl}}");
552  DrawInPad(p1, 2, dndetaEmp, "nostack", (sys == 2 ? kLogy : 0),
553  "\\mathrm{d}N_{\\mathrm{ch}}/\\mathrm{d}\\eta|"
554  "_{\\mathrm{prim}}");
555  DrawInPad(p1, 3, leftRight, "nostack", 0, "Left/Right");
556  p1->GetPad(1)->SetGridx();
557  p1->GetPad(2)->SetGridx();
558  p1->GetPad(3)->SetGridx();
559 
560  PrintCanvas(Form("%s results - stacks", base.Data()));
561  }
562 
563  Int_t cLow = centAxis ? centAxis->GetXmin() : 0;
564  Int_t cHigh = centAxis ? -centAxis->GetXmax() : -100;
565  DrawCentRes(c, base, cLow, cHigh);
566  if (onlyMB || !centAxis) {
567  // Info("", "Returning dndeta for MB");
568  dndeta = MakeMBStack(c, base);
569  return dndeta;
570  }
571 
572  for (Int_t i = 1; i <= centAxis->GetNbins(); i++)
573  DrawCentRes(c, base, centAxis->GetBinLowEdge(i),
574  centAxis->GetBinUpEdge(i));
575 
576  return dndeta;
577  }
578  //____________________________________________________________________
579  THStack* MakeMBStack(const TCollection* sums, const TString& base)
580  {
581  TString title("results");
582  TCollection* c = GetCentCollection(sums, base, 0, -1, title);
583  if (!c) return 0;
584 
585  TH1* dndeta = GetH1(c, Form("dndeta%s",base.Data()));
586  if (!dndeta) return 0;
587 
588  THStack* ret = new THStack("dndetaMB", title);
589  ret->Add(dndeta);
590 
591  if (base.EqualTo("MCTruth")) {
592  dndeta = GetH1(c, "dndetaTruth");
593  if (dndeta) ret->Add(dndeta);
594  }
595  return ret;
596  }
597 
598  //____________________________________________________________________
599  void DrawCentRes(const TCollection* sums, const TString& base,
600  Float_t cLow, Float_t cHigh)
601  {
602  // Info("DrawCentRes", "Drawing centrality results [%d,%d] in %s (%s)",
603  // cLow, cHigh, sums->GetName(), base.Data());
604  TString title("results");
605  TCollection* c = GetCentCollection(sums, base, cLow, cHigh, title);
606  if (!c) return;
607 
608 
609  TH1* trig = GetH1(c, "triggers");
610  TH1* norm = GetH1(c, Form("norm%s",base.Data()));
611  TH1* dndeta = GetH1(c, Form("dndeta%s",base.Data()));
612  TH1* dndetaEmp = GetH1(c, Form("dndeta%sEmp",base.Data()));
613  TH2* d2ndetadphi = GetH2(c, Form("d2Ndetadphi%s", base.Data()));
614  if (!trig || !norm || !dndeta || !d2ndetadphi) return;
615  if (norm->GetEntries() <= 0) return;
616 
617  norm->SetFillColor(kGreen+1);
618  norm->SetFillStyle(3001);
619 
620 
621  fBody->Divide(2, (dndetaEmp ? 4 : 3), 0.05, 0);
622 
623  Int_t trP = 1;
624  TVirtualPad* p = fBody->GetPad(trP);
625  p->SetBottomMargin(0.15);
626  p->SetLeftMargin(0.15);
627  if (trP > 2) p->SetTopMargin(0.05);
628 
629  DrawInPad(fBody, trP, trig, "HIST TEXT");
630  DrawInPad(fBody, 2, d2ndetadphi, "colz", 0,
631  "d^{2}#it{N}/d#it{#eta}d#it{phi}|_{incl}");
632  DrawInPad(fBody, 4, norm, "", 0, "Normalization");
633  DrawInPad(fBody, 6, dndeta, "", 0,
634  "d#it{N}_{ch}/d#it{#eta}|_{incl}");
635 
636  fBody->GetPad(2)->SetGridx();
637  fBody->GetPad(4)->SetGridx();
638  fBody->GetPad(6)->SetGridx();
639  fBody->GetPad(2)->SetLeftMargin(0.15);
640  fBody->GetPad(4)->SetLeftMargin(0.15);
641  fBody->GetPad(6)->SetLeftMargin(0.15);
642  fBody->GetPad(4)->SetRightMargin(0.15);
643  fBody->GetPad(6)->SetRightMargin(0.15);
644  if (dndetaEmp) {
645  DrawInPad(fBody, 8, dndetaEmp, "", 0,
646  "d#it{N}_{ch}/d#it{#eta}|_{prim}");
647  fBody->GetPad(8)->SetGridx();
648  fBody->GetPad(8)->SetLeftMargin(0.15);
649  fBody->GetPad(8)->SetRightMargin(0.15);
650  }
651 
652  TObject* normCalc = GetObject(c, "normCalc");
653  TString calc = normCalc ? normCalc->GetTitle() : "?";
654 
655  // Beautify the text
656  calc.ReplaceAll("beta", "#beta");
657  calc.ReplaceAll("eps", "#varepsilon");
658  const char* sufs[] = { "all", "acc", "trg", "vtx", "B", "A", "C", "E",
659  "V", "T", 0 };
660  const char** suf = sufs;
661  while (*suf) {
662  calc.ReplaceAll(Form("_%s", *suf), Form("_{%s}", *suf));
663  suf++;
664  }
665 
666  p = fBody->cd(3);
667  p->SetPad(p->GetXlowNDC(), 0,
668  p->GetXlowNDC()+p->GetWNDC(), p->GetYlowNDC()+p->GetHNDC());
669  fBody->GetPad(5)->Delete();
670  if (dndetaEmp) fBody->GetPad(7)->Delete();
671  TObjArray* lines = calc.Tokenize("\n");
672  // TPaveText* disp = new TPaveText(.1,.1,.9,.9, "NDC");
673  TIter next(lines);
674  TObjString* sline = 0;
675  Double_t y = .95;
676  Double_t xSave = fParName->GetX();
677  Int_t aSave = fParName->GetTextAlign();
678  Double_t tSave = fParVal->GetTextSize();
679  fParName->SetTextAlign(33);
680  fParName->SetX(fParVal->GetX()-.05);
681  while ((sline = static_cast<TObjString*>(next()))) {
682  // disp->AddText(line->GetName());
683  TString& line = sline->String();
684  Ssiz_t eq = line.Last('=');
685  if (eq == kNPOS) {
686  DrawParameter(y, line, "", .6*tSave);
687  continue;
688  }
689  TString name = line(0, eq);
690  TString val = line(eq+1,line.Length()-eq-1);
691  DrawParameter(y, name.Strip(TString::kBoth),
692  val.Strip(TString::kBoth),
693  .6*tSave);
694 
695  }
696  fParName->SetTextAlign(aSave);
697  fParName->SetX(xSave);
698  // disp->SetBorderSize(0);
699  // disp->SetBorderSize(0);
700  // disp->SetFillStyle(0);
701  // DrawInPad(fBody, 3, disp);
702  // fBody->cd();
703 
704  PrintCanvas(title);
705 
706  DrawCentResDetails(c, title);
707  }
708  //____________________________________________________________________
709  void DrawCentResDetails(const TCollection* sums, const TString& base)
710  {
711  TString title = TString::Format("%s - details: ", base.Data());
712 
713  TCollection* c = GetCollection(sums, "partial");
714  if (!c) {
715  Warning("", "Collection partical not found in %s", sums->GetName());
716  sums->ls();
717  return;
718  }
719 
720  fBody->Divide(3, 1, 0.05, 0);
721 
722  const char* typs[] = { "", "0", "All" };
723  const char* tits[] = { "Non-zero events", "Zero events", "Weighted sum" };
724  for (Int_t i = 1; i <= 3; i++) {
725  const char* suf = typs[i-1];
726  TVirtualPad* p = fBody->cd(i);
727  p->SetTopMargin(0.10);
728 
729  TLatex* ltx = new TLatex(0.5, .99, tits[i-1]);
730  ltx->SetNDC();
731  ltx->SetTextAlign(23);
732  ltx->SetTextSize(0.05);
733  ltx->Draw();
734 
735  TH1* sum = GetH2(c, Form("sum%s", suf));
736  TH1* norm = GetH1(c, Form("norm%s", suf));
737  TH1* phi = GetH1(c, Form("phi%s", suf));
738 
739  norm->SetFillColor(kGreen+1);
740  norm->SetFillStyle(3002);
741  phi->SetFillColor(kBlue+1);
742  phi->SetFillStyle(3001);
743 
744  p->Divide(1, 3, 0, 0);
745  DrawInPad(p, 1, sum, sum->Integral()>0 ? "col" : "", 0,
746  "d^{2}#it{N}_{ch}/d#it{#varphi}d#it{#eta}");
747  DrawInPad(p, 2, GetH1(c, Form("average%s", suf)), "", 0,
748  "d#it{N}_{ch}/d#it{#eta}");
749  DrawInPad(p, 3, norm, "", 0, "#eta-coverage/#varphi-acceptance");
750  DrawInPad(p, 3, phi, "same", kLegend);
751  }
752  PrintCanvas(title);
753  }
754 
755  //____________________________________________________________________
756  THStack* CleanStack(const THStack* stack, TLegend* l, const TAxis* axis)
757  {
758  if (!stack) return 0;
759  THStack* ret = new THStack(stack->GetName(), stack->GetTitle());
760  TList* hists = stack->GetHists();
761  TIter next(hists);
762  TH1* h = 0;
763  Int_t j = 0;
764  Bool_t ok = false;
765  while ((h = static_cast<TH1*>(next()))) {
766  TString name(h->GetTitle());
767  TString nme(h->GetName());
768  if (nme.Contains("_mirror", TString::kIgnoreCase)) {
769  // Printf("Ignore %s/%s in stack", nme.Data(), name.Data());
770  continue;
771  }
772  if (l && !ok) {
773  j++;
774  if (axis) {
775  Int_t bin = j; // axis ? TMath::Min(j, axis->GetNbins()) : 1;
776  if (j > axis->GetNbins())
777  name = "0% - 100%";
778  else
779  name.Form("%3d%% - %3d%%",
780  Int_t(axis->GetBinLowEdge(bin)),
781  Int_t(axis->GetBinUpEdge(bin)));
782  ok = axis->GetBinUpEdge(bin) > 100;
783  }
784  else {
785  name.ReplaceAll("ALICE", "");
786  name.ReplaceAll("dNdeta", " - work in progress");
787  }
788  // Printf("Adding entry %d: %s/%s", j, nme.Data(), name.Data());
789  TLegendEntry* e = l->AddEntry("dummy", name, "f");
790  e->SetFillStyle(1001);
791  e->SetFillColor(h->GetMarkerColor());
792  }
793  ret->Add(h);
794  }
795  return ret;
796  }
797 };
798 //
799 // EOF
800 //
const char * filename
Definition: TestFCM.C:1
THStack * DrawRes(TDirectory *top, const TString &base, Bool_t onlyMB)
TLatex * fParName
double Double_t
Definition: External.C:58
Base class for classes to draw summaries.
const char * title
Definition: MakeQAPdf.C:26
static TH1 * GetH1(const TObject *parent, const TString &name, Bool_t verb=true)
TAxis * GetCentAxis(const TCollection *parent, Bool_t verbose=false)
static TH2 * GetH2(const TObject *parent, const TString &name, Bool_t verb=true)
static void TriggerString(ULong_t trigger, TString &str)
TCanvas * c
Definition: TestFitELoss.C:172
TCollection * GetCentCollection(const TCollection *sums, const TString &base, Float_t cLow, Float_t cHigh, TString &title)
void DrawCentResDetails(const TCollection *sums, const TString &base)
void DrawBoth(THStack *rC, THStack *rF, THStack *rM, TFile *file, Bool_t onlyMB)
void Run(const char *fname="forward_dndeta.root", UShort_t flags=kNormal)
static void SysString(UShort_t sys, TString &str)
Int_t cH
Definition: Combine.C:26
void DrawSums(TDirectory *top, const TString &base, bool onlyMB)
int Int_t
Definition: External.C:63
float Float_t
Definition: External.C:68
const char * ColName(const char *prefix, bool results=false)
static Bool_t GetParameter(const TObject *c, const TString &name, Short_t &value, Bool_t verb=true)
static TObject * GetObject(const TObject *parent, const TString &name, Bool_t verb=true)
unsigned long ULong_t
Definition: External.C:38
const char * pdfName
Definition: DrawAnaELoss.C:30
void CreateCanvas(const TString &pname, Bool_t landscape=false, Bool_t pdf=true, Bool_t useTop=true)
void DrawParameter(Double_t &y, const TString &name, const TString &value, Double_t size=0)
void DrawResTitle(TCollection *c, Double_t &y, Bool_t onlyMB)
void PrintCanvas(const TString &title, Float_t size=.7)
static void SNNString(UShort_t sNN, TString &str)
void CloseCanvas()
THStack * MakeMBStack(const TCollection *sums, const TString &base)
TLatex * fParVal
Definition: External.C:220
void DrawTitlePage(TFile *file, Bool_t mc, Bool_t onlyMB)
TFile * file
unsigned short UShort_t
Definition: External.C:28
TObject * DrawInPad(TVirtualPad *c, Int_t padNo, TObject *h, Option_t *opts="", UInt_t flags=0x0, const char *title="")
bool Bool_t
Definition: External.C:53
static TCollection * GetCollection(const TObject *parent, const TString &name, Bool_t verb=true)
THStack * CleanStack(const THStack *stack, TLegend *l, const TAxis *axis)
Definition: External.C:196
void DrawCentRes(const TCollection *sums, const TString &base, Float_t cLow, Float_t cHigh)
static THStack * GetStack(const TObject *parent, const TString &name, const char *sub=0, Bool_t verb=true)
void DrawCentSum(const TCollection *sums, const TString &base, Float_t cLow, Float_t cHigh)