]>
Commit | Line | Data |
---|---|---|
e1f47419 | 1 | /** |
ffca499d | 2 | * @file DrawdNdeta.C |
0ccdab7b | 3 | * @author Christian Holm Christensen <cholm@nbi.dk> |
ffca499d | 4 | * @date Wed Mar 23 14:07:10 2011 |
5 | * | |
5bb5d1f6 | 6 | * @brief Script to visualise the dN/deta for pp and PbPb |
e1f47419 | 7 | * |
ffca499d | 8 | * This script is independent of any AliROOT code - and should stay |
9 | * that way. | |
10 | * | |
970b1a8a | 11 | * The script is <i>very</i> long - sigh - the joy of drawing |
5bb5d1f6 | 12 | * things nicely in ROOT |
ffca499d | 13 | * |
bd6f5206 | 14 | * @ingroup pwglf_forward_dndeta |
e1f47419 | 15 | */ |
b2e7f2d6 | 16 | #include <TH1.h> |
f709afcf | 17 | #include <TColor.h> |
b2e7f2d6 | 18 | #include <THStack.h> |
19 | #include <TGraphErrors.h> | |
20 | #include <TGraphAsymmErrors.h> | |
21 | #include <TMultiGraph.h> | |
22 | #include <TFile.h> | |
23 | #include <TList.h> | |
24 | #include <TString.h> | |
25 | #include <TError.h> | |
26 | #include <TSystem.h> | |
27 | #include <TROOT.h> | |
28 | #include <TMath.h> | |
29 | #include <TCanvas.h> | |
30 | #include <TPad.h> | |
31 | #include <TStyle.h> | |
32 | #include <TLegend.h> | |
33 | #include <TLegendEntry.h> | |
34 | #include <TLatex.h> | |
35 | #include <TImage.h> | |
33ab4a11 | 36 | #include <TRandom.h> |
bfab35d9 | 37 | #include <TParameter.h> |
7095962e | 38 | #include <TGClient.h> |
33ab4a11 | 39 | #include <fstream> |
298a7264 | 40 | #include <iostream> |
290052e7 | 41 | /** Systematic error color */ |
7095962e CHC |
42 | // #define SYSERR_COLOR kGray; |
43 | // #define SYSERR_COLOR kBlue-10 | |
44 | // #define SYSERR_COLOR kCyan-10 | |
45 | #define SYSERR_COLOR TColor::GetColor(220, 220, 255) | |
290052e7 | 46 | /** Systematic error style */ |
37079f20 | 47 | #define SYSERR_STYLE 1001 |
b2e7f2d6 | 48 | |
797161e8 | 49 | Double_t myFunc(Double_t* xp, Double_t* pp); |
50 | ||
fe0e93e6 | 51 | |
2f92ff0e | 52 | /** |
53 | * Class to draw dN/deta results | |
54 | * | |
bd6f5206 | 55 | * @ingroup pwglf_forward_tasks_dndeta |
56 | * @ingroup pwglf_forward_dndeta | |
2f92ff0e | 57 | */ |
b2e7f2d6 | 58 | struct dNdetaDrawer |
59 | { | |
7095962e | 60 | enum EFlags { |
fe0e93e6 | 61 | kShowRatios = 0x00001, |
62 | kShowLeftRight = 0x00002, | |
63 | kShowSysError = 0x00004, | |
64 | kShowRings = 0x00008, | |
65 | kCutEdges = 0x00010, | |
66 | kRemoveOuters = 0x00020, | |
67 | kUseFinalMC = 0x00040, | |
68 | kUseEmpirical = 0x00080, | |
69 | kForceMB = 0x00100, | |
70 | kMirror = 0x00200, | |
71 | kExport = 0x00400, | |
72 | kAddExec = 0x00800, | |
73 | kOldFormat = 0x01000, | |
74 | kVerbose = 0x02000, | |
75 | kHiRes = 0x04000, | |
76 | kExtraWhite = 0x08000, | |
77 | kLogo = 0x10000, | |
78 | kDefaultOptions = 0x1CE07 | |
7095962e CHC |
79 | }; |
80 | enum EOutFormat { | |
fe0e93e6 | 81 | kPNG = 0x1, |
82 | kPDF = 0x2, | |
83 | kROOT = 0x4, | |
84 | kScript = 0x8, | |
7095962e CHC |
85 | kAllFormats = 0xF |
86 | }; | |
fe0e93e6 | 87 | struct MarkerUtil |
88 | { | |
89 | /** | |
90 | * Marker styles | |
91 | */ | |
92 | enum { | |
93 | kSolid = 0x000, | |
94 | kHollow = 0x001, | |
95 | kCircle = 0x002, | |
96 | kSquare = 0x004, | |
97 | kUpTriangle = 0x006, | |
98 | kDownTriangle = 0x008, | |
99 | kDiamond = 0x00a, | |
100 | kCross = 0x00c, | |
101 | kStar = 0x00e | |
102 | }; | |
103 | /** | |
104 | * Get the marker style from option bits | |
105 | * | |
106 | * @param bits Option bits | |
107 | * | |
108 | * @return Marker style | |
109 | */ | |
110 | static Int_t GetMarkerStyle(UShort_t bits) | |
111 | { | |
112 | Int_t base = bits & (0xFE); | |
113 | Bool_t hollow = bits & kHollow; | |
114 | switch (base) { | |
115 | case kCircle: return (hollow ? 24 : 20); | |
116 | case kSquare: return (hollow ? 25 : 21); | |
117 | case kUpTriangle: return (hollow ? 26 : 22); | |
118 | case kDownTriangle: return (hollow ? 32 : 23); | |
119 | case kDiamond: return (hollow ? 27 : 33); | |
120 | case kCross: return (hollow ? 28 : 34); | |
121 | case kStar: return (hollow ? 30 : 29); | |
122 | } | |
123 | return 1; | |
124 | } | |
125 | /** | |
126 | * Get the marker option bits from a style | |
127 | * | |
128 | * @param style Style | |
129 | * | |
130 | * @return option bits | |
131 | */ | |
132 | static UShort_t GetMarkerBits(Int_t style) | |
133 | { | |
134 | UShort_t bits = 0; | |
135 | switch (style) { | |
136 | case 24: case 25: case 26: case 27: case 28: case 30: case 32: | |
137 | bits |= kHollow; break; | |
138 | } | |
139 | switch (style) { | |
140 | case 20: case 24: bits |= kCircle; break; | |
141 | case 21: case 25: bits |= kSquare; break; | |
142 | case 22: case 26: bits |= kUpTriangle; break; | |
143 | case 23: case 32: bits |= kDownTriangle; break; | |
144 | case 27: case 33: bits |= kDiamond; break; | |
145 | case 28: case 34: bits |= kCross; break; | |
146 | case 29: case 30: bits |= kStar; break; | |
147 | } | |
148 | return bits; | |
149 | } | |
150 | /** | |
151 | * Flip an option bit | |
152 | * | |
153 | * @param style Style parameter | |
154 | * | |
155 | * @return New style | |
156 | */ | |
157 | static Int_t FlipHollowStyle(Int_t style) | |
158 | { | |
159 | UShort_t bits = GetMarkerBits(style); | |
160 | Int_t ret = GetMarkerStyle(bits ^ kHollow); | |
161 | return ret; | |
162 | } | |
163 | }; | |
7095962e | 164 | |
2f92ff0e | 165 | /** |
166 | * POD of data for range zooming | |
167 | */ | |
b2e7f2d6 | 168 | struct RangeParam |
169 | { | |
170 | TAxis* fMasterAxis; // Master axis | |
171 | TAxis* fSlave1Axis; // First slave axis | |
172 | TVirtualPad* fSlave1Pad; // First slave pad | |
173 | TAxis* fSlave2Axis; // Second slave axis | |
174 | TVirtualPad* fSlave2Pad; // Second slave pad | |
175 | }; | |
176 | //__________________________________________________________________ | |
2f92ff0e | 177 | /** |
178 | * Constructor | |
179 | * | |
180 | */ | |
b2e7f2d6 | 181 | dNdetaDrawer() |
7095962e CHC |
182 | : fOptions(kDefaultOptions), |
183 | fFormats(kAllFormats), | |
184 | fShowOthers(0), // Show other data | |
185 | // Settings | |
186 | fRebin(0), // Rebinning factor | |
187 | fFwdSysErr(0.076), // Systematic error in forward range | |
188 | fCenSysErr(0), // Systematic error in central range | |
189 | fTitle(""), // Title on plot | |
190 | fBase(""), // Optional base name of output files | |
191 | fClusterScale(""), // Scaling of clusters to tracklets | |
192 | // Read (or set) information | |
193 | fTrigString(0), // Trigger string (read, or set) | |
194 | fNormString(0), // Normalisation string (read, or set) | |
195 | fSNNString(0), // Energy string (read, or set) | |
196 | fSysString(0), // Collision system string (read or set) | |
197 | fVtxAxis(0), // Vertex cuts (read or set) | |
198 | fCentAxis(0), // Centrality axis | |
199 | fTriggerEff(1), // Trigger efficency | |
200 | fExtTriggerEff(false), // True if fTriggerEff was read | |
f71c22f5 CHC |
201 | fCentMin(0), // Least centrality to plot |
202 | fCentMax(100), // Largest centrality to plot | |
7095962e CHC |
203 | // Resulting plots |
204 | fResults(0), // Stack of results | |
205 | fRatios(0), // Stack of ratios | |
206 | fLeftRight(0), // Left-right asymmetry | |
207 | fOthers(0), // Older data | |
208 | fTriggers(0), // Number of triggers | |
209 | fTruth(0), // Pointer to truth | |
210 | fRangeParam(0) // Parameter object for range zoom | |
b2e7f2d6 | 211 | { |
212 | fRangeParam = new RangeParam; | |
213 | fRangeParam->fMasterAxis = 0; | |
214 | fRangeParam->fSlave1Axis = 0; | |
215 | fRangeParam->fSlave1Pad = 0; | |
216 | fRangeParam->fSlave2Axis = 0; | |
217 | fRangeParam->fSlave2Pad = 0; | |
7095962e CHC |
218 | |
219 | TColor* sysErr = gROOT->GetColor(kSysErrColor); | |
220 | sysErr->SetAlpha(0.7); | |
b2e7f2d6 | 221 | } |
290052e7 | 222 | /** |
223 | * Cpoy constructor | |
224 | */ | |
33ab4a11 | 225 | dNdetaDrawer(const dNdetaDrawer&) {} |
290052e7 | 226 | /** |
227 | * Assignment operator | |
228 | * | |
229 | * | |
230 | * @return Reference to this object | |
231 | */ | |
33ab4a11 | 232 | dNdetaDrawer& operator=(const dNdetaDrawer&) { return *this; } |
233 | ||
ffca499d | 234 | //__________________________________________________________________ |
290052e7 | 235 | /** |
236 | * Destructor | |
237 | */ | |
ffca499d | 238 | virtual ~dNdetaDrawer() |
239 | { | |
240 | if (fRatios && fRatios->GetHists()) fRatios->GetHists()->Delete(); | |
241 | if (fResults && fResults->GetHists()) fResults->GetHists()->Delete(); | |
242 | ||
243 | if (fTrigString) { delete fTrigString; fTrigString = 0; } | |
244 | if (fSNNString) { delete fSNNString; fSNNString = 0; } | |
245 | if (fSysString) { delete fSysString; fSysString = 0; } | |
246 | if (fVtxAxis) { delete fVtxAxis; fVtxAxis = 0; } | |
247 | if (fCentAxis) { delete fCentAxis; fCentAxis = 0; } | |
248 | if (fResults) { delete fResults; fResults = 0; } | |
249 | if (fRatios) { delete fRatios; fRatios = 0; } | |
250 | if (fOthers) { delete fOthers; fOthers = 0; } | |
251 | if (fTriggers) { delete fTriggers; fTriggers = 0; } | |
252 | fRangeParam = 0; | |
253 | } | |
254 | ||
2f92ff0e | 255 | //================================================================== |
7095962e | 256 | void SetShowOthers(UInt_t others) { fShowOthers = others; } |
2f92ff0e | 257 | /** |
258 | * Set the rebinning factor | |
259 | * | |
260 | * @param x Rebinning factor (must be a divisor in the number of bins) | |
261 | */ | |
b2e7f2d6 | 262 | void SetRebin(UShort_t x) { fRebin = x; } |
2f92ff0e | 263 | //__________________________________________________________________ |
2f92ff0e | 264 | /** |
265 | * Set the title of the plot | |
266 | * | |
267 | * @param x Title | |
268 | */ | |
b2e7f2d6 | 269 | void SetTitle(TString x) { fTitle = x; } |
797161e8 | 270 | //__________________________________________________________________ |
d43c6cd0 | 271 | /** |
272 | * Set the base name of the output files | |
273 | * | |
274 | * @param x Base name | |
275 | */ | |
276 | void SetBase(TString x) { fBase = x; } | |
d43c6cd0 | 277 | //__________________________________________________________________ |
797161e8 | 278 | /** |
279 | * Set the systematic error in the forward region | |
280 | * | |
281 | * @param e Systematic error in the forward region | |
282 | */ | |
283 | void SetForwardSysError(Double_t e=0) { fFwdSysErr = e; } | |
284 | //__________________________________________________________________ | |
285 | /** | |
286 | * Set the systematic error in the forward region | |
287 | * | |
288 | * @param e Systematic error in the forward region | |
289 | */ | |
290 | void SetCentralSysError(Double_t e=0) { fCenSysErr = e; } | |
73b32206 | 291 | /** |
292 | * Set the 'Final MC' correction file. This is needed if the | |
293 | * secondary maps where produced using the old code | |
294 | * | |
295 | * @param file Filename | |
296 | */ | |
b615e211 | 297 | void SetFinalMC(const TString& file) { fFinalMC = file; } |
73b32206 | 298 | /** |
299 | * Set the file that contains the empirical correction. This is | |
300 | * needed when the secondary maps was generated with an in-accurate | |
301 | * geometry, and when we're analysing data at nominal interaction | |
302 | * points | |
303 | * | |
304 | * @param file Filename | |
305 | */ | |
b615e211 | 306 | void SetEmpirical(const TString& file) { fEmpirical = file; } |
2f92ff0e | 307 | /* @} */ |
308 | //================================================================== | |
309 | /** | |
310 | * @{ | |
311 | * @name Override settings from input | |
312 | */ | |
313 | /** | |
314 | * Override setting from file | |
315 | * | |
316 | * @param sNN Center of mass energy per nucleon pair (GeV) | |
317 | */ | |
318 | void SetSNN(UShort_t sNN) | |
319 | { | |
cdb9e20f | 320 | fSNNString = new TNamed("sNN", Form("%04dGeV", sNN)); |
2f92ff0e | 321 | fSNNString->SetUniqueID(sNN); |
322 | } | |
323 | //__________________________________________________________________ | |
324 | /** | |
325 | * Set the collision system | |
326 | * - 1: pp | |
327 | * - 2: PbPb | |
328 | * | |
329 | * @param sys collision system | |
330 | */ | |
331 | void SetSys(UShort_t sys) | |
332 | { | |
fe52e455 | 333 | fSysString = new TNamed("sys", (sys == 1 ? "pp" : |
4d97b4bb | 334 | sys == 2 ? "PbPb" : |
335 | sys == 3 ? "pPb" : | |
336 | "unknown")); | |
fe52e455 | 337 | fSysString->SetUniqueID(sys); |
2f92ff0e | 338 | } |
339 | //__________________________________________________________________ | |
340 | /** | |
341 | * Set the vertex range in centimeters | |
342 | * | |
343 | * @param vzMin Min @f$ v_z@f$ | |
344 | * @param vzMax Max @f$ v_z@f$ | |
345 | */ | |
346 | void SetVertexRange(Double_t vzMin, Double_t vzMax) | |
347 | { | |
348 | fVtxAxis = new TAxis(10, vzMin, vzMax); | |
349 | fVtxAxis->SetName("vtxAxis"); | |
350 | fVtxAxis->SetTitle(Form("v_{z}#in[%+5.1f,%+5.1f]cm", vzMin, vzMax)); | |
351 | } | |
b2e7f2d6 | 352 | //__________________________________________________________________ |
f71c22f5 CHC |
353 | /** |
354 | * Set the centrality range in centimeters | |
355 | * | |
356 | * @param centMin Min @f$ v_z@f$ | |
357 | * @param centMax Max @f$ v_z@f$ | |
358 | */ | |
359 | void SetCentralityRange(UShort_t centMin, UShort_t centMax) | |
360 | { | |
361 | fCentMin = centMin; | |
362 | fCentMax = centMax; | |
363 | } | |
364 | //__________________________________________________________________ | |
bfab35d9 | 365 | /** |
366 | * Set the trigger mask (overrides what's in the file) | |
367 | * | |
368 | * @param trig Trigger mask (0x1: INEL, 0x2: INEL>0, 0x4: NSD) | |
369 | */ | |
2f92ff0e | 370 | void SetTrigger(UShort_t trig) |
371 | { | |
372 | fTrigString = new TNamed("trigString", (trig & 0x1 ? "INEL" : | |
373 | trig & 0x2 ? "INEL>0" : | |
374 | trig & 0x4 ? "NSD" : | |
375 | "unknown")); | |
376 | fTrigString->SetUniqueID(trig); | |
377 | } | |
bfab35d9 | 378 | //__________________________________________________________________ |
379 | /** | |
380 | * Set the trigger efficiency - if set, then scale result histograms | |
381 | * by this factor | |
382 | * | |
383 | * @param eff @f$\varepsilon_{T}@f$ | |
384 | */ | |
d43c6cd0 | 385 | void SetTriggerEfficiency(Float_t eff) |
386 | { | |
387 | fTriggerEff = eff; | |
388 | fExtTriggerEff = false; | |
389 | } | |
2f92ff0e | 390 | |
391 | //================================================================== | |
392 | /** | |
393 | * @{ | |
394 | * @name Main steering functions | |
395 | */ | |
396 | /** | |
397 | * Run the code to produce the final result. | |
398 | * | |
399 | * @param filename File containing the data | |
400 | */ | |
7095962e CHC |
401 | void Run(const char* filename="forward_dndeta.root", |
402 | UInt_t flags=kDefaultOptions, UInt_t formats=kAllFormats) | |
b2e7f2d6 | 403 | { |
7095962e CHC |
404 | fOptions = flags; |
405 | fFormats = formats; | |
406 | SetForwardSysError(flags & kShowSysError ? 0.076 : 0); | |
407 | SetFinalMC (flags & kUseFinalMC ? "forward_dndetamc.root" : ""); | |
408 | SetEmpirical (flags & kUseEmpirical ? "EmpiricalCorrection.root" : ""); | |
409 | // SetBase(base); | |
410 | ||
e308a636 | 411 | Double_t max = 0, rmax=0, amax=0; |
b2e7f2d6 | 412 | |
e308a636 | 413 | gStyle->SetPalette(1); |
e1f47419 | 414 | |
e308a636 | 415 | // --- Open input file ------------------------------------------- |
416 | TFile* file = TFile::Open(filename, "READ"); | |
417 | if (!file) { | |
73b32206 | 418 | Error("Run", "Cannot open %s", filename); |
e308a636 | 419 | return; |
420 | } | |
73b32206 | 421 | Info("Run", "Drawing results from %s", file->GetName()); |
29befb1a | 422 | |
e308a636 | 423 | // --- Get forward list ------------------------------------------ |
c8b1a7db | 424 | TList* forward = static_cast<TList*>(file->Get("ForwarddNdetaResults")); |
e308a636 | 425 | if (!forward) { |
c8b1a7db | 426 | Error("Run", "Couldn't find list ForwarddNdetaResults"); |
e308a636 | 427 | return; |
428 | } | |
c8b1a7db | 429 | TList* sums = static_cast<TList*>(file->Get("ForwarddNdetaSums")); |
bfab35d9 | 430 | if (!sums) { |
c8b1a7db | 431 | Error("Run", "Couldn't find list ForwarddNdetaSums"); |
bfab35d9 | 432 | return; |
433 | } | |
434 | TParameter<bool>* p = | |
435 | static_cast<TParameter<bool>*>(sums->FindObject("empirical")); | |
436 | if (p && p->GetVal() && !fEmpirical.IsNull()) { | |
437 | Warning("Run", "Empirical correction already applied"); | |
438 | fEmpirical = "__task__"; | |
439 | } | |
e308a636 | 440 | // --- Get information on the run -------------------------------- |
441 | FetchInformation(forward); | |
73b32206 | 442 | |
443 | // --- Print settings -------------------------------------------- | |
444 | Info("Run", "Settings for the drawer:\n" | |
445 | " Show ratios: %5s\n" | |
446 | " Show Left/right: %5s\n" | |
447 | " Show rings: %5s\n" | |
448 | " Export to file: %5s\n" | |
449 | " Cut edges when rebinning: %5s\n" | |
450 | " Remove outer rings: %5s\n" | |
451 | " Mirror to un-covered regions: %5s\n" | |
452 | " Force minimum bias: %5s\n" | |
453 | " Show other results: 0x%03x\n" | |
454 | " Rebinning factor: %5d\n" | |
455 | " Forward systematic error: %5.1f%%\n" | |
456 | " Central systematic error: %5.1f%%\n" | |
457 | " Title on plot: %s\n" | |
458 | " Scaling of clusters to tracklets: %s\n" | |
459 | " Final MC correction file: %s\n" | |
d43c6cd0 | 460 | " Empirical correction file: %s", |
7095962e CHC |
461 | ((fOptions & kShowRatios) ? "yes" : "no"), |
462 | ((fOptions & kShowLeftRight) ? "yes" : "no"), | |
463 | ((fOptions & kShowRings) ? "yes" : "no"), | |
464 | ((fOptions & kExport) ? "yes" : "no"), | |
465 | ((fOptions & kCutEdges) ? "yes" : "no"), | |
466 | ((fOptions & kRemoveOuters) ? "yes" : "no"), | |
467 | ((fOptions & kMirror) ? "yes" : "no"), | |
468 | ((fOptions & kForceMB) ? "yes" : "no"), | |
73b32206 | 469 | fShowOthers, fRebin, (100*fFwdSysErr), (100*fCenSysErr), |
470 | fTitle.Data(), fClusterScale.Data(), fFinalMC.Data(), | |
d43c6cd0 | 471 | fEmpirical.Data()); |
73b32206 | 472 | |
e308a636 | 473 | // --- Set the macro pathand load other data script -------------- |
d43c6cd0 | 474 | TString savPath(gROOT->GetMacroPath()); |
bd6f5206 | 475 | gROOT->SetMacroPath(Form("%s:$(ALICE_ROOT)/PWGLF/FORWARD/analysis2", |
e308a636 | 476 | gROOT->GetMacroPath())); |
bfab35d9 | 477 | // Always recompile |
d43c6cd0 | 478 | if (!gROOT->GetClass("RefData")) |
a19faec0 | 479 | gROOT->LoadMacro("OtherData.C+"); |
d43c6cd0 | 480 | gROOT->SetMacroPath(savPath); |
e308a636 | 481 | |
482 | // --- Get the central results ----------------------------------- | |
c8b1a7db | 483 | TList* clusters = static_cast<TList*>(file->Get("CentraldNdetaResults")); |
484 | if (!clusters) Warning("Run", "Couldn't find list CentraldNdetaResults"); | |
b2e7f2d6 | 485 | |
9ecab72f | 486 | // --- Get the central results ----------------------------------- |
c8b1a7db | 487 | TList* mcTruth = static_cast<TList*>(file->Get("MCTruthdNdetaResults")); |
488 | if (!mcTruth) Warning("Run", "Couldn't find list MCTruthdNdetaResults"); | |
9ecab72f | 489 | |
e308a636 | 490 | // --- Make our containtes --------------------------------------- |
491 | fResults = new THStack("results", "Results"); | |
492 | fRatios = new THStack("ratios", "Ratios"); | |
493 | fLeftRight = new THStack("asymmetry", "Left-right asymmetry"); | |
494 | fOthers = new TMultiGraph(); | |
b615e211 | 495 | |
496 | // --- Try to open the final MC file, and find relevant lists ---- | |
497 | TList* forwardMC = 0; | |
29befb1a | 498 | // TList* centralMC = 0; |
b615e211 | 499 | if (!fFinalMC.IsNull()) { |
500 | TFile* finalMC = TFile::Open(fFinalMC, "READ"); | |
501 | if (!finalMC) { | |
502 | Warning("Run", "Failed to open file %s for final MC corrections", | |
503 | fFinalMC.Data()); | |
504 | } | |
505 | else { | |
c8b1a7db | 506 | forwardMC = static_cast<TList*>(finalMC->Get("ForwarddNdetaResults")); |
b615e211 | 507 | if (!forwardMC) |
c8b1a7db | 508 | Warning("Run","Couldn't find list ForwarddNdetaResults for final MC"); |
29befb1a | 509 | #if 0 |
c8b1a7db | 510 | centralMC = static_cast<TList*>(finalMC->Get("CentradNdetalResults")); |
29befb1a | 511 | if (!centralMC) |
c8b1a7db | 512 | Warning("Run","Couldn't find list CentraldNdetaResults for final MC"); |
29befb1a | 513 | #endif |
b615e211 | 514 | } |
515 | } | |
516 | if (!forwardMC) fFinalMC = ""; | |
517 | ||
518 | // --- Try to get the emperical correction ----------------------- | |
519 | TGraphErrors* empCorr = 0; | |
bfab35d9 | 520 | if (!fEmpirical.IsNull() && !fEmpirical.EqualTo("__task__")) { |
5934a3e3 | 521 | if (gSystem->AccessPathName(fEmpirical.Data())) { // Not found here |
522 | fEmpirical = | |
bfab35d9 | 523 | gSystem->ExpandPathName(Form("$ALICE_ROOT/PWGLF/FORWARD/" |
524 | "corrections/Empirical/%s", | |
525 | fEmpirical.Data())); | |
5934a3e3 | 526 | if (gSystem->AccessPathName(fEmpirical.Data())) { // Not found here |
73b32206 | 527 | Warning("Run", "Couldn't get empirical correction file"); |
5934a3e3 | 528 | fEmpirical = ""; |
529 | } | |
b615e211 | 530 | } |
5934a3e3 | 531 | if (!fEmpirical.IsNull()) { |
532 | TFile* empirical = TFile::Open(fEmpirical, "READ"); | |
533 | if (!empirical) { | |
73b32206 | 534 | Warning("Run", "couldn't open empirical correction file: %s", |
5934a3e3 | 535 | fEmpirical.Data()); |
536 | fEmpirical = ""; | |
537 | } | |
538 | const char* empPath = "fmdfull/average"; | |
539 | empCorr = static_cast<TGraphErrors*>(empirical->Get(empPath)); | |
540 | if (!empCorr) { | |
73b32206 | 541 | Warning("Run", "Didn't find the graph %s in %s", |
5934a3e3 | 542 | empPath, fEmpirical.Data()); |
543 | fEmpirical = ""; | |
544 | } | |
b615e211 | 545 | } |
546 | } | |
bfab35d9 | 547 | if (!empCorr && !fEmpirical.EqualTo("__task__")) fEmpirical = ""; |
b615e211 | 548 | |
e308a636 | 549 | // --- Loop over input data -------------------------------------- |
cd3b3781 | 550 | TObjArray truths; |
7095962e CHC |
551 | FetchTopResults(mcTruth, 0, 0, "MCTruth", max, rmax, amax,truths); |
552 | TObjArray* fwdA = FetchTopResults(forward, forwardMC, empCorr, "Forward", | |
b615e211 | 553 | max, rmax, amax,truths); |
7095962e | 554 | TObjArray* cenA = FetchTopResults(clusters, 0, 0, "Central", |
b615e211 | 555 | max, rmax, amax,truths); |
b2e7f2d6 | 556 | |
e308a636 | 557 | // --- Get trigger information ----------------------------------- |
bfab35d9 | 558 | // TList* sums = static_cast<TList*>(file->Get("ForwardSums")); |
e308a636 | 559 | if (sums) { |
7095962e CHC |
560 | TList* all = (fOptions & kOldFormat ? sums : |
561 | static_cast<TList*>(sums->FindObject("all"))); | |
e308a636 | 562 | if (all) { |
7095962e | 563 | fTriggers = FetchHistogram(all, "triggers"); |
e308a636 | 564 | if (!fTriggers) all->ls(); |
565 | } | |
566 | else { | |
567 | Warning("Run", "List all not found in ForwardSums"); | |
568 | sums->ls(); | |
569 | } | |
570 | } | |
571 | else { | |
572 | Warning("Run", "No ForwardSums directory found in %s", file->GetName()); | |
573 | file->ls(); | |
574 | } | |
575 | ||
576 | // --- Check our stacks ------------------------------------------ | |
577 | if (!fResults->GetHists() || | |
578 | fResults->GetHists()->GetEntries() <= 0) { | |
579 | Error("Run", "No histograms in result stack!"); | |
580 | return; | |
581 | } | |
582 | if (!fOthers->GetListOfGraphs() || | |
583 | fOthers->GetListOfGraphs()->GetEntries() <= 0) { | |
584 | Warning("Run", "No other data found - disabling that"); | |
797161e8 | 585 | fShowOthers = 0; |
e308a636 | 586 | } |
587 | if (!fRatios->GetHists() || | |
588 | fRatios->GetHists()->GetEntries() <= 0) { | |
589 | Warning("Run", "No ratio data found - disabling that"); | |
590 | // fRatios->ls(); | |
7095962e | 591 | fOptions &= ~kShowRatios; |
e308a636 | 592 | } |
593 | if (!fLeftRight->GetHists() || | |
594 | fLeftRight->GetHists()->GetEntries() <= 0) { | |
595 | Warning("Run", "No left/right data found - disabling that"); | |
596 | // fLeftRight->ls(); | |
7095962e | 597 | fOptions &= ~kShowLeftRight; |
e308a636 | 598 | } |
797161e8 | 599 | if (fFwdSysErr > 0) { |
600 | if (fCenSysErr <= 0) fCenSysErr = fFwdSysErr; | |
601 | for (Int_t i = 0; i < fwdA->GetEntriesFast(); i++) { | |
602 | TH1* fwd = static_cast<TH1*>(fwdA->At(i)); | |
603 | TH1* cen = static_cast<TH1*>(cenA->At(i)); | |
604 | CorrectForward(fwd); | |
605 | CorrectCentral(cen); | |
606 | Double_t low, high; | |
607 | TH1* tmp = Merge(cen, fwd, low, high); | |
608 | TF1* f = FitMerged(tmp, low, high); | |
609 | MakeSysError(tmp, cen, fwd, f); | |
610 | delete f; | |
7095962e CHC |
611 | if (fOptions & kVerbose) |
612 | Info("", "Adding systematic error histogram %s", tmp->GetName()); | |
797161e8 | 613 | fResults->GetHists()->AddFirst(tmp, "e5"); |
298a7264 | 614 | |
7095962e | 615 | if (!(fOptions & kMirror)) continue; |
298a7264 | 616 | |
797161e8 | 617 | TH1* tmp2 = Symmetrice(tmp); |
618 | tmp2->SetFillColor(tmp->GetFillColor()); | |
619 | tmp2->SetFillStyle(tmp->GetFillStyle()); | |
620 | tmp2->SetMarkerStyle(tmp->GetMarkerStyle()); | |
621 | tmp2->SetLineWidth(tmp->GetLineWidth()); | |
622 | fResults->GetHists()->AddFirst(tmp2, "e5"); | |
623 | fResults->Modified(); | |
624 | } | |
625 | } | |
626 | delete fwdA; | |
627 | delete cenA; | |
e308a636 | 628 | |
629 | // --- Close the input file -------------------------------------- | |
630 | file->Close(); | |
631 | ||
632 | ||
b2e7f2d6 | 633 | |
e308a636 | 634 | // --- Plot results ---------------------------------------------- |
635 | Plot(max, rmax, amax); | |
b2e7f2d6 | 636 | } |
e1f47419 | 637 | |
b2e7f2d6 | 638 | //__________________________________________________________________ |
2f92ff0e | 639 | /** |
e1f47419 | 640 | * Fetch the information on the run from the results list |
2f92ff0e | 641 | * |
e1f47419 | 642 | * @param results Results list |
2f92ff0e | 643 | */ |
e1f47419 | 644 | void FetchInformation(const TList* results) |
b2e7f2d6 | 645 | { |
2f92ff0e | 646 | if (!fTrigString) |
0be6c8cd | 647 | fTrigString = static_cast<TNamed*>(results->FindObject("trigger")); |
648 | if (!fNormString) | |
649 | fNormString = static_cast<TNamed*>(results->FindObject("scheme")); | |
2f92ff0e | 650 | if (!fSNNString) |
651 | fSNNString = static_cast<TNamed*>(results->FindObject("sNN")); | |
652 | if (!fSysString) | |
653 | fSysString = static_cast<TNamed*>(results->FindObject("sys")); | |
654 | if (!fVtxAxis) | |
655 | fVtxAxis = static_cast<TAxis*>(results->FindObject("vtxAxis")); | |
e308a636 | 656 | if (!fCentAxis) |
657 | fCentAxis = static_cast<TAxis*>(results->FindObject("centAxis")); | |
d43c6cd0 | 658 | if (fTriggerEff <= 0 || (1-fTriggerEff)<1e-6) { |
659 | TParameter<double>* eff = | |
660 | static_cast<TParameter<double>*>(results->FindObject("triggerEff")); | |
661 | if (eff) { | |
662 | fTriggerEff = eff->GetVal(); | |
663 | fExtTriggerEff = true; | |
664 | } | |
665 | if (fTriggerEff <= 0) SetTriggerEfficiency(1); | |
666 | } | |
c25b5e1b | 667 | |
668 | TNamed* options = static_cast<TAxis*>(results->FindObject("options")); | |
0be6c8cd | 669 | if (!fTrigString) fTrigString = new TNamed("trigger", "unknown"); |
670 | if (!fNormString) fNormString = new TNamed("scheme", "unknown"); | |
bad9a3c1 | 671 | if (!fSNNString) fSNNString = new TNamed("sNN", "unknown"); |
672 | if (!fSysString) fSysString = new TNamed("sys", "unknown"); | |
673 | if (!fVtxAxis) { | |
674 | fVtxAxis = new TAxis(1,0,0); | |
675 | fVtxAxis->SetName("vtxAxis"); | |
676 | fVtxAxis->SetTitle("v_{z} range unspecified"); | |
677 | } | |
f71c22f5 CHC |
678 | if (fCentAxis) { |
679 | TArrayD bins(fCentAxis->GetNbins()+1); | |
680 | Int_t nBins = 0; | |
681 | Double_t high = -1; | |
682 | for (Int_t i = 1; i <= fCentAxis->GetNbins(); i++) { | |
683 | Double_t binLow = fCentAxis->GetBinLowEdge(i); | |
684 | Double_t binHigh = fCentAxis->GetBinUpEdge(i); | |
685 | if (binLow < fCentMin-.5) continue; | |
686 | if (binHigh > fCentMax+.5) continue; | |
687 | high = binHigh; | |
688 | bins[nBins] = binLow; | |
689 | nBins++; | |
690 | } | |
691 | bins[nBins] = high; | |
692 | fCentAxis->Set(nBins, bins.GetArray()); | |
693 | } | |
694 | ||
b2e7f2d6 | 695 | |
7095962e CHC |
696 | if (fOptions & kVerbose) { |
697 | TString centTxt("none"); | |
698 | if (fCentAxis) { | |
699 | Int_t nCent = fCentAxis->GetNbins(); | |
700 | centTxt = Form("%d bins", nCent); | |
701 | for (Int_t i = 0; i <= nCent; i++) | |
702 | centTxt.Append(Form("%c%d", i == 0 ? ' ' : '-', | |
703 | int(fCentAxis->GetXbins()->At(i)))); | |
704 | } | |
705 | Info("FetchInformation", | |
706 | "Initialized for\n" | |
707 | " Trigger: %-30s (0x%x)\n" | |
708 | " Efficiency: %-6.4f\n" | |
709 | " sqrt(sNN): %-30s (%dGeV)\n" | |
710 | " System: %-30s (%d)\n" | |
711 | " Vz range: %-30s (%f,%f)\n" | |
712 | " Normalization: %-30s (%d)\n" | |
713 | " Centrality: %s\n" | |
714 | " Options: %s", | |
715 | fTrigString->GetTitle(), fTrigString->GetUniqueID(), | |
716 | fTriggerEff, | |
717 | fSNNString->GetTitle(), fSNNString->GetUniqueID(), | |
718 | fSysString->GetTitle(), fSysString->GetUniqueID(), | |
719 | fVtxAxis->GetTitle(), fVtxAxis->GetXmin(), fVtxAxis->GetXmax(), | |
720 | fNormString->GetTitle(), fNormString->GetUniqueID(), | |
721 | centTxt.Data(), (options ? options->GetTitle() : "none")); | |
722 | } | |
29befb1a | 723 | if (fSysString->GetUniqueID() == 3) { |
7095962e | 724 | Info("FetchTopResults", "Left/Right assymmetry, mirror, and systematic " |
29befb1a | 725 | "errors explicitly disabled for pPb"); |
7095962e CHC |
726 | fOptions &= ~kShowLeftRight; |
727 | fOptions &= ~kMirror; | |
728 | fFwdSysErr = 0; | |
729 | fCenSysErr = 0; | |
29befb1a | 730 | } |
e1f47419 | 731 | } |
e308a636 | 732 | //__________________________________________________________________ |
733 | TMultiGraph* FetchOthers(UShort_t centLow, UShort_t centHigh) | |
734 | { | |
735 | TMultiGraph* thisOther = 0; | |
797161e8 | 736 | if (fShowOthers == 0) return 0; |
e308a636 | 737 | |
e308a636 | 738 | UShort_t sys = (fSysString ? fSysString->GetUniqueID() : 0); |
739 | UShort_t trg = (fTrigString ? fTrigString->GetUniqueID() : 0); | |
740 | UShort_t snn = (fSNNString ? fSNNString->GetUniqueID() : 0); | |
d43c6cd0 | 741 | Long_t ret = |
742 | gROOT->ProcessLine(Form("RefData::GetData(%d,%d,%d,%d,%d,%d);", | |
743 | sys,snn,trg,centLow,centHigh,fShowOthers)); | |
7095962e CHC |
744 | if (!ret) { |
745 | Warning("FetchOthers", | |
746 | "No other data for %s %s %s %3d%%-%3d%% central (0x%x)", | |
747 | fSysString ? fSysString->GetTitle() : "unknown", | |
748 | fTrigString ? fTrigString->GetTitle() : "unknown", | |
749 | fSNNString ? fSNNString->GetTitle() : "unknown", | |
750 | centLow, centHigh, fShowOthers); | |
751 | return 0; | |
752 | } | |
797161e8 | 753 | |
754 | thisOther = reinterpret_cast<TMultiGraph*>(ret); | |
e308a636 | 755 | return thisOther; |
756 | } | |
e1f47419 | 757 | //__________________________________________________________________ |
758 | /** | |
7095962e | 759 | * Get the results from the top-level list (MC, SPD, FMD) |
e1f47419 | 760 | * |
66cf95f2 | 761 | * @param list List |
762 | * @param mcList List of histograms from MC | |
763 | * @param empCorr Emperical correction if any | |
764 | * @param name name | |
765 | * @param max On return, maximum of data | |
766 | * @param rmax On return, maximum of ratios | |
767 | * @param amax On return, maximum of left-right comparisons | |
768 | * @param truths List of MC truths to compare to. | |
290052e7 | 769 | * |
770 | * @return Array of results | |
e1f47419 | 771 | */ |
797161e8 | 772 | TObjArray* |
7095962e CHC |
773 | FetchTopResults(const TList* list, |
774 | const TList* mcList, | |
775 | TGraphErrors* empCorr, | |
776 | const char* name, | |
777 | Double_t& max, | |
778 | Double_t& rmax, | |
779 | Double_t& amax, | |
780 | TObjArray& truths) | |
e1f47419 | 781 | { |
3d532276 | 782 | if (!list) return 0; |
0ccdab7b | 783 | UShort_t n = HasCent() ? fCentAxis->GetNbins() : 1; |
784 | // // Info("FetchTopResults","got %d centrality bins", n); | |
785 | // if (n == 0) { | |
786 | // TH1* h = FetchOne(list, mcList, empCorr, name, "all", | |
787 | // FetchOthers(0,0), -1000, 0, | |
788 | // max, rmax, amax, truths); | |
789 | // if (!h) return 0; | |
790 | // TObjArray* a = new TObjArray; | |
791 | // // Info("FetchTopResults", "Adding %s to result stack", h->GetName()); | |
792 | // a->AddAt(h, 0); | |
793 | // return a; | |
794 | // } | |
e1f47419 | 795 | |
797161e8 | 796 | TObjArray* a = new TObjArray; |
cd3b3781 | 797 | truths.Expand(n); |
e308a636 | 798 | for (UShort_t i = 0; i < n; i++) { |
0ccdab7b | 799 | UShort_t centLow = 0; |
800 | UShort_t centHigh = 0; | |
801 | TString lname = "all"; | |
802 | Int_t col = -1000; | |
803 | TString centTxt = ""; | |
804 | if (HasCent()) { | |
805 | centLow = fCentAxis->GetBinLowEdge(i+1); | |
806 | centHigh = fCentAxis->GetBinUpEdge(i+1); | |
807 | lname = Form("cent%03d_%03d", centLow, centHigh); | |
808 | col = GetCentralityColor(i+1); | |
809 | centTxt = Form("%3d%%-%3d%% central", centLow, centHigh); | |
810 | } | |
cd3b3781 | 811 | TH1* tt = static_cast<TH1*>(truths.At(i)); |
812 | TH1* ot = tt; | |
b615e211 | 813 | TH1* h = FetchOne(list, mcList, empCorr, name, lname, |
814 | FetchOthers(centLow,centHigh), col, | |
0ccdab7b | 815 | centTxt.Data(), max, rmax, amax, tt); |
b615e211 | 816 | if (!h) continue; |
0ccdab7b | 817 | |
818 | if (tt != ot) { | |
cd3b3781 | 819 | truths.AddAt(tt, i); |
820 | } | |
7095962e | 821 | // Info("FetchTopResults", "Adding %p to result stack", h); |
797161e8 | 822 | a->AddAt(h, i); |
e1f47419 | 823 | } |
0ccdab7b | 824 | if (a->GetEntries() <= 0) { |
825 | delete a; | |
826 | a = 0; | |
827 | } | |
797161e8 | 828 | return a; |
e308a636 | 829 | } |
830 | //__________________________________________________________________ | |
290052e7 | 831 | /** |
7095962e | 832 | * Steer retrieval one centrality bin results |
290052e7 | 833 | * |
7095962e CHC |
834 | * @param list Input list |
835 | * @param mcList Possible MC list | |
836 | * @param empCorr Possible empirical correction | |
837 | * @param name Name of bing | |
838 | * @param folderName What sub-folder to get | |
839 | * @param others What else to plot | |
840 | * @param col Color | |
841 | * @param txt Centrality text | |
842 | * @param max Current maximum, on return new maximum | |
843 | * @param rmax Current range maximum, on return new maximum | |
844 | * @param amax Current A maximum, on return new maximum | |
845 | * @param truth Possible MC truth histogram | |
290052e7 | 846 | * |
7095962e | 847 | * @return |
290052e7 | 848 | */ |
b615e211 | 849 | TH1* FetchOne(const TList* list, |
850 | const TList* mcList, | |
851 | TGraphErrors* empCorr, | |
852 | const char* name, | |
853 | const char* folderName, | |
854 | TMultiGraph* others, | |
855 | Int_t col, | |
7095962e | 856 | const char* /* txt */, |
b615e211 | 857 | Double_t& max, |
858 | Double_t& rmax, | |
859 | Double_t& amax, | |
860 | TH1*& truth) | |
861 | { | |
7095962e CHC |
862 | TList* folder = (fOptions & kOldFormat ? const_cast<TList*>(list) : |
863 | static_cast<TList*>(list->FindObject(folderName))); | |
b615e211 | 864 | if (!folder) { |
865 | Error("FetchResults", "Couldn't find list '%s' in %s", | |
866 | folderName, list->GetName()); | |
867 | return 0; | |
868 | } | |
869 | TList* mcFolder = 0; | |
870 | if (mcList) { | |
871 | mcFolder = static_cast<TList*>(mcList->FindObject(folderName)); | |
872 | if (!mcFolder) | |
873 | Warning("FetchResults", | |
874 | "Didn't find the list '%s' in %s for final MC correction", | |
875 | folderName, mcList->GetName()); | |
876 | } | |
7095962e CHC |
877 | if (fOptions & kVerbose) { |
878 | TObject* normCalc = folder->FindObject("normCalc"); | |
879 | if (normCalc) Info("FetchOne", "%s:\n%s", | |
880 | folderName, normCalc->GetTitle()); | |
881 | } | |
882 | TH1* h = FetchCentResults(folder, mcFolder, empCorr, name, | |
883 | others, col, folderName, max, rmax, amax, truth); | |
b615e211 | 884 | return h; |
885 | } | |
e308a636 | 886 | //__________________________________________________________________ |
887 | /** | |
888 | * Fetch results for a particular centrality bin | |
889 | * | |
890 | * @param list List | |
66cf95f2 | 891 | * @param mcList List of MC results |
892 | * @param empCorr Emperical correction if any | |
e308a636 | 893 | * @param name Name |
ffca499d | 894 | * @param thisOther Other graphs |
895 | * @param color Color | |
896 | * @param centTxt Centrality text | |
e308a636 | 897 | * @param max On return, data maximum |
898 | * @param rmax On return, ratio maximum | |
899 | * @param amax On return, left-right maximum | |
66cf95f2 | 900 | * @param truth MC truth to compare to or possibly update |
290052e7 | 901 | * |
902 | * @return Histogram of results | |
e308a636 | 903 | */ |
7095962e CHC |
904 | TH1* FetchCentResults(const TList* list, |
905 | const TList* mcList, | |
906 | TGraphErrors* empCorr, | |
907 | const char* name, | |
908 | TMultiGraph* thisOther, | |
909 | Int_t color, | |
910 | const char* centTxt, | |
911 | Double_t& max, | |
912 | Double_t& rmax, | |
913 | Double_t& amax, | |
914 | TH1*& truth) | |
b2e7f2d6 | 915 | { |
cd3b3781 | 916 | |
7095962e CHC |
917 | TH1* dndeta = FetchHistogram(list, Form("dndeta%s", name)); |
918 | TH1* dndetaMC = FetchHistogram(list, Form("dndeta%sMC", name)); | |
919 | TH1* dndetaTruth = FetchHistogram(list, "dndetaTruth"); | |
0ccdab7b | 920 | Info("", "dN/deta truth from %s: %p", list->GetName(), dndetaTruth); |
921 | Info("", "dN/deta truth from external: %p", truth); | |
b615e211 | 922 | |
7095962e CHC |
923 | if (mcList && FetchHistogram(mcList, "finalMCCorr")) |
924 | Warning("FetchCentResults", "dNdeta already corrected for final MC"); | |
b615e211 | 925 | else |
926 | CorrectFinalMC(dndeta, mcList); | |
927 | ||
928 | CorrectEmpirical(dndeta, empCorr); | |
bfab35d9 | 929 | CorrectTriggerEff(dndeta); |
930 | CorrectTriggerEff(dndetaMC); | |
b615e211 | 931 | |
e1f47419 | 932 | TH1* dndetaSym = 0; |
933 | TH1* dndetaMCSym = 0; | |
e308a636 | 934 | SetAttributes(dndeta, color); |
298a7264 | 935 | SetAttributes(dndetaMC, HasCent() ? color : color+2); |
e308a636 | 936 | SetAttributes(dndetaTruth,color); |
937 | SetAttributes(dndetaSym, color); | |
298a7264 | 938 | SetAttributes(dndetaMCSym,HasCent() ? color : color+2); |
939 | if (dndetaMC && HasCent()) | |
5bb5d1f6 | 940 | dndetaMC->SetMarkerStyle(dndetaMC->GetMarkerStyle()+2); |
298a7264 | 941 | if (dndetaMCSym && HasCent()) |
5bb5d1f6 | 942 | dndetaMCSym->SetMarkerStyle(dndetaMCSym->GetMarkerStyle()+2); |
298a7264 | 943 | if (dndetaTruth && HasCent()) { |
5bb5d1f6 | 944 | dndetaTruth->SetMarkerStyle(34); |
9ecab72f | 945 | dndetaTruth->SetMarkerColor(kYellow-1); |
946 | } | |
298a7264 | 947 | if (dndetaTruth) { |
948 | dndetaTruth->SetLineColor(kBlack); | |
cd3b3781 | 949 | dndetaTruth->SetFillColor(kBlack); |
950 | dndetaTruth->SetFillStyle(3002); | |
951 | // dndetaTruth->SetLineColor(kBlack); | |
298a7264 | 952 | } |
e308a636 | 953 | ModifyTitle(dndeta, centTxt); |
954 | ModifyTitle(dndetaMC, centTxt); | |
955 | ModifyTitle(dndetaTruth,centTxt); | |
956 | ModifyTitle(dndetaSym, centTxt); | |
957 | ModifyTitle(dndetaMCSym,centTxt); | |
b615e211 | 958 | |
e1f47419 | 959 | |
cd3b3781 | 960 | max = TMath::Max(max, AddHistogram(fResults, dndetaTruth, "e5")); |
e308a636 | 961 | max = TMath::Max(max, AddHistogram(fResults, dndetaMC, dndetaMCSym)); |
962 | max = TMath::Max(max, AddHistogram(fResults, dndeta, dndetaSym)); | |
5bb5d1f6 | 963 | |
9ecab72f | 964 | if (dndetaTruth) { |
cd3b3781 | 965 | truth = dndetaTruth; |
5bb5d1f6 | 966 | } |
9ecab72f | 967 | else { |
7095962e | 968 | if ((fOptions & kShowRings)) { |
9ecab72f | 969 | THStack* rings = static_cast<THStack*>(list->FindObject("dndetaRings")); |
970 | if (rings) { | |
971 | TIter next(rings->GetHists()); | |
972 | TH1* hist = 0; | |
973 | while ((hist = static_cast<TH1*>(next()))) | |
974 | max = TMath::Max(max, AddHistogram(fResults, hist)); | |
ffca499d | 975 | } |
0ccdab7b | 976 | } // If show rings |
7095962e | 977 | // Info("FetchCentResults", "Got %p, %p, %p from %s with name %s, max=%f", |
9ecab72f | 978 | // dndeta, dndetaMC, dndetaTruth, list->GetName(), name, max); |
979 | ||
7095962e | 980 | if ((fOptions & kShowLeftRight)) { |
9ecab72f | 981 | fLeftRight->Add(Asymmetry(dndeta, amax)); |
982 | fLeftRight->Add(Asymmetry(dndetaMC, amax)); | |
9ecab72f | 983 | } |
984 | ||
985 | if (thisOther) { | |
986 | TIter next(thisOther->GetListOfGraphs()); | |
987 | TGraph* g = 0; | |
988 | while ((g = static_cast<TGraph*>(next()))) { | |
989 | fRatios->Add(Ratio(dndeta, g, rmax)); | |
990 | fRatios->Add(Ratio(dndetaSym, g, rmax)); | |
991 | SetAttributes(g, color); | |
992 | ModifyTitle(g, centTxt); | |
993 | if (!fOthers->GetListOfGraphs() || | |
994 | !fOthers->GetListOfGraphs()->FindObject(g->GetName())) { | |
995 | max = TMath::Max(max,TMath::MaxElement(g->GetN(), g->GetY())); | |
996 | fOthers->Add(g); | |
997 | } | |
998 | } | |
999 | // fOthers->Add(thisOther); | |
0ccdab7b | 1000 | } // if others for this |
1001 | } // if not truth | |
5bb5d1f6 | 1002 | if (dndetaMC) { |
1003 | fRatios->Add(Ratio(dndeta, dndetaMC, rmax)); | |
1004 | fRatios->Add(Ratio(dndetaSym, dndetaMCSym, rmax)); | |
1005 | } | |
cd3b3781 | 1006 | if (truth) { |
0ccdab7b | 1007 | Info("", "Forming ratio to truth:\n\t%s\n\t%s", |
1008 | dndeta->GetName(), | |
1009 | truth->GetName()); | |
cd3b3781 | 1010 | fRatios->Add(Ratio(dndeta, truth, rmax)); |
1011 | fRatios->Add(Ratio(dndetaSym, truth, rmax)); | |
b2e7f2d6 | 1012 | } |
797161e8 | 1013 | return dndeta; |
b2e7f2d6 | 1014 | } |
b615e211 | 1015 | //__________________________________________________________________ |
1016 | void CorrectFinalMC(TH1* dndeta, const TList* mcList) | |
1017 | { | |
1018 | if (!dndeta) return; | |
1019 | if (!mcList) return; | |
1020 | ||
7095962e CHC |
1021 | TH1* dndetaMC = FetchHistogram(mcList, dndeta->GetName()); |
1022 | TH1* dndetaTruth = FetchHistogram(mcList, "dndetaTruth"); | |
b615e211 | 1023 | if (!dndetaMC || !dndetaTruth) return; |
1024 | ||
1025 | TH1* corr = static_cast<TH1*>(dndetaMC->Clone("finalMCCorr")); | |
1026 | corr->Divide(dndetaTruth); | |
1027 | ||
1028 | Info("CorrectFinalMC", "Correcting dN/deta with final MC correction"); | |
1029 | dndeta->Divide(corr); | |
1030 | } | |
1031 | //__________________________________________________________________ | |
1032 | void CorrectEmpirical(TH1* dndeta, const TGraphErrors* empCorr) | |
1033 | { | |
1034 | if (!dndeta) return; | |
1035 | if (!empCorr) return; | |
1036 | ||
1037 | Info("CorrectEmpirical", "Doing empirical correction of dN/deta"); | |
1038 | TAxis* xAxis = dndeta->GetXaxis(); | |
1039 | for (Int_t i = 1; i <= xAxis->GetNbins(); i++) { | |
1040 | Double_t x = xAxis->GetBinCenter(i); | |
1041 | Double_t y = dndeta->GetBinContent(i); | |
1042 | Double_t c = empCorr->Eval(x); | |
1043 | dndeta->SetBinContent(i, y / c); | |
1044 | } | |
1045 | } | |
bfab35d9 | 1046 | //__________________________________________________________________ |
1047 | void CorrectTriggerEff(TH1* dndeta) | |
1048 | { | |
d43c6cd0 | 1049 | if (fExtTriggerEff) return; |
bfab35d9 | 1050 | if (!dndeta) return; |
1051 | if (fTriggerEff <= 0 || fTriggerEff >= 1) return; | |
1052 | dndeta->Scale(fTriggerEff); | |
1053 | } | |
b2e7f2d6 | 1054 | //__________________________________________________________________ |
2f92ff0e | 1055 | /** |
1056 | * Plot the results | |
2f92ff0e | 1057 | * @param max Max value |
2f92ff0e | 1058 | * @param rmax Maximum diviation from 1 of ratios |
2f92ff0e | 1059 | * @param amax Maximum diviation from 1 of asymmetries |
1060 | */ | |
e308a636 | 1061 | void Plot(Double_t max, |
b2e7f2d6 | 1062 | Double_t rmax, |
b2e7f2d6 | 1063 | Double_t amax) |
1064 | { | |
1065 | gStyle->SetOptTitle(0); | |
d43c6cd0 | 1066 | gStyle->SetTitleFont(kFont, "xyz"); |
1067 | gStyle->SetLabelFont(kFont, "xyz"); | |
b2e7f2d6 | 1068 | |
7095962e CHC |
1069 | Int_t h = (gROOT->IsBatch() ? |
1070 | ((fOptions & kHiRes) ? 1600 : 900) : | |
1071 | gClient->GetDisplayHeight()); | |
1072 | Int_t w = h; // h / TMath::Sqrt(2); | |
e308a636 | 1073 | Double_t y1 = 0; |
1074 | Double_t y2 = 0; | |
1075 | Double_t y3 = 0; | |
7095962e | 1076 | if (!(fOptions & kShowRatios)) w *= 1.3; |
e308a636 | 1077 | else y1 = 0.3; |
7095962e | 1078 | if (!(fOptions & kShowLeftRight)) w *= 1.3; |
e308a636 | 1079 | else { |
1080 | Double_t y11 = y1; | |
1081 | y1 = (y11 > 0.0001 ? 0.4 : 0.2); | |
ffca499d | 1082 | y2 = (y11 > 0.0001 ? 0.2 : 0.3); |
e308a636 | 1083 | } |
b2e7f2d6 | 1084 | TCanvas* c = new TCanvas("Results", "Results", w, h); |
1085 | c->SetFillColor(0); | |
1086 | c->SetBorderSize(0); | |
1087 | c->SetBorderMode(0); | |
1088 | ||
e308a636 | 1089 | PlotResults(max, y1); |
b2e7f2d6 | 1090 | c->cd(); |
1091 | ||
e308a636 | 1092 | PlotRatios(rmax, y2, y1); |
b2e7f2d6 | 1093 | c->cd( ); |
1094 | ||
e308a636 | 1095 | PlotLeftRight(amax, y3, y2); |
b2e7f2d6 | 1096 | c->cd(); |
bad9a3c1 | 1097 | |
1098 | ||
1099 | Int_t vMin = fVtxAxis->GetXmin(); | |
1100 | Int_t vMax = fVtxAxis->GetXmax(); | |
1101 | TString trg(fTrigString->GetTitle()); | |
e308a636 | 1102 | Int_t nev = 0; |
ffca499d | 1103 | if (fTriggers) nev = fTriggers->GetBinContent(1); |
7095962e | 1104 | if (HasCent()) trg = "CENT"; |
bad9a3c1 | 1105 | trg = trg.Strip(TString::kBoth); |
f709afcf | 1106 | trg.ReplaceAll(" ", "_"); |
1107 | trg.ReplaceAll(">", "Gt"); | |
29befb1a | 1108 | trg.ReplaceAll("&", "AND"); |
1109 | trg.ReplaceAll("|", "OR"); | |
d43c6cd0 | 1110 | if (fBase.IsNull()) |
1111 | fBase = "dndeta_<sys>_<snn>_<trig>_<ipmin><ipmax>cm_<nev>ev"; | |
1112 | fBase.ReplaceAll("<sys>", fSysString->GetTitle()); | |
1113 | fBase.ReplaceAll("<snn>", fSNNString->GetTitle()); | |
1114 | fBase.ReplaceAll("<trig>", trg.Data()); | |
1115 | fBase.ReplaceAll("<ipmin>", Form("%c%02d",vMin<0?'m':'p',TMath::Abs(vMin))); | |
1116 | fBase.ReplaceAll("<ipmax>", Form("%c%02d",vMax<0?'m':'p',TMath::Abs(vMax))); | |
1117 | fBase.ReplaceAll("<nev>", Form("%09d", nev)); | |
7095962e CHC |
1118 | if ((fFormats & kPNG)) c->SaveAs(Form("%s.png", fBase.Data())); |
1119 | if ((fFormats & kROOT)) c->SaveAs(Form("%s.root", fBase.Data())); | |
1120 | if ((fFormats & kScript))c->SaveAs(Form("%s.C", fBase.Data())); | |
1121 | if ((fFormats & kPDF)) c->SaveAs(Form("%s.pdf", fBase.Data())); | |
1122 | if (fOptions & kExport) { | |
1123 | TString exp(fBase); | |
1124 | exp.ReplaceAll("dndeta", "export"); | |
1125 | exp.ReplaceAll("dNdeta", "export"); | |
1126 | Export(exp); | |
1127 | } | |
b2e7f2d6 | 1128 | } |
fe0e93e6 | 1129 | //__________________________________________________________________ |
1130 | /** | |
1131 | * Plot the title on a pad | |
1132 | * | |
1133 | * @param p Pad to draw in | |
1134 | * @param bottom Bottom or top of pad | |
1135 | */ | |
1136 | void PlotTitle(TVirtualPad* p, Double_t yd, Bool_t bottom=true) | |
1137 | { | |
1138 | // Put a title on top | |
1139 | p->cd(); | |
1140 | fTitle.ReplaceAll("@", " "); | |
1141 | Double_t s = 1/yd/1.2; | |
1142 | TLatex* tit = new TLatex((bottom ? kRightMargin : p->GetLeftMargin()), | |
1143 | (bottom ? 0.01 : .99), fTitle.Data()); | |
1144 | tit->SetNDC(); | |
1145 | tit->SetTextFont(kFont); | |
1146 | tit->SetTextAlign(bottom ? 11 : 13); | |
1147 | tit->SetTextSize(s*0.045); | |
1148 | tit->SetTextColor(kAlicePurple); | |
1149 | tit->Draw(); | |
1150 | } | |
1151 | ||
b2e7f2d6 | 1152 | //__________________________________________________________________ |
e308a636 | 1153 | /** |
1154 | * Build main legend | |
1155 | * | |
66cf95f2 | 1156 | * @param stack Stack to include |
1157 | * @param mg (optional) Multi graph to include | |
1158 | * @param x1 Lower X coordinate in the range [0,1] | |
1159 | * @param y1 Lower Y coordinate in the range [0,1] | |
1160 | * @param x2 Upper X coordinate in the range [0,1] | |
1161 | * @param y2 Upper Y coordinate in the range [0,1] | |
1162 | * @param forceCol If non-zero, force this many columns | |
e308a636 | 1163 | */ |
ffca499d | 1164 | void BuildLegend(THStack* stack, TMultiGraph* mg, |
cd3b3781 | 1165 | Double_t x1, Double_t y1, Double_t x2, Double_t y2, |
1166 | Int_t forceCol=0) | |
e308a636 | 1167 | { |
1168 | TLegend* l = new TLegend(x1,y1,x2,y2); | |
cd3b3781 | 1169 | Int_t nCol = forceCol; |
1170 | if (nCol <= 0) nCol = HasCent() ? 1 : 2; | |
1171 | l->SetNColumns(nCol); | |
e308a636 | 1172 | l->SetFillColor(0); |
1173 | l->SetFillStyle(0); | |
1174 | l->SetBorderSize(0); | |
d43c6cd0 | 1175 | l->SetTextFont(kFont); |
1176 | l->SetTextColor(kAliceBlue); | |
e308a636 | 1177 | |
797161e8 | 1178 | // Loop over items in stack and get unique items, while ignoring |
1179 | // mirrored data and systematic error bands | |
ffca499d | 1180 | TIter next(stack->GetHists()); |
5bb5d1f6 | 1181 | TH1* hist = 0; |
1182 | TObjArray unique; | |
1183 | unique.SetOwner(); | |
797161e8 | 1184 | Bool_t sysErrSeen = false; |
fe0e93e6 | 1185 | Bool_t mirrorSeen = false; |
5bb5d1f6 | 1186 | while ((hist = static_cast<TH1*>(next()))) { |
797161e8 | 1187 | TString t(hist->GetTitle()); |
1188 | TString n(hist->GetName()); | |
1189 | n.ToLower(); | |
fe0e93e6 | 1190 | if (t.Contains("mirrored")) { mirrorSeen = true; continue; } |
797161e8 | 1191 | if (n.Contains("syserror")) { sysErrSeen = true; continue; } |
1192 | if (unique.FindObject(t.Data())) continue; | |
298a7264 | 1193 | TObjString* s1 = new TObjString(hist->GetTitle()); |
1194 | s1->SetUniqueID(((hist->GetMarkerStyle() & 0xFFFF) << 16) | | |
1195 | ((hist->GetMarkerColor() & 0xFFFF) << 0)); | |
1196 | unique.Add(s1); | |
5bb5d1f6 | 1197 | // l->AddEntry(hist, hist->GetTitle(), "pl"); |
e308a636 | 1198 | } |
ffca499d | 1199 | if (mg) { |
797161e8 | 1200 | // If we have other stuff, scan for unique names |
ffca499d | 1201 | TIter nexto(mg->GetListOfGraphs()); |
5bb5d1f6 | 1202 | TGraph* g = 0; |
1203 | while ((g = static_cast<TGraph*>(nexto()))) { | |
1204 | TString n(g->GetTitle()); | |
ffca499d | 1205 | if (n.Contains("mirrored")) continue; |
5bb5d1f6 | 1206 | if (unique.FindObject(n.Data())) continue; |
298a7264 | 1207 | TObjString* s2 = new TObjString(n); |
1208 | s2->SetUniqueID(((g->GetMarkerStyle() & 0xFFFF) << 16) | | |
1209 | ((g->GetMarkerColor() & 0xFFFF) << 0)); | |
1210 | unique.Add(s2); | |
5bb5d1f6 | 1211 | // l->AddEntry(hist, hist->GetTitle(), "pl"); |
ffca499d | 1212 | } |
e308a636 | 1213 | } |
797161e8 | 1214 | |
1215 | // Add legend entries for unique items only | |
5bb5d1f6 | 1216 | TIter nextu(&unique); |
1217 | TObject* s = 0; | |
1218 | Int_t i = 0; | |
1219 | while ((s = nextu())) { | |
1220 | TLegendEntry* dd = l->AddEntry(Form("data%2d", i++), | |
1221 | s->GetName(), "lp"); | |
1222 | Int_t style = (s->GetUniqueID() >> 16) & 0xFFFF; | |
1223 | Int_t color = (s->GetUniqueID() >> 0) & 0xFFFF; | |
1224 | dd->SetLineColor(kBlack); | |
298a7264 | 1225 | if (HasCent()) dd->SetMarkerColor(kBlack); |
5bb5d1f6 | 1226 | else dd->SetMarkerColor(color); |
1227 | dd->SetMarkerStyle(style); | |
1228 | } | |
797161e8 | 1229 | if (sysErrSeen) { |
1230 | // Add entry for systematic errors | |
1231 | TLegendEntry* d0 = l->AddEntry("d0", Form("%4.1f%% Systematic error", | |
1232 | 100*fFwdSysErr), "f"); | |
7095962e CHC |
1233 | d0->SetLineColor(kSysErrColor); |
1234 | d0->SetMarkerColor(kSysErrColor); | |
1235 | d0->SetFillColor(kSysErrColor); | |
e42bc740 | 1236 | d0->SetFillStyle(SYSERR_STYLE); |
797161e8 | 1237 | d0->SetMarkerStyle(0); |
1238 | d0->SetLineWidth(0); | |
1239 | i++; | |
1240 | } | |
cd3b3781 | 1241 | if (nCol == 2 && i % 2 == 1) { |
797161e8 | 1242 | // To make sure the 'data' and 'mirrored' entries are on a line |
1243 | // by themselves | |
1244 | TLegendEntry* dd = l->AddEntry("dd", " ", ""); | |
1245 | dd->SetTextSize(0); | |
1246 | dd->SetFillStyle(0); | |
1247 | dd->SetFillColor(0); | |
1248 | dd->SetLineWidth(0); | |
1249 | dd->SetLineColor(0); | |
1250 | dd->SetMarkerSize(0); | |
1251 | } | |
fe0e93e6 | 1252 | if (mirrorSeen /* (fOptions & kMirror) */) { |
cd3b3781 | 1253 | // Add entry for 'data' |
1254 | TLegendEntry* d1 = l->AddEntry("d1", "Data", "lp"); | |
1255 | d1->SetLineColor(kBlack); | |
1256 | d1->SetMarkerColor(kBlack); | |
1257 | d1->SetMarkerStyle(20); | |
1258 | ||
1259 | // Add entry for 'mirrored data' | |
298a7264 | 1260 | TLegendEntry* d2 = l->AddEntry("d2", "Mirrored data", "lp"); |
1261 | d2->SetLineColor(kBlack); | |
1262 | d2->SetMarkerColor(kBlack); | |
1263 | d2->SetMarkerStyle(24); | |
1264 | } | |
e308a636 | 1265 | |
1266 | l->Draw(); | |
1267 | } | |
1268 | //__________________________________________________________________ | |
5bb5d1f6 | 1269 | /** |
1270 | * Build centrality legend | |
1271 | * | |
5bb5d1f6 | 1272 | * @param x1 Lower X coordinate in the range [0,1] |
1273 | * @param y1 Lower Y coordinate in the range [0,1] | |
1274 | * @param x2 Upper X coordinate in the range [0,1] | |
1275 | * @param y2 Upper Y coordinate in the range [0,1] | |
1276 | */ | |
1277 | void BuildCentLegend(Double_t x1, Double_t y1, Double_t x2, Double_t y2) | |
1278 | { | |
298a7264 | 1279 | if (!HasCent()) return; |
5bb5d1f6 | 1280 | |
fe0e93e6 | 1281 | if (fCentAxis->GetNbins() <= 4) y1 += .15; |
5bb5d1f6 | 1282 | TLegend* l = new TLegend(x1,y1,x2,y2); |
1283 | l->SetNColumns(1); | |
1284 | l->SetFillColor(0); | |
1285 | l->SetFillStyle(0); | |
1286 | l->SetBorderSize(0); | |
d43c6cd0 | 1287 | l->SetTextFont(kFont); |
1288 | l->SetTextColor(kAliceBlue); | |
5bb5d1f6 | 1289 | |
1290 | Int_t n = fCentAxis->GetNbins(); | |
1291 | for (Int_t i = 1; i <= n; i++) { | |
1292 | Double_t low = fCentAxis->GetBinLowEdge(i); | |
1293 | Double_t upp = fCentAxis->GetBinUpEdge(i); | |
1294 | TLegendEntry* e = l->AddEntry(Form("dummy%02d", i), | |
1295 | Form("%3d%% - %3d%%", | |
1296 | int(low), int(upp)), "pl"); | |
1297 | e->SetMarkerColor(GetCentralityColor(i)); | |
1298 | } | |
1299 | l->Draw(); | |
1300 | } | |
1301 | //__________________________________________________________________ | |
fe0e93e6 | 1302 | void AttachExec(TVirtualPad* p, const char* plot, UShort_t id, |
1303 | Bool_t isBottom) | |
1304 | { | |
1305 | if (!(fOptions & kAddExec)) return; | |
1306 | ||
1307 | if (isBottom) { | |
1308 | fRangeParam->fMasterAxis = FindXAxis(p, plot); | |
1309 | p->AddExec("range", Form("RangeExec((dNdetaDrawer::RangeParam*)%p)", | |
1310 | fRangeParam)); | |
1311 | } | |
1312 | else { | |
1313 | if (id == 1) { | |
1314 | fRangeParam->fSlave1Axis = FindXAxis(p, plot); | |
1315 | fRangeParam->fSlave1Pad = p; | |
1316 | } | |
1317 | else if (id == 2) { | |
1318 | fRangeParam->fSlave2Axis = FindXAxis(p, plot); | |
1319 | fRangeParam->fSlave2Pad = p; | |
1320 | } | |
1321 | } | |
1322 | } | |
1323 | //__________________________________________________________________ | |
2f92ff0e | 1324 | /** |
1325 | * Plot the results | |
1326 | * | |
2f92ff0e | 1327 | * @param max Maximum |
1328 | * @param yd Bottom position of pad | |
1329 | */ | |
e308a636 | 1330 | void PlotResults(Double_t max, Double_t yd) |
b2e7f2d6 | 1331 | { |
fe0e93e6 | 1332 | // --- Make a sub-pad for the result itself ---------------------- |
b2e7f2d6 | 1333 | TPad* p1 = new TPad("p1", "p1", 0, yd, 1.0, 1.0, 0, 0, 0); |
fe0e93e6 | 1334 | p1->SetTopMargin(kRightMargin); |
b2e7f2d6 | 1335 | p1->SetBorderSize(0); |
1336 | p1->SetBorderMode(0); | |
1337 | p1->SetBottomMargin(yd > 0.001 ? 0.001 : 0.1); | |
d43c6cd0 | 1338 | p1->SetRightMargin(kRightMargin); |
7095962e | 1339 | if ((fOptions & kShowLeftRight) || (fOptions & kShowRatios)) p1->SetGridx(); |
b2e7f2d6 | 1340 | p1->SetTicks(1,1); |
1341 | p1->SetNumber(1); | |
1342 | p1->Draw(); | |
1343 | p1->cd(); | |
33ab4a11 | 1344 | |
fe0e93e6 | 1345 | // --- Figure out min/max ---------------------------------------- |
7095962e CHC |
1346 | // Info("PlotResults", "Plotting results with max=%f", max); |
1347 | fResults->SetMaximum((fOptions & kExtraWhite ? 1.4 : 1.15)*max); | |
cd3b3781 | 1348 | fResults->SetMinimum(yd > 0.00001 ? -0.02*max : 0); |
73b32206 | 1349 | // fResults->SetMinimum(yd > 0.00001 ? -0.02*max : 0); |
b2e7f2d6 | 1350 | |
fe0e93e6 | 1351 | // --- Fix up axis ----------------------------------------------- |
1352 | Double_t yyd = (1-yd)*(yd > .001 ? 1 : .9 / 1.2); | |
1353 | FixAxis(fResults, yyd, | |
1354 | "1/#it{N}#kern[.1]{d#it{N}_{ch}/d#it{#eta}}" | |
1355 | //"#frac{1}{#it{N}}#kern[.1]{#frac{d#it{N}_{ch}}{d#it{#eta}}}" | |
1356 | ); | |
b2e7f2d6 | 1357 | |
fe0e93e6 | 1358 | // --- Clear pad and re-draw ------------------------------------ |
b2e7f2d6 | 1359 | p1->Clear(); |
e308a636 | 1360 | fResults->DrawClone("nostack e1"); |
b2e7f2d6 | 1361 | |
fe0e93e6 | 1362 | // --- Draw other data ------------------------------------------- |
797161e8 | 1363 | if (fShowOthers != 0) { |
b2e7f2d6 | 1364 | TGraphAsymmErrors* o = 0; |
e308a636 | 1365 | TIter nextG(fOthers->GetListOfGraphs()); |
b2e7f2d6 | 1366 | while ((o = static_cast<TGraphAsymmErrors*>(nextG()))) |
1367 | o->DrawClone("same p"); | |
1368 | } | |
1369 | ||
fe0e93e6 | 1370 | // --- Make a legend in the result pad --------------------------- |
1371 | Double_t x1 = p1->GetLeftMargin()+.08; | |
1372 | Double_t x2 = 1-p1->GetRightMargin()-.08; | |
7095962e CHC |
1373 | Double_t y1 = p1->GetBottomMargin()+.01; |
1374 | Double_t y2 = .35; | |
fe0e93e6 | 1375 | Int_t nC = 2; |
7095962e | 1376 | if (HasCent()) { |
fe0e93e6 | 1377 | if (fCentAxis->GetNbins() <= 4) { |
1378 | x1 = p1->GetLeftMargin()+.15; | |
1379 | x2 = 1-p1->GetRightMargin()-.15; | |
1380 | y2 = .2; | |
7095962e CHC |
1381 | } |
1382 | else { | |
fe0e93e6 | 1383 | nC = 1; |
1384 | if (fOptions & kExtraWhite) { | |
1385 | x1 = .40; | |
1386 | x2 = .60; | |
1387 | y1 = .70; | |
1388 | y2 = 1-p1->GetTopMargin()-0.06; | |
1389 | } | |
1390 | else { | |
1391 | x1 = .70; | |
1392 | x2 = 1-p1->GetRightMargin()-.01; | |
1393 | y1 = .50; | |
1394 | y2 = 1-p1->GetTopMargin()-0.16; | |
1395 | } | |
7095962e | 1396 | } |
fe0e93e6 | 1397 | } |
1398 | BuildLegend(fResults, fOthers, x1, y1, x2, y2, nC); | |
b2e7f2d6 | 1399 | |
fe0e93e6 | 1400 | // --- Parameters for stuff on the right ------------------------- |
1401 | Double_t yTop = 1-p1->GetTopMargin()-.02; | |
1402 | Double_t xR = .95; | |
1403 | Double_t yR = yTop; | |
b2e7f2d6 | 1404 | |
fe0e93e6 | 1405 | // --- Put a nice label in the plot ------------------------------ |
b2e7f2d6 | 1406 | TString eS; |
bad9a3c1 | 1407 | UShort_t snn = fSNNString->GetUniqueID(); |
1408 | const char* sys = fSysString->GetTitle(); | |
5bb5d1f6 | 1409 | if (snn == 2750) snn = 2760; |
b2e7f2d6 | 1410 | if (snn > 1000) eS = Form("%4.2fTeV", float(snn)/1000); |
1411 | else eS = Form("%3dGeV", snn); | |
fe0e93e6 | 1412 | TLatex* tt = new TLatex(xR, yR, Form("%s #sqrt{s%s}=%s, %s", |
1413 | sys, | |
1414 | (HasCent() ? "_{NN}" : ""), | |
1415 | eS.Data(), | |
1416 | HasCent() ? "by centrality" : | |
1417 | fTrigString->GetTitle())); | |
d43c6cd0 | 1418 | tt->SetTextColor(kAliceBlue); |
b2e7f2d6 | 1419 | tt->SetNDC(); |
d43c6cd0 | 1420 | tt->SetTextFont(kFont); |
b2e7f2d6 | 1421 | tt->SetTextAlign(33); |
1422 | tt->Draw(); | |
fe0e93e6 | 1423 | yR -= tt->GetTextSize() + .01; |
bfab35d9 | 1424 | |
fe0e93e6 | 1425 | // --- Put number of accepted events on the plot ----------------- |
e308a636 | 1426 | Int_t nev = 0; |
ffca499d | 1427 | if (fTriggers) nev = fTriggers->GetBinContent(1); |
fe0e93e6 | 1428 | TLatex* et = new TLatex(xR, yR, Form("%d events", nev)); |
d43c6cd0 | 1429 | et->SetTextColor(kAliceBlue); |
b2e7f2d6 | 1430 | et->SetNDC(); |
d43c6cd0 | 1431 | et->SetTextFont(kFont); |
b2e7f2d6 | 1432 | et->SetTextAlign(33); |
1433 | et->Draw(); | |
fe0e93e6 | 1434 | yR -= et->GetTextSize() + .01; |
b2e7f2d6 | 1435 | |
fe0e93e6 | 1436 | // --- Put vertex axis on the plot ------------------------------- |
b2e7f2d6 | 1437 | if (fVtxAxis) { |
fe0e93e6 | 1438 | TLatex* vt = new TLatex(xR, yR, fVtxAxis->GetTitle()); |
b2e7f2d6 | 1439 | vt->SetNDC(); |
d43c6cd0 | 1440 | vt->SetTextFont(kFont); |
b2e7f2d6 | 1441 | vt->SetTextAlign(33); |
d43c6cd0 | 1442 | vt->SetTextColor(kAliceBlue); |
b2e7f2d6 | 1443 | vt->Draw(); |
fe0e93e6 | 1444 | yR -= vt->GetTextSize() + .01; |
b2e7f2d6 | 1445 | } |
1446 | // results->Draw("nostack e1 same"); | |
1447 | ||
fe0e93e6 | 1448 | // --- Put statement on corrections used on the plot ------------- |
b615e211 | 1449 | TString corrs; |
1450 | if (!fEmpirical.IsNull()) corrs.Append("Emperical"); | |
1451 | if (!fFinalMC.IsNull()) { | |
1452 | if (!corrs.IsNull()) corrs.Append("+"); | |
1453 | corrs.Append("Final MC"); | |
1454 | } | |
1455 | ||
1456 | if (!corrs.IsNull()) { | |
1457 | corrs.Append(" correction"); | |
1458 | if (corrs.Index("+") != kNPOS) corrs.Append("s"); | |
fe0e93e6 | 1459 | TLatex* em = new TLatex(xR, yR, corrs); |
b615e211 | 1460 | em->SetNDC(); |
d43c6cd0 | 1461 | em->SetTextFont(kFont); |
b615e211 | 1462 | em->SetTextAlign(33); |
d43c6cd0 | 1463 | em->SetTextColor(kAliceBlue); |
b615e211 | 1464 | em->Draw(); |
fe0e93e6 | 1465 | yR -= em->GetTextSize() + .01; |
b615e211 | 1466 | } |
fe0e93e6 | 1467 | |
1468 | // --- Put trigger efficiency on the plot (only pp and if != 1) -- | |
bfab35d9 | 1469 | if (fTriggerEff > 0 && fTriggerEff <= 1 && !HasCent()) { |
fe0e93e6 | 1470 | TLatex* ef = new TLatex(xR, yR, Form("#varepsilon_{%s} = %5.3f", |
bfab35d9 | 1471 | fTrigString->GetTitle(), |
1472 | fTriggerEff)); | |
1473 | ef->SetNDC(); | |
d43c6cd0 | 1474 | ef->SetTextFont(kFont); |
bfab35d9 | 1475 | ef->SetTextAlign(33); |
d43c6cd0 | 1476 | ef->SetTextColor(kAliceBlue); |
bfab35d9 | 1477 | ef->Draw(); |
fe0e93e6 | 1478 | yR -= ef->GetTextSize() + .01; |
1479 | } | |
1480 | ||
1481 | // --- Put logo on the plot -------------------------------------- | |
1482 | if (fOptions & kLogo) { | |
1483 | TString savPath(gROOT->GetMacroPath()); | |
1484 | TString fwd("$(ALICE_ROOT)/PWGLF/FORWARD/analysis2"); | |
1485 | gROOT->SetMacroPath(Form("%s:%s/scripts", | |
1486 | gROOT->GetMacroPath(), fwd.Data())); | |
1487 | // Always recompile | |
1488 | if (!gROOT->GetClass("AliceLogo")) | |
a19faec0 | 1489 | gROOT->LoadMacro("AliceLogo.C+"); |
fe0e93e6 | 1490 | gROOT->SetMacroPath(savPath); |
1491 | ||
1492 | if (gROOT->GetClass("AliceLogo")) { | |
1493 | yR -= .22; | |
1494 | p1->cd(); | |
1495 | p1->Range(0,0,1,1); | |
1496 | gROOT->ProcessLine("AliceLogo* al = new AliceLogo();"); | |
1497 | gROOT->ProcessLine(Form("al->Draw(0,.88,%f,.2, 0, 0);", yR)); | |
1498 | yR -= .01; | |
1499 | } | |
bfab35d9 | 1500 | } |
b2e7f2d6 | 1501 | |
fe0e93e6 | 1502 | |
1503 | // --- Parameter for stuff on the left --------------------------- | |
1504 | Double_t yL = yTop; | |
1505 | Double_t xL = .12; | |
b2e7f2d6 | 1506 | |
fe0e93e6 | 1507 | // --- Mark as work in progress ---------------------------------- |
1508 | TLatex* pt = new TLatex(xL, yL, "Work in progress"); | |
b2e7f2d6 | 1509 | pt->SetNDC(); |
d43c6cd0 | 1510 | pt->SetTextFont(62); |
f709afcf | 1511 | // pt->SetTextSize(); |
d43c6cd0 | 1512 | pt->SetTextColor(kAliceRed); |
b2e7f2d6 | 1513 | pt->SetTextAlign(13); |
1514 | pt->Draw(); | |
fe0e93e6 | 1515 | yL -= pt->GetTextSize()+.01; |
1516 | ||
1517 | TDatime now; | |
1518 | TLatex* dt = new TLatex(xL, yL, now.AsSQLString()); | |
d43c6cd0 | 1519 | dt->SetNDC(); |
1520 | dt->SetTextFont(42); | |
1521 | dt->SetTextSize(0.04); | |
1522 | dt->SetTextColor(kAliceBlue); // kAliceRed); | |
1523 | dt->SetTextAlign(13); | |
1524 | dt->Draw(); | |
fe0e93e6 | 1525 | yL -= dt->GetTextSize()+.01; |
d43c6cd0 | 1526 | |
fe0e93e6 | 1527 | // --- Possible centrality legend -------------------------------- |
1528 | BuildCentLegend(xL, yL-.4, xL+.23, yL); | |
1529 | ||
1530 | // --- Attach Zoom executor -------------------------------------- | |
1531 | AttachExec(p1, fResults->GetName(), 1, false); | |
1532 | ||
1533 | // --- Possibly add title ---------------------------------------- | |
1534 | if (yd < 0.0001) PlotTitle(p1, yyd, true); | |
b2e7f2d6 | 1535 | |
b2e7f2d6 | 1536 | p1->cd(); |
1537 | } | |
1538 | //__________________________________________________________________ | |
2f92ff0e | 1539 | /** |
1540 | * Plot the ratios | |
1541 | * | |
2f92ff0e | 1542 | * @param max Maximum diviation from 1 |
1543 | * @param y1 Lower y coordinate of pad | |
1544 | * @param y2 Upper y coordinate of pad | |
1545 | */ | |
e308a636 | 1546 | void PlotRatios(Double_t max, Double_t y1, Double_t y2) |
b2e7f2d6 | 1547 | { |
7095962e | 1548 | if ((fOptions & kShowRatios) == 0) return; |
e308a636 | 1549 | |
b2e7f2d6 | 1550 | bool isBottom = (y1 < 0.0001); |
1551 | Double_t yd = y2 - y1; | |
fe0e93e6 | 1552 | |
1553 | // --- Make a sub-pad for the ratios ----------------------------- | |
b2e7f2d6 | 1554 | TPad* p2 = new TPad("p2", "p2", 0, y1, 1.0, y2, 0, 0, 0); |
1555 | p2->SetTopMargin(0.001); | |
d43c6cd0 | 1556 | p2->SetRightMargin(kRightMargin); |
b2e7f2d6 | 1557 | p2->SetBottomMargin(isBottom ? 1/yd * 0.07 : 0.0001); |
1558 | p2->SetGridx(); | |
1559 | p2->SetTicks(1,1); | |
1560 | p2->SetNumber(2); | |
1561 | p2->Draw(); | |
1562 | p2->cd(); | |
1563 | ||
fe0e93e6 | 1564 | // --- Fix up axis ----------------------------------------------- |
797161e8 | 1565 | FixAxis(fRatios, yd, "Ratios", 7); |
b2e7f2d6 | 1566 | |
fe0e93e6 | 1567 | // --- Set up min/max -------------------------------------------- |
e308a636 | 1568 | fRatios->SetMaximum(1+TMath::Max(.22,1.05*max)); |
1569 | fRatios->SetMinimum(1-TMath::Max(.32,1.05*max)); | |
fe0e93e6 | 1570 | |
1571 | // --- Clear pad and draw ---------------------------------------- | |
b2e7f2d6 | 1572 | p2->Clear(); |
e308a636 | 1573 | fRatios->DrawClone("nostack e1"); |
b2e7f2d6 | 1574 | |
1575 | ||
fe0e93e6 | 1576 | // --- Make a legend --------------------------------------------- |
ffca499d | 1577 | BuildLegend(fRatios, 0, .15,p2->GetBottomMargin()+.01,.9, |
cd3b3781 | 1578 | isBottom ? .6 : .4, 2); |
fe0e93e6 | 1579 | |
1580 | // --- Make a nice band from 0.9 to 1.1 -------------------------- | |
b2e7f2d6 | 1581 | TGraphErrors* band = new TGraphErrors(2); |
e1f47419 | 1582 | band->SetPoint(0, fResults->GetXaxis()->GetXmin(), 1); |
1583 | band->SetPoint(1, fResults->GetXaxis()->GetXmax(), 1); | |
b2e7f2d6 | 1584 | band->SetPointError(0, 0, .1); |
1585 | band->SetPointError(1, 0, .1); | |
1586 | band->SetFillColor(kYellow+2); | |
1587 | band->SetFillStyle(3002); | |
1588 | band->SetLineStyle(2); | |
1589 | band->SetLineWidth(1); | |
1590 | band->Draw("3 same"); | |
1591 | band->DrawClone("X L same"); | |
1592 | ||
fe0e93e6 | 1593 | // --- Replot the ratios on top ---------------------------------- |
e308a636 | 1594 | fRatios->DrawClone("nostack e1 same"); |
b2e7f2d6 | 1595 | |
fe0e93e6 | 1596 | // --- Some more stuff ------------------------------------------- |
1597 | AttachExec(p2, fRatios->GetName(), 2, isBottom); | |
1598 | if (isBottom) PlotTitle(p2, yd, true); | |
b2e7f2d6 | 1599 | } |
1600 | //__________________________________________________________________ | |
2f92ff0e | 1601 | /** |
1602 | * Plot the asymmetries | |
1603 | * | |
2f92ff0e | 1604 | * @param max Maximum diviation from 1 |
1605 | * @param y1 Lower y coordinate of pad | |
1606 | * @param y2 Upper y coordinate of pad | |
1607 | */ | |
e308a636 | 1608 | void PlotLeftRight(Double_t max, Double_t y1, Double_t y2) |
b2e7f2d6 | 1609 | { |
7095962e | 1610 | if (!(fOptions & kShowLeftRight)) return; |
e308a636 | 1611 | |
b2e7f2d6 | 1612 | bool isBottom = (y1 < 0.0001); |
1613 | Double_t yd = y2 - y1; | |
fe0e93e6 | 1614 | // --- Make a sub-pad for the asymmetry -------------------------- |
b2e7f2d6 | 1615 | TPad* p3 = new TPad("p3", "p3", 0, y1, 1.0, y2, 0, 0, 0); |
1616 | p3->SetTopMargin(0.001); | |
d43c6cd0 | 1617 | p3->SetRightMargin(kRightMargin); |
b2e7f2d6 | 1618 | p3->SetBottomMargin(isBottom ? 1/yd * 0.07 : 0.0001); |
1619 | p3->SetGridx(); | |
1620 | p3->SetTicks(1,1); | |
1621 | p3->SetNumber(2); | |
1622 | p3->Draw(); | |
1623 | p3->cd(); | |
1624 | ||
fe0e93e6 | 1625 | // --- Fix up axis ----------------------------------------------- |
797161e8 | 1626 | FixAxis(fLeftRight, yd, "Right/Left", 4); |
b2e7f2d6 | 1627 | |
fe0e93e6 | 1628 | // ---- Setup min/max -------------------------------------------- |
e308a636 | 1629 | fLeftRight->SetMaximum(1+TMath::Max(.12,1.05*max)); |
1630 | fLeftRight->SetMinimum(1-TMath::Max(.15,1.05*max)); | |
fe0e93e6 | 1631 | |
1632 | // --- Clear pad and redraw -------------------------------------- | |
b2e7f2d6 | 1633 | p3->Clear(); |
e308a636 | 1634 | fLeftRight->DrawClone("nostack e1"); |
b2e7f2d6 | 1635 | |
1636 | ||
fe0e93e6 | 1637 | // --- Make a legend --------------------------------------------- |
298a7264 | 1638 | Double_t xx1 = (HasCent() ? .7 : .15); |
1639 | Double_t xx2 = (HasCent() ? 1-p3->GetRightMargin()-.01 : .90); | |
5bb5d1f6 | 1640 | Double_t yy1 = p3->GetBottomMargin()+.01; |
298a7264 | 1641 | Double_t yy2 = (HasCent() ? 1-p3->GetTopMargin()-.01-.15 : .5); |
5bb5d1f6 | 1642 | BuildLegend(fLeftRight, 0, xx1, yy1, xx2, yy2); |
fe0e93e6 | 1643 | |
1644 | // --- Make a nice band from 0.9 to 1.1 -------------------------- | |
b2e7f2d6 | 1645 | TGraphErrors* band = new TGraphErrors(2); |
e1f47419 | 1646 | band->SetPoint(0, fResults->GetXaxis()->GetXmin(), 1); |
1647 | band->SetPoint(1, fResults->GetXaxis()->GetXmax(), 1); | |
b2e7f2d6 | 1648 | band->SetPointError(0, 0, .05); |
1649 | band->SetPointError(1, 0, .05); | |
1650 | band->SetFillColor(kYellow+2); | |
1651 | band->SetFillStyle(3002); | |
1652 | band->SetLineStyle(2); | |
1653 | band->SetLineWidth(1); | |
1654 | band->Draw("3 same"); | |
1655 | band->DrawClone("X L same"); | |
1656 | ||
fe0e93e6 | 1657 | // --- Re-draw over ---------------------------------------------- |
e308a636 | 1658 | fLeftRight->DrawClone("nostack e1 same"); |
fe0e93e6 | 1659 | |
1660 | // --- Misc stuff ------------------------------------------------ | |
1661 | AttachExec(p3, fLeftRight->GetName(), 2, isBottom); | |
1662 | if (isBottom) PlotTitle(p3, yd, true); | |
b2e7f2d6 | 1663 | } |
2f92ff0e | 1664 | /** @} */ |
1665 | //================================================================== | |
1666 | /** | |
1667 | * @{ | |
1668 | * @name Data utility functions | |
1669 | */ | |
7095962e CHC |
1670 | //__________________________________________________________________ |
1671 | /** | |
1672 | * Get the color for a centrality bin | |
1673 | * | |
1674 | * @param bin Centrality bin | |
1675 | * | |
1676 | * @return Color | |
1677 | */ | |
1678 | Int_t GetCentralityColor(Int_t bin) const | |
1679 | { | |
1680 | if (fCentAxis->GetNbins() < 6) { | |
1681 | switch (bin) { | |
1682 | case 1: return kRed+2; | |
1683 | case 2: return kGreen+2; | |
1684 | case 3: return kBlue+1; | |
1685 | case 4: return kCyan+1; | |
1686 | case 5: return kMagenta+1; | |
1687 | case 6: return kYellow+2; | |
1688 | } | |
1689 | } | |
1690 | UShort_t centLow = fCentAxis->GetBinLowEdge(bin); | |
1691 | UShort_t centHigh = fCentAxis->GetBinUpEdge(bin); | |
1692 | Float_t fc = (centLow+double(centHigh-centLow)/2) / 100; | |
1693 | Int_t nCol = gStyle->GetNumberOfColors(); | |
1694 | Int_t icol = TMath::Min(nCol-1,int(fc * nCol + .5)); | |
1695 | Int_t col = gStyle->GetColorPalette(icol); | |
1696 | //Info("GetCentralityColor","%3d: %3d-%3d -> %3d",bin,centLow,centHigh,col); | |
1697 | return col; | |
1698 | } | |
1699 | //__________________________________________________________________ | |
1700 | /** | |
1701 | * Set attributed on a histogram | |
1702 | * | |
1703 | * @param h Histogram | |
1704 | * @param color Color | |
1705 | */ | |
1706 | void SetAttributes(TH1* h, Int_t color) | |
1707 | { | |
1708 | if (!h) return; | |
1709 | if (color < 0) return; | |
1710 | // h->SetLineColor(color); | |
1711 | h->SetMarkerColor(color); | |
1712 | // h->SetFillColor(color); | |
1713 | } | |
1714 | //__________________________________________________________________ | |
1715 | /** | |
1716 | * Set attributed on a graph | |
1717 | * | |
1718 | * @param g Graph | |
1719 | * @param color Color | |
1720 | */ | |
1721 | void SetAttributes(TGraph* g, Int_t color) | |
1722 | { | |
1723 | if (!g) return; | |
1724 | if (color < 0) return; | |
1725 | // g->SetLineColor(color); | |
1726 | g->SetMarkerColor(color); | |
1727 | // g->SetFillColor(color); | |
1728 | } | |
1729 | //__________________________________________________________________ | |
1730 | /** | |
1731 | * Modify the title | |
1732 | * | |
1733 | */ | |
1734 | void ModifyTitle(TNamed* h, const char* centTxt) | |
1735 | { | |
1736 | if (!h) return; | |
1737 | ||
1738 | TString title(h->GetTitle()); | |
1739 | title.ReplaceAll("ALICE ",""); | |
1740 | if (title.Contains("Central")) | |
1741 | title.ReplaceAll("CentraldNdeta", "SPD clusters"); | |
1742 | if (title.Contains("Forward")) | |
1743 | title.ReplaceAll("ForwarddNdeta", "FMD"); | |
1744 | h->SetTitle(title); | |
1745 | ||
1746 | if (centTxt && centTxt[0] != '\0') { | |
1747 | TString name(h->GetName()); | |
1748 | name.Append(Form("_%s", centTxt)); | |
1749 | h->SetName(name); | |
1750 | } | |
1751 | ||
1752 | return; | |
1753 | // if (!centTxt || !h) return; | |
1754 | // h->SetTitle(Form("%s, %s", h->GetTitle(), centTxt)); | |
1755 | } | |
2f92ff0e | 1756 | /** |
1757 | * Get a result from the passed list | |
1758 | * | |
1759 | * @param list List to search | |
1760 | * @param name Object name to search for | |
1761 | * | |
290052e7 | 1762 | * @return Histogram |
2f92ff0e | 1763 | */ |
7095962e | 1764 | TH1* FetchHistogram(const TList* list, const char* name) const |
2f92ff0e | 1765 | { |
1766 | if (!list) return 0; | |
1767 | ||
e308a636 | 1768 | TH1* ret = static_cast<TH1*>(list->FindObject(name)); |
5bb5d1f6 | 1769 | #if 0 |
e1f47419 | 1770 | if (!ret) { |
1771 | // all->ls(); | |
1772 | Warning("GetResult", "Histogram %s not found", name); | |
2f92ff0e | 1773 | } |
5bb5d1f6 | 1774 | #endif |
b2e7f2d6 | 1775 | |
e1f47419 | 1776 | return ret; |
2f92ff0e | 1777 | } |
b2e7f2d6 | 1778 | //__________________________________________________________________ |
2f92ff0e | 1779 | /** |
1780 | * Add a histogram to the stack after possibly rebinning it | |
1781 | * | |
1782 | * @param stack Stack to add to | |
1783 | * @param hist histogram | |
1784 | * @param option Draw options | |
1785 | * | |
1786 | * @return Maximum of histogram | |
1787 | */ | |
e308a636 | 1788 | Double_t AddHistogram(THStack* stack, TH1* hist, Option_t* option="") const |
b2e7f2d6 | 1789 | { |
1790 | // Check if we have input | |
1791 | if (!hist) return 0; | |
1792 | ||
1793 | // Rebin if needed | |
1794 | Rebin(hist); | |
1795 | ||
1796 | stack->Add(hist, option); | |
1797 | return hist->GetMaximum(); | |
1798 | } | |
1799 | //__________________________________________________________________ | |
2f92ff0e | 1800 | /** |
1801 | * Add a histogram to the stack after possibly rebinning it | |
1802 | * | |
1803 | * @param stack Stack to add to | |
1804 | * @param hist histogram | |
1805 | * @param option Draw options | |
1806 | * @param sym On return, the data symmetriced (added to stack) | |
1807 | * | |
1808 | * @return Maximum of histogram | |
1809 | */ | |
e308a636 | 1810 | Double_t AddHistogram(THStack* stack, TH1* hist, TH1*& sym, |
1811 | Option_t* option="") const | |
b2e7f2d6 | 1812 | { |
1813 | // Check if we have input | |
1814 | if (!hist) return 0; | |
1815 | ||
1816 | // Rebin if needed | |
1817 | Rebin(hist); | |
1818 | stack->Add(hist, option); | |
1819 | ||
1820 | // Now symmetrice the histogram | |
7095962e | 1821 | if ((fOptions & kMirror)) { |
298a7264 | 1822 | sym = Symmetrice(hist); |
1823 | stack->Add(sym, option); | |
1824 | } | |
b2e7f2d6 | 1825 | |
1826 | return hist->GetMaximum(); | |
1827 | } | |
1828 | ||
1829 | //__________________________________________________________________ | |
1830 | /** | |
1831 | * Rebin a histogram | |
1832 | * | |
1833 | * @param h Histogram to rebin | |
b2e7f2d6 | 1834 | */ |
1835 | virtual void Rebin(TH1* h) const | |
1836 | { | |
1837 | if (fRebin <= 1) return; | |
1838 | ||
1839 | Int_t nBins = h->GetNbinsX(); | |
1840 | if(nBins % fRebin != 0) { | |
1841 | Warning("Rebin", "Rebin factor %d is not a devisor of current number " | |
1842 | "of bins %d in the histogram %s", fRebin, nBins, h->GetName()); | |
1843 | return; | |
1844 | } | |
1845 | ||
1846 | // Make a copy | |
1847 | TH1* tmp = static_cast<TH1*>(h->Clone("tmp")); | |
1848 | tmp->Rebin(fRebin); | |
1849 | tmp->SetDirectory(0); | |
e28f5fc5 | 1850 | tmp->Reset(); |
b2e7f2d6 | 1851 | // The new number of bins |
1852 | Int_t nBinsNew = nBins / fRebin; | |
1853 | for(Int_t i = 1;i<= nBinsNew; i++) { | |
1854 | Double_t content = 0; | |
1855 | Double_t sumw = 0; | |
1856 | Double_t wsum = 0; | |
1857 | Int_t nbins = 0; | |
1858 | for(Int_t j = 1; j<=fRebin;j++) { | |
1859 | Int_t bin = (i-1)*fRebin + j; | |
1860 | Double_t c = h->GetBinContent(bin); | |
1861 | ||
1862 | if (c <= 0) continue; | |
1863 | ||
7095962e | 1864 | if ((fOptions & kCutEdges)) { |
b2e7f2d6 | 1865 | if (h->GetBinContent(bin+1)<=0 || |
73b32206 | 1866 | h->GetBinContent(bin-1)<=0) { |
b2e7f2d6 | 1867 | Warning("Rebin", "removing bin %d=%f of %s (%d=%f,%d=%f)", |
1868 | bin, c, h->GetName(), | |
1869 | bin+1, h->GetBinContent(bin+1), | |
1870 | bin-1, h->GetBinContent(bin-1)); | |
1871 | continue; | |
1872 | } | |
1873 | } | |
1874 | Double_t e = h->GetBinError(bin); | |
1875 | Double_t w = 1 / (e*e); // 1/c/c | |
1876 | content += c; | |
1877 | sumw += w; | |
1878 | wsum += w * c; | |
1879 | nbins++; | |
1880 | } | |
1881 | ||
e28f5fc5 | 1882 | if(content > 0 && nbins > 1 ) { |
b2e7f2d6 | 1883 | tmp->SetBinContent(i, wsum / sumw); |
1884 | tmp->SetBinError(i,1./TMath::Sqrt(sumw)); | |
1885 | } | |
1886 | } | |
1887 | ||
1888 | // Finally, rebin the histogram, and set new content | |
1889 | h->Rebin(fRebin); | |
1890 | h->Reset(); | |
1891 | for(Int_t i = 1; i<= nBinsNew; i++) { | |
1892 | h->SetBinContent(i,tmp->GetBinContent(i)); | |
1893 | h->SetBinError(i, tmp->GetBinError(i)); | |
1894 | } | |
1895 | ||
1896 | delete tmp; | |
1897 | } | |
1898 | //__________________________________________________________________ | |
1899 | /** | |
1900 | * Make an extension of @a h to make it symmetric about 0 | |
1901 | * | |
1902 | * @param h Histogram to symmertrice | |
1903 | * | |
1904 | * @return Symmetric extension of @a h | |
1905 | */ | |
1906 | TH1* Symmetrice(const TH1* h) const | |
1907 | { | |
1908 | Int_t nBins = h->GetNbinsX(); | |
1909 | TH1* s = static_cast<TH1*>(h->Clone(Form("%s_mirror", h->GetName()))); | |
1910 | s->SetTitle(Form("%s (mirrored)", h->GetTitle())); | |
1911 | s->Reset(); | |
1912 | s->SetBins(nBins, -h->GetXaxis()->GetXmax(), -h->GetXaxis()->GetXmin()); | |
1913 | s->SetMarkerColor(h->GetMarkerColor()); | |
1914 | s->SetMarkerSize(h->GetMarkerSize()); | |
1915 | s->SetMarkerStyle(h->GetMarkerStyle()+4); | |
1916 | s->SetFillColor(h->GetFillColor()); | |
1917 | s->SetFillStyle(h->GetFillStyle()); | |
1918 | s->SetDirectory(0); | |
1919 | ||
1920 | // Find the first and last bin with data | |
1921 | Int_t first = nBins+1; | |
1922 | Int_t last = 0; | |
1923 | for (Int_t i = 1; i <= nBins; i++) { | |
1924 | if (h->GetBinContent(i) <= 0) continue; | |
1925 | first = TMath::Min(first, i); | |
1926 | last = TMath::Max(last, i); | |
1927 | } | |
1928 | ||
1929 | Double_t xfirst = h->GetBinCenter(first-1); | |
1930 | Int_t f1 = h->GetXaxis()->FindBin(-xfirst); | |
1931 | Int_t l2 = s->GetXaxis()->FindBin(xfirst); | |
1932 | for (Int_t i = f1, j=l2; i <= last; i++,j--) { | |
1933 | s->SetBinContent(j, h->GetBinContent(i)); | |
1934 | s->SetBinError(j, h->GetBinError(i)); | |
1935 | } | |
1936 | // Fill in overlap bin | |
1937 | s->SetBinContent(l2+1, h->GetBinContent(first)); | |
1938 | s->SetBinError(l2+1, h->GetBinError(first)); | |
1939 | return s; | |
1940 | } | |
1941 | //__________________________________________________________________ | |
2f92ff0e | 1942 | /** |
1943 | * Calculate the left-right asymmetry of input histogram | |
1944 | * | |
1945 | * @param h Input histogram | |
1946 | * @param max On return, the maximum distance from 1 of the histogram | |
1947 | * | |
1948 | * @return Asymmetry | |
1949 | */ | |
1950 | TH1* Asymmetry(TH1* h, Double_t& max) | |
1951 | { | |
1952 | if (!h) return 0; | |
1953 | ||
1954 | TH1* ret = static_cast<TH1*>(h->Clone(Form("%s_leftright", h->GetName()))); | |
1955 | // Int_t oBins = h->GetNbinsX(); | |
1956 | // Double_t high = h->GetXaxis()->GetXmax(); | |
1957 | // Double_t low = h->GetXaxis()->GetXmin(); | |
1958 | // Double_t dBin = (high - low) / oBins; | |
1959 | // Int_t tBins = Int_t(2*high/dBin+.5); | |
1960 | // ret->SetBins(tBins, -high, high); | |
e308a636 | 1961 | ret->SetDirectory(0); |
2f92ff0e | 1962 | ret->Reset(); |
1963 | ret->SetTitle(Form("%s (+/-)", h->GetTitle())); | |
1964 | ret->SetYTitle("Right/Left"); | |
1965 | Int_t nBins = h->GetNbinsX(); | |
1966 | for (Int_t i = 1; i <= nBins; i++) { | |
1967 | Double_t x = h->GetBinCenter(i); | |
1968 | if (x > 0) break; | |
1969 | ||
1970 | Double_t c1 = h->GetBinContent(i); | |
1971 | Double_t e1 = h->GetBinError(i); | |
1972 | if (c1 <= 0) continue; | |
1973 | ||
1974 | Int_t j = h->FindBin(-x); | |
1975 | if (j <= 0 || j > nBins) continue; | |
1976 | ||
1977 | Double_t c2 = h->GetBinContent(j); | |
1978 | Double_t e2 = h->GetBinError(j); | |
1979 | ||
1980 | Double_t c12 = c1*c1; | |
1981 | Double_t e = TMath::Sqrt((e2*e2*c1*c1+e1*e1*c2*c2)/(c12*c12)); | |
1982 | ||
1983 | Int_t k = ret->FindBin(x); | |
1984 | ret->SetBinContent(k, c2/c1); | |
1985 | ret->SetBinError(k, e); | |
1986 | } | |
1987 | max = TMath::Max(max, RatioMax(ret)); | |
1988 | ||
1989 | return ret; | |
1990 | } | |
1991 | //__________________________________________________________________ | |
1992 | /** | |
1993 | * Transform a graph into a histogram | |
1994 | * | |
290052e7 | 1995 | * @param g Graph |
2f92ff0e | 1996 | * |
290052e7 | 1997 | * @return Histogram |
2f92ff0e | 1998 | */ |
1999 | TH1* Graph2Hist(const TGraphAsymmErrors* g) const | |
2000 | { | |
2001 | Int_t nBins = g->GetN(); | |
2002 | TArrayF bins(nBins+1); | |
2003 | Double_t dx = 0; | |
2004 | for (Int_t i = 0; i < nBins; i++) { | |
2005 | Double_t x = g->GetX()[i]; | |
2006 | Double_t exl = g->GetEXlow()[i]; | |
2007 | Double_t exh = g->GetEXhigh()[i]; | |
2008 | bins.fArray[i] = x-exl; | |
2009 | bins.fArray[i+1] = x+exh; | |
2010 | Double_t dxi = exh+exl; | |
2011 | if (i == 0) dx = dxi; | |
2012 | else if (dxi != dx) dx = 0; | |
2013 | } | |
2014 | TString name(g->GetName()); | |
2015 | TString title(g->GetTitle()); | |
2016 | TH1D* h = 0; | |
2017 | if (dx != 0) { | |
2018 | h = new TH1D(name.Data(), title.Data(), nBins, bins[0], bins[nBins]); | |
2019 | } | |
2020 | else { | |
2021 | h = new TH1D(name.Data(), title.Data(), nBins, bins.fArray); | |
2022 | } | |
2023 | h->SetMarkerStyle(g->GetMarkerStyle()); | |
2024 | h->SetMarkerColor(g->GetMarkerColor()); | |
2025 | h->SetMarkerSize(g->GetMarkerSize()); | |
2026 | ||
2027 | return h; | |
2028 | } | |
2029 | /* @} */ | |
2030 | //================================================================== | |
2031 | /** | |
2032 | * @{ | |
2033 | * @name Ratio utility functions | |
2034 | */ | |
2035 | /** | |
2036 | * Get the maximum diviation from 1 in the passed ratio | |
2037 | * | |
2038 | * @param h Ratio histogram | |
2039 | * | |
2040 | * @return Max diviation from 1 | |
2041 | */ | |
b2e7f2d6 | 2042 | Double_t RatioMax(TH1* h) const |
2043 | { | |
2044 | Int_t nBins = h->GetNbinsX(); | |
2045 | Double_t ret = 0; | |
2046 | for (Int_t i = 1; i <= nBins; i++) { | |
2047 | Double_t c = h->GetBinContent(i); | |
2048 | if (c == 0) continue; | |
2049 | Double_t e = h->GetBinError(i); | |
2050 | Double_t d = TMath::Abs(1-c-e); | |
2051 | ret = TMath::Max(d, ret); | |
2052 | } | |
2053 | return ret; | |
2054 | } | |
2055 | //__________________________________________________________________ | |
e1f47419 | 2056 | /** |
2057 | * Make the ratio of h1 to h2 | |
2058 | * | |
ffca499d | 2059 | * @param o1 First object (numerator) |
2060 | * @param o2 Second object (denominator) | |
2061 | * @param max Maximum diviation from 1 | |
e1f47419 | 2062 | * |
ffca499d | 2063 | * @return o1 / o2 |
e1f47419 | 2064 | */ |
2065 | TH1* Ratio(const TObject* o1, const TObject* o2, Double_t& max) const | |
2066 | { | |
b30dee70 | 2067 | if (!o1 || !o2) return 0; |
2068 | ||
fe0e93e6 | 2069 | bool mirror = false; |
2070 | TString n1(o1->GetName()); | |
2071 | if (n1.Contains("mirror")) mirror = true; | |
2072 | ||
5bb5d1f6 | 2073 | TH1* r = 0; |
2074 | const TAttMarker* m1 = 0; | |
2075 | const TAttMarker* m2 = 0; | |
e1f47419 | 2076 | const TH1* h1 = dynamic_cast<const TH1*>(o1); |
2077 | if (h1) { | |
5bb5d1f6 | 2078 | m1 = h1; |
e1f47419 | 2079 | const TH1* h2 = dynamic_cast<const TH1*>(o2); |
ffca499d | 2080 | if (h2) { |
5bb5d1f6 | 2081 | m2 = h2; |
2082 | r = RatioHH(h1,h2); | |
ffca499d | 2083 | } |
2084 | else { | |
2085 | const TGraph* g2 = dynamic_cast<const TGraph*>(o2); | |
2086 | if (g2) { | |
5bb5d1f6 | 2087 | m2 = g2; |
2088 | r = RatioHG(h1,g2); | |
ffca499d | 2089 | } |
2090 | } | |
e1f47419 | 2091 | } |
ffca499d | 2092 | else { |
2093 | const TGraphAsymmErrors* g1 = dynamic_cast<const TGraphAsymmErrors*>(o1); | |
2094 | if (g1) { | |
5bb5d1f6 | 2095 | m1 = g1; |
ffca499d | 2096 | const TGraphAsymmErrors* g2 = |
2097 | dynamic_cast<const TGraphAsymmErrors*>(o2); | |
2098 | if (g2) { | |
5bb5d1f6 | 2099 | m2 = g2; |
2100 | r = RatioGG(g1, g2); | |
ffca499d | 2101 | } |
2102 | } | |
e1f47419 | 2103 | } |
ffca499d | 2104 | if (!r) { |
c8b1a7db | 2105 | // Warning("Ratio", "Don't know how to divide a %s (%s) with a %s (%s)", |
2106 | // o1->ClassName(),o1->GetName(),o2->ClassName(),o2->GetName()); | |
ffca499d | 2107 | return 0; |
2108 | } | |
2109 | // Check that the histogram isn't empty | |
2110 | if (r->GetEntries() <= 0) { | |
2111 | delete r; | |
2112 | r = 0; | |
2113 | } | |
5bb5d1f6 | 2114 | if (r) { |
fe0e93e6 | 2115 | UShort_t m2bits = MarkerUtil::GetMarkerBits(m2->GetMarkerStyle()); |
2116 | if (mirror) m2bits |= MarkerUtil::kHollow; | |
2117 | else m2bits &= ~MarkerUtil::kHollow; | |
2118 | r->SetMarkerStyle(MarkerUtil::GetMarkerStyle(m2bits)); | |
9ecab72f | 2119 | r->SetMarkerColor(m1->GetMarkerColor()); |
cd3b3781 | 2120 | if (TString(o2->GetName()).Contains("truth", TString::kIgnoreCase)) |
2121 | r->SetMarkerStyle(m1->GetMarkerStyle()); | |
5bb5d1f6 | 2122 | r->SetMarkerSize(0.9*m1->GetMarkerSize()); |
2123 | r->SetName(Form("%s_over_%s", o1->GetName(), o2->GetName())); | |
2124 | r->SetTitle(Form("%s / %s", o1->GetTitle(), o2->GetTitle())); | |
2125 | r->SetDirectory(0); | |
2126 | max = TMath::Max(RatioMax(r), max); | |
2127 | } | |
2128 | ||
ffca499d | 2129 | return r; |
e1f47419 | 2130 | } |
2131 | //__________________________________________________________________ | |
b2e7f2d6 | 2132 | /** |
2133 | * Compute the ratio of @a h to @a g. @a g is evaluated at the bin | |
2134 | * centers of @a h | |
2135 | * | |
2136 | * @param h Numerator | |
2137 | * @param g Divisor | |
2138 | * | |
2139 | * @return h/g | |
2140 | */ | |
5bb5d1f6 | 2141 | TH1* RatioHG(const TH1* h, const TGraph* g) const |
b2e7f2d6 | 2142 | { |
2143 | if (!h || !g) return 0; | |
2144 | ||
2145 | TH1* ret = static_cast<TH1*>(h->Clone("tmp")); | |
b2e7f2d6 | 2146 | ret->Reset(); |
b2e7f2d6 | 2147 | Double_t xlow = g->GetX()[0]; |
2148 | Double_t xhigh = g->GetX()[g->GetN()-1]; | |
2149 | if (xlow > xhigh) { Double_t t = xhigh; xhigh = xlow; xlow = t; } | |
2150 | ||
2151 | for (Int_t i = 1; i <= h->GetNbinsX(); i++) { | |
2152 | Double_t c = h->GetBinContent(i); | |
2153 | if (c <= 0) continue; | |
2154 | ||
2155 | Double_t x = h->GetBinCenter(i); | |
2156 | if (x < xlow || x > xhigh) continue; | |
2157 | ||
2158 | Double_t f = g->Eval(x); | |
2159 | if (f <= 0) continue; | |
2160 | ||
2161 | ret->SetBinContent(i, c / f); | |
2162 | ret->SetBinError(i, h->GetBinError(i) / f); | |
2163 | } | |
b2e7f2d6 | 2164 | return ret; |
2165 | } | |
2166 | //__________________________________________________________________ | |
2167 | /** | |
2168 | * Make the ratio of h1 to h2 | |
2169 | * | |
2170 | * @param h1 First histogram (numerator) | |
2171 | * @param h2 Second histogram (denominator) | |
2172 | * | |
2173 | * @return h1 / h2 | |
2174 | */ | |
5bb5d1f6 | 2175 | TH1* RatioHH(const TH1* h1, const TH1* h2) const |
b2e7f2d6 | 2176 | { |
2177 | if (!h1 || !h2) return 0; | |
c8b1a7db | 2178 | Bool_t bad = false; |
2179 | if (h1->GetNbinsX() != h2->GetNbinsX()) { | |
2180 | Error("RatioHH", "They have differnet number of bins"); | |
2181 | bad = true; | |
2182 | } | |
2183 | for (Int_t i = 1; i <= h1->GetNbinsX(); i++) { | |
2184 | if (h1->GetXaxis()->GetBinLowEdge(i) != | |
2185 | h2->GetXaxis()->GetBinLowEdge(i)) { | |
2186 | // Error("RatioHH", "They have incompatible variable bins"); | |
2187 | bad = true; | |
2188 | break; | |
2189 | } | |
2190 | } | |
2191 | if (bad) return 0; | |
2192 | ||
5bb5d1f6 | 2193 | TH1* t1 = static_cast<TH1*>(h1->Clone("tmp")); |
c8b1a7db | 2194 | // Printf("Dividing %s with %s", h1->GetName(), h2->GetName()); |
b2e7f2d6 | 2195 | t1->Divide(h2); |
b2e7f2d6 | 2196 | return t1; |
2197 | } | |
2198 | //__________________________________________________________________ | |
2199 | /** | |
2200 | * Calculate the ratio of two graphs - g1 / g2 | |
2201 | * | |
2202 | * @param g1 Numerator | |
2203 | * @param g2 Denominator | |
2204 | * | |
2205 | * @return g1 / g2 in a histogram | |
2206 | */ | |
ffca499d | 2207 | TH1* RatioGG(const TGraphAsymmErrors* g1, |
5bb5d1f6 | 2208 | const TGraphAsymmErrors* g2) const |
b2e7f2d6 | 2209 | { |
2210 | Int_t nBins = g1->GetN(); | |
2211 | TArrayF bins(nBins+1); | |
2212 | Double_t dx = 0; | |
2213 | for (Int_t i = 0; i < nBins; i++) { | |
2214 | Double_t x = g1->GetX()[i]; | |
2215 | Double_t exl = g1->GetEXlow()[i]; | |
2216 | Double_t exh = g1->GetEXhigh()[i]; | |
2217 | bins.fArray[i] = x-exl; | |
2218 | bins.fArray[i+1] = x+exh; | |
2219 | Double_t dxi = exh+exl; | |
2220 | if (i == 0) dx = dxi; | |
2221 | else if (dxi != dx) dx = 0; | |
2222 | } | |
b2e7f2d6 | 2223 | TH1* h = 0; |
2224 | if (dx != 0) { | |
5bb5d1f6 | 2225 | h = new TH1F("tmp", "tmp", nBins, bins[0], bins[nBins]); |
b2e7f2d6 | 2226 | } |
2227 | else { | |
5bb5d1f6 | 2228 | h = new TH1F("tmp", "tmp", nBins, bins.fArray); |
b2e7f2d6 | 2229 | } |
b2e7f2d6 | 2230 | |
2231 | Double_t low = g2->GetX()[0]; | |
2232 | Double_t high = g2->GetX()[g2->GetN()-1]; | |
2233 | if (low > high) { Double_t t = low; low = high; high = t; } | |
2234 | for (Int_t i = 0; i < nBins; i++) { | |
2235 | Double_t x = g1->GetX()[i]; | |
2236 | if (x < low-dx || x > high+dx) continue; | |
2237 | Double_t c1 = g1->GetY()[i]; | |
2238 | Double_t e1 = g1->GetErrorY(i); | |
2239 | Double_t c2 = g2->Eval(x); | |
2240 | ||
2241 | h->SetBinContent(i+1, c1 / c2); | |
2242 | h->SetBinError(i+1, e1 / c2); | |
2243 | } | |
b2e7f2d6 | 2244 | return h; |
2f92ff0e | 2245 | } |
2246 | /* @} */ | |
2247 | //================================================================== | |
2248 | /** | |
2249 | * @{ | |
2250 | * @name Graphics utility functions | |
2251 | */ | |
2252 | /** | |
2253 | * Find an X axis in a pad | |
2254 | * | |
2255 | * @param p Pad | |
2256 | * @param name Histogram to find axis for | |
2257 | * | |
2258 | * @return Found axis or null | |
2259 | */ | |
2260 | TAxis* FindXAxis(TVirtualPad* p, const char* name) | |
b2e7f2d6 | 2261 | { |
2f92ff0e | 2262 | TObject* o = p->GetListOfPrimitives()->FindObject(name); |
2263 | if (!o) { | |
2264 | Warning("FindXAxis", "%s not found in pad", name); | |
2265 | return 0; | |
b2e7f2d6 | 2266 | } |
2f92ff0e | 2267 | THStack* stack = dynamic_cast<THStack*>(o); |
2268 | if (!stack) { | |
2269 | Warning("FindXAxis", "%s is not a THStack", name); | |
2270 | return 0; | |
b2e7f2d6 | 2271 | } |
2f92ff0e | 2272 | if (!stack->GetHistogram()) { |
2273 | Warning("FindXAxis", "%s has no histogram", name); | |
2274 | return 0; | |
b2e7f2d6 | 2275 | } |
2f92ff0e | 2276 | TAxis* ret = stack->GetHistogram()->GetXaxis(); |
2277 | return ret; | |
b2e7f2d6 | 2278 | } |
2f92ff0e | 2279 | |
b2e7f2d6 | 2280 | //__________________________________________________________________ |
2f92ff0e | 2281 | /** |
2282 | * Fix the apperance of the axis in a stack | |
2283 | * | |
2284 | * @param stack stack of histogram | |
c6115ede | 2285 | * @param yd How the canvas is cut |
2f92ff0e | 2286 | * @param ytitle Y axis title |
2f92ff0e | 2287 | * @param ynDiv Divisions on Y axis |
c6115ede | 2288 | * @param force Whether to draw the stack first or not |
2f92ff0e | 2289 | */ |
797161e8 | 2290 | void FixAxis(THStack* stack, Double_t yd, const char* ytitle, |
2f92ff0e | 2291 | Int_t ynDiv=210, Bool_t force=true) |
b2e7f2d6 | 2292 | { |
e308a636 | 2293 | if (!stack) { |
2294 | Warning("FixAxis", "No stack passed for %s!", ytitle); | |
2295 | return; | |
2296 | } | |
2f92ff0e | 2297 | if (force) stack->Draw("nostack e1"); |
b2e7f2d6 | 2298 | |
2f92ff0e | 2299 | TH1* h = stack->GetHistogram(); |
e308a636 | 2300 | if (!h) { |
2301 | Warning("FixAxis", "Stack %s has no histogram", stack->GetName()); | |
2302 | return; | |
2303 | } | |
797161e8 | 2304 | Double_t s = 1/yd/1.2; |
2305 | // Info("FixAxis", "for %s, s=1/%f=%f", stack->GetName(), yd, s); | |
2306 | ||
d43c6cd0 | 2307 | h->SetXTitle("#it{#eta}"); |
2f92ff0e | 2308 | h->SetYTitle(ytitle); |
2309 | TAxis* xa = h->GetXaxis(); | |
2310 | TAxis* ya = h->GetYaxis(); | |
797161e8 | 2311 | |
2312 | // Int_t npixels = 20 | |
2313 | // Float_t dy = gPad->AbsPixeltoY(0) - gPad->AbsPixeltoY(npixels); | |
2314 | // Float_t ts = dy/(gPad->GetY2() - gPad->GetY1()); | |
2315 | ||
2f92ff0e | 2316 | if (xa) { |
797161e8 | 2317 | // xa->SetTitle(h->GetXTitle()); |
2f92ff0e | 2318 | // xa->SetTicks("+-"); |
2319 | xa->SetTitleSize(s*xa->GetTitleSize()); | |
2320 | xa->SetLabelSize(s*xa->GetLabelSize()); | |
2321 | xa->SetTickLength(s*xa->GetTickLength()); | |
797161e8 | 2322 | // xa->SetTitleOffset(xa->GetTitleOffset()/s); |
e308a636 | 2323 | |
2324 | if (stack != fResults) { | |
2325 | TAxis* rxa = fResults->GetXaxis(); | |
2326 | xa->Set(rxa->GetNbins(), rxa->GetXmin(), rxa->GetXmax()); | |
2327 | } | |
2f92ff0e | 2328 | } |
2329 | if (ya) { | |
797161e8 | 2330 | // ya->SetTitle(h->GetYTitle()); |
2f92ff0e | 2331 | ya->SetDecimals(); |
2332 | // ya->SetTicks("+-"); | |
2333 | ya->SetNdivisions(ynDiv); | |
2334 | ya->SetTitleSize(s*ya->GetTitleSize()); | |
fe0e93e6 | 2335 | ya->SetTitleOffset(/*1.15*/1.4*ya->GetTitleOffset()/s); |
2f92ff0e | 2336 | ya->SetLabelSize(s*ya->GetLabelSize()); |
b2e7f2d6 | 2337 | } |
b2e7f2d6 | 2338 | } |
797161e8 | 2339 | //__________________________________________________________________ |
2340 | /** | |
2341 | * Merge two histograms into one | |
2342 | * | |
2343 | * @param cen Central part | |
2344 | * @param fwd Forward part | |
2345 | * @param xlow On return, lower eta bound | |
2346 | * @param xhigh On return, upper eta bound | |
2347 | * | |
2348 | * @return Newly allocated histogram or null | |
2349 | */ | |
2350 | TH1* | |
2351 | Merge(const TH1* cen, const TH1* fwd, Double_t& xlow, Double_t& xhigh) | |
2352 | { | |
2353 | TH1* tmp = static_cast<TH1*>(fwd->Clone("tmp")); | |
2354 | TString name(fwd->GetName()); | |
2355 | name.ReplaceAll("Forward", "Merged"); | |
2356 | tmp->SetName(name); | |
2357 | ||
2358 | // tmp->SetMarkerStyle(28); | |
2359 | // tmp->SetMarkerColor(kBlack); | |
2360 | tmp->SetDirectory(0); | |
2361 | xlow = 100; | |
2362 | xhigh = -100; | |
2363 | for (Int_t i = 1; i <= tmp->GetNbinsX(); i++) { | |
2364 | Double_t cc = cen->GetBinContent(i); | |
2365 | Double_t cf = fwd->GetBinContent(i); | |
2366 | Double_t ec = cen->GetBinError(i); | |
2367 | Double_t ef = fwd->GetBinError(i); | |
2368 | Double_t nc = cf; | |
2369 | Double_t ne = ef; | |
2370 | if (cc < 0.001 && cf < 0.01) continue; | |
2371 | xlow = TMath::Min(tmp->GetXaxis()->GetBinLowEdge(i),xlow); | |
2372 | xhigh = TMath::Max(tmp->GetXaxis()->GetBinUpEdge(i),xhigh); | |
2373 | if (cc > 0.001) { | |
2374 | nc = cc; | |
2375 | ne = ec; | |
2376 | if (cf > 0.001) { | |
2377 | nc = (cf + cc) / 2; | |
2378 | ne = TMath::Sqrt(ec*ec + ef*ef); | |
2379 | } | |
2380 | } | |
2381 | tmp->SetBinContent(i, nc); | |
2382 | tmp->SetBinError(i, ne); | |
2383 | } | |
2384 | return tmp; | |
2385 | } | |
2386 | //____________________________________________________________________ | |
2387 | /** | |
2388 | * Fit @f$g(x;A_1,A_2,\sigma_1,\sigma_2)@f$ to histogram data | |
2389 | * | |
2390 | * @param tmp Histogram | |
2391 | * @param xlow Lower x bound | |
2392 | * @param xhigh Upper x bound | |
2393 | * | |
2394 | * @return Fitted function | |
2395 | */ | |
2396 | TF1* | |
2397 | FitMerged(TH1* tmp, Double_t xlow, Double_t xhigh) | |
2398 | { | |
2399 | TF1* tmpf = new TF1("tmpf", "gaus", xlow, xhigh); | |
2400 | tmp->Fit(tmpf, "NQ", ""); | |
2401 | tmp->SetDirectory(0); | |
2402 | ||
2403 | TF1* fit = new TF1("f", myFunc, xlow, xhigh, 4); | |
2404 | fit->SetParNames("a_{1}", "a_{2}", "#sigma_{1}", "#sigma_{2}"); | |
2405 | fit->SetParameters(tmpf->GetParameter(0), | |
2406 | .2, | |
2407 | tmpf->GetParameter(2), | |
2408 | tmpf->GetParameter(2)/4); | |
2409 | fit->SetParLimits(3, 0, 100); | |
2410 | fit->SetParLimits(4, 0, 100); | |
7095962e CHC |
2411 | const char* fitOpts = (fOptions & kVerbose ? "0W" : "Q0W"); |
2412 | tmp->Fit(fit,fitOpts,""); | |
797161e8 | 2413 | |
2414 | delete tmpf; | |
2415 | return fit; | |
2416 | } | |
2417 | //____________________________________________________________________ | |
2418 | /** | |
2419 | * Make band of systematic errors | |
2420 | * | |
2421 | * @param tmp Histogram | |
c6115ede | 2422 | * @param cen Central |
2423 | * @param fwd Forward | |
797161e8 | 2424 | * @param fit Fit |
2425 | */ | |
2426 | void | |
2427 | MakeSysError(TH1* tmp, TH1* cen, TH1* fwd, TF1* fit) | |
2428 | { | |
2429 | for (Int_t i = 1; i <= tmp->GetNbinsX(); i++) { | |
2430 | Double_t tc = tmp->GetBinContent(i); | |
2431 | if (tc < 0.01) continue; | |
2432 | Double_t fc = fwd->GetBinContent(i); | |
2433 | Double_t cc = cen->GetBinContent(i); | |
2434 | Double_t sysErr = fFwdSysErr; | |
2435 | if (cc > .01 && fc > 0.01) | |
2436 | sysErr = (fFwdSysErr+fCenSysErr) / 2; | |
2437 | else if (cc > .01) | |
2438 | sysErr = fCenSysErr; | |
2439 | Double_t x = tmp->GetXaxis()->GetBinCenter(i); | |
2440 | Double_t y = fit->Eval(x); | |
2441 | tmp->SetBinContent(i, y); | |
2442 | tmp->SetBinError(i,sysErr*y); | |
2443 | } | |
2444 | TString name(tmp->GetName()); | |
2445 | name.ReplaceAll("Merged", "SysError"); | |
2446 | tmp->SetName(name); | |
7095962e CHC |
2447 | tmp->SetMarkerColor(kSysErrColor); |
2448 | tmp->SetLineColor(kSysErrColor); | |
2449 | tmp->SetFillColor(kSysErrColor); | |
e42bc740 | 2450 | tmp->SetFillStyle(SYSERR_STYLE); |
797161e8 | 2451 | tmp->SetMarkerStyle(0); |
2452 | tmp->SetLineWidth(0); | |
2453 | } | |
2454 | void CorrectForward(TH1* h) const | |
2455 | { | |
7095962e | 2456 | if (!(fOptions & kRemoveOuters)) return; |
797161e8 | 2457 | |
2458 | for (Int_t i = 1; i <= h->GetNbinsX(); i++) { | |
2459 | Double_t eta = h->GetBinCenter(i); | |
2460 | if (TMath::Abs(eta) < 2.3) { | |
2461 | h->SetBinContent(i, 0); | |
2462 | h->SetBinError(i, 0); | |
2463 | } | |
2464 | } | |
2465 | } | |
2466 | void CorrectCentral(TH1* h) const | |
2467 | { | |
2468 | if (fClusterScale.IsNull()) return; | |
2469 | TString t(h->GetTitle()); | |
2470 | Info("CorrectCentral", "Replacing Central with Tracklets in %s", t.Data()); | |
2471 | t.ReplaceAll("Central", "Tracklets"); | |
2472 | h->SetTitle(t); | |
2473 | ||
2474 | TF1* cf = new TF1("clusterScale", fClusterScale); | |
2475 | for (Int_t i = 1; i <= h->GetNbinsX(); i++) { | |
2476 | Double_t eta = h->GetBinCenter(i); | |
2477 | Double_t f = cf->Eval(eta); | |
2478 | Double_t c = h->GetBinContent(i); | |
2479 | if (f < .1) f = 1; | |
2480 | h->SetBinContent(i, c / f); | |
2481 | } | |
2482 | delete cf; | |
2483 | } | |
33ab4a11 | 2484 | //____________________________________________________________________ |
2485 | void Export(const char* basename) | |
2486 | { | |
2487 | TString bname(basename); | |
2488 | bname.ReplaceAll(" ", "_"); | |
2489 | bname.ReplaceAll("-", "_"); | |
2490 | TString fname(Form("%s.C", bname.Data())); | |
2491 | ||
298a7264 | 2492 | std::ofstream outf(fname.Data()); |
2493 | if (!outf) { | |
33ab4a11 | 2494 | Error("Export", "Failed to open output file %s", fname.Data()); |
2495 | return; | |
2496 | } | |
7095962e | 2497 | Info("Export", "Exporting data to %s", fname.Data()); |
298a7264 | 2498 | outf << "// Create by dNdetaDrawer\n" |
2499 | << "void " << bname << "(THStack* stack, TLegend* l, Int_t m)\n" | |
2500 | << "{" | |
2501 | << " Int_t ma[] = { 24, 25, 26, 32,\n" | |
2502 | << " 20, 21, 22, 33,\n" | |
2503 | << " 34, 30, 29, 0, \n" | |
2504 | << " 23, 27 };\n" | |
2505 | << " Int_t mm = ((m < 20 || m > 34) ? 0 : ma[m-20]);\n\n"; | |
33ab4a11 | 2506 | TList* hists = fResults->GetHists(); |
2507 | TIter next(hists); | |
2508 | TH1* hist = 0; | |
2509 | while ((hist = static_cast<TH1*>(next()))) { | |
2510 | TString hname = hist->GetName(); | |
2511 | hname.Append(Form("_%04x", (gRandom->Integer(0xffff) & 0xffff))); | |
2512 | hist->SetName(hname); | |
2513 | hist->GetListOfFunctions()->Clear(); | |
298a7264 | 2514 | hist->SavePrimitive(outf, "nodraw"); |
33ab4a11 | 2515 | bool mirror = hname.Contains("mirror"); |
2516 | bool syserr = hname.Contains("SysError"); | |
2517 | if (!syserr) | |
298a7264 | 2518 | outf << " " << hname << "->SetMarkerStyle(" |
2519 | << (mirror ? "mm" : "m") << ");\n"; | |
33ab4a11 | 2520 | else |
298a7264 | 2521 | outf << " " << hname << "->SetMarkerStyle(1);\n"; |
2522 | outf << " stack->Add(" << hname | |
2523 | << (syserr ? ",\"e5\"" : "") << ");\n\n"; | |
33ab4a11 | 2524 | } |
2525 | UShort_t snn = fSNNString->GetUniqueID(); | |
2526 | // const char* sys = fSysString->GetTitle(); | |
2527 | TString eS; | |
2528 | if (snn == 2750) snn = 2760; | |
2529 | if (snn < 1000) eS = Form("%3dGeV", snn); | |
2530 | else if (snn % 1000 == 0) eS = Form("%dTeV", snn/1000); | |
2531 | else eS = Form("%4.2fTeV", float(snn)/1000); | |
298a7264 | 2532 | outf << " if (l) {\n" |
2533 | << " TLegendEntry* e = l->AddEntry(\"\",\"" << eS << "\",\"pl\");\n" | |
2534 | << " e->SetMarkerStyle(m);\n" | |
2535 | << " e->SetMarkerColor(kBlack);\n" | |
2536 | << " }\n" | |
2537 | << "}\n" << std::endl; | |
33ab4a11 | 2538 | } |
09d8cf0a | 2539 | /* @} */ |
73b32206 | 2540 | /** |
2541 | * Check if we have centrality dependent information, and we're not | |
2542 | * forcing to use minimum bias | |
2543 | * | |
2544 | * | |
2545 | * @return True if we should do centrality dependent ploting | |
2546 | */ | |
c8b1a7db | 2547 | Bool_t HasCent() const |
2548 | { | |
7095962e | 2549 | return fCentAxis && fCentAxis->GetNbins() > 0 && !(fOptions & kForceMB); |
c8b1a7db | 2550 | } |
2f92ff0e | 2551 | |
b2e7f2d6 | 2552 | |
2553 | ||
2554 | //__________________________________________________________________ | |
33ab4a11 | 2555 | /** |
2556 | * @{ | |
2557 | * @name Options | |
2558 | */ | |
7095962e CHC |
2559 | UInt_t fOptions; // Options |
2560 | UInt_t fFormats; // Output formats | |
33ab4a11 | 2561 | UShort_t fShowOthers; // Show other data |
2562 | /* @} */ | |
2563 | /** | |
2564 | * @{ | |
2565 | * @name Settings | |
2566 | */ | |
2567 | UShort_t fRebin; // Rebinning factor | |
2568 | Double_t fFwdSysErr; // Systematic error in forward range | |
2569 | Double_t fCenSysErr; // Systematic error in central range | |
e308a636 | 2570 | TString fTitle; // Title on plot |
d43c6cd0 | 2571 | TString fBase; // Base name of output |
33ab4a11 | 2572 | TString fClusterScale; // Scaling of clusters to tracklets |
b615e211 | 2573 | TString fFinalMC; // Final MC correction file name |
2574 | TString fEmpirical; // Empirical correction file name | |
33ab4a11 | 2575 | /* @} */ |
2576 | /** | |
2577 | * @{ | |
2578 | * @name Read (or set) information | |
2579 | */ | |
e308a636 | 2580 | TNamed* fTrigString; // Trigger string (read, or set) |
0be6c8cd | 2581 | TNamed* fNormString; // Normalisation string (read, or set) |
e308a636 | 2582 | TNamed* fSNNString; // Energy string (read, or set) |
2583 | TNamed* fSysString; // Collision system string (read or set) | |
2584 | TAxis* fVtxAxis; // Vertex cuts (read or set) | |
2585 | TAxis* fCentAxis; // Centrality axis | |
bfab35d9 | 2586 | Float_t fTriggerEff; // Trigger efficiency |
d43c6cd0 | 2587 | Bool_t fExtTriggerEff;// True if read externally |
f71c22f5 CHC |
2588 | UShort_t fCentMin; // Least centrality to plot |
2589 | UShort_t fCentMax; // Largest centrality to plot | |
33ab4a11 | 2590 | /* @} */ |
2591 | /** | |
2592 | * @{ | |
2593 | * @name Resulting plots | |
2594 | */ | |
e308a636 | 2595 | THStack* fResults; // Stack of results |
2596 | THStack* fRatios; // Stack of ratios | |
2597 | THStack* fLeftRight; // Left-right asymmetry | |
2598 | TMultiGraph* fOthers; // Older data | |
2599 | TH1* fTriggers; // Number of triggers | |
9ecab72f | 2600 | TH1* fTruth; // Pointer to truth |
33ab4a11 | 2601 | /* @} */ |
e308a636 | 2602 | RangeParam* fRangeParam; // Parameter object for range zoom |
d43c6cd0 | 2603 | |
2604 | static const Float_t kRightMargin; | |
2605 | static const Int_t kFont; | |
2606 | static const Int_t kAliceBlue; | |
2607 | static const Int_t kAliceRed; | |
2608 | static const Int_t kAlicePurple; | |
2609 | static const Int_t kAliceYellow; | |
7095962e | 2610 | static const Int_t kSysErrColor; |
b2e7f2d6 | 2611 | }; |
2612 | ||
fe0e93e6 | 2613 | const Float_t dNdetaDrawer::kRightMargin = 0.01; |
d43c6cd0 | 2614 | const Int_t dNdetaDrawer::kFont = 42; // 132 for serif |
2615 | const Int_t dNdetaDrawer::kAliceBlue = TColor::GetColor(40, 58, 68); | |
2616 | const Int_t dNdetaDrawer::kAliceRed = TColor::GetColor(226, 0, 26); | |
2617 | const Int_t dNdetaDrawer::kAlicePurple = TColor::GetColor(202, 71, 67); | |
2618 | const Int_t dNdetaDrawer::kAliceYellow = TColor::GetColor(238, 125, 17); | |
7095962e | 2619 | const Int_t dNdetaDrawer::kSysErrColor = SYSERR_COLOR; |
d43c6cd0 | 2620 | |
2621 | ||
797161e8 | 2622 | //____________________________________________________________________ |
2623 | /** | |
2624 | * Function to calculate | |
2625 | * @f[ | |
2626 | * g(x;A_1,A_2,\sigma_1,\sigma_2) = | |
2627 | * A_1\left(\frac{1}{2\pi\sigma_1}e^{(x/\sigma_1)^2} - | |
2628 | * A_2\frac{1}{2\pi\sigma_2}e^{(x/\sigma_2)^2}\right) | |
2629 | * @f] | |
2630 | * | |
2631 | * @param xp Pointer to x array | |
2632 | * @param pp Pointer to parameter array | |
2633 | * | |
2634 | * @return @f$g(x;A_1,A_2,\sigma_1,\sigma_2)@f$ | |
2635 | */ | |
2636 | Double_t myFunc(Double_t* xp, Double_t* pp) | |
2637 | { | |
2638 | Double_t x = xp[0]; | |
2639 | Double_t a1 = pp[0]; | |
2640 | Double_t a2 = pp[1]; | |
2641 | Double_t s1 = pp[2]; | |
2642 | Double_t s2 = pp[3]; | |
2643 | return a1*(TMath::Gaus(x, 0, s1) - a2 * TMath::Gaus(x, 0, s2)); | |
2644 | } | |
2645 | ||
b2e7f2d6 | 2646 | //=== Stuff for auto zooming ========================================= |
290052e7 | 2647 | /** |
2648 | * Update canvas range | |
2649 | * | |
2650 | * @param p Parameter | |
2651 | */ | |
b2e7f2d6 | 2652 | void UpdateRange(dNdetaDrawer::RangeParam* p) |
2653 | { | |
2654 | if (!p) { | |
2655 | Warning("UpdateRange", "No parameters %p", p); | |
2656 | return; | |
2657 | } | |
2658 | if (!p->fMasterAxis) { | |
2659 | Warning("UpdateRange", "No master axis %p", p->fMasterAxis); | |
2660 | return; | |
2661 | } | |
2662 | Int_t first = p->fMasterAxis->GetFirst(); | |
2663 | Int_t last = p->fMasterAxis->GetLast(); | |
2664 | Double_t x1 = p->fMasterAxis->GetBinCenter(first); | |
2665 | Double_t x2 = p->fMasterAxis->GetBinCenter(last); | |
b2e7f2d6 | 2666 | |
2667 | if (p->fSlave1Axis) { | |
2668 | Int_t i1 = p->fSlave1Axis->FindBin(x1); | |
2669 | Int_t i2 = p->fSlave1Axis->FindBin(x2); | |
2670 | p->fSlave1Axis->SetRange(i1, i2); | |
2671 | p->fSlave1Pad->Modified(); | |
2672 | p->fSlave1Pad->Update(); | |
2673 | } | |
2674 | if (p->fSlave2Axis) { | |
2675 | Int_t i1 = p->fSlave2Axis->FindBin(x1); | |
2676 | Int_t i2 = p->fSlave2Axis->FindBin(x2); | |
2677 | p->fSlave2Axis->SetRange(i1, i2); | |
2678 | p->fSlave2Pad->Modified(); | |
2679 | p->fSlave2Pad->Update(); | |
2680 | } | |
2681 | TCanvas* c = gPad->GetCanvas(); | |
2682 | c->cd(); | |
2683 | } | |
2684 | ||
2685 | //____________________________________________________________________ | |
290052e7 | 2686 | /** |
2687 | * Called when user changes X range | |
2688 | * | |
2689 | * @param p Parameter | |
2690 | */ | |
b2e7f2d6 | 2691 | void RangeExec(dNdetaDrawer::RangeParam* p) |
2692 | { | |
2693 | // Event types: | |
2694 | // 51: Mouse move | |
2695 | // 53: | |
2696 | // 1: Button down | |
2697 | // 21: Mouse drag | |
2698 | // 11: Mouse release | |
2699 | // dNdetaDrawer::RangeParam* p = | |
2700 | // reinterpret_cast<dNdetaDrawer::RangeParam*>(addr); | |
2701 | Int_t event = gPad->GetEvent(); | |
2702 | TObject *select = gPad->GetSelected(); | |
2703 | if (event == 53) { | |
2704 | UpdateRange(p); | |
2705 | return; | |
2706 | } | |
2707 | if (event != 11 || !select || select != p->fMasterAxis) return; | |
2708 | UpdateRange(p); | |
2709 | } | |
2710 | ||
33ab4a11 | 2711 | //=== Steering functions |
2712 | //============================================== | |
290052e7 | 2713 | /** |
2714 | * Display usage information | |
2715 | * | |
2716 | */ | |
33ab4a11 | 2717 | void |
2718 | Usage() | |
2719 | { | |
7095962e CHC |
2720 | std::ostream& o = std::cout; |
2721 | o << "Usage: DrawdNdeta(FILE,TITLE,REBIN,OTHERS,FLAGS," | |
2722 | << "SNN,SYS,TRIG,IPZMIN,IPZMAX,BASE,FMT)\n" | |
2723 | << " const char* FILE File name to open (\"forward_dndeta.root\")\n" | |
2724 | << " const char* TITLE Title to put on plot (\"\")\n" | |
2725 | << " UShort_t REBIN Rebinning factor (1)\n" | |
2726 | << " UShort_t OTHERS Other data to draw - more below (0x7)\n" | |
2727 | << " UShort_t FLAGS Visualisation flags - more below (0x7)\n" | |
2728 | << " UShort_t SYS (optional) 1:pp, 2:PbPb, 3:pPb\n" | |
2729 | << " UShort_t SNN (optional) sqrt(s_NN) in GeV\n" | |
2730 | << " UShort_t TRIG (optional) 1: INEL, 2: INEL>0, 4: NSD, ...\n" | |
2731 | << " Float_t EFF (optional) Trigger efficiency\n" | |
2732 | << " Float_t IPZMIN (optional) Least z coordinate of IP\n" | |
2733 | << " Float_t IPZMAX (optional) Largest z coordinate of IP\n" | |
2734 | << " const char* BASE (optional) base name of output files\n" | |
2735 | << " UShort_t FMT (optional) Output formats\n" | |
2736 | << "\n"; | |
2737 | o << " OTHERS is a bit mask of\n" | |
2738 | << " 0x1 Show UA5 data (INEL,NSD, ppbar, 900GeV)\n" | |
2739 | << " 0x2 Show CMS data (NSD, pp)\n" | |
2740 | << " 0x4 Show published ALICE data (INEL,INEL>0,NSD, pp)\n" | |
2741 | << " 0x8 Show event genertor data\n" | |
2742 | << "\n"; | |
2743 | o << " FMT is a bit mask of\n" | |
2744 | << " 0x1 Make PNG output\n" | |
2745 | << " 0x2 Make PDF output\n" | |
2746 | << " 0x4 Make ROOT file output\n" | |
2747 | << " 0x8 Make ROOT script output\n" | |
2748 | << "\n"; | |
2749 | o << " FLAGS is a bit mask of\n" | |
2750 | << " 0x1 Show ratios of data to other data and possibly MC\n" | |
2751 | << " 0x2 Show left-right asymmetry\n" | |
2752 | << " 0x4 Show systematic error band\n" | |
2753 | << " 0x8 Show individual ring results (INEL only)\n" | |
2754 | << " 0x10 Cut edges when rebinning\n" | |
2755 | << " 0x20 Remove FMDxO points\n" | |
2756 | << " 0x40 Apply `final MC' correction\n" | |
2757 | << " 0x80 Apply `Emperical' correction\n" | |
2758 | << " 0x100 Force use of MB\n" | |
2759 | << " 0x200 Mirror data\n" | |
2760 | << " 0x400 Export results to script\n" | |
2761 | << " 0x800 Add code to do combined zooms on eta axis\n" | |
2762 | << " 0x1000 Assume old-style input\n" | |
2763 | << " 0x2000 Be verbose\n" | |
2764 | << " 0x4000 Hi-res batch output\n" | |
2765 | << " 0x8000 Add aditional white-space above results\n" | |
2766 | << "\n"; | |
2767 | o << "0x200 requires the file forward_dndetamc.root\n" | |
2768 | << "0x400 requires the file EmpiricalCorrection.root\n" | |
2769 | << "To specify that you want ratios, force MB, apply empirical " | |
2770 | << "correction, and export to script, set flags to\n\n" | |
2771 | << " 0x1|0x100|0x80|0x400=0x581\n" | |
2772 | << std::endl; | |
2773 | ||
33ab4a11 | 2774 | } |
2775 | ||
2776 | //____________________________________________________________________ | |
ffca499d | 2777 | /** |
2778 | * Draw @f$ dN/d\eta@f$ | |
2779 | * | |
2780 | * @param filename File name | |
ffca499d | 2781 | * @param title Title |
2782 | * @param rebin Rebinning factor | |
c6115ede | 2783 | * @param others What other data to show |
2784 | * @param flags Flags | |
ffca499d | 2785 | * @param sNN (optional) Collision energy [GeV] |
2786 | * @param sys (optional) Collision system (1: pp, 2: PbPb) | |
2787 | * @param trg (optional) Trigger (1: INEL, 2: INEL>0, 4: NSD) | |
c8b1a7db | 2788 | * @param eff (optional) Trigger efficiency |
ffca499d | 2789 | * @param vzMin Least @f$ v_z@f$ |
2790 | * @param vzMax Largest @f$ v_z@f$ | |
c8b1a7db | 2791 | * @param base Base name |
2792 | * @param outflg Output flags | |
ffca499d | 2793 | * |
bd6f5206 | 2794 | * @ingroup pwglf_forward_dndeta |
ffca499d | 2795 | */ |
b2e7f2d6 | 2796 | void |
2797 | DrawdNdeta(const char* filename="forward_dndeta.root", | |
bad9a3c1 | 2798 | const char* title="", |
b2e7f2d6 | 2799 | UShort_t rebin=5, |
797161e8 | 2800 | UShort_t others=0x7, |
fe0e93e6 | 2801 | UInt_t flags=dNdetaDrawer::kDefaultOptions, |
fe52e455 | 2802 | UShort_t sNN=0, |
2803 | UShort_t sys=0, | |
ffca499d | 2804 | UShort_t trg=0, |
bfab35d9 | 2805 | Float_t eff=0, |
f71c22f5 CHC |
2806 | UShort_t centMin=0, |
2807 | UShort_t centMax=100, | |
fe52e455 | 2808 | Float_t vzMin=999, |
d43c6cd0 | 2809 | Float_t vzMax=-999, |
2810 | const char* base="", | |
7095962e | 2811 | UShort_t outflg=dNdetaDrawer::kAllFormats) |
b2e7f2d6 | 2812 | { |
33ab4a11 | 2813 | TString fname(filename); |
2814 | fname.ToLower(); | |
2815 | if (fname.CompareTo("help") == 0 || | |
2816 | fname.CompareTo("--help") == 0) { | |
2817 | Usage(); | |
2818 | return; | |
2819 | } | |
ffca499d | 2820 | dNdetaDrawer* pd = new dNdetaDrawer; |
2821 | dNdetaDrawer& d = *pd; | |
b2e7f2d6 | 2822 | d.SetRebin(rebin); |
b2e7f2d6 | 2823 | d.SetTitle(title); |
797161e8 | 2824 | d.SetShowOthers(others); |
d43c6cd0 | 2825 | d.SetBase(base); |
797161e8 | 2826 | // d.fClusterScale = "1.06 -0.003*x +0.0119*x*x"; |
2f92ff0e | 2827 | // Do the below if your input data does not contain these settings |
fe52e455 | 2828 | if (sNN > 0) d.SetSNN(sNN); // Collision energy per nucleon pair (GeV) |
2829 | if (sys > 0) d.SetSys(sys); // Collision system (1:pp, 2:PbPB) | |
2830 | if (trg > 0) d.SetTrigger(trg); // Collision trigger (1:INEL, 2:INEL>0, 4:NSD) | |
bfab35d9 | 2831 | if (eff > 0) d.SetTriggerEfficiency(eff); // Trigger efficiency |
fe52e455 | 2832 | if (vzMin < 999 && vzMax > -999) |
2833 | d.SetVertexRange(vzMin,vzMax); // Collision vertex range (cm) | |
f71c22f5 | 2834 | d.SetCentralityRange(centMin,centMax); // Collision vertex range (cm) |
7095962e | 2835 | d.Run(filename, flags, outflg); |
b2e7f2d6 | 2836 | } |
2f92ff0e | 2837 | //____________________________________________________________________ |
2838 | // | |
2839 | // EOF | |
2840 | // | |
2841 |