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

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