source: tags/Mars-V0.8.4/mhbase/MH3.cc

Last change on this file was 2737, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 19.2 KB
Line 
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-2002
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 TH1F, TH2F or TH3F 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// The axis binning is retrieved from the parameter list, too. Create a
44// MBinning with the name "Binning" plus the name of your MH3 container
45// plus the axis name ("X", "Y" or "Z") and add it to the parameter list.
46//
47// If you want to use a different unit for histogramming use SetScaleX,
48// SetScaleY and SetScaleZ.
49//
50// For example:
51// MH3 myhist("MHillas.fLength");
52// myhist.SetName("MyHist");
53// myhist.SetScaleX(geomcam.GetConvMm2Deg()); //convert length to degree
54// MBinning bins("BinningMyHistX");
55// bins.SetEdges(10, 0, 150);
56// plist.AddToList(&myhist);
57// plist.AddToList(&bins);
58//
59/////////////////////////////////////////////////////////////////////////////
60#include "MH3.h"
61
62#include <ctype.h> // tolower
63#include <fstream>
64
65#include <TPad.h>
66#include <TStyle.h>
67#include <TCanvas.h>
68
69#include <TH2.h>
70#include <TH3.h>
71#include <TProfile.h>
72#include <TProfile2D.h>
73
74#include "MLog.h"
75#include "MLogManip.h"
76
77#include "MParList.h"
78#include "MBinning.h"
79#include "MDataChain.h"
80
81ClassImp(MH3);
82
83using namespace std;
84
85static const TString gsDefName = "MH3";
86static const TString gsDefTitle = "Container for a %dD Mars Histogram";
87
88// --------------------------------------------------------------------------
89//
90// Default constructor.
91//
92MH3::MH3(const unsigned int dim)
93 : fDimension(dim>3?3:dim), fHist(NULL)
94{
95 switch (fDimension)
96 {
97 case 1:
98 fHist = new TH1F;
99 fHist->SetYTitle("Counts");
100 break;
101 case 2:
102 fHist = new TH2F;
103 fHist->SetZTitle("Counts");
104 break;
105 case 3:
106 fHist = new TH3F;
107 break;
108 }
109
110 fData[0] = NULL;
111 fData[1] = NULL;
112 fData[2] = NULL;
113
114 fName = gsDefName;
115 fTitle = Form(gsDefTitle.Data(), 1);
116
117 if (fHist)
118 {
119 fHist->SetDirectory(NULL);
120 fHist->UseCurrentStyle();
121 }
122
123 fScale[0] = 1;
124 fScale[1] = 1;
125 fScale[2] = 1;
126}
127
128// --------------------------------------------------------------------------
129//
130// Creates an TH1F. memberx is filled into the X-bins. For a more detailed
131// description see the class description above.
132//
133MH3::MH3(const char *memberx)
134 : fDimension(1)
135{
136 fHist = new TH1F;
137
138 fData[0] = new MDataChain(memberx);
139 fData[1] = NULL;
140 fData[2] = NULL;
141
142 fName = gsDefName;
143 fTitle = Form(gsDefTitle.Data(), 1);
144
145 fHist->UseCurrentStyle();
146 fHist->SetDirectory(NULL);
147 fHist->SetYTitle("Counts");
148
149 fScale[0] = 1;
150 fScale[1] = 1;
151 fScale[2] = 1;
152}
153
154// --------------------------------------------------------------------------
155//
156// Creates an TH2F. memberx is filled into the X-bins. membery is filled
157// into the Y-bins. For a more detailed description see the class
158// description above.
159//
160MH3::MH3(const char *memberx, const char *membery)
161 : fDimension(2)
162{
163 fHist = new TH2F;
164
165 fData[0] = new MDataChain(memberx);
166 fData[1] = new MDataChain(membery);
167 fData[2] = NULL;
168
169 fName = gsDefName;
170 fTitle = Form(gsDefTitle.Data(), 2);
171
172 fHist->UseCurrentStyle();
173 fHist->SetDirectory(NULL);
174 fHist->SetZTitle("Counts");
175
176 fScale[0] = 1;
177 fScale[1] = 1;
178 fScale[2] = 1;
179}
180
181// --------------------------------------------------------------------------
182//
183// Creates an TH3F. memberx is filled into the X-bins. membery is filled
184// into the Y-bins. membery is filled into the Z-bins. For a more detailed
185// description see the class description above.
186//
187MH3::MH3(const char *memberx, const char *membery, const char *memberz)
188 : fDimension(3)
189{
190 fHist = new TH3F;
191
192 fData[0] = new MDataChain(memberx);
193 fData[1] = new MDataChain(membery);
194 fData[2] = new MDataChain(memberz);
195
196 fName = gsDefName;
197 fTitle = Form(gsDefTitle.Data(), 3);
198
199 fHist->UseCurrentStyle();
200 fHist->SetDirectory(NULL);
201
202 fScale[0] = 1;
203 fScale[1] = 1;
204 fScale[2] = 1;
205}
206
207// --------------------------------------------------------------------------
208//
209// Deletes the histogram
210//
211MH3::~MH3()
212{
213 delete fHist;
214
215 for (int i=0; i<3; i++)
216 if (fData[i])
217 delete fData[i];
218}
219
220// --------------------------------------------------------------------------
221//
222// Return the data members used by the data chain to be used in
223// MTask::AddBranchToList
224//
225TString MH3::GetDataMember() const
226{
227 TString str=fData[0]->GetDataMember();
228 if (fData[1])
229 {
230 str += ";";
231 str += fData[1]->GetDataMember();
232 }
233 if (fData[2])
234 {
235 str += ";";
236 str += fData[2]->GetDataMember();
237 }
238 return str;
239}
240
241// --------------------------------------------------------------------------
242//
243// Setup the Binning for the histograms automatically if the correct
244// instances of MBinning are found in the parameter list
245// For a more detailed description see class description above.
246//
247Bool_t MH3::SetupFill(const MParList *plist)
248{
249 // reset histogram (necessary if the same eventloop is run more than once)
250 fHist->Reset();
251
252 TString bname("Binning");
253 bname += fName;
254
255 MBinning *binsx = NULL;
256 MBinning *binsy = NULL;
257 MBinning *binsz = NULL;
258
259 switch (fDimension)
260 {
261 case 3:
262 binsz = (MBinning*)plist->FindObject(bname+"Z", "MBinning");
263 if (!binsz)
264 {
265 *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
266 return kFALSE;
267 }
268 if (binsz->IsLogarithmic())
269 fHist->SetBit(kIsLogz);
270 if (fData[2]) fHist->SetZTitle(fData[2]->GetTitle());
271 if (fData[2] && !fData[2]->PreProcess(plist))
272 return kFALSE;
273 case 2:
274 binsy = (MBinning*)plist->FindObject(bname+"Y", "MBinning");
275 if (!binsy)
276 {
277 *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl;
278 return kFALSE;
279 }
280 if (binsy->IsLogarithmic())
281 fHist->SetBit(kIsLogy);
282 if (fData[1]) fHist->SetYTitle(fData[1]->GetTitle());
283 if (fData[1] && !fData[1]->PreProcess(plist))
284 return kFALSE;
285 case 1:
286 binsx = (MBinning*)plist->FindObject(bname+"X", "MBinning");
287 if (!binsx)
288 {
289 if (fDimension==1)
290 binsx = (MBinning*)plist->FindObject(bname, "MBinning");
291
292 if (!binsx)
293 {
294 *fLog << err << dbginf << "Neither MBinning '" << bname << "X' nor '" << bname << "' found... aborting." << endl;
295 return kFALSE;
296 }
297
298 }
299 if (binsx->IsLogarithmic())
300 fHist->SetBit(kIsLogx);
301
302 if (fData[0]!=NULL) fHist->SetXTitle(fData[0]->GetTitle());
303 if (fData[0] && !fData[0]->PreProcess(plist))
304 return kFALSE;
305 }
306
307 fHist->SetName(fName);
308
309 TString title("Histogram for ");
310 title += fName;
311
312 switch (fDimension)
313 {
314 case 1:
315 fHist->SetTitle(title+" (1D)");
316 SetBinning(fHist, binsx);
317 return kTRUE;
318 case 2:
319 fHist->SetTitle(title+" (2D)");
320 SetBinning((TH2*)fHist, binsx, binsy);
321 return kTRUE;
322 case 3:
323 fHist->SetTitle(title+" (3D)");
324 SetBinning((TH3*)fHist, binsx, binsy, binsz);
325 return kTRUE;
326 }
327 cout << "Still alive...?" << endl;
328 return kTRUE;
329}
330
331// --------------------------------------------------------------------------
332//
333// Set the name of the histogram ant the MH3 container
334//
335void MH3::SetName(const char *name)
336{
337 fHist->SetName(name);
338 MParContainer::SetName(name);
339}
340
341// --------------------------------------------------------------------------
342//
343// Set the title of the histogram ant the MH3 container
344//
345void MH3::SetTitle(const char *title)
346{
347 fHist->SetTitle(title);
348 MParContainer::SetTitle(title);
349}
350
351// --------------------------------------------------------------------------
352//
353// Fills the one, two or three data members into our histogram
354//
355Bool_t MH3::Fill(const MParContainer *par, const Stat_t w)
356{
357 Double_t x=0;
358 Double_t y=0;
359 Double_t z=0;
360
361 switch (fDimension)
362 {
363 case 3:
364 z = fData[2]->GetValue()*fScale[2];
365 case 2:
366 y = fData[1]->GetValue()*fScale[1];
367 case 1:
368 x = fData[0]->GetValue()*fScale[0];
369 }
370
371 switch (fDimension)
372 {
373 case 3:
374 ((TH3*)fHist)->Fill(x, y, z, w);
375 return kTRUE;
376 case 2:
377 ((TH2*)fHist)->Fill(x, y, w);
378 return kTRUE;
379 case 1:
380 fHist->Fill(x, w);
381 return kTRUE;
382 }
383
384 return kFALSE;
385}
386/*
387// --------------------------------------------------------------------------
388//
389// Set the palette you wanna use:
390// - you could set the root "Pretty Palette Violet->Red" by
391// gStyle->SetPalette(1, 0), but in some cases this may look
392// confusing
393// - The maximum colors root allowes us to set by ourself
394// is 50 (idx: 51-100). This colors are set to a grayscaled
395// palette
396// - the number of contours must be two less than the number
397// of palette entries
398//
399void MHStarMap::PrepareDrawing() const
400{
401 const Int_t numg = 32; // number of gray scaled colors
402 const Int_t numw = 32; // number of white
403
404 Int_t palette[numg+numw];
405
406 //
407 // The first half of the colors are white.
408 // This is some kind of optical background supression
409 //
410 gROOT->GetColor(51)->SetRGB(1, 1, 1);
411
412 Int_t i;
413 for (i=0; i<numw; i++)
414 palette[i] = 51;
415
416 //
417 // now the (gray) scaled part is coming
418 //
419 for (;i<numw+numg; i++)
420 {
421 const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
422
423 gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
424 palette[i] = 52+i;
425 }
426
427 //
428 // Set the palette and the number of contour levels
429 //
430 gStyle->SetPalette(numg+numw, palette);
431 fStarMap->SetContour(numg+numw-2);
432}
433*/
434// --------------------------------------------------------------------------
435//
436// Setup a inversed deep blue sea palette for the fCenter histogram.
437//
438void MH3::SetColors() const
439{
440 // FIXME: This must be redone each time the canvas is repainted....
441 gStyle->SetPalette(51, NULL);
442 Int_t c[50];
443 for (int i=0; i<50; i++)
444 c[49-i] = gStyle->GetColorPalette(i);
445 gStyle->SetPalette(50, c);
446}
447
448// --------------------------------------------------------------------------
449//
450// Draw clone of histogram. So that the object can be deleted
451//
452// Possible options are:
453// PROFX: Draw a x-profile into the histogram (for 2D histograms only)
454// PROFY: Draw a y-profile into the histogram (for 2D histograms only)
455// ONLY: Draw the profile histogram only (for 2D histograms only)
456//
457// If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
458// eg this is set when applying a logarithmic MBinning
459//
460// and the histogram is still visible in the canvas.
461// The cloned object are deleted together with the canvas if the canvas is
462// destroyed. If you want to handle destroying the canvas you can get a
463// pointer to it from this function
464//
465/*
466TObject *MH3::DrawClone(Option_t *opt) const
467{
468 TString str(opt);
469
470 TVirtualPad *c = gPad;
471 if (!str.Contains("nonew", TString::kIgnoreCase))
472 {
473 c = MH::MakeDefCanvas(fHist);
474
475 //
476 // This is necessary to get the expected bahviour of DrawClone
477 //
478 gROOT->SetSelectedPad(NULL);
479 }
480
481 if (str.Contains("COL", TString::kIgnoreCase))
482 SetColors();
483
484 Bool_t only = str.Contains("ONLY", TString::kIgnoreCase) && fDimension==2;
485 if (!only)
486 fHist->DrawCopy(opt);
487
488 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
489 {
490 TProfile *p = ((TH2*)fHist)->ProfileX("_pfx", -1, 9999, "s");
491 p->SetLineColor(kBlue);
492 p->Draw(only?"":"same");
493 p->SetBit(kCanDelete);
494 p->SetDirectory(NULL);
495 }
496 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
497 {
498 TProfile *p = ((TH2*)fHist)->ProfileY("_pfy", -1, 9999, "s");
499 p->SetLineColor(kBlue);
500 p->Draw(only?"":"same");
501 p->SetBit(kCanDelete);
502 p->SetDirectory(NULL);
503 }
504
505 if (fHist->TestBit(kIsLogx)) c->SetLogx();
506 if (fHist->TestBit(kIsLogy)) c->SetLogy();
507 if (fHist->TestBit(kIsLogz)) c->SetLogz();
508
509 c->Modified();
510 c->Update();
511
512 return c;
513}
514*/
515
516// --------------------------------------------------------------------------
517//
518// Creates a new canvas and draws the histogram into it.
519//
520// Possible options are:
521// PROFX: Draw a x-profile into the histogram (for 2D histograms only)
522// PROFY: Draw a y-profile into the histogram (for 2D histograms only)
523// ONLY: Draw the profile histogram only (for 2D histograms only)
524//
525// If the kIsLog?-Bit is set the axis is displayed lkogarithmically.
526// eg this is set when applying a logarithmic MBinning
527//
528// Be careful: The histogram belongs to this object and won't get deleted
529// together with the canvas.
530//
531void MH3::Draw(Option_t *opt)
532{
533 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(fHist);
534 pad->SetBorderMode(0);
535
536 AppendPad("");
537
538 TString str(opt);
539
540 // FIXME: Do it in Paint()
541 if (str.Contains("COL", TString::kIgnoreCase))
542 SetColors();
543
544 fHist->SetFillStyle(4000);
545
546 Bool_t only = str.Contains("ONLY", TString::kIgnoreCase) && fDimension==2;
547 if (!only)
548 fHist->Draw(opt);
549
550 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
551 {
552 TProfile *p = ((TH2*)fHist)->ProfileX("_pfx", -1, 9999, "s");
553 p->SetLineColor(kBlue);
554 p->Draw(only?"":"same");
555 p->SetBit(kCanDelete);
556 p->SetDirectory(NULL);
557 }
558 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
559 {
560 TProfile *p = ((TH2*)fHist)->ProfileY("_pfy", -1, 9999, "s");
561 p->SetLineColor(kBlue);
562 p->Draw(only?"":"same");
563 p->SetBit(kCanDelete);
564 p->SetDirectory(NULL);
565 }
566
567 if (fHist->TestBit(kIsLogx)) pad->SetLogx();
568 if (fHist->TestBit(kIsLogy)) pad->SetLogy();
569 if (fHist->TestBit(kIsLogz)) pad->SetLogz();
570
571 pad->Modified();
572 pad->Update();
573}
574
575// --------------------------------------------------------------------------
576//
577// Implementation of SavePrimitive. Used to write the call to a constructor
578// to a macro. In the original root implementation it is used to write
579// gui elements to a macro-file.
580//
581void MH3::StreamPrimitive(ofstream &out) const
582{
583 TString name = GetUniqueName();
584
585 out << " MH3 " << name << "(\"";
586 out << fData[0]->GetRule() << "\"";
587 if (fDimension>1)
588 out << ", \"" << fData[1]->GetRule() << "\"";
589 if (fDimension>2)
590 out << ", \"" << fData[2]->GetRule() << "\"";
591
592 out << ");" << endl;
593
594 if (fName!=gsDefName)
595 out << " " << name << ".SetName(\"" << fName << "\");" << endl;
596
597 if (fTitle!=Form(gsDefTitle.Data(), fDimension))
598 out << " " << name << ".SetTitle(\"" << fTitle << "\");" << endl;
599
600 switch (fDimension)
601 {
602 case 3:
603 if (fScale[2]!=1)
604 out << " " << name << ".SetScaleZ(" << fScale[2] << ");" << endl;
605 case 2:
606 if (fScale[1]!=1)
607 out << " " << name << ".SetScaleY(" << fScale[1] << ");" << endl;
608 case 1:
609 if (fScale[0]!=1)
610 out << " " << name << ".SetScaleX(" << fScale[0] << ");" << endl;
611 }
612}
613
614// --------------------------------------------------------------------------
615//
616// Used to rebuild a MH3 object of the same type (data members,
617// dimension, ...)
618//
619MParContainer *MH3::New() const
620{
621 MH3 *h = NULL;
622
623 if (fData[0] == NULL)
624 h=new MH3(fDimension);
625 else
626 switch (fDimension)
627 {
628 case 1:
629 h=new MH3(fData[0]->GetRule());
630 break;
631 case 2:
632 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule());
633 break;
634 case 3:
635 h=new MH3(fData[0]->GetRule(), fData[1]->GetRule(), fData[2]->GetRule());
636 break;
637 }
638 switch (fDimension)
639 {
640 case 3:
641 h->SetScaleZ(fScale[2]);
642 case 2:
643 h->SetScaleY(fScale[1]);
644 case 1:
645 h->SetScaleX(fScale[0]);
646 }
647 return h;
648}
649
650TString MH3::GetRule(const Char_t axis) const
651{
652 switch (tolower(axis))
653 {
654 case 'x':
655 return fData[0] ? fData[0]->GetRule() : TString("");
656 case 'y':
657 return fData[1] ? fData[1]->GetRule() : TString("");
658 case 'z':
659 return fData[2] ? fData[2]->GetRule() : TString("");
660 default:
661 return "<n/a>";
662 }
663}
664
665// --------------------------------------------------------------------------
666//
667// Returns the total number of bins in a histogram (excluding under- and
668// overflow bins)
669//
670Int_t MH3::GetNbins() const
671{
672 Int_t num = 1;
673
674 switch (fDimension)
675 {
676 case 3:
677 num *= fHist->GetNbinsZ()+2;
678 case 2:
679 num *= fHist->GetNbinsY()+2;
680 case 1:
681 num *= fHist->GetNbinsX()+2;
682 }
683
684 return num;
685}
686
687// --------------------------------------------------------------------------
688//
689// Returns the total number of bins in a histogram (excluding under- and
690// overflow bins) Return -1 if bin is underflow or overflow
691//
692Int_t MH3::FindFixBin(Double_t x, Double_t y, Double_t z) const
693{
694 const TAxis &axex = *fHist->GetXaxis();
695 const TAxis &axey = *fHist->GetYaxis();
696 const TAxis &axez = *fHist->GetZaxis();
697
698 Int_t binz = 0;
699 Int_t biny = 0;
700 Int_t binx = 0;
701
702 switch (fDimension)
703 {
704 case 3:
705 binz = axez.FindFixBin(z);
706 if (binz>axez.GetLast() || binz<axez.GetFirst())
707 return -1;
708 case 2:
709 biny = axey.FindFixBin(y);
710 if (biny>axey.GetLast() || biny<axey.GetFirst())
711 return -1;
712 case 1:
713 binx = axex.FindFixBin(x);
714 if (binx<axex.GetFirst() || binx>axex.GetLast())
715 return -1;
716 }
717
718 const Int_t nx = fHist->GetNbinsX()+2;
719 const Int_t ny = fHist->GetNbinsY()+2;
720
721 return binx + nx*(biny +ny*binz);
722}
Note: See TracBrowser for help on using the repository browser.