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

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