]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWG2/FORWARD/analysis2/qa/DrawELossPoisson.C
27fdbe494781d14048e09bd10d3ad23a1509aae7
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / qa / DrawELossPoisson.C
1 /**
2  * @file   DrawELossPoisson.C
3  * @author Christian Holm Christensen <cholm@dalsgaard.hehi.nbi.dk>
4  * @date   Thu Jul  7 10:24:58 2011
5  * 
6  * @brief  A script to draw the Poisson vs Energy Loss correlation 
7  * 
8  * 
9  */
10
11 /** 
12  * Draw the poisson @f$N_{ch}@f$ estimate against the @f$\Delta@f$
13  * @f$N_{ch}@f$ estimate and do a regression line between the two 
14  * 
15  * @param p            List 
16  * @param d            Detector
17  * @param r            Ring
18  * @param xmin         Minimum
19  * @param xmax         Maximum
20  * 
21  * @return The regression coefficient 
22  */
23 Double_t
24 DrawRingELossPoisson(TList* p, UShort_t d, Char_t r, 
25                      Double_t xmin, Double_t xmax)
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("DrawELossPoisson", "List FMD%d%c not found in %s",d,r,p->GetName());
32     return;
33   }
34   
35   TH2* corr = static_cast<TH2D*>(ring->FindObject("elossVsPoisson"));
36   if (!corr) { 
37     Error("DrawRingELossPoisson", "Histogram esdEloss not found in FMD%d%c",
38           d, r);
39     return;
40   }
41   TPad* pad = static_cast<TPad*>(gPad);
42   pad->SetGridy();
43   pad->SetGridx();
44   pad->SetLogz();
45   if (d == 3) { 
46     pad->SetPad(pad->GetXlowNDC(), pad->GetYlowNDC(), .99, 
47                  pad->GetYlowNDC()+pad->GetHNDC());
48     pad->SetRightMargin(0.15);
49   }
50   pad->SetFillColor(0);
51   if (xmax < 0) xmax = corr->GetXaxis()->GetXmax();
52   corr->GetXaxis()->SetRangeUser(xmin,xmax);
53   corr->GetYaxis()->SetRangeUser(xmin,xmax);
54   corr->SetTitle(Form("FMD%d%c",d,r));
55   Info("", "Entries: %d, integral: %f", corr->GetEntries(), corr->Integral()); 
56   corr->Draw("colz");
57
58   // Calculate the linear regression 
59   Double_t dx    = (xmax-xmin);
60   Double_t rxy   = corr->GetCorrelationFactor();
61   Double_t sx    = corr->GetRMS(1);
62   Double_t sy    = corr->GetRMS(2);
63   Double_t sx2   = sx*sx;
64   Double_t sy2   = sy*sy;
65   Double_t cxy   = corr->GetCovariance();
66   Double_t mx    = corr->GetMean(1);
67   Double_t my    = corr->GetMean(2);
68   Double_t delta = 1;
69 #if 0
70   Double_t beta  = rxy * sy / sx;
71 #else
72   Double_t beta  = ((sy2 - delta*sx2 + 
73                      TMath::Sqrt(TMath::Power(sy2-delta*sx2,2) + 
74                                  4*delta*cxy*cxy))
75                     / 2 / cxy);
76 #endif
77   Double_t alpha = my - beta * mx;
78   TF1* f = new TF1("f", "[0]+[1]*x", xmin, xmax);
79   f->SetParameters(alpha, beta);
80   f->SetLineWidth(1);
81   f->SetLineStyle(1);
82   f->Draw("same");
83
84   TLine* l = new TLine(xmin,xmin,xmax,xmax);
85   l->SetLineWidth(1);
86   l->SetLineStyle(2);
87   l->SetLineColor(kBlack);
88   l->Draw();
89   // corr->GetXaxis()->SetRangeUser(0, 2);
90
91   Info("", "FMD%d%c correlation coefficient: %9.5f%% "
92        "line y = %f + %f * x", d, r, 100*rxy, alpha, beta);
93
94   Double_t x = pad->GetLeftMargin()+.01;
95   Double_t y = 1-pad->GetTopMargin()-.01;
96   TLatex* ltx = new TLatex(x, y, Form("FMD%d%c", d, r));
97   ltx->SetNDC();
98   ltx->SetTextAlign(13);
99   ltx->SetTextSize(0.08);
100   ltx->Draw();
101   y -= 0.12;
102   ltx->SetTextSize(0.06);
103   ltx->DrawLatex(x,y,"Deming regression: y=#alpha+#beta x");
104   x += .02;
105   y -= .06;
106   ltx->SetTextSize(0.05);
107   ltx->DrawLatex(x, y, Form("#alpha=%5.3f", alpha));
108   y -= .06;
109   ltx->DrawLatex(x, y, Form("#beta= %5.3f", beta));
110
111   TLinearFitter* fitter = new TLinearFitter(1);
112   fitter->SetFormula("1 ++ x");
113   for (Int_t i = 1; i <= corr->GetNbinsX(); i++) { 
114     x = corr->GetXaxis()->GetBinCenter(i);
115     if (x < xmin || x > xmax) continue;
116     for (Int_t j = 1; j <= corr->GetNbinsY(); j++) {
117       y = corr->GetYaxis()->GetBinCenter(j);
118       if (y < xmin || y > xmax) continue;
119       Double_t c = corr->GetBinContent(i,j);
120       if (c < .1) continue;
121       fitter->AddPoint(&x, y, c);
122     }
123   }
124   Bool_t robust = false;
125   if (robust)
126     fitter->EvalRobust();
127   else 
128     fitter->Eval();
129   for (Int_t i = 0; i < 2; i++) { 
130     std::cout << i << "\t" 
131               << fitter->GetParName(i) << "\t"
132               << fitter->GetParameter(i) << "\t";
133     if (!robust) 
134       std::cout << fitter->GetParError(i); 
135     std::cout << std::endl;
136   }
137   Double_t chi2 = fitter->GetChisquare();
138   Int_t    ndf  = (fitter->GetNpoints() - 
139                    fitter->GetNumberFreeParameters() );
140   std::cout << "chi2/ndf: " << chi2 << '/' << ndf 
141             << '=' << chi2 / ndf << std::endl;
142
143   // corr->Scale(1. / corr->GetMaximum());
144   pad->cd();
145   return corr->GetCorrelationFactor();
146 }
147
148 /** 
149  * Draw the correlation between the Poisson @f$N_{ch}@f$ estimate
150  * and the @f$\Delta@f$ @f$N_{ch}@f$ estimate and do a regression
151  * line between the two for each ring
152  * 
153  * @param filename File to read
154  * @param xmax     Minimum X
155  * @param xmin     Maximum X 
156  */
157 void
158 DrawELossPoisson(const char* filename="forward.root", 
159                  Double_t xmax=-1,
160                  Double_t xmin=-1)
161 {
162   gStyle->SetPalette(1);
163   gStyle->SetOptFit(0);
164   gStyle->SetOptStat(0);
165   gStyle->SetOptTitle(0);
166   gStyle->SetTitleW(.4);
167   gStyle->SetTitleH(.1);
168   gStyle->SetTitleX(.4);
169   // gStyle->SetTitleY(.1);
170   gStyle->SetTitleColor(0);
171   gStyle->SetTitleStyle(0);
172   gStyle->SetTitleBorderSize(0);
173   
174   TFile* file = TFile::Open(filename, "READ");
175   if (!file) { 
176     Error("DrawELossPoisson", "failed to open %s", filename);
177     return;
178   }
179
180   TList* forward = static_cast<TList*>(file->Get("Forward"));
181   if (!forward) { 
182     Error("DrawELossPoisson", "List Forward not found in %s", filename);
183     return;
184   }
185
186   TList* dc = static_cast<TList*>(forward->FindObject("fmdDensityCalculator"));
187   if (!dc) { 
188     Error("DrawELossPoisson", "List fmdDensityCalculator not found in Forward");
189     return;
190   }
191   
192   TCanvas* c = new TCanvas("elossVsPoisson", 
193                            "N_ch from ELoss versus from Poisson", 900, 700);
194   c->SetFillColor(0);
195   c->SetBorderSize(0);
196   c->SetBorderMode(0);
197   c->SetHighLightColor(0);
198   c->Divide(3, 2, 0, 0);
199
200   Double_t corrs[5];
201   c->cd(1); corrs[0] = DrawRingELossPoisson(dc, 1, 'I', xmin, xmax);
202   c->cd(2); corrs[1] = DrawRingELossPoisson(dc, 2, 'I', xmin, xmax);
203   c->cd(5); corrs[2] = DrawRingELossPoisson(dc, 2, 'O', xmin, xmax);
204   c->cd(3); corrs[3] = DrawRingELossPoisson(dc, 3, 'I', xmin, xmax);
205   c->cd(6); corrs[4] = DrawRingELossPoisson(dc, 3, 'O', xmin, xmax);
206
207   TVirtualPad* p = c->cd(4);
208   p->SetTopMargin(0.05);
209   p->SetRightMargin(0.10);
210   p->SetLeftMargin(0.15);
211   p->SetBottomMargin(0.15);
212   p->SetFillColor(0);
213
214   TH1D* hc = new TH1D("corrs", "Correlation factors", 5, .5, 5.5);
215   hc->SetFillColor(kRed+1);
216   hc->SetFillStyle(3001);
217   hc->SetMinimum(0.0);
218   hc->SetMaximum(1.3);
219   hc->GetXaxis()->SetBinLabel(1,"FMD1i"); hc->SetBinContent(1,corrs[0]);
220   hc->GetXaxis()->SetBinLabel(2,"FMD2i"); hc->SetBinContent(2,corrs[1]);
221   hc->GetXaxis()->SetBinLabel(3,"FMD2o"); hc->SetBinContent(3,corrs[2]);
222   hc->GetXaxis()->SetBinLabel(4,"FMD3i"); hc->SetBinContent(4,corrs[3]);
223   hc->GetXaxis()->SetBinLabel(5,"FMD3o"); hc->SetBinContent(5,corrs[4]);
224   hc->GetXaxis()->SetLabelSize(0.08);
225   hc->GetYaxis()->SetTitle("r");
226   hc->SetMarkerSize(1.5);
227   hc->Draw("text hist");
228
229   // TH2D* highCuts = static_cast<TH2D*>(dc->FindObject("highCuts"));
230   // if (highCuts) highCuts->Draw("colz");
231   c->cd();
232   c->SaveAs("elossVsPoisson.png");
233 }
234 //
235 // EOF
236 //