1 | /* ======================================================================== *\
|
---|
2 | !
|
---|
3 | ! *
|
---|
4 | ! * This file is part of MARS, the MAGIC Analysis and Reconstruction
|
---|
5 | ! * Software. It is distributed to you in the hope that it can be a useful
|
---|
6 | ! * and timesaving tool in analysing Data of imaging Cerenkov telescopes.
|
---|
7 | ! * It is distributed WITHOUT ANY WARRANTY.
|
---|
8 | ! *
|
---|
9 | ! * Permission to use, copy, modify and distribute this software and its
|
---|
10 | ! * documentation for any purpose is hereby granted without fee,
|
---|
11 | ! * provided that the above copyright notice appear in all copies and
|
---|
12 | ! * that both that copyright notice and this permission notice appear
|
---|
13 | ! * in supporting documentation. It is provided "as is" without express
|
---|
14 | ! * or implied warranty.
|
---|
15 | ! *
|
---|
16 | !
|
---|
17 | !
|
---|
18 | ! Author(s): Thomas Bretz 2002 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | !
|
---|
20 | ! Copyright: MAGIC Software Development, 2000-2010
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | /////////////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MH3
|
---|
28 | //
|
---|
29 | // With this histogram you can fill a histogram with up to three
|
---|
30 | // variables from Mars parameter containers in an eventloop.
|
---|
31 | //
|
---|
32 | // In the constructor you can give up to three variables which should be
|
---|
33 | // filled in the histogram. Dependend on the number of given variables
|
---|
34 | // (data members) a TH1D, TH2D or TH3D is created.
|
---|
35 | // Specify the data mamber like the following:
|
---|
36 | // "MHillas.fLength"
|
---|
37 | // Where MHillas is the name of the parameter container in the parameter
|
---|
38 | // list and fLength is the name of the data member which should be filled
|
---|
39 | // in the histogram. Assuming that your MHillas container has a different
|
---|
40 | // name (MyHillas) the name to give would be:
|
---|
41 | // "MyHillas.fLength"
|
---|
42 | //
|
---|
43 | // If you want to use a different unit for histogramming use SetScaleX,
|
---|
44 | // SetScaleY and SetScaleZ.
|
---|
45 | //
|
---|
46 | //
|
---|
47 | // Binning/Binning name
|
---|
48 | // =====================
|
---|
49 | //
|
---|
50 | // The axis binning is retrieved from the parameter list, too. Create a
|
---|
51 | // MBinning with the name "Binning" plus the name of your MH3 container
|
---|
52 | // plus the axis name ("X", "Y" or "Z") and add it to the parameter list.
|
---|
53 | //
|
---|
54 | // If the binning should have a different name than the histogram name
|
---|
55 | // the binning name can be added to the name, eg.:
|
---|
56 | // SetName("MyHistName;MyXBins;MyYBins")
|
---|
57 | // Instead of BinningMyHistName[XYZ] the parameter list will be searched
|
---|
58 | // for BinningMyXBinning, BinningMyYBins and BinningMyHistNameZ
|
---|
59 | //
|
---|
60 | // If you don't want to use a MBinning object from the parameter list
|
---|
61 | // you can also set one directly, for example
|
---|
62 | // MBinning bins(10, 0, 1);
|
---|
63 | // mh3.SetBinningX(&bins);
|
---|
64 | // You must not delete the MBinning object before the class has been
|
---|
65 | // PreProcessed.
|
---|
66 | //
|
---|
67 | //
|
---|
68 | // Axis titles
|
---|
69 | // ===========
|
---|
70 | //
|
---|
71 | // 1) If no other title is given the rule for this axis is used.
|
---|
72 | // 2) If the MBinning used for this axis has a non-default Title
|
---|
73 | // (MBinning::HasTitle) this title is used for the corresponding axis
|
---|
74 | // 3) If the MH3 has a non-default title (MH3::SetTitle called)
|
---|
75 | // this title is set as the histogram title. It can be used to overwrite
|
---|
76 | // the axis titles. For more information see TH1::SetTitle, eg.
|
---|
77 | // SetTitle("MyHist;x[mm];y[cm];Counts");
|
---|
78 | //
|
---|
79 | //
|
---|
80 | // Labels
|
---|
81 | // ======
|
---|
82 | //
|
---|
83 | // To use labels at an axis you have to initialize this for the axis
|
---|
84 | // by either calling InitLabels(Labels_t) or setiting a DefaultLabel.
|
---|
85 | // For the axis for which the labels have been initialized the
|
---|
86 | // number returned by the given corresponding phrase is converted
|
---|
87 | // to int with TMath::Nint and used as a label. If you want to replace
|
---|
88 | // this id by a named label you can call DefineLabel to do that.
|
---|
89 | // Several ids can be replaced by the same label. If you define
|
---|
90 | // named labels for every label which was not defined the default
|
---|
91 | // is used, if any, otherwise an unnamed label is created.
|
---|
92 | //
|
---|
93 | // In the case of an axis with labels the axis-title cannot be
|
---|
94 | // set via a MBinning, because the MBinning is not evaluated.
|
---|
95 | //
|
---|
96 | // Please note that for some reason not all combinations of
|
---|
97 | // labels, dimensions and weights are available in the root-
|
---|
98 | // histogram classes. Please check the MH3::Fill function to see
|
---|
99 | // whether your combination is supported.
|
---|
100 | //
|
---|
101 | //
|
---|
102 | // Examples:
|
---|
103 | // =========
|
---|
104 | //
|
---|
105 | // 1) MH3 myhist("MHillas.fLength");
|
---|
106 | // myhist.SetName("MyHist");
|
---|
107 | // myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree
|
---|
108 | // MBinning bins("BinningMyHistX", "Title for my x-axis [Hz]");
|
---|
109 | // bins.SetEdges(10, 0, 150);
|
---|
110 | // plist.AddToList(&bins);
|
---|
111 | //
|
---|
112 | // 2) MH3 myhist("MHillas.fLength");
|
---|
113 | // myhist.SetName("MyHist;MyX");
|
---|
114 | // myhist.SetTitle("Histogram Title;X-Title [mm];Counts");
|
---|
115 | // MBinning bins("BinningMyX");
|
---|
116 | // bins.SetEdges(10, 0, 150);
|
---|
117 | // plist.AddToList(&bins);
|
---|
118 | //
|
---|
119 | // 3) MH3 myhist("MTriggerPatter.GetUnprescaled");
|
---|
120 | // myhist.SetWeight("1./MRawRunHeader.GetRunLength");
|
---|
121 | // myhist.SetTitle("Rate of the trigger pattern [Hz];Run Number;Trigger Pattern;Rate [Hz]");
|
---|
122 | // myhist.InitLabels(MH3::kLabelsXY);
|
---|
123 | // myhist.DefaultLabelY("UNKNOWN"); // Lvl1
|
---|
124 | // myhist.DefineLabelY( 1, "Trig"); // Lvl1
|
---|
125 | // myhist.DefineLabelY( 2, "Cal"); // Cal
|
---|
126 | // myhist.DefineLabelY( 4, "Trig"); // Lvl2
|
---|
127 | // myhist.DefineLabelY( 8, "Ped"); // Ped
|
---|
128 | //
|
---|
129 | //
|
---|
130 | // Class Version 1:
|
---|
131 | // ----------------
|
---|
132 | // - MData *fData[3];
|
---|
133 | // + MData *fData[4];
|
---|
134 | //
|
---|
135 | // Class Version 2:
|
---|
136 | // ----------------
|
---|
137 | // - MDataChain *fData[3]; // Object from which the data is filled
|
---|
138 | // + MData *fData[3]; // Object from which the data is filled
|
---|
139 | //
|
---|
140 | // Class Version 3:
|
---|
141 | // ----------------
|
---|
142 | // - Byte_t fStyleBits
|
---|
143 | // + MBinning fBins[3]
|
---|
144 | //
|
---|
145 | // Class Version 5:
|
---|
146 | // ----------------
|
---|
147 | // + TFormula *fConversion
|
---|
148 | //
|
---|
149 | // Class Version 6:
|
---|
150 | // ----------------
|
---|
151 | // + MData *fWeight;
|
---|
152 | //
|
---|
153 | // Class Version 7:
|
---|
154 | // ----------------
|
---|
155 | // + MData *fData[4]
|
---|
156 | // + Double_t fScale[4]
|
---|
157 | // - MData *fData[3]
|
---|
158 | // - Double_t fScale[3]
|
---|
159 | //
|
---|
160 | /////////////////////////////////////////////////////////////////////////////
|
---|
161 | #include "MH3.h"
|
---|
162 |
|
---|
163 | #include <ctype.h> // tolower
|
---|
164 | #include <stdlib.h> // atoi (Ubuntu 8.10)
|
---|
165 | #include <fstream>
|
---|
166 |
|
---|
167 | #include <TMath.h>
|
---|
168 | #include <TFormula.h>
|
---|
169 |
|
---|
170 | #include <THashList.h>
|
---|
171 | #include <TObjString.h>
|
---|
172 |
|
---|
173 | //#include <TPad.h>
|
---|
174 | #include <TStyle.h>
|
---|
175 | #include <TCanvas.h>
|
---|
176 |
|
---|
177 | #include <TH2.h>
|
---|
178 | #include <TH3.h>
|
---|
179 | #include <TProfile.h>
|
---|
180 | #include <TProfile2D.h>
|
---|
181 | #include <TProfile3D.h>
|
---|
182 |
|
---|
183 | #include "MLog.h"
|
---|
184 | #include "MLogManip.h"
|
---|
185 |
|
---|
186 | #include "MString.h"
|
---|
187 |
|
---|
188 | #include "MParList.h"
|
---|
189 | #include "MBinning.h"
|
---|
190 | #include "MDataPhrase.h"
|
---|
191 |
|
---|
192 | ClassImp(MH3);
|
---|
193 |
|
---|
194 | using namespace std;
|
---|
195 |
|
---|
196 | const TString MH3::gsDefName = "MH3";
|
---|
197 | const TString MH3::gsDefTitle = "Container for a n-D Mars Histogram";
|
---|
198 |
|
---|
199 | // --------------------------------------------------------------------------
|
---|
200 | //
|
---|
201 | // Set fStyleBits to 0, Reset fBins to NULL, fScale to 1, set name and title
|
---|
202 | // to gsDefName and gsDefTitle and if fHist!=NULL UseCurrentStyle and
|
---|
203 | // SetDirectory(0)
|
---|
204 | //
|
---|
205 | void MH3::Init()
|
---|
206 | {
|
---|
207 | fStyleBits = 0;
|
---|
208 |
|
---|
209 | fWeight = NULL;
|
---|
210 |
|
---|
211 | fData[0] = NULL;
|
---|
212 | fData[1] = NULL;
|
---|
213 | fData[2] = NULL;
|
---|
214 | fData[3] = NULL;
|
---|
215 |
|
---|
216 | fBins[0] = NULL;
|
---|
217 | fBins[1] = NULL;
|
---|
218 | fBins[2] = NULL;
|
---|
219 |
|
---|
220 | fScale[0] = 1;
|
---|
221 | fScale[1] = 1;
|
---|
222 | fScale[2] = 1;
|
---|
223 | fScale[3] = 1;
|
---|
224 |
|
---|
225 | fConversion = NULL;
|
---|
226 |
|
---|
227 | fName = gsDefName;
|
---|
228 | fTitle = gsDefTitle;
|
---|
229 |
|
---|
230 | if (!fHist)
|
---|
231 | return;
|
---|
232 |
|
---|
233 | fHist->UseCurrentStyle();
|
---|
234 | fHist->SetDirectory(NULL);
|
---|
235 | }
|
---|
236 |
|
---|
237 | // --------------------------------------------------------------------------
|
---|
238 | //
|
---|
239 | // Default constructor.
|
---|
240 | //
|
---|
241 | MH3::MH3(const Int_t dim, Type_t type) : fDimension(dim), fHist(NULL)
|
---|
242 | {
|
---|
243 | // FIXME?
|
---|
244 | switch (type)
|
---|
245 | {
|
---|
246 | case kHistogram:
|
---|
247 | if (fDimension>3)
|
---|
248 | fDimension=3;
|
---|
249 | break;
|
---|
250 | case kProfile:
|
---|
251 | case kProfileSpread:
|
---|
252 | fDimension = -TMath::Abs(fDimension);
|
---|
253 | if (fDimension<-2)
|
---|
254 | fDimension = -2;
|
---|
255 | break;
|
---|
256 | }
|
---|
257 |
|
---|
258 | switch (fDimension)
|
---|
259 | {
|
---|
260 | case 1:
|
---|
261 | fHist = new TH1D;
|
---|
262 | fHist->SetYTitle("Counts");
|
---|
263 | break;
|
---|
264 | case -1:
|
---|
265 | fHist = new TProfile;
|
---|
266 | fHist->SetYTitle("Average");
|
---|
267 | if (type==kProfileSpread)
|
---|
268 | static_cast<TProfile*>(fHist)->SetErrorOption("s");
|
---|
269 | break;
|
---|
270 | case 2:
|
---|
271 | fHist = new TH2D;
|
---|
272 | fHist->SetZTitle("Counts");
|
---|
273 | break;
|
---|
274 | case -2:
|
---|
275 | fHist = new TProfile2D;
|
---|
276 | fHist->SetZTitle("Average");
|
---|
277 | static_cast<TProfile2D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
|
---|
278 | break;
|
---|
279 | case 3:
|
---|
280 | fHist = new TH3D;
|
---|
281 | break;
|
---|
282 | case -3:
|
---|
283 | fHist = new TProfile2D;
|
---|
284 | fHist->SetZTitle("Average");
|
---|
285 | static_cast<TProfile2D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
|
---|
286 | break;
|
---|
287 | }
|
---|
288 |
|
---|
289 | Init();
|
---|
290 | }
|
---|
291 |
|
---|
292 | // --------------------------------------------------------------------------
|
---|
293 | //
|
---|
294 | // Creates an TH1D. memberx is filled into the X-bins. For a more detailed
|
---|
295 | // description see the class description above.
|
---|
296 | //
|
---|
297 | MH3::MH3(const char *memberx, Type_t type) : fDimension(1)
|
---|
298 | {
|
---|
299 | fHist = new TH1D;
|
---|
300 | fHist->SetYTitle("Counts");
|
---|
301 |
|
---|
302 | Init();
|
---|
303 |
|
---|
304 | fData[0] = new MDataPhrase(memberx);
|
---|
305 | }
|
---|
306 |
|
---|
307 | // --------------------------------------------------------------------------
|
---|
308 | //
|
---|
309 | // Adapt a given histogram
|
---|
310 | //
|
---|
311 | MH3::MH3(const TH1 &h1) : fDimension(1)
|
---|
312 | {
|
---|
313 | if (h1.InheritsFrom(TH3::Class()))
|
---|
314 | fDimension = 3;
|
---|
315 | if (h1.InheritsFrom(TH2::Class()))
|
---|
316 | fDimension = 2;
|
---|
317 |
|
---|
318 | if (h1.InheritsFrom(TProfile3D::Class()))
|
---|
319 | {
|
---|
320 | fDimension = -3;
|
---|
321 | *fLog << warn << "WARNING - MH3::MH3(TH1&) does not support TProfile3D." << endl;
|
---|
322 | }
|
---|
323 | if (h1.InheritsFrom(TProfile2D::Class()))
|
---|
324 | fDimension = -2;
|
---|
325 | if (h1.InheritsFrom(TProfile::Class()))
|
---|
326 | fDimension = -1;
|
---|
327 |
|
---|
328 | fHist = (TH1*)h1.Clone();
|
---|
329 |
|
---|
330 | Init(); // Before without SeUseCurrentStyle!
|
---|
331 |
|
---|
332 | switch (fDimension)
|
---|
333 | {
|
---|
334 | case 3:
|
---|
335 | case -2:
|
---|
336 | fData[2] = new MDataPhrase(h1.GetZaxis()->GetTitle());
|
---|
337 | case 2:
|
---|
338 | case -1:
|
---|
339 | fData[1] = new MDataPhrase(h1.GetYaxis()->GetTitle());
|
---|
340 | case 1:
|
---|
341 | fData[0] = new MDataPhrase(h1.GetXaxis()->GetTitle());
|
---|
342 | }
|
---|
343 | }
|
---|
344 |
|
---|
345 | // --------------------------------------------------------------------------
|
---|
346 | //
|
---|
347 | // Creates an TH2D. memberx is filled into the X-bins. membery is filled
|
---|
348 | // into the Y-bins. For a more detailed description see the class
|
---|
349 | // description above.
|
---|
350 | //
|
---|
351 | MH3::MH3(const char *memberx, const char *membery, Type_t type)
|
---|
352 | : fDimension(type==kHistogram?2:-1)
|
---|
353 | {
|
---|
354 |
|
---|
355 | switch (fDimension)
|
---|
356 | {
|
---|
357 | case 2:
|
---|
358 | fHist = static_cast<TH1*>(new TH2D);
|
---|
359 | break;
|
---|
360 | case -1:
|
---|
361 | fHist = static_cast<TH1*>(new TProfile);
|
---|
362 | if (type==kProfileSpread)
|
---|
363 | static_cast<TProfile*>(fHist)->SetErrorOption("s");
|
---|
364 |
|
---|
365 | break;
|
---|
366 | }
|
---|
367 |
|
---|
368 | fHist->SetZTitle(fDimension>0?"Counts":"Average");
|
---|
369 |
|
---|
370 | Init();
|
---|
371 |
|
---|
372 | fData[0] = new MDataPhrase(memberx);
|
---|
373 | fData[1] = new MDataPhrase(membery);
|
---|
374 | }
|
---|
375 |
|
---|
376 | // --------------------------------------------------------------------------
|
---|
377 | //
|
---|
378 | // Creates an TH3D. memberx is filled into the X-bins. membery is filled
|
---|
379 | // into the Y-bins. membery is filled into the Z-bins. For a more detailed
|
---|
380 | // description see the class description above.
|
---|
381 | //
|
---|
382 | MH3::MH3(const char *memberx, const char *membery, const char *memberz, Type_t type)
|
---|
383 | : fDimension(type==kHistogram?3:-2)
|
---|
384 | {
|
---|
385 | switch (fDimension)
|
---|
386 | {
|
---|
387 | case 3:
|
---|
388 | fHist = static_cast<TH1*>(new TH3D);
|
---|
389 | break;
|
---|
390 | case -2:
|
---|
391 | fHist = static_cast<TH1*>(new TProfile2D);
|
---|
392 | static_cast<TProfile2D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
|
---|
393 | }
|
---|
394 |
|
---|
395 | Init();
|
---|
396 |
|
---|
397 | fData[0] = new MDataPhrase(memberx);
|
---|
398 | fData[1] = new MDataPhrase(membery);
|
---|
399 | fData[2] = new MDataPhrase(memberz);
|
---|
400 | }
|
---|
401 |
|
---|
402 | // --------------------------------------------------------------------------
|
---|
403 | //
|
---|
404 | // Creates an TH3D. memberx is filled into the X-bins. membery is filled
|
---|
405 | // into the Y-bins. membery is filled into the Z-bins. Weight is used as a
|
---|
406 | // weight for the profile histogram. For a more detailed description see the
|
---|
407 | // class description above.
|
---|
408 | //
|
---|
409 | MH3::MH3(const char *memberx, const char *membery, const char *memberz, const char *weight, Type_t type)
|
---|
410 | : fDimension(-3)
|
---|
411 | {
|
---|
412 | fHist = static_cast<TH1*>(new TProfile3D);
|
---|
413 | static_cast<TProfile3D*>(fHist)->BuildOptions(0, 0, type==kProfileSpread?"s":"");
|
---|
414 |
|
---|
415 | Init();
|
---|
416 |
|
---|
417 | fData[0] = new MDataPhrase(memberx);
|
---|
418 | fData[1] = new MDataPhrase(membery);
|
---|
419 | fData[2] = new MDataPhrase(memberz);
|
---|
420 | fData[3] = new MDataPhrase(weight);
|
---|
421 | }
|
---|
422 |
|
---|
423 | // --------------------------------------------------------------------------
|
---|
424 | //
|
---|
425 | // Deletes the histogram
|
---|
426 | //
|
---|
427 | MH3::~MH3()
|
---|
428 | {
|
---|
429 | if (fHist)
|
---|
430 | delete fHist;
|
---|
431 |
|
---|
432 | if (fConversion)
|
---|
433 | delete fConversion;
|
---|
434 |
|
---|
435 | if (fWeight)
|
---|
436 | delete fWeight;
|
---|
437 |
|
---|
438 | for (int i=0; i<4; i++)
|
---|
439 | if (fData[i])
|
---|
440 | delete fData[i];
|
---|
441 |
|
---|
442 | for (int i=0; i<3; i++)
|
---|
443 | if (fLabels[i].GetDefault())
|
---|
444 | delete fLabels[i].GetDefault();
|
---|
445 | }
|
---|
446 |
|
---|
447 | // --------------------------------------------------------------------------
|
---|
448 | //
|
---|
449 | // You can set a weight as a phrase additionally to the one given
|
---|
450 | // as an argument to Fill (most likely from MFillH). The two weights
|
---|
451 | // are multiplied together.
|
---|
452 | //
|
---|
453 | void MH3::SetWeight(const char *phrase)
|
---|
454 | {
|
---|
455 | if (fWeight)
|
---|
456 | delete fWeight;
|
---|
457 | fWeight = new MDataPhrase(phrase);
|
---|
458 | }
|
---|
459 |
|
---|
460 | // --------------------------------------------------------------------------
|
---|
461 | //
|
---|
462 | // Set a function which is applied to the histogram before it is displayed.
|
---|
463 | // Note, that it only effects the displayed histogram.
|
---|
464 | //
|
---|
465 | // e.g. SetConversion("sqrt(x)");
|
---|
466 | //
|
---|
467 | Bool_t MH3::SetConversion(const char *func)
|
---|
468 | {
|
---|
469 | if (TString(func).IsNull())
|
---|
470 | {
|
---|
471 | delete fConversion;
|
---|
472 | fConversion = 0;
|
---|
473 | return kTRUE;
|
---|
474 | }
|
---|
475 |
|
---|
476 | fConversion = new TFormula;
|
---|
477 |
|
---|
478 | // Must have a name otherwise all axis labels disappear like a miracle
|
---|
479 | fConversion->SetName("ConversionFunction");
|
---|
480 | if (fConversion->Compile(func))
|
---|
481 | {
|
---|
482 | *fLog << err << dbginf << "Syntax Error: TFormula::Compile failed for " << func << endl;
|
---|
483 | delete fConversion;
|
---|
484 | fConversion = 0;
|
---|
485 | return kFALSE;
|
---|
486 | }
|
---|
487 |
|
---|
488 | gROOT->GetListOfFunctions()->Remove(fConversion);
|
---|
489 |
|
---|
490 | return kTRUE;
|
---|
491 | }
|
---|
492 |
|
---|
493 | // --------------------------------------------------------------------------
|
---|
494 | //
|
---|
495 | // The axis label is centered and the labeling of the axis is initialized.
|
---|
496 | //
|
---|
497 | // This function must not be called after any label has been created!
|
---|
498 | //
|
---|
499 | void MH3::InitLabels(TAxis &x) const
|
---|
500 | {
|
---|
501 | x.CenterTitle();
|
---|
502 | x.SetBinLabel(1, "");
|
---|
503 | x.LabelsOption("h"); // FIXME: Is "a" thread safe? (Paint and Fill?)
|
---|
504 | x.GetLabels()->Delete();
|
---|
505 | }
|
---|
506 |
|
---|
507 | // --------------------------------------------------------------------------
|
---|
508 | //
|
---|
509 | // Depending on the bits set the InitLabels(TAxis&) function for
|
---|
510 | // the corresponding axes are called. In any case the kCanRebin bit
|
---|
511 | // is set.
|
---|
512 | //
|
---|
513 | // This function must not be called after any label has been created!
|
---|
514 | //
|
---|
515 | void MH3::InitLabels(Labels_t type) const
|
---|
516 | {
|
---|
517 | if (!fHist)
|
---|
518 | return;
|
---|
519 |
|
---|
520 | if (type&kLabelsX && fHist->GetXaxis())
|
---|
521 | InitLabels(*fHist->GetXaxis());
|
---|
522 |
|
---|
523 | if (type&kLabelsY && fHist->GetYaxis())
|
---|
524 | InitLabels(*fHist->GetYaxis());
|
---|
525 |
|
---|
526 | if (type&kLabelsZ && fHist->GetZaxis())
|
---|
527 | InitLabels(*fHist->GetZaxis());
|
---|
528 |
|
---|
529 | #if ROOT_VERSION_CODE < ROOT_VERSION(6,00,00)
|
---|
530 | if (type&kLabelsXYZ)
|
---|
531 | fHist->SetBit(TH1::kCanRebin);
|
---|
532 | #else
|
---|
533 | // https://root.cern.ch/content/main-histogram-changes-root-6
|
---|
534 | if (type&kLabelsXYZ)
|
---|
535 | fHist->SetCanExtend(TH1::kAllAxes);
|
---|
536 | #endif
|
---|
537 | }
|
---|
538 |
|
---|
539 | // --------------------------------------------------------------------------
|
---|
540 | //
|
---|
541 | // Return the corresponding Labels_t describing for which axis
|
---|
542 | // axis-labels are switched on.
|
---|
543 | //
|
---|
544 | MH3::Labels_t MH3::GetLabels() const
|
---|
545 | {
|
---|
546 | UInt_t type = kNoLabels;
|
---|
547 | if (fHist->GetXaxis() && fHist->GetXaxis()->GetLabels())
|
---|
548 | type |= kLabelsX;
|
---|
549 | if (fHist->GetYaxis() && fHist->GetYaxis()->GetLabels())
|
---|
550 | type |= kLabelsY;
|
---|
551 | if (fHist->GetZaxis() && fHist->GetZaxis()->GetLabels())
|
---|
552 | type |= kLabelsZ;
|
---|
553 | return (Labels_t)type;
|
---|
554 | }
|
---|
555 |
|
---|
556 | // --------------------------------------------------------------------------
|
---|
557 | //
|
---|
558 | // Calls the LabelsDeflate from the histogram for all axes.
|
---|
559 | // LabelsDeflate will just do nothing if the axis has no labels
|
---|
560 | // initialized.
|
---|
561 | //
|
---|
562 | void MH3::DeflateLabels() const
|
---|
563 | {
|
---|
564 | fHist->LabelsDeflate("X");
|
---|
565 | fHist->LabelsDeflate("Y");
|
---|
566 | fHist->LabelsDeflate("Z");
|
---|
567 | }
|
---|
568 |
|
---|
569 | // --------------------------------------------------------------------------
|
---|
570 | //
|
---|
571 | // Returns the named label corresponding to the given value
|
---|
572 | // and the given axis. The names are defined with the
|
---|
573 | // DefineLabel-functions. if no name is defined the value
|
---|
574 | // is converted to a string with %d and TMath::Nint.
|
---|
575 | // If names are defined, but not for the given value, the default
|
---|
576 | // label is returned instead. If no default is defined the
|
---|
577 | // %d-converted string is returned.
|
---|
578 | //
|
---|
579 | const char *MH3::GetLabel(Int_t axe, Double_t val) const
|
---|
580 | {
|
---|
581 | const Int_t v = TMath::Nint(val);
|
---|
582 |
|
---|
583 | if (fLabels[axe].GetSize())
|
---|
584 | {
|
---|
585 | const char *l = fLabels[axe].GetObjName(v);
|
---|
586 | if (l)
|
---|
587 | return l;
|
---|
588 | }
|
---|
589 |
|
---|
590 | return Form("%d", v);
|
---|
591 | }
|
---|
592 |
|
---|
593 | // --------------------------------------------------------------------------
|
---|
594 | //
|
---|
595 | // Return the data members used by the data chain to be used in
|
---|
596 | // MTask::AddBranchToList
|
---|
597 | //
|
---|
598 | TString MH3::GetDataMember() const
|
---|
599 | {
|
---|
600 | TString str=fData[0]->GetDataMember();
|
---|
601 |
|
---|
602 | for (int i=1; i<4; i++)
|
---|
603 | if (fData[i])
|
---|
604 | {
|
---|
605 | str += ";";
|
---|
606 | str += fData[i]->GetDataMember();
|
---|
607 | }
|
---|
608 |
|
---|
609 | return str;
|
---|
610 | }
|
---|
611 |
|
---|
612 | // --------------------------------------------------------------------------
|
---|
613 | //
|
---|
614 | // Setup the Binning for the histograms automatically if the correct
|
---|
615 | // instances of MBinning are found in the parameter list
|
---|
616 | // For a more detailed description see class description above.
|
---|
617 | //
|
---|
618 | Bool_t MH3::SetupFill(const MParList *plist)
|
---|
619 | {
|
---|
620 | // reset histogram (necessary if the same eventloop is run more than once)
|
---|
621 | if (!TestBit(kDoNotReset))
|
---|
622 | fHist->Reset();
|
---|
623 |
|
---|
624 | // Tokenize name into name and binnings names
|
---|
625 | TObjArray *tok = fName.Tokenize(";");
|
---|
626 |
|
---|
627 | const TString name = (*tok)[0] ? (*tok)[0]->GetName() : fName.Data();
|
---|
628 |
|
---|
629 | TString bx = (*tok)[1] ? (*tok)[1]->GetName() : Form("%sX", name.Data());
|
---|
630 | TString by = (*tok)[2] ? (*tok)[2]->GetName() : Form("%sY", name.Data());
|
---|
631 | TString bz = (*tok)[3] ? (*tok)[3]->GetName() : Form("%sZ", name.Data());
|
---|
632 |
|
---|
633 | bx.Prepend("Binning");
|
---|
634 | by.Prepend("Binning");
|
---|
635 | bz.Prepend("Binning");
|
---|
636 |
|
---|
637 | delete tok;
|
---|
638 |
|
---|
639 | MBinning *binsx = NULL;
|
---|
640 | MBinning *binsy = NULL;
|
---|
641 | MBinning *binsz = NULL;
|
---|
642 |
|
---|
643 | const Labels_t labels = GetLabels();
|
---|
644 |
|
---|
645 | switch (TMath::Abs(fDimension))
|
---|
646 | {
|
---|
647 | case 3:
|
---|
648 | if (fData[2])
|
---|
649 | fHist->SetZTitle(fData[2]->GetTitle());
|
---|
650 | if (!(labels&kLabelsZ))
|
---|
651 | {
|
---|
652 | binsz = fBins[2] ? fBins[2] : (MBinning*)plist->FindObject(bz, "MBinning");
|
---|
653 | if (!binsz)
|
---|
654 | {
|
---|
655 | *fLog << err << dbginf << "MBinning '" << bz << "' not found... aborting." << endl;
|
---|
656 | return kFALSE;
|
---|
657 | }
|
---|
658 | if (binsz->HasTitle())
|
---|
659 | fHist->SetZTitle(binsz->GetTitle());
|
---|
660 | if (binsz->IsLogarithmic())
|
---|
661 | fHist->SetBit(kIsLogz);
|
---|
662 | }
|
---|
663 | case 2:
|
---|
664 | if (fData[1])
|
---|
665 | fHist->SetYTitle(fData[1]->GetTitle());
|
---|
666 | if (!(labels&kLabelsY))
|
---|
667 | {
|
---|
668 | binsy = fBins[1] ? fBins[1] : (MBinning*)plist->FindObject(by, "MBinning");
|
---|
669 | if (!binsy)
|
---|
670 | {
|
---|
671 | *fLog << err << dbginf << "MBinning '" << by << "' not found... aborting." << endl;
|
---|
672 | return kFALSE;
|
---|
673 | }
|
---|
674 | if (binsy->HasTitle())
|
---|
675 | fHist->SetYTitle(binsy->GetTitle());
|
---|
676 | if (binsy->IsLogarithmic())
|
---|
677 | fHist->SetBit(kIsLogy);
|
---|
678 | }
|
---|
679 | case 1:
|
---|
680 | if (fData[0]!=NULL)
|
---|
681 | fHist->SetXTitle(fData[0]->GetTitle());
|
---|
682 | if (!(labels&kLabelsX))
|
---|
683 | {
|
---|
684 | binsx = fBins[0] ? fBins[0] : (MBinning*)plist->FindObject(bx, "MBinning");
|
---|
685 | if (!binsx)
|
---|
686 | {
|
---|
687 | if (fDimension==1)
|
---|
688 | binsx = (MBinning*)plist->FindObject("Binning"+fName, "MBinning");
|
---|
689 |
|
---|
690 | if (!binsx)
|
---|
691 | {
|
---|
692 | *fLog << err << dbginf << "Neither '" << bx << "' nor 'Binning" << fName << "' found... aborting." << endl;
|
---|
693 | return kFALSE;
|
---|
694 | }
|
---|
695 | }
|
---|
696 | if (binsx->HasTitle())
|
---|
697 | fHist->SetXTitle(binsx->GetTitle());
|
---|
698 | if (binsx->IsLogarithmic())
|
---|
699 | fHist->SetBit(kIsLogx);
|
---|
700 | }
|
---|
701 | }
|
---|
702 |
|
---|
703 | // PreProcess existing fData members
|
---|
704 | for (int i=0; i<4; i++)
|
---|
705 | if (fData[i] && !fData[i]->PreProcess(plist))
|
---|
706 | return kFALSE;
|
---|
707 |
|
---|
708 | if (fWeight && !fWeight->PreProcess(plist))
|
---|
709 | return kFALSE;
|
---|
710 |
|
---|
711 | TString title(fDimension>0?"Histogram":"Profile");
|
---|
712 | title += " for ";
|
---|
713 | title += name;
|
---|
714 | title += Form(" (%dD)", TMath::Abs(fDimension));
|
---|
715 |
|
---|
716 | fHist->SetName(name);
|
---|
717 | fHist->SetTitle(fTitle==gsDefTitle ? title : fTitle);
|
---|
718 | fHist->SetDirectory(0);
|
---|
719 |
|
---|
720 | // This is for the case we have set lables
|
---|
721 | const MBinning def(1, 0, 1);
|
---|
722 | if (!binsx)
|
---|
723 | binsx = const_cast<MBinning*>(&def);
|
---|
724 | if (!binsy)
|
---|
725 | binsy = const_cast<MBinning*>(&def);
|
---|
726 | if (!binsz)
|
---|
727 | binsz = const_cast<MBinning*>(&def);
|
---|
728 |
|
---|
729 | // set binning
|
---|
730 | switch (TMath::Abs(fDimension))
|
---|
731 | {
|
---|
732 | case 1:
|
---|
733 | SetBinning(*fHist, *binsx);
|
---|
734 | return kTRUE;
|
---|
735 | case 2:
|
---|
736 | SetBinning(static_cast<TH2&>(*fHist), *binsx, *binsy);
|
---|
737 | return kTRUE;
|
---|
738 | case 3:
|
---|
739 | SetBinning(static_cast<TH3&>(*fHist), *binsx, *binsy, *binsz);
|
---|
740 | return kTRUE;
|
---|
741 | }
|
---|
742 |
|
---|
743 | *fLog << err << "ERROR - MH3 has " << TMath::Abs(fDimension) << " dimensions!" << endl;
|
---|
744 | return kFALSE;
|
---|
745 | }
|
---|
746 |
|
---|
747 | // --------------------------------------------------------------------------
|
---|
748 | //
|
---|
749 | // Set the name of the histogram ant the MH3 container
|
---|
750 | //
|
---|
751 | void MH3::SetName(const char *name)
|
---|
752 | {
|
---|
753 | if (fHist)
|
---|
754 | {
|
---|
755 | if (gPad)
|
---|
756 | {
|
---|
757 | const TString pfx(MString::Format("%sProfX", fHist->GetName()));
|
---|
758 | const TString pfy(MString::Format("%sProfY", fHist->GetName()));
|
---|
759 |
|
---|
760 | TProfile *p = 0;
|
---|
761 | if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
|
---|
762 | p->SetName(MString::Format("%sProfX", name));
|
---|
763 | if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
|
---|
764 | p->SetName(MString::Format("%sProfY", name));
|
---|
765 | }
|
---|
766 |
|
---|
767 | fHist->SetName(name);
|
---|
768 | fHist->SetDirectory(0);
|
---|
769 |
|
---|
770 | }
|
---|
771 | MParContainer::SetName(name);
|
---|
772 | }
|
---|
773 |
|
---|
774 | // --------------------------------------------------------------------------
|
---|
775 | //
|
---|
776 | // Set the title of the histogram ant the MH3 container
|
---|
777 | //
|
---|
778 | void MH3::SetTitle(const char *title)
|
---|
779 | {
|
---|
780 | if (fHist)
|
---|
781 | fHist->SetTitle(title);
|
---|
782 | MParContainer::SetTitle(title);
|
---|
783 | }
|
---|
784 |
|
---|
785 | // --------------------------------------------------------------------------
|
---|
786 | //
|
---|
787 | // Fills the one, two or three data members into our histogram
|
---|
788 | //
|
---|
789 | Int_t MH3::Fill(const MParContainer *par, const Stat_t ww)
|
---|
790 | {
|
---|
791 | // Get Information about labels (UInt_t, to supress warning about
|
---|
792 | // unhandeled cases in switch)
|
---|
793 | const UInt_t type = GetLabels();
|
---|
794 |
|
---|
795 | // Get values for axis
|
---|
796 | Double_t x=0;
|
---|
797 | Double_t y=0;
|
---|
798 | Double_t z=0;
|
---|
799 | Double_t t=0;
|
---|
800 | Double_t w=ww;
|
---|
801 |
|
---|
802 | switch (fDimension)
|
---|
803 | {
|
---|
804 | case -3:
|
---|
805 | t = fData[3]->GetValue()*fScale[3];
|
---|
806 | case -2:
|
---|
807 | case 3:
|
---|
808 | z = fData[2]->GetValue()*fScale[2];
|
---|
809 | case -1:
|
---|
810 | case 2:
|
---|
811 | y = fData[1]->GetValue()*fScale[1];
|
---|
812 | case 1:
|
---|
813 | x = fData[0]->GetValue()*fScale[0];
|
---|
814 | }
|
---|
815 |
|
---|
816 | if (fWeight)
|
---|
817 | w *= fWeight->GetValue();
|
---|
818 |
|
---|
819 | // If label option is set, convert value to label
|
---|
820 | TString labelx, labely, labelz;
|
---|
821 | if (type&kLabelsX)
|
---|
822 | labelx = GetLabel(0, x);
|
---|
823 | if (type&kLabelsY)
|
---|
824 | labely = GetLabel(1, y);
|
---|
825 | if (type&kLabelsZ)
|
---|
826 | labelz = GetLabel(2, z);
|
---|
827 |
|
---|
828 | // Fill histogram
|
---|
829 | switch (fDimension)
|
---|
830 | {
|
---|
831 | case 3:
|
---|
832 | switch (type)
|
---|
833 | {
|
---|
834 | case kNoLabels:
|
---|
835 | static_cast<TH3*>(fHist)->Fill(x, y, z, w);
|
---|
836 | return kTRUE;
|
---|
837 | case kLabelsX:
|
---|
838 | #if ROOT_VERSION_CODE < ROOT_VERSION(5,32,01)
|
---|
839 | static_cast<TH3*>(fHist)->Fill(labelx, y, z);
|
---|
840 | #endif
|
---|
841 | return kTRUE;
|
---|
842 | case kLabelsY:
|
---|
843 | static_cast<TH3*>(fHist)->Fill(x, labely, z, w);
|
---|
844 | return kTRUE;
|
---|
845 | case kLabelsZ:
|
---|
846 | static_cast<TH3*>(fHist)->Fill(x, y, labelz, w);
|
---|
847 | return kTRUE;
|
---|
848 | case kLabelsXY:
|
---|
849 | static_cast<TH3*>(fHist)->Fill(labelx, labely, z, w);
|
---|
850 | return kTRUE;
|
---|
851 | case kLabelsXZ:
|
---|
852 | static_cast<TH3*>(fHist)->Fill(labelx, y, labelz, w);
|
---|
853 | return kTRUE;
|
---|
854 | case kLabelsYZ:
|
---|
855 | static_cast<TH3*>(fHist)->Fill(x, labely, labelz, w);
|
---|
856 | return kTRUE;
|
---|
857 | case kLabelsXYZ:
|
---|
858 | static_cast<TH3*>(fHist)->Fill(labelx, labely, labelz, w);
|
---|
859 | return kTRUE;
|
---|
860 | }
|
---|
861 | break;
|
---|
862 | case 2:
|
---|
863 | switch (type)
|
---|
864 | {
|
---|
865 | case kNoLabels:
|
---|
866 | static_cast<TH2*>(fHist)->Fill(x, y, w);
|
---|
867 | return kTRUE;
|
---|
868 | case kLabelsX:
|
---|
869 | static_cast<TH2*>(fHist)->Fill(x, labely, w);
|
---|
870 | return kTRUE;
|
---|
871 | case kLabelsY:
|
---|
872 | static_cast<TH2*>(fHist)->Fill(labelx, y, w);
|
---|
873 | return kTRUE;
|
---|
874 | case kLabelsXY:
|
---|
875 | static_cast<TH2*>(fHist)->Fill(labelx, labely, w);
|
---|
876 | return kTRUE;
|
---|
877 | }
|
---|
878 | break;
|
---|
879 | case 1:
|
---|
880 | switch (type)
|
---|
881 | {
|
---|
882 | case kNoLabels:
|
---|
883 | fHist->Fill(x, w);
|
---|
884 | return kTRUE;
|
---|
885 | case kLabelsX:
|
---|
886 | fHist->Fill(labelx, w);
|
---|
887 | return kTRUE;
|
---|
888 | }
|
---|
889 | break;
|
---|
890 | case -1:
|
---|
891 | switch (type)
|
---|
892 | {
|
---|
893 | case kNoLabels:
|
---|
894 | static_cast<TProfile*>(fHist)->Fill(x, y, w);
|
---|
895 | return kTRUE;
|
---|
896 | case kLabelsX:
|
---|
897 | static_cast<TProfile*>(fHist)->Fill(labelx, y, w);
|
---|
898 | return kTRUE;
|
---|
899 | }
|
---|
900 | break;
|
---|
901 | case -2:
|
---|
902 | switch (type)
|
---|
903 | {
|
---|
904 | case kNoLabels:
|
---|
905 | static_cast<TProfile2D*>(fHist)->Fill(x, y, z, w);
|
---|
906 | return kTRUE;
|
---|
907 | case kLabelsX:
|
---|
908 | static_cast<TProfile2D*>(fHist)->Fill(labelx, y, z);
|
---|
909 | return kTRUE;
|
---|
910 | case kLabelsY:
|
---|
911 | static_cast<TProfile2D*>(fHist)->Fill(x, labely, z);
|
---|
912 | return kTRUE;
|
---|
913 | case kLabelsXY:
|
---|
914 | static_cast<TProfile2D*>(fHist)->Fill(labelx, labely, z);
|
---|
915 | return kTRUE;
|
---|
916 | }
|
---|
917 | break;
|
---|
918 | case -3:
|
---|
919 | switch (type)
|
---|
920 | {
|
---|
921 | case kNoLabels:
|
---|
922 | static_cast<TProfile3D*>(fHist)->Fill(x, y, z, t, w);
|
---|
923 | return kTRUE;
|
---|
924 | default:
|
---|
925 | *fLog << err << "ERROR - Labels not supported in TProfile3D." << endl;
|
---|
926 | }
|
---|
927 | break;
|
---|
928 | }
|
---|
929 |
|
---|
930 | *fLog << err << "MH3::Fill: ERROR - A fatal error occured." << endl;
|
---|
931 | return kERROR;
|
---|
932 | }
|
---|
933 |
|
---|
934 | // --------------------------------------------------------------------------
|
---|
935 | //
|
---|
936 | // If an auto range bit is set the histogram range of the corresponding
|
---|
937 | // axis is set to show only the non-empty bins (with a single empty bin
|
---|
938 | // on both sides)
|
---|
939 | //
|
---|
940 | Bool_t MH3::Finalize()
|
---|
941 | {
|
---|
942 | DeflateLabels();
|
---|
943 |
|
---|
944 | Bool_t autorangex=TESTBIT(fStyleBits, 0);
|
---|
945 | Bool_t autorangey=TESTBIT(fStyleBits, 1);
|
---|
946 | //Bool_t autorangez=TESTBIT(fStyleBits, 2);
|
---|
947 |
|
---|
948 | Int_t lo, hi;
|
---|
949 | if (autorangex)
|
---|
950 | {
|
---|
951 | GetRangeX(*fHist, lo, hi);
|
---|
952 | fHist->GetXaxis()->SetRange(lo-2, hi+1);
|
---|
953 | }
|
---|
954 | if (autorangey)
|
---|
955 | {
|
---|
956 | GetRangeY(*fHist, lo, hi);
|
---|
957 | fHist->GetYaxis()->SetRange(lo-2, hi+1);
|
---|
958 | }
|
---|
959 | /*
|
---|
960 | if (autorangez)
|
---|
961 | {
|
---|
962 | GetRangeZ(*fHist, lo, hi);
|
---|
963 | fHist->GetZaxis()->SetRange(lo-2, hi+1);
|
---|
964 | }
|
---|
965 | */
|
---|
966 | return kTRUE;
|
---|
967 | }
|
---|
968 |
|
---|
969 | // --------------------------------------------------------------------------
|
---|
970 | //
|
---|
971 | // Apply the conversion function to the contents stored in fHist and
|
---|
972 | // store the result in h.
|
---|
973 | //
|
---|
974 | // In the case of a TProfile we keep it a TProfile to keep the mean and
|
---|
975 | // rms in y displayed. To get the error correctly we have to reverse
|
---|
976 | // the calculation done in TProfile::GetBinError of course.
|
---|
977 | //
|
---|
978 | void MH3::Convert(TH1 &h) const
|
---|
979 | {
|
---|
980 | const Bool_t prof = h.InheritsFrom(TProfile::Class()) || h.InheritsFrom(TProfile2D::Class()) || h.InheritsFrom(TProfile3D::Class());
|
---|
981 |
|
---|
982 | for (Int_t z=0; z<=h.GetNbinsZ()+1; z++)
|
---|
983 | for (Int_t y=0; y<=h.GetNbinsY()+1; y++)
|
---|
984 | for (Int_t x=0; x<=h.GetNbinsX()+1; x++)
|
---|
985 | {
|
---|
986 | h.SetBinContent(x, y, z, fConversion->Eval(fHist->GetBinContent(x, y, z)));
|
---|
987 |
|
---|
988 | if (prof)
|
---|
989 | h.SetBinError(x, y, z, TMath::Hypot(fConversion->Eval(fHist->GetBinContent(x, y, z)),
|
---|
990 | fConversion->Eval(fHist->GetBinError( x, y, z))));
|
---|
991 | else
|
---|
992 | h.SetBinError(x, y, z, fConversion->Eval(fHist->GetBinError(x, y, z)));
|
---|
993 | }
|
---|
994 |
|
---|
995 | TProfile *p1 = dynamic_cast<TProfile*>(fHist);
|
---|
996 | if (p1)
|
---|
997 | for (Int_t i=0; i<p1->GetSize(); i++)
|
---|
998 | static_cast<TProfile&>(h).SetBinEntries(i, p1->GetBinEntries(i)>0 ? 1 : 0);
|
---|
999 |
|
---|
1000 | TProfile2D *p2 = dynamic_cast<TProfile2D*>(fHist);
|
---|
1001 | if (p2)
|
---|
1002 | for (Int_t i=0; i<p2->GetSize(); i++)
|
---|
1003 | static_cast<TProfile2D&>(h).SetBinEntries(i, p2->GetBinEntries(i)>0 ? 1 : 0);
|
---|
1004 |
|
---|
1005 | TProfile3D *p3 = dynamic_cast<TProfile3D*>(fHist);
|
---|
1006 | if (p3)
|
---|
1007 | for (Int_t i=0; i<p3->GetSize(); i++)
|
---|
1008 | static_cast<TProfile3D&>(h).SetBinEntries(i, p3->GetBinEntries(i)>0 ? 1 : 0);
|
---|
1009 | }
|
---|
1010 |
|
---|
1011 | // --------------------------------------------------------------------------
|
---|
1012 | //
|
---|
1013 | // FIXME
|
---|
1014 | //
|
---|
1015 | void MH3::Paint(Option_t *o)
|
---|
1016 | {
|
---|
1017 | TProfile *p=0;
|
---|
1018 |
|
---|
1019 | if (TMath::Abs(fDimension)==2)
|
---|
1020 | MH::SetPalette("pretty");
|
---|
1021 |
|
---|
1022 | if (fConversion)
|
---|
1023 | {
|
---|
1024 | TH1 *h = 0;
|
---|
1025 | if ((h=dynamic_cast<TH1*>(gPad->FindObject(fHist->GetName()))))
|
---|
1026 | Convert(*h);
|
---|
1027 | }
|
---|
1028 |
|
---|
1029 | const TString pfx(MString::Format("%sProfX", fHist->GetName()));
|
---|
1030 | if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
|
---|
1031 | {
|
---|
1032 | Int_t col = p->GetLineColor();
|
---|
1033 | p = ((TH2*)fHist)->ProfileX(pfx, -1, -1, "s");
|
---|
1034 | p->SetLineColor(col);
|
---|
1035 | }
|
---|
1036 |
|
---|
1037 | const TString pfy(MString::Format("%sProfY", fHist->GetName()));
|
---|
1038 | if ((p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
|
---|
1039 | {
|
---|
1040 | Int_t col = p->GetLineColor();
|
---|
1041 | p = ((TH2*)fHist)->ProfileY(pfy, -1, -1, "s");
|
---|
1042 | p->SetLineColor(col);
|
---|
1043 | }
|
---|
1044 | /*
|
---|
1045 | if (fHist->TestBit(kIsLogx) && fHist->GetEntries()>0) gPad->SetLogx();
|
---|
1046 | if (fHist->TestBit(kIsLogy) && fHist->GetEntries()>0) gPad->SetLogy();
|
---|
1047 | if (fHist->TestBit(kIsLogz) && fHist->GetEntries()>0) gPad->SetLogz();
|
---|
1048 | */
|
---|
1049 | }
|
---|
1050 |
|
---|
1051 | // --------------------------------------------------------------------------
|
---|
1052 | //
|
---|
1053 | // If Xmax is < 3000*Xmin SetMoreLogLabels is called. If Xmax<5000
|
---|
1054 | // the exponent is switched off (SetNoExponent)
|
---|
1055 | //
|
---|
1056 | void MH3::HandleLogAxis(TAxis &axe) const
|
---|
1057 | {
|
---|
1058 | if (axe.GetXmax()>3000*axe.GetXmin())
|
---|
1059 | return;
|
---|
1060 |
|
---|
1061 | axe.SetMoreLogLabels();
|
---|
1062 | if (axe.GetXmax()<5000)
|
---|
1063 | axe.SetNoExponent();
|
---|
1064 | }
|
---|
1065 |
|
---|
1066 | // --------------------------------------------------------------------------
|
---|
1067 | //
|
---|
1068 | // Creates a new canvas and draws the histogram into it.
|
---|
1069 | //
|
---|
1070 | // Possible options are:
|
---|
1071 | // PROFX: Draw a x-profile into the histogram (for 2D histograms only)
|
---|
1072 | // PROFY: Draw a y-profile into the histogram (for 2D histograms only)
|
---|
1073 | // ONLY: Draw the profile histogram only (for 2D histograms only)
|
---|
1074 | // BLUE: Draw the profile in blue color instead of the histograms
|
---|
1075 | // line color
|
---|
1076 | //
|
---|
1077 | // If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
|
---|
1078 | // eg this is set when applying a logarithmic MBinning
|
---|
1079 | //
|
---|
1080 | // Be careful: The histogram belongs to this object and won't get deleted
|
---|
1081 | // together with the canvas.
|
---|
1082 | //
|
---|
1083 | void MH3::Draw(Option_t *opt)
|
---|
1084 | {
|
---|
1085 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist);
|
---|
1086 | pad->SetBorderMode(0);
|
---|
1087 | pad->SetGridx();
|
---|
1088 | pad->SetGridy();
|
---|
1089 |
|
---|
1090 | AppendPad();
|
---|
1091 |
|
---|
1092 | if (fHist->TestBit(kIsLogx))
|
---|
1093 | {
|
---|
1094 | pad->SetLogx();
|
---|
1095 | HandleLogAxis(*fHist->GetXaxis());
|
---|
1096 | }
|
---|
1097 | if (fHist->TestBit(kIsLogy))
|
---|
1098 | {
|
---|
1099 | pad->SetLogy();
|
---|
1100 | HandleLogAxis(*fHist->GetYaxis());
|
---|
1101 | }
|
---|
1102 | if (fHist->TestBit(kIsLogz))
|
---|
1103 | {
|
---|
1104 | pad->SetLogz();
|
---|
1105 | HandleLogAxis(*fHist->GetZaxis());
|
---|
1106 | }
|
---|
1107 |
|
---|
1108 | fHist->SetFillStyle(4000);
|
---|
1109 |
|
---|
1110 | TString str(opt);
|
---|
1111 | str.ToLower();
|
---|
1112 | str.ReplaceAll(" ", "");
|
---|
1113 |
|
---|
1114 | if (GetLabels())
|
---|
1115 | {
|
---|
1116 | if (str.IsNull() && fDimension==2)
|
---|
1117 | str = "colz";
|
---|
1118 |
|
---|
1119 | if (str.Contains("box", TString::kIgnoreCase) && fDimension==2)
|
---|
1120 | fHist->SetLineColor(kBlue);
|
---|
1121 | }
|
---|
1122 |
|
---|
1123 | const Bool_t only = str.Contains("only") && TMath::Abs(fDimension)==2;
|
---|
1124 | const Bool_t same = str.Contains("same") && TMath::Abs(fDimension)<3;
|
---|
1125 | const Bool_t blue = str.Contains("blue") && TMath::Abs(fDimension)==2;
|
---|
1126 | const Bool_t profx = str.Contains("profx") && TMath::Abs(fDimension)==2;
|
---|
1127 | const Bool_t profy = str.Contains("profy") && TMath::Abs(fDimension)==2;
|
---|
1128 |
|
---|
1129 | str.ReplaceAll("only", "");
|
---|
1130 | str.ReplaceAll("blue", "");
|
---|
1131 | str.ReplaceAll("profx", "");
|
---|
1132 | str.ReplaceAll("profy", "");
|
---|
1133 |
|
---|
1134 | if (same && TMath::Abs(fDimension)==1)
|
---|
1135 | {
|
---|
1136 | fHist->SetLineColor(kBlue);
|
---|
1137 | fHist->SetMarkerColor(kBlue);
|
---|
1138 | }
|
---|
1139 |
|
---|
1140 |
|
---|
1141 | TH1 *h = fHist;
|
---|
1142 |
|
---|
1143 | if (fConversion)
|
---|
1144 | {
|
---|
1145 | h = static_cast<TH1*>(fHist->Clone());
|
---|
1146 |
|
---|
1147 | h->SetDirectory(0);
|
---|
1148 | h->SetBit(kCanDelete);
|
---|
1149 |
|
---|
1150 | Convert(*h);
|
---|
1151 | }
|
---|
1152 |
|
---|
1153 | // FIXME: We may have to remove all our own options from str!
|
---|
1154 | if (!only)
|
---|
1155 | h->Draw(str);
|
---|
1156 |
|
---|
1157 | TProfile *p=0;
|
---|
1158 | if (profx)
|
---|
1159 | {
|
---|
1160 | const TString pfx(MString::Format("%sProfX", h->GetName()));
|
---|
1161 |
|
---|
1162 | if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfx))))
|
---|
1163 | *fLog << warn << "TProfile " << pfx << " already in pad." << endl;
|
---|
1164 |
|
---|
1165 | p = ((TH2*)h)->ProfileX(pfx, -1, -1, "s");
|
---|
1166 | p->UseCurrentStyle();
|
---|
1167 | p->SetLineColor(blue ? kBlue : h->GetLineColor());
|
---|
1168 | p->SetBit(kCanDelete);
|
---|
1169 | p->SetDirectory(NULL);
|
---|
1170 | p->SetXTitle(h->GetXaxis()->GetTitle());
|
---|
1171 | p->SetYTitle(h->GetYaxis()->GetTitle());
|
---|
1172 | p->Draw(only&&!same?"":"same");
|
---|
1173 | }
|
---|
1174 | if (profy)
|
---|
1175 | {
|
---|
1176 | const TString pfy(MString::Format("%sProfY", h->GetName()));
|
---|
1177 |
|
---|
1178 | if (same && (p=dynamic_cast<TProfile*>(gPad->FindObject(pfy))))
|
---|
1179 | *fLog << warn << "TProfile " << pfy << " already in pad." << endl;
|
---|
1180 |
|
---|
1181 | p = ((TH2*)h)->ProfileY(pfy, -1, -1, "s");
|
---|
1182 | p->UseCurrentStyle();
|
---|
1183 | p->SetLineColor(blue ? kBlue : h->GetLineColor());
|
---|
1184 | p->SetBit(kCanDelete);
|
---|
1185 | p->SetDirectory(NULL);
|
---|
1186 | p->SetYTitle(h->GetXaxis()->GetTitle());
|
---|
1187 | p->SetXTitle(h->GetYaxis()->GetTitle());
|
---|
1188 | p->Draw(only&&!same?"":"same");
|
---|
1189 | }
|
---|
1190 |
|
---|
1191 | //AppendPad("log");
|
---|
1192 | }
|
---|
1193 |
|
---|
1194 | // --------------------------------------------------------------------------
|
---|
1195 | //
|
---|
1196 | // Implementation of SavePrimitive. Used to write the call to a constructor
|
---|
1197 | // to a macro. In the original root implementation it is used to write
|
---|
1198 | // gui elements to a macro-file.
|
---|
1199 | //
|
---|
1200 | void MH3::StreamPrimitive(ostream &out) const
|
---|
1201 | {
|
---|
1202 | TString name = GetUniqueName();
|
---|
1203 |
|
---|
1204 | out << " MH3 " << name << "(\"";
|
---|
1205 | out << fData[0]->GetRule() << "\"";
|
---|
1206 | if (fDimension>1 || fDimension<0)
|
---|
1207 | out << ", \"" << fData[1]->GetRule() << "\"";
|
---|
1208 | if (fDimension>2 || fDimension<-1)
|
---|
1209 | out << ", \"" << fData[2]->GetRule() << "\"";
|
---|
1210 |
|
---|
1211 | out << ");" << endl;
|
---|
1212 |
|
---|
1213 | if (fName!=gsDefName)
|
---|
1214 | out << " " << name << ".SetName(\"" << fName << "\");" << endl;
|
---|
1215 |
|
---|
1216 | if (fTitle!=gsDefTitle)
|
---|
1217 | out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl;
|
---|
1218 |
|
---|
1219 | if (fWeight)
|
---|
1220 | out << " " << name << ".SetWeight(\"" << fWeight->GetRule() << "\");" << endl;
|
---|
1221 |
|
---|
1222 | switch (fDimension)
|
---|
1223 | {
|
---|
1224 | case -3:
|
---|
1225 | if (fScale[3]!=1)
|
---|
1226 | out << " " << name << ".SetScaleT(" << fScale[3] << ");" << endl;
|
---|
1227 | case -2:
|
---|
1228 | case 3:
|
---|
1229 | if (fScale[2]!=1)
|
---|
1230 | out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl;
|
---|
1231 | case -1:
|
---|
1232 | case 2:
|
---|
1233 | if (fScale[1]!=1)
|
---|
1234 | out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl;
|
---|
1235 | case 1:
|
---|
1236 | if (fScale[0]!=1)
|
---|
1237 | out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl;
|
---|
1238 | }
|
---|
1239 | }
|
---|
1240 |
|
---|
1241 | MH3::Type_t MH3::GetType() const
|
---|
1242 | {
|
---|
1243 | switch (fDimension)
|
---|
1244 | {
|
---|
1245 | case -1:
|
---|
1246 | return TString(static_cast<TProfile*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
|
---|
1247 | case -2:
|
---|
1248 | return TString(static_cast<TProfile2D*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
|
---|
1249 | case -3:
|
---|
1250 | return TString(static_cast<TProfile3D*>(fHist)->GetErrorOption())=="s" ? kProfileSpread : kProfile;
|
---|
1251 | }
|
---|
1252 | return kHistogram;
|
---|
1253 | }
|
---|
1254 |
|
---|
1255 | // --------------------------------------------------------------------------
|
---|
1256 | //
|
---|
1257 | // Used to rebuild a MH3 object of the same type (data members,
|
---|
1258 | // dimension, ...)
|
---|
1259 | //
|
---|
1260 | MParContainer *MH3::New() const
|
---|
1261 | {
|
---|
1262 | // FIXME: TREAT THE NEW OPTIONS CORRECTLY (PROFILE, LABELS)
|
---|
1263 |
|
---|
1264 | MH3 *h = NULL;
|
---|
1265 |
|
---|
1266 | if (fData[0] == NULL)
|
---|
1267 | h=new MH3(fDimension);
|
---|
1268 | else
|
---|
1269 | switch (fDimension)
|
---|
1270 | {
|
---|
1271 | case 1:
|
---|
1272 | h=new MH3(fData[0]->GetRule());
|
---|
1273 | break;
|
---|
1274 | case 2:
|
---|
1275 | case -1:
|
---|
1276 | h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), GetType());
|
---|
1277 | break;
|
---|
1278 | case 3:
|
---|
1279 | case -2:
|
---|
1280 | h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule(), GetType());
|
---|
1281 | break;
|
---|
1282 | case -3:
|
---|
1283 | h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule(), fData[3]->GetRule());
|
---|
1284 | break;
|
---|
1285 | }
|
---|
1286 |
|
---|
1287 | switch (fDimension)
|
---|
1288 | {
|
---|
1289 | case -3:
|
---|
1290 | h->SetScaleZ(fScale[3]);
|
---|
1291 | case -2:
|
---|
1292 | case 3:
|
---|
1293 | h->SetScaleZ(fScale[2]);
|
---|
1294 | case 2:
|
---|
1295 | case -1:
|
---|
1296 | h->SetScaleY(fScale[1]);
|
---|
1297 | case 1:
|
---|
1298 | h->SetScaleX(fScale[0]);
|
---|
1299 | }
|
---|
1300 |
|
---|
1301 | if (fWeight)
|
---|
1302 | h->SetWeight(fWeight->GetRule());
|
---|
1303 |
|
---|
1304 | return h;
|
---|
1305 | }
|
---|
1306 |
|
---|
1307 | // --------------------------------------------------------------------------
|
---|
1308 | //
|
---|
1309 | // FIXME
|
---|
1310 | //
|
---|
1311 | TString MH3::GetRule(const Char_t axis) const
|
---|
1312 | {
|
---|
1313 | switch (tolower(axis))
|
---|
1314 | {
|
---|
1315 | case 'x':
|
---|
1316 | case 'X':
|
---|
1317 | return fData[0] ? fData[0]->GetRule() : TString("");
|
---|
1318 | case 'y':
|
---|
1319 | case 'Y':
|
---|
1320 | return fData[1] ? fData[1]->GetRule() : TString("");
|
---|
1321 | case 'z':
|
---|
1322 | case 'Z':
|
---|
1323 | return fData[2] ? fData[2]->GetRule() : TString("");
|
---|
1324 | case 't':
|
---|
1325 | case 'T':
|
---|
1326 | return fData[3] ? fData[3]->GetRule() : TString("");
|
---|
1327 | case 'w':
|
---|
1328 | case 'W':
|
---|
1329 | return fWeight ? fWeight->GetRule() : TString("");
|
---|
1330 | default:
|
---|
1331 | return "<n/a>";
|
---|
1332 | }
|
---|
1333 | }
|
---|
1334 |
|
---|
1335 | // --------------------------------------------------------------------------
|
---|
1336 | //
|
---|
1337 | // Returns the total number of bins in a histogram (excluding under- and
|
---|
1338 | // overflow bins)
|
---|
1339 | //
|
---|
1340 | Int_t MH3::GetNbins() const
|
---|
1341 | {
|
---|
1342 | Int_t num = 1;
|
---|
1343 |
|
---|
1344 | switch (TMath::Abs(fDimension))
|
---|
1345 | {
|
---|
1346 | case 3:
|
---|
1347 | num *= fHist->GetNbinsZ()+2;
|
---|
1348 | case 2:
|
---|
1349 | num *= fHist->GetNbinsY()+2;
|
---|
1350 | case 1:
|
---|
1351 | num *= fHist->GetNbinsX()+2;
|
---|
1352 | }
|
---|
1353 |
|
---|
1354 | return num;
|
---|
1355 | }
|
---|
1356 |
|
---|
1357 | // --------------------------------------------------------------------------
|
---|
1358 | //
|
---|
1359 | // Returns the total number of bins in a histogram (excluding under- and
|
---|
1360 | // overflow bins) Return -1 if bin is underflow or overflow
|
---|
1361 | //
|
---|
1362 | Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const
|
---|
1363 | {
|
---|
1364 | const TAxis &axex = *fHist->GetXaxis();
|
---|
1365 | const TAxis &axey = *fHist->GetYaxis();
|
---|
1366 | const TAxis &axez = *fHist->GetZaxis();
|
---|
1367 |
|
---|
1368 | Int_t binz = 0;
|
---|
1369 | Int_t biny = 0;
|
---|
1370 | Int_t binx = 0;
|
---|
1371 |
|
---|
1372 | switch (fDimension)
|
---|
1373 | {
|
---|
1374 | case 3:
|
---|
1375 | case -2:
|
---|
1376 | binz = axez.FindFixBin(z);
|
---|
1377 | if (binz>axez.GetLast() || binz<axez.GetFirst())
|
---|
1378 | return -1;
|
---|
1379 | case 2:
|
---|
1380 | case -1:
|
---|
1381 | biny = axey.FindFixBin(y);
|
---|
1382 | if (biny>axey.GetLast() || biny<axey.GetFirst())
|
---|
1383 | return -1;
|
---|
1384 | case 1:
|
---|
1385 | binx = axex.FindFixBin(x);
|
---|
1386 | if (binx<axex.GetFirst() || binx>axex.GetLast())
|
---|
1387 | return -1;
|
---|
1388 | }
|
---|
1389 |
|
---|
1390 | const Int_t nx = fHist->GetNbinsX()+2;
|
---|
1391 | const Int_t ny = fHist->GetNbinsY()+2;
|
---|
1392 |
|
---|
1393 | return binx + nx*(biny +ny*binz);
|
---|
1394 | }
|
---|
1395 |
|
---|
1396 | // --------------------------------------------------------------------------
|
---|
1397 | //
|
---|
1398 | // Return the MObjLookup corresponding to the axis/character.
|
---|
1399 | // Note that only lower-case charecters (x, y, z) are supported.
|
---|
1400 | // If for the axis no labels were set, the corresponding
|
---|
1401 | // InitLabels is called.
|
---|
1402 | //
|
---|
1403 | MObjLookup *MH3::GetLabels(char axe)
|
---|
1404 | {
|
---|
1405 | if (!fHist)
|
---|
1406 | return 0;
|
---|
1407 |
|
---|
1408 | TAxis *x = 0;
|
---|
1409 |
|
---|
1410 | switch (axe)
|
---|
1411 | {
|
---|
1412 | case 'x':
|
---|
1413 | x = fHist->GetXaxis();
|
---|
1414 | break;
|
---|
1415 | case 'y':
|
---|
1416 | x = fHist->GetYaxis();
|
---|
1417 | break;
|
---|
1418 | case 'z':
|
---|
1419 | x = fHist->GetZaxis();
|
---|
1420 | break;
|
---|
1421 | }
|
---|
1422 |
|
---|
1423 | if (!x)
|
---|
1424 | return 0;
|
---|
1425 |
|
---|
1426 | const Int_t idx = axe-'x';
|
---|
1427 |
|
---|
1428 | if (!x->GetLabels())
|
---|
1429 | switch (idx)
|
---|
1430 | {
|
---|
1431 | case 0:
|
---|
1432 | InitLabels(kLabelsX);
|
---|
1433 | break;
|
---|
1434 | case 1:
|
---|
1435 | InitLabels(kLabelsY);
|
---|
1436 | break;
|
---|
1437 | case 2:
|
---|
1438 | InitLabels(kLabelsZ);
|
---|
1439 | break;
|
---|
1440 | }
|
---|
1441 |
|
---|
1442 | return &fLabels[idx];
|
---|
1443 | }
|
---|
1444 |
|
---|
1445 | // --------------------------------------------------------------------------
|
---|
1446 | //
|
---|
1447 | // Set a default label which is used if no other is found in the list
|
---|
1448 | // of labels. if a default was set already it is overwritten. If the
|
---|
1449 | // axis has not yet been initialized to use labels it it now.
|
---|
1450 | //
|
---|
1451 | void MH3::DefaultLabel(char axe, const char *name)
|
---|
1452 | {
|
---|
1453 | MObjLookup *arr = GetLabels(axe);
|
---|
1454 | if (!arr)
|
---|
1455 | return;
|
---|
1456 |
|
---|
1457 | if (arr->GetDefault())
|
---|
1458 | {
|
---|
1459 | delete arr->GetDefault();
|
---|
1460 | arr->SetDefault(0);
|
---|
1461 | }
|
---|
1462 |
|
---|
1463 | if (name)
|
---|
1464 | arr->SetDefault(new TObjString(name));
|
---|
1465 | }
|
---|
1466 |
|
---|
1467 | // --------------------------------------------------------------------------
|
---|
1468 | //
|
---|
1469 | // Define a name for a label. More than one label can have the same
|
---|
1470 | // name. If the axis has not yet been initialized to use labels
|
---|
1471 | // it it now.
|
---|
1472 | //
|
---|
1473 | void MH3::DefineLabel(char axe, Int_t label, const char *name)
|
---|
1474 | {
|
---|
1475 | MObjLookup *arr = GetLabels(axe);
|
---|
1476 |
|
---|
1477 | if (!arr || !name)
|
---|
1478 | return;
|
---|
1479 |
|
---|
1480 | if (arr->GetObj(label)!=arr->GetDefault())
|
---|
1481 | return;
|
---|
1482 |
|
---|
1483 | arr->Add(label, name);
|
---|
1484 | }
|
---|
1485 |
|
---|
1486 | // --------------------------------------------------------------------------
|
---|
1487 | //
|
---|
1488 | // Define names for labels, like
|
---|
1489 | // 1=Trig;2=Cal;4=Ped;8=Lvl2
|
---|
1490 | // More than one label can have the same name. If the axis has not
|
---|
1491 | // yet been initialized to use labels it it now.
|
---|
1492 | //
|
---|
1493 | // A default cannot be set here. Use DefaultLabel instead.
|
---|
1494 | //
|
---|
1495 | void MH3::DefineLabels(char axe, const TString &labels)
|
---|
1496 | {
|
---|
1497 | TObjArray *arr = labels.Tokenize(';');
|
---|
1498 |
|
---|
1499 | for (int i=0; i<arr->GetEntries(); i++)
|
---|
1500 | {
|
---|
1501 | const char *s = (*arr)[0]->GetName();
|
---|
1502 | const char *v = strchr(s, '=');
|
---|
1503 |
|
---|
1504 | if (v)
|
---|
1505 | DefineLabel(axe, atoi(s), v+1);
|
---|
1506 | }
|
---|
1507 |
|
---|
1508 | delete arr;
|
---|
1509 | }
|
---|
1510 |
|
---|
1511 | void MH3::RecursiveRemove(TObject *obj)
|
---|
1512 | {
|
---|
1513 | if (obj==fHist)
|
---|
1514 | fHist = 0;
|
---|
1515 |
|
---|
1516 | MH::RecursiveRemove(obj);
|
---|
1517 | }
|
---|