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

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