source: trunk/MagicSoft/Mars/mhbase/MH3.cc@ 7900

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