]> git.uio.no Git - u/mrichter/AliRoot.git/blame - PWG2/FORWARD/analysis2/scripts/DrawCorrELoss.C
Fix up some drawing stuff
[u/mrichter/AliRoot.git] / PWG2 / FORWARD / analysis2 / scripts / DrawCorrELoss.C
CommitLineData
1c762251 1/**
2 * @file
3 *
4 * Scripts to draw energy loss fits from correction object file
5 *
6 * @ingroup pwg2_forward_analysis_scripts
7 */
8/**
9 * Clear canvas
10 *
11 * @param c Canvas to clear
12 *
13 * @ingroup pwg2_forward_analysis_scripts
14 */
cbd6464b 15void
16ClearCanvas(TCanvas* c)
17{
18 c->SetLeftMargin(.1);
19 c->SetRightMargin(.05);
20 c->SetBottomMargin(.1);
21 c->SetTopMargin(.05);
22 c->Clear();
23}
24
1c762251 25/**
26 * Draw energy loss fits to a multi-page PDF.
27 *
28 * @par Input:
29 * The input file is expected to contain a AliFMDCorrELossFit object
0f84fefb 30 * named @i elossfits in the top level directory.
1c762251 31 *
32 * @para Output:
33 * A multi-page PDF. Note, that the PDF generated by ROOT in this way
34 * is broken (cannot be read by Acrobat Reader on Windows and MacOSX)
35 * and one should pass it through a filter to correct these problems.
36 *
37 * @param fname File name
38 * @param option Drawing options
39 *
40 * @ingroup pwg2_forward_analysis_scripts
41 */
cbd6464b 42void
7dacdf84 43DrawCorrELoss(const char* fname, const char* option="summary error")
cbd6464b 44{
1c762251 45 //__________________________________________________________________
46 // Load libraries and object
cbd6464b 47 gROOT->Macro("$ALICE_ROOT/PWG2/FORWARD/analysis2/scripts/LoadLibs.C");
48
49 TFile* file = TFile::Open(fname, "READ");
50 if (!file) {
b0e36b4a 51 Error("DrawCorrELoss", "Failed to open %s", fname);
cbd6464b 52 return;
53 }
7dacdf84 54 TString pname(gSystem->BaseName(fname));
cbd6464b 55 pname.ReplaceAll(".root", ".pdf");
56
57 AliFMDCorrELossFit* fits =
58 static_cast<AliFMDCorrELossFit*>(file->Get("elossfits"));
59 if (!fits) {
b0e36b4a 60 Error("DrawCorrELoss", "Object 'elossfits' not found in %s", fname);
cbd6464b 61 return;
62 }
63
7dacdf84 64 //__________________________________________________________________
65 TString opts(option);
66 opts.ToLower();
67 Bool_t summary = opts.Contains("summary");
68 if (summary) opts.ReplaceAll("summary", "");
69
1c762251 70 //__________________________________________________________________
71 // Create a canvas
7dacdf84 72 Int_t h = 800;
73 Int_t w = h / TMath::Sqrt(2);
74 if (summary) w = h;
75 TCanvas* c = new TCanvas("c", "c", w, h);
cbd6464b 76 c->SetFillColor(0);
77 c->SetBorderSize(0);
78 c->SetBorderMode(0);
79 c->Print(Form("%s[", pname.Data()));
80
81 gStyle->SetOptStat(0);
7dacdf84 82 gStyle->SetTitleFillColor(0);
cbd6464b 83 gStyle->SetTitleStyle(0);
84 gStyle->SetTitleBorderSize(0);
85 gStyle->SetTitleX(.7);
86 gStyle->SetTitleY(1);
87 gStyle->SetTitleW(.3);
88 gStyle->SetTitleH(.1);
89 gStyle->SetFrameFillColor(kWhite);
90 gStyle->SetFrameBorderSize(1);
91 gStyle->SetFrameBorderMode(1);
92
93 ClearCanvas(c);
1c762251 94 //__________________________________________________________________
7dacdf84 95 if (!summary) {
96 // Create a title page
97 TLatex* ll = new TLatex(.5,.8, fname);
98 ll->SetTextAlign(22);
99 ll->SetTextSize(0.05);
100 ll->SetNDC();
101 ll->Draw();
102
103 TLatex* l = new TLatex(.5,.8, fname);
104 l->SetNDC();
105 l->SetTextSize(0.03);
106 l->SetTextFont(132);
107 l->SetTextAlign(12);
108 l->DrawLatex(0.2, 0.70, "1^{st} page is a summary of fit parameters");
109 l->DrawLatex(0.2, 0.67, "2^{nd} page is a summary of relative errors");
110 l->DrawLatex(0.2, 0.64, "Subsequent pages shows the fitted functions");
111 l->DrawLatex(0.3, 0.60, "Black line is the full fitted function");
112 l->DrawLatex(0.3, 0.57, "Coloured lines are the individual N-mip comp.");
113 l->DrawLatex(0.3, 0.54, "Full drawn lines correspond to used components");
114 l->DrawLatex(0.3, 0.51, "Dashed lines correspond to ignored components");
115 l->DrawLatex(0.2, 0.47, "Each component has the form");
116 l->DrawLatex(0.3, 0.42, "f_{n}(x; #Delta, #xi, #sigma') = "
117 "#int_{-#infty}^{+#infty}d#Delta' "
118 "landau(x; #Delta, #xi)gaus(x; #Delta', #sigma')");
119 l->DrawLatex(0.2, 0.37, "The full function is given by");
120 l->DrawLatex(0.3, 0.32, "f_{N}(x; #Delta, #xi, #sigma', #bf{a}) = "
121 "#sum_{i=1}^{N} a_{i} "
122 "f_{i}(x; #Delta_{i}, #xi_{i}, #sigma_{i}')");
123 l->DrawLatex(0.3, 0.26, "#Delta_{i} = i (#Delta_{1} + #xi_{1} log(i))");
124 l->DrawLatex(0.3, 0.23, "#xi_{i} = i #xi_{1}");
125 l->DrawLatex(0.3, 0.20, "#sigma_{i} = #sqrt{i} #sigma_{1}");
126 l->DrawLatex(0.3, 0.17, "#sigma_{n} #dot{=} 0");
127 l->DrawLatex(0.3, 0.14, "#sigma'^{2} = #sigma^{2}_{n} + #sigma^{2}");
128 l->DrawLatex(0.3, 0.11, "a_{1} = 1");
129 c->Print(pname.Data(), "Title:Title page");
130
131 ClearCanvas(c);
132 }
cbd6464b 133
1c762251 134 //__________________________________________________________________
135 // Draw overview page
7dacdf84 136 fits->Draw(opts.Data());
cbd6464b 137 c->Print(pname.Data(), "Title:Fit overview");
138
7dacdf84 139 if (summary) {
140 c->Print(Form("%s]", pname.Data()));
141 TString pngName(pname.Data());
142 pngName.ReplaceAll(".pdf", ".png");
143 c->Print(pngName.Data());
144 return;
145 }
cbd6464b 146 ClearCanvas(c);
7dacdf84 147
1c762251 148 //__________________________________________________________________
149 // Draw relative parameter errors
7dacdf84 150 fits->Draw("relative");
cbd6464b 151 c->Print(pname.Data(), "Title:Relative parameter errors");
152
1c762251 153 //__________________________________________________________________
154 // Draw all fits individually
cbd6464b 155 Int_t nPad = 6;
156 for (UShort_t d=1; d<=3; d++) {
157 UShort_t nQ = (d == 1 ? 1 : 2);
158 for (UShort_t q = 0; q < nQ; q++) {
159 Char_t r = (q == 0 ? 'I' : 'O');
160
161 TObjArray* ra = fits->GetRingArray(d, r);
162 if (!ra) continue;
163
164 AliFMDCorrELossFit::ELossFit* fit = 0;
165 TIter next(ra);
166 Int_t i = 0;
167 while ((fit = static_cast<AliFMDCorrELossFit::ELossFit*>(next()))) {
168 if ((i % nPad) == 0) {
169 ClearCanvas(c);
170 c->Divide(2,nPad/2,0,0);
171 }
172 TVirtualPad* p = c->cd((i % nPad) + 1);
173 p->SetLogy();
174 p->SetGridx();
175 p->SetGridy();
176 p->SetFillColor(kWhite);
177 fit->Draw("comp");
178
179 if ((i % nPad) == (nPad-1))
180 c->Print(pname.Data(),
181 Form("Title:FMD%d%c page %d", d, r, (i/nPad)+1));
182 i++;
183 }
184 if (i % nPad != 0)
185 c->Print(pname.Data(),
186 Form("Title:FMD%d%c page %d", d, r, (i/nPad)+1));
187 }
188 }
189
1c762251 190 //__________________________________________________________________
191 // Close output file
cbd6464b 192 c->Print(Form("%s]", pname.Data()));
193}
1c762251 194//
195// EOF
196//