source: trunk/MagicSoft/Mars/mhist/MH3.cc@ 2124

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