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