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