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

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