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

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