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

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