]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/qa/DrawRecAnaEloss.C
Fixed references from PWG2 -> PWGLF - very efficiently done using ETags.
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / qa / DrawRecAnaEloss.C
1 /**
2  * @file   DrawRecAnaEloss.C
3  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
4  * @date   Thu Jul  7 10:58:50 2011
5  * 
6  * @brief  Draw energ-loss before/after merging and used in the
7  * density calculations 
8  *
9  * @deprecated Use QATrender instead
10  * @ingroup pwglf_forward_scripts_qa
11  */
12
13 #ifndef __CINT__
14 #include <TList.h>
15 #include <TH1.h>
16 #include <TH2.h>
17 #include <TCanvas.h>
18 #include <TLatex.h>
19 #include <TLine.h>
20 #include <TFile.h>
21 #include <TError.h>
22 #include <TParameter.h>
23 #include <TStyle.h>
24 #else
25 class TLatex;
26 #endif
27
28 /** 
29  * Draw some text
30  * 
31  * @param l    LaTeX object
32  * @param x    x coordinate
33  * @param y    y coordinate (incremented on return)
34  * @param c1   First string
35  * @param c2   Second string 
36  *
37  * @deprecated Use QATrender instead
38  * @ingroup pwglf_forward_scripts_qa
39  */
40 void 
41 DrawText(TLatex* l, Double_t x, Double_t& y, const char* c1, const char* c2)
42 {
43   y -= 1.2*l->GetTextSize();
44   l->DrawLatex(x,    y, c1);
45   l->DrawLatex(x+.4, y, c2);
46 }
47 /** 
48  * Draw the energy loss before/after mergin for a single ring
49  * 
50  * @param p      List 1 
51  * @param p2     List 2
52  * @param lowCut Low cut
53  * @param d      Detector
54  * @param r      Ring 
55  *
56  * @deprecated Use QATrender instead
57  * @ingroup pwglf_forward_scripts_qa
58  */
59 void
60 DrawRingRecAnaEloss(TList* p, TList* p2, Double_t lowCut, UShort_t d, Char_t r)
61 {
62   if (!p) return;
63
64   TList* ring = static_cast<TList*>(p->FindObject(Form("FMD%d%c",d,r)));
65   if (!ring) { 
66     Error("DrawRecAnaEloss", "List FMD%d%c not found in %s",d,r,p->GetName());
67     return;
68   }
69   TList* ring2 = static_cast<TList*>(p2->FindObject(Form("FMD%d%c",d,r)));
70   if (!ring2){
71     Error("DrawRecAnaEloss","List FMD%d%c not found in %s",d,r,p2->GetName());
72     return;
73   }
74
75   TH1* before = static_cast<TH1D*>(ring->FindObject("esdEloss"));
76   if (!before) { 
77     Error("DrawRingRecAnaEloss", "Histogram esdEloss not found in FMD%d%c",
78           d, r);
79     return;
80   }
81   TH1* after = static_cast<TH1D*>(ring->FindObject("anaEloss"));
82   if (!after) { 
83     Error("DrawRingRecAnaEloss", "Histogram anaEloss not found in FMD%d%c",
84           d, r);
85     return;
86   }
87   TH1* presented = static_cast<TH1D*>(ring2->FindObject("eloss"));
88   if (!presented) {
89     Error("DrawRingRecAnaEloss", "Histogram eloss not found in FMD%d%c",
90           d, r);
91     return;
92   }
93   TH1* used = static_cast<TH1D*>(ring2->FindObject("elossUsed"));
94   if (!used) {
95     Error("DrawRingRecAnaEloss", "Histogram elossUsed not found in FMD%d%c",
96           d, r);
97     return;
98   }
99
100
101   Int_t low = before->GetXaxis()->FindBin(lowCut);
102   Int_t ib  = Int_t(before->Integral(low,before->GetNbinsX()));
103   Int_t ia  = Int_t(after->Integral(low,after->GetNbinsX()));
104   Int_t ip  = Int_t(presented->Integral(low,presented->GetNbinsX()));
105   Int_t iu  = Int_t(used->Integral(low,used->GetNbinsX()));
106   // before->Scale(1. / ib);
107   // after->Scale(1. / ib);
108   // presented->Scale(1. / ib);
109   // used->Scale(1. / ib);
110
111   gPad->SetLogy();
112   gPad->SetFillColor(0);
113   before->SetTitle(Form("FMD%d%c",d,r));
114   before->Draw("");
115   after->Draw("same");
116   presented->Draw("same");
117   used->Draw("same");
118   
119   // ib           = before->Integral(low,before->GetNbinsX());
120   // ia           = after->Integral(low,after->GetNbinsX());
121   // ip           = presented->Integral(low,presented->GetNbinsX());
122   // iu           = used->Integral(low,used->GetNbinsX());
123   Double_t ts  = 0.03;
124   Double_t x   = gPad->GetLeftMargin() + .25;
125   Double_t y   = 1-gPad->GetTopMargin()-gStyle->GetTitleH()+ts;
126   TLatex*  ltx = new TLatex(x, y, Form("FMD%d%c", d, r));
127   ltx->SetNDC();
128   ltx->SetTextAlign(13);
129   ltx->SetTextSize(ts);
130   // ltx->Draw();
131   // ltx->SetTextSize(.05);
132   TString inte(Form("Integral [%4.2f,#infty]", lowCut));
133   DrawText(ltx, x, y, Form("%s before:", inte.Data()), Form("%9d", ib));
134   DrawText(ltx, x, y, Form("%s after:",  inte.Data()), Form("%9d", ia));
135   DrawText(ltx, x, y, Form("%s input:",  inte.Data()), Form("%9d", ip));
136   DrawText(ltx, x, y, Form("%s user:",   inte.Data()), Form("%9d", iu));
137   TLine* l = new TLine;
138   l->SetLineWidth(1);
139   l->DrawLineNDC(x, y-0.9*ts, 1-gPad->GetRightMargin()-0.01, y-0.9*ts);
140   if (ib != 0 && ia != 0) {
141     DrawText(ltx, x, y, "Change (merging)", Form("%5.1f%%", (100.*(ia-ib))/ib));
142     DrawText(ltx, x, y, "Change (input)",   Form("%5.1f%% (%5.1f%%)", 
143                                                  (100.*(ip-ia))/ia,
144                                                  (100.*(ip-ib))/ib));
145     DrawText(ltx, x, y, "Change (use)",     Form("%5.1f%% (%5.1f%%)", 
146                                                  (100.*(iu-ip))/ip,
147                                                  (100.*(iu-ib))/ib));
148   }
149   before->GetXaxis()->SetRangeUser(0, 4);
150   gPad->cd();
151 }
152
153 /** 
154  * Draw energy loss before/after merging 
155  * 
156  * @param filename Input file name
157  * @param folder   Input folder name in file 
158  *
159  * @deprecated Use QATrender instead
160  * @ingroup pwglf_forward_scripts_qa
161  */
162 void
163 DrawRecAnaEloss(const char* filename="forward.root",
164                 const char* folder="ForwardResults")
165 {
166   gStyle->SetPalette(1);
167   gStyle->SetOptFit(0);
168   gStyle->SetOptStat(0);
169   gStyle->SetOptTitle(1);
170   gStyle->SetTitleW(.4);
171   gStyle->SetTitleH(.1);
172   gStyle->SetTitleColor(0);
173   gStyle->SetTitleStyle(0);
174   gStyle->SetTitleBorderSize(0);
175   gStyle->SetTitleX(.6);
176   
177   TFile* file = TFile::Open(filename, "READ");
178   if (!file) { 
179     Error("DrawRecAnaEloss", "failed to open %s", filename);
180     return;
181   }
182
183   TList* forward = static_cast<TList*>(file->Get(folder));
184   if (!forward) { 
185     Error("DrawRecAnaEloss", "List %s not found in %s", folder, filename);
186     return;
187   }
188
189   TList* sf = static_cast<TList*>(forward->FindObject("fmdSharingFilter"));
190   if (!sf) { 
191     Error("DrawRecAnaEloss", "List fmdSharingFilter not found in Forward");
192     return;
193   }
194   TList* dc = static_cast<TList*>(forward->FindObject("fmdDensityCalculator"));
195   if (!dc) {
196     Error("DrawRecAnaEloss","List fmdDensityCalculator not found in Forward");
197     return;
198   }
199
200   TParameter<double>* lowCut = 
201     static_cast<TParameter<double>*>(sf->FindObject("lowCut"));
202   Double_t low = (lowCut ? lowCut->GetVal() : 0.15);
203   if (!lowCut)
204     Warning("DrawRecAnaEloss", "Low cut not found in %s, assuming %f",
205             sf->GetName(), low);
206   TCanvas* c = new TCanvas("recAnaELoss", 
207                            "Reconstructed and Analysed energy loss", 900, 700);
208   c->SetFillColor(0);
209   c->SetBorderSize(0);
210   c->SetLeftMargin(0.15);
211   c->SetRightMargin(0.02);
212   c->SetTopMargin(0.02);
213   c->Divide(3, 2, 0, 0);
214   
215   c->cd(1); DrawRingRecAnaEloss(sf, dc, low, 1, 'I');
216   c->cd(2); DrawRingRecAnaEloss(sf, dc, low, 2, 'I');
217   c->cd(5); DrawRingRecAnaEloss(sf, dc, low, 2, 'O');
218   c->cd(3); DrawRingRecAnaEloss(sf, dc, low, 3, 'I');
219   c->cd(6); DrawRingRecAnaEloss(sf, dc, low, 3, 'O');
220   TVirtualPad* p = c->cd(4);
221   // p->SetTopMargin(0.05);
222   p->SetRightMargin(0.15);
223   p->SetFillColor(0);
224   TH2D* highCuts = static_cast<TH2D*>(sf->FindObject("highCuts"));
225   if (highCuts) highCuts->Draw("colz");
226   c->cd();
227   c->SaveAs("recAnaELoss.png");
228 }
229  
230 //
231 // EOF
232 //