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