]> git.uio.no Git - u/mrichter/AliRoot.git/blob - PWGLF/FORWARD/analysis2/corrs/DrawCorrAcc.C
Various fixes to make the script easier to use
[u/mrichter/AliRoot.git] / PWGLF / FORWARD / analysis2 / corrs / DrawCorrAcc.C
1 /**
2  * @file 
3  * 
4  * Scripts to draw energy loss fits from correction object file 
5  *
6  * @ingroup pwglf_forward_scripts_corr
7  */
8 #ifndef __CINT__
9 # include <TH1.h>
10 # include <TH2.h>
11 # include <THStack.h>
12 # include <TObjArray.h>
13 # include <TList.h>
14 # include <TFile.h>
15 # include <TError.h>
16 # include <TCanvas.h>
17 # include <TLatex.h>
18 # include <TMath.h>
19 # include <TStyle.h>
20 # include "AliForwardCorrectionManager.h"
21 # include "AliFMDCorrAcceptance.h"
22 # include "AliForwardUtil.h"
23 #else
24 class TCanvas;
25 class TObjArray;
26 #endif
27 /** 
28  * Clear canvas 
29  * 
30  * @param c Canvas to clear 
31  *
32  * @ingroup pwglf_forward_scripts_corr
33  */
34 void
35 ClearCanvas(TCanvas* c)
36 {
37   c->SetLeftMargin(.1);
38   c->SetRightMargin(.05);
39   c->SetBottomMargin(.1);
40   c->SetTopMargin(.05);
41   c->Clear();
42 }
43
44 Bool_t
45 DrawSummary(TObjArray* stacks, TCanvas* c)
46 {
47   if (!stacks) return false;
48   // --- Make summary page -------------------------------------------
49   ClearCanvas(c);
50   Int_t nVtx = 10; // stacks->GetEntries();
51
52   c->Divide(3, (nVtx+2)/3, 0, 0);
53   Int_t ipad = 0;
54   for (UShort_t v = 1; v <= nVtx; v++) {
55     ipad++;
56     
57     if (ipad == 1 || ipad == 12) ipad++;
58
59     THStack*     stack = static_cast<THStack*>(stacks->At(v-1));
60     if (!stack) { 
61       Error("", "No stack at v=%d", v-1);
62       continue;
63     }
64     TVirtualPad* pad   = c->cd(ipad);
65     if (!pad) { 
66       Error("", "No pad at %d", ipad);
67       continue;
68     }
69     pad->SetFillColor(kWhite);
70     
71     stack->SetMaximum(1.2);
72     stack->Draw("nostack hist");
73   }
74   c->Modified();
75   c->Update();
76   c->cd();
77   return true;
78 }
79
80 /** 
81  * Draw energy loss fits to a multi-page PDF. 
82  *
83  * @par Input: 
84  * The input file is expected to contain a AliFMDCorrELossFit object
85  * named @c elossfits in the top level directory.
86  * 
87  * @par Output: 
88  * A multi-page PDF.  Note, that the PDF generated by ROOT in this way
89  * is broken (cannot be read by Acrobat Reader on Windows and MacOSX)
90  * and one should pass it through a filter to correct these problems.
91  * 
92  * @param fname   File name 
93  * @param option  Drawing options 
94  *
95  * @ingroup pwglf_forward_scripts_corr
96  */
97 void
98 DrawCorrAcc(const char* fname, const char* option="colz")
99 {
100   //__________________________________________________________________
101   // Load libraries and object 
102 #ifdef __CINT__
103   gROOT->Macro("$ALICE_ROOT/PWGLF/FORWARD/analysis2/scripts/LoadLibs.C");
104 #endif
105
106   TFile* file = TFile::Open(fname, "READ");
107   if (!file) { 
108     Error("DrawCorrAcc", "Failed to open %s", fname);
109     return;
110   }
111   TString pname(fname);
112   pname.ReplaceAll(".root", ".pdf");
113
114   const char* objName = 
115     AliForwardCorrectionManager::Instance()
116     .GetObjectName(AliForwardCorrectionManager::kAcceptance);
117   AliFMDCorrAcceptance* corr = 
118     static_cast<AliFMDCorrAcceptance*>(file->Get(objName));
119   if (!corr) { 
120     Error("DrawCorrAcc", "Object '%s' not found in %s", objName, fname);
121     return;
122   }
123
124   //__________________________________________________________________
125   // Create a canvas
126   TCanvas* c = new TCanvas("c", "c", 800 / TMath::Sqrt(2), 800);
127   c->SetFillColor(0);
128   c->SetBorderSize(0);
129   c->SetBorderMode(0);
130   c->Print(Form("%s[", pname.Data()));
131   
132   gStyle->SetOptStat(0);
133   gStyle->SetTitleColor(0);
134   gStyle->SetTitleStyle(0);
135   gStyle->SetTitleBorderSize(0);
136   gStyle->SetTitleX(.5);
137   gStyle->SetTitleY(1);
138   gStyle->SetTitleW(.8);
139   gStyle->SetTitleH(.09);
140   gStyle->SetFrameFillColor(kWhite);
141   gStyle->SetFrameBorderSize(1);
142   gStyle->SetFrameBorderMode(1);
143   gStyle->SetPalette(1);
144
145   ClearCanvas(c);
146   //__________________________________________________________________
147   // Create a title page 
148   TLatex* ll = new TLatex(.5,.8, fname);
149   ll->SetTextAlign(22);
150   ll->SetTextSize(0.03);
151   ll->SetNDC();
152   ll->Draw();
153
154   TLatex* l = new TLatex(.5,.8, fname);
155   l->SetNDC();
156   l->SetTextSize(0.03);
157   l->SetTextFont(132);
158   l->SetTextAlign(12);
159   l->DrawLatex(0.2, 0.70, "Acceptance due to dead channels");
160   l->SetTextAlign(22);
161   l->DrawLatex(0.5, 0.60, "c_{v,r}(#eta,#phi)=#frac{"
162                "#sum active strips#in(#eta,#phi)}{"
163                "#sum strips#in(#eta,#phi)}");
164   
165   c->Print(pname.Data(), "Title:Title page");
166
167   ClearCanvas(c);
168
169   //__________________________________________________________________
170   // Draw all corrections
171   const TAxis& vtxAxis = corr->GetVertexAxis();
172   Int_t        nVtx    = vtxAxis.GetNbins();
173
174   // --- Create stacks for summaries ---------------------------------
175   TObjArray* stacks  = new TObjArray(nVtx);
176   TObjArray* stacks2 = (corr->HasOverflow() ? new TObjArray(nVtx) : 0);
177   for (UShort_t v = 1; v <= nVtx; v++) { 
178     THStack* stack = new THStack(Form("vtx%02d", v),
179                                  Form("%+5.1f<v_{z}<%+5.1f",
180                                       vtxAxis.GetBinLowEdge(v),
181                                       vtxAxis.GetBinUpEdge(v)));
182     stacks->AddAt(stack, v-1);
183     if (!stacks2) continue;
184     stacks2->AddAt(stack->Clone(), v-1);
185   }
186
187   // --- Loop over detectors -----------------------------------------
188   for (UShort_t d = 1; d <= 3; d++) {
189     UShort_t     nQ = (d == 1 ? 1 : 2);
190     for (UShort_t q = 0; q < nQ; q++) { 
191       Char_t r = (q == 0 ? 'I' : 'O');
192
193       ClearCanvas(c);
194       c->Divide(2, (nVtx+1)/2);
195       for (UShort_t v=1; v <= nVtx; v++) { 
196         TVirtualPad* p = c->cd(v);
197         p->SetFillColor(kWhite);
198       
199         TH2* h2 = corr->GetCorrection(d, r, v);
200         if (!h2) { 
201           Warning("DrawCorrAcc", "No correction for r=%c, v=%d", r, v);
202           continue;
203         }
204         h2->Draw(option);
205
206         Int_t nY = h2->GetNbinsY();
207         TH1* hh = h2->ProjectionX(Form("FMD%d%c", d, r), 1, nY);
208         hh->Scale(1. / nY);
209         hh->SetDirectory(0);
210         hh->SetMarkerColor(AliForwardUtil::RingColor(d, r));
211         hh->SetLineColor(AliForwardUtil::RingColor(d, r));
212         hh->SetFillColor(AliForwardUtil::RingColor(d, r));
213         hh->SetFillStyle(3001);
214
215         THStack* stack = static_cast<THStack*>(stacks->At(v-1));
216         if (!stack) { 
217           Error("", "No stack at v=%d", v-1);
218           continue;
219         }
220         stack->Add(hh);
221
222         if (!stacks2) {
223           Warning("", "No phi acceptance defined");
224           continue;
225         }
226         stack = static_cast<THStack*>(stacks2->At(v-1));
227         if (!stack) { 
228           Error("", "No stack at v=%d", v-1);
229           continue;
230         }
231         TH1* hp = corr->GetPhiAcceptance(d, r, v);
232         if (!hp) { 
233           Error("", "No phi acceptance at v=%d", v-1);
234           continue;
235         }
236         hp->SetDirectory(0);
237         hp->SetMarkerColor(AliForwardUtil::RingColor(d, r));
238         hp->SetLineColor(AliForwardUtil::RingColor(d, r));
239         hp->SetFillColor(AliForwardUtil::RingColor(d, r));
240         hp->SetFillStyle(3001);
241         Info("", "Adding phi acceptance plot %d", hp->GetEntries());
242         stack->Add(hp);
243
244       }
245       c->Print(pname.Data(), Form("Title:FMD%d%c", d, r));
246     }
247   }
248
249   if (DrawSummary(stacks2, c))  
250     c->Print(pname.Data(), "Title:Summary2");
251   if (DrawSummary(stacks, c))  
252     c->Print(pname.Data(), "Title:Summary");
253
254   //__________________________________________________________________
255   // Close output file 
256   c->Print(Form("%s]", pname.Data()));
257 }
258 //
259 // EOF
260 //