]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/qa/DrawMCResult.C
Fixed references from PWG2 -> PWGLF - very efficiently done using ETags.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / qa / DrawMCResult.C
1 /**
2  * @file   DrawMCResult.C
3  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
4  * @date   Thu Jul  7 10:57:01 2011
5  * 
6  * @brief  Script to draw steps (deprecated version - use DrawSteps.C)
7  * 
8  * 
9  * @deprecated Use DrawSteps.C instead
10  * @ingroup pwglf_forward_scripts_qa
11  */
12 /** 
13  * 
14  * 
15  * @param forward Forward folder
16  * @param sub     Sub-folder to get
17  * @param name    Name of stack
18  * 
19  * @return Stack, or null
20  *
21  * @deprecated Use DrawSteps.C instead
22  * @ingroup pwglf_forward_scripts_qa
23  */
24 THStack*
25 GetStack(const TList& forward,  const char* sub, const char* name)
26 {
27   TList* lsub = static_cast<TList*>(forward.FindObject(sub));
28   if (!lsub) { 
29     Warning("GetStack", "Sub list %s not found in %s", sub, forward.GetName());
30     return 0;
31   }
32   THStack* ret = static_cast<THStack*>(lsub->FindObject(name));
33   if (!ret) 
34     Warning("GetStack" "Stack %s not found in %s", name, sub);
35   return ret;
36 }
37
38 /** 
39  * Rebin a histogram
40  * 
41  * @param h      Histogram
42  * @param rebin  Rebinning factor
43  * 
44  * @return Pointer to histogram
45  *
46  * @deprecated Use DrawSteps.C instead
47  * @ingroup pwglf_forward_scripts_qa
48  */
49 TH1* 
50 Rebin(TH1* h, Int_t rebin)
51 {
52   if (rebin <= 1) return h;
53   h->Rebin(rebin);
54   h->Scale(1. / rebin);
55   return h;
56 }
57
58 /** 
59  * Take ratio of two histograms
60  * 
61  * @param h1  Numerator
62  * @param h2  Denominator
63  * 
64  * @return Newly allocated histogram containg ratio
65  *
66  * @deprecated Use DrawSteps.C instead
67  * @ingroup pwglf_forward_scripts_qa
68  */
69 TH1*
70 Ratio(const TH1* h1, const TH1* h2)
71 {
72   if (!h1) return;
73   if (!h2) return;
74   
75   TH1* copy = static_cast<TH1*>(h2->Clone("tmp"));
76   copy->SetName(Form("%s_%s", h2->GetName(), h1->GetName()));
77   copy->SetTitle(Form("%s/%s", h2->GetTitle(), h1->GetTitle()));
78   copy->SetDirectory(0);
79   copy->Divide(h1);
80
81   return copy;
82 }
83 /** 
84  * Take ratio of histograms in stacks
85  * 
86  * @param r   Return stack of ratios
87  * @param h1  Numerators
88  * @param h2  Denominators
89  * 
90  * @return How many histograms in the return stack
91  *
92  * @deprecated Use DrawSteps.C instead
93  * @ingroup pwglf_forward_scripts_qa
94  */
95 Int_t 
96 Ratio(THStack* r, const THStack* h1, const THStack* h2)
97 {
98   if (!h1) return 0;
99   if (!h2) return 0;
100
101   int n1 = h1->GetHists()->GetEntries();
102   int n2 = h2->GetHists()->GetEntries();
103   int nH = 0;
104   for (int i = 0; i < n1 && i < n2; i++) { 
105     TH1* hh1 = static_cast<TH1*>(h1->GetHists()->At(i));
106     TH1* hh2 = static_cast<TH1*>(h2->GetHists()->At(i));
107     TH1* h   = Ratio(hh1, hh2);
108     if (!h) continue;
109     nH++;
110     r->Add(h);
111   }
112   return nH;
113 }
114 /** 
115  * Draw MC results
116  * 
117  * @param filename  Input file name
118  * @param rebin     Rebinning factor
119  * @param ratios    Whether to show ratios
120  *
121  * @deprecated Use DrawSteps.C instead
122  * @ingroup pwglf_forward_scripts_qa
123  */
124 void
125 DrawMCResult(const char* filename="forward.root", Int_t rebin=1,
126              Bool_t ratios=true)
127 {
128   gStyle->SetPalette(1);
129   gStyle->SetOptFit(0);
130   gStyle->SetOptStat(0);
131
132   TFile* file = TFile::Open(filename, "READ");
133   if (!file) { 
134     Error("DrawMCResult", "failed to open %s", filename);
135     return;
136   }
137
138   TList* forward = static_cast<TList*>(file->Get("Forward"));
139   if (!forward) { 
140     Error("DrawMCResult", "List Forward not found in %s", filename);
141     return;
142   }
143   THStack* res    = GetStack(*forward, "ringResults", "all");
144   THStack* mcRes  = GetStack(*forward, "mcRingResults", "all");
145   THStack* deltas = GetStack(*forward, "fmdSharingFilter", "sums");
146   THStack* nchs   = GetStack(*forward, "fmdDensityCalculator", "sums");
147   THStack* prims  = GetStack(*forward, "fmdCorrector", "sums");
148   
149   TH1* sumEta = static_cast<TH1*>(forward->FindObject("mcSumEta"));
150   if (!sumEta) { 
151     Warning("DrawMCResults", "mcSumEta not found in Forward");
152   }
153
154
155   gStyle->SetTitleBorderSize(0);
156   gStyle->SetTitleFillColor(0);
157   gStyle->SetTitleStyle(0);
158   // gStyle->SetTitleColor(kBlack);
159
160
161   TCanvas* c = new TCanvas("c", "C", 900, 700);
162   c->SetFillColor(0);
163   c->SetBorderSize(0);
164   c->SetTopMargin(0.05);
165   c->SetRightMargin(0.05);
166
167   Double_t y1 = (ratios ? .3 : 0);
168   TPad* p1 = new TPad("p1", "p1", 0, y1, 1, 1, 0, 0, 0); 
169   p1->SetBottomMargin(ratios ? 0.01 : .10);
170   p1->SetFillColor(0);
171   p1->SetBorderSize(0);
172   p1->SetTopMargin(0.05);
173   p1->SetRightMargin(0.05);
174   p1->Draw();
175   p1->cd();
176
177   THStack* all = new THStack("all", "Analysis steps");
178   if (res) {
179     res->SetTitle("dN_{ch}/d#eta");
180     TIter next(res->GetHists());
181     TH1* h = 0;
182     while ((h = static_cast<TH1*>(next()))) all->Add(Rebin(h,rebin));
183   }
184   if (mcRes) {
185     mcRes->SetTitle("Track-Refs");
186     TIter next(mcRes->GetHists());
187     TH1* h = 0;
188     while ((h = static_cast<TH1*>(next()))) all->Add(Rebin(h,rebin));
189   }
190   if (deltas) {
191     deltas->SetTitle("#sum_{} #Delta/#Delta_{mip}");
192     TIter next(deltas->GetHists());
193     TH1* h = 0;
194     while ((h = static_cast<TH1*>(next()))) { 
195       h->SetMarkerStyle(25);
196       all->Add(Rebin(h,rebin));
197     }
198   }
199   if (nchs) {
200     nchs->SetTitle("#sum_{} N_{ch,incl}");
201     TIter next(nchs->GetHists());
202     TH1* h = 0;
203     while ((h = static_cast<TH1*>(next()))) { 
204       h->SetMarkerStyle(21);
205       all->Add(Rebin(h,rebin));
206     }
207   }
208   if (prims) {
209     prims->SetTitle("#sum_{} N_{ch,primary}");
210     TIter next(prims->GetHists());
211     TH1* h = 0;
212     while ((h = static_cast<TH1*>(next()))) { 
213       h->SetMarkerStyle(22);
214       all->Add(Rebin(h,rebin));
215     }
216   }
217   if (sumEta) all->Add(sumEta);
218   all->Draw("nostack");
219   all->GetHistogram()->SetXTitle("#eta");
220   all->GetHistogram()->SetYTitle("signal");
221   all->GetHistogram()->GetXaxis()->SetLabelFont(132);
222   all->GetHistogram()->GetXaxis()->SetTitleFont(132);
223   all->GetHistogram()->GetYaxis()->SetLabelFont(132);
224   all->GetHistogram()->GetYaxis()->SetTitleFont(132);
225   c->SetGridx();
226
227   TLegend* l = new TLegend(.33, .2, .53, .9);
228   TLegendEntry* e = 0;
229   l->SetFillColor(0);
230   l->SetFillStyle(0);
231   l->SetBorderSize(0);
232   l->SetNColumns(1);
233   l->SetTextFont(132);
234   Int_t i = 0;
235   if (res) { 
236     TIter next(res->GetHists());
237     TH1* h = 0;
238     while ((h = static_cast<TH1*>(next()))) {
239       e = l->AddEntry(Form("dummy%02d", i++),h->GetTitle(),"pl");
240       e->SetMarkerStyle(20);
241       e->SetMarkerColor(h->GetMarkerColor());
242     }
243     e = l->AddEntry(Form("dummy%02d", i++),res->GetTitle(),"pl");
244     e->SetMarkerStyle(20);
245     e->SetMarkerColor(kBlack);
246   }
247   if (deltas) { 
248     e = l->AddEntry(Form("dummy%02d", i++), deltas->GetTitle(),"pl");
249     TH1* h = static_cast<TH1*>(deltas->GetHists()->At(0));
250     e->SetMarkerStyle(h->GetMarkerStyle());
251     e->SetMarkerColor(kBlack);
252   }
253   if (nchs) { 
254     e = l->AddEntry(Form("dummy%02d",i++),nchs->GetTitle(),"pl");
255     TH1* h = static_cast<TH1*>(nchs->GetHists()->At(0));
256     e->SetMarkerStyle(h->GetMarkerStyle());
257     e->SetMarkerColor(kBlack);
258   }
259   if (prims) { 
260     e = l->AddEntry(Form("dummy%02d", i++), prims->GetTitle(),"pl");
261     TH1* h = static_cast<TH1*>(prims->GetHists()->At(0));
262     e->SetMarkerStyle(h->GetMarkerStyle());
263     e->SetMarkerColor(kBlack);
264   }
265
266   if (mcRes) { 
267     e = l->AddEntry(Form("dummy%02d", i++), mcRes->GetTitle(), "pl");
268     TH1* h = static_cast<TH1*>(mcRes->GetHists()->At(0));
269     e->SetMarkerStyle(h->GetMarkerStyle());
270     e->SetMarkerColor(kBlack);
271   }
272
273   if (sumEta) l->AddEntry(sumEta);
274   l->Draw();
275
276
277   if (!ratios) return;
278
279   c->cd();
280   TPad* p2 = new TPad("p2", "p2", 0, 0, 1, y1, 0, 0, 0); 
281   p2->SetTopMargin(0);
282   p2->SetFillColor(0);
283   p2->SetBorderSize(0);
284   p2->SetRightMargin(0.05);
285   p2->Draw();
286   p2->cd();
287
288   THStack* rs = new THStack("ratios", "Ratios");
289   Int_t nDN = Ratio(rs, deltas, nchs);
290   Int_t nNR = Ratio(rs, nchs, res);
291   Int_t nRP = Ratio(rs, res, prims);
292   rs->Draw("nostack");
293
294   TLegend* ll = new TLegend(.38, .2, .48, .9);
295   ll->SetFillColor(0);
296   ll->SetFillStyle(0);
297   ll->SetBorderSize(0);
298   ll->SetNColumns(1);
299   ll->SetTextFont(132);
300   if (nDN) {
301     e = ll->AddEntry("d1", Form("#frac{%s}{%s}", 
302                                               nchs->GetTitle(), 
303                                               deltas->GetTitle()), "lp");
304     e->SetMarkerStyle(21);
305   }
306   if (nNR) {
307     e = ll->AddEntry("d2", Form("#frac{%s}{%s}", res->GetTitle(), 
308                                 nchs->GetTitle()), "lp");
309     e->SetMarkerStyle(20);
310   }
311   if (nRP) {
312     e = ll->AddEntry("d3", Form("#frac{%s}{%s}", prims->GetTitle(), 
313                             res->GetTitle()), "lp");
314     e->SetMarkerStyle(22);
315   }
316   ll->Draw();
317   
318
319 }
320  
321 //
322 // EOF
323 //