source: trunk/MagicSoft/Mars/mhist/MHHillas.cc@ 1210

Last change on this file since 1210 was 1210, checked in by tbretz, 23 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 4.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 2001 <mailto:tbretz@uni-sw.gwdg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHHillas
28//
29// This class contains histograms for every Hillas parameter
30//
31/////////////////////////////////////////////////////////////////////////////
32
33#include "MHHillas.h"
34
35#include <math.h>
36
37#include <TH1.h>
38#include <TPad.h>
39#include <TCanvas.h>
40
41#include "MHillas.h"
42#include "MParList.h"
43
44ClassImp(MHHillas);
45
46// --------------------------------------------------------------------------
47//
48// Setup four histograms for Width, Length
49//
50MHHillas::MHHillas(const char *name, const char *title)
51{
52 //
53 // set the name and title of this object
54 //
55 fName = name ? name : "MHHillas";
56 fTitle = title ? title : "Container for Hillas histograms";
57
58 //
59 // loop over all Pixels and create two histograms
60 // one for the Low and one for the High gain
61 // connect all the histogram with the container fHist
62 //
63 fWidth = new TH1F("Width [mm]", "Width of Hillas", 100, 0, 300);
64 fLength = new TH1F("Length [mm]", "Length of Hillas", 100, 0, 300);
65
66 fLength->SetDirectory(NULL);
67 fWidth->SetDirectory(NULL);
68
69 fLength->GetXaxis()->SetTitle("Length [mm]");
70 fWidth->GetXaxis()->SetTitle("Width [mm]");
71
72 fLength->GetYaxis()->SetTitle("Counts");
73 fWidth->GetYaxis()->SetTitle("Counts");
74}
75
76// --------------------------------------------------------------------------
77//
78// Delete the four histograms
79//
80MHHillas::~MHHillas()
81{
82 delete fWidth;
83 delete fLength;
84}
85
86// --------------------------------------------------------------------------
87//
88// Setup the Binning for the histograms automatically if the correct
89// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
90// are found in the parameter list
91//
92Bool_t MHHillas::SetupFill(const MParList *plist)
93{
94 const MBinning* binsw = (MBinning*)plist->FindObject("BinningWidth");
95 const MBinning* binsl = (MBinning*)plist->FindObject("BinningLength");
96
97 if (binsw)
98 SetBinning(fWidth, binsw);
99
100 if (binsl)
101 SetBinning(fLength, binsl);
102
103 return kTRUE;
104}
105
106// --------------------------------------------------------------------------
107//
108// Fill the four histograms with data from a MHillas-Container.
109// Be careful: Only call this with an object of type MHillas
110//
111Bool_t MHHillas::Fill(const MParContainer *par)
112{
113 const MHillas &h = *(MHillas*)par;
114
115 fWidth ->Fill(h.GetWidth());
116 fLength->Fill(h.GetLength());
117
118 return kTRUE;
119}
120
121// --------------------------------------------------------------------------
122//
123// Draw clones of all four histograms. So that the object can be deleted
124// and the histograms are still visible in the canvas.
125// The cloned object are deleted together with the canvas if the canvas is
126// destroyed. If you want to handle dostroying the canvas you can get a
127// pointer to it from this function
128//
129TObject *MHHillas::DrawClone(Option_t *opt) const
130{
131 TCanvas *c = MakeDefCanvas("Hillas", "Histograms of Hillas Parameters",
132 350, 500);
133 c->Divide(1, 2);
134
135 gROOT->SetSelectedPad(NULL);
136
137 //
138 // This is necessary to get the expected bahviour of DrawClone
139 //
140 c->cd(1);
141 fLength->DrawCopy();
142
143 c->cd(2);
144 fWidth->DrawCopy();
145
146 c->Modified();
147 c->Update();
148
149 return c;
150}
151
152// --------------------------------------------------------------------------
153//
154// Creates a new canvas and draws the four histograms into it.
155// Be careful: The histograms belongs to this object and won't get deleted
156// together with the canvas.
157//
158void MHHillas::Draw(Option_t *)
159{
160 if (!gPad)
161 MakeDefCanvas("Hillas", "Histograms of Hillas Parameters", 350, 500);
162
163 gPad->Divide(1, 2);
164
165 gPad->cd(1);
166 fLength->Draw();
167
168 gPad->cd(2);
169 fWidth->Draw();
170
171 gPad->Modified();
172 gPad->Update();
173}
Note: See TracBrowser for help on using the repository browser.