source: trunk/MagicSoft/Mars/mhist/MH.cc@ 1290

Last change on this file since 1290 was 1283, checked in by tbretz, 24 years ago
*** empty log message ***
File size: 9.2 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 07/2001 <mailto:tbretz@uni-sw.gwdg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26// //
27// MH //
28// //
29// This is a base tasks for mars histograms. It defines a common interface //
30// for filling the histograms with events (MH::Fill) which is used by a //
31// common 'filler' And a SetupFill member function which may be used //
32// by MFillH. The idea is: //
33// 1) If your Histogram can become filled by one single container //
34// (like MHHillas) you overload MH::Fill and it gets called with //
35// a pointer to the container with which it should be filled. //
36// //
37// 2) You histogram needs several containers to get filled. Than you //
38// have to overload MH::SetupFill and get the necessary objects from //
39// the parameter list. Use this objects in Fill to fill your //
40// histogram. //
41// //
42// If you want to create your own histogram class the new class must be //
43// derived from MH (instead of the base MParContainer) and you must //
44// the fill function of MH. This is the function which is called to fill //
45// the histogram(s) by the data of a corresponding parameter container. //
46// //
47//////////////////////////////////////////////////////////////////////////////
48
49#include "MH.h"
50
51#include <TH1.h>
52#include <TCanvas.h>
53
54#include "MBinning.h"
55
56ClassImp(MH);
57
58// --------------------------------------------------------------------------
59//
60// Default Constructor. It sets name and title only. Typically you won't
61// need to change this.
62//
63MH::MH(const char *name, const char *title)
64{
65 //
66 // set the name and title of this object
67 //
68 fName = name ? name : "MH";
69 fTitle = title ? title : "Base class for Mars histograms";
70}
71
72// --------------------------------------------------------------------------
73//
74// This is a function which should replace the creation of default
75// canvases like root does. Because this is inconvinient in some aspects.
76// need to change this.
77// You can specify a name for the default canvas and a title. Also
78// width and height can be given.
79// MakeDefCanvas looks for a canvas with the given name. If now name is
80// given the DefCanvasName of root is used. If no such canvas is existing
81// it is created and returned. If such a canvas already exists a new canvas
82// with a name plus anumber is created (the number is calculated by the
83// number of all existing canvases plus one)
84//
85TCanvas *MH::MakeDefCanvas(TString name, const char *title,
86 const UInt_t w, const UInt_t h)
87{
88 const TList *list = (TList*)gROOT->GetListOfCanvases();
89
90 if (name.IsNull())
91 name = gROOT->GetDefCanvasName();
92
93 if (list->FindObject(name))
94 name += Form(" <%d>", list->GetSize()+1);
95
96 return new TCanvas(name, title, w, h);
97}
98
99// --------------------------------------------------------------------------
100//
101// This function works like MakeDefCanvas(name, title, w, h) but name
102// and title are retrieved from the given TObject.
103//
104TCanvas *MH::MakeDefCanvas(const TObject *obj,
105 const UInt_t w, const UInt_t h)
106{
107 return MakeDefCanvas(obj->GetName(), obj->GetTitle(), w, h);
108}
109
110void MH::SetBinning(TH1 *h, const MBinning *binsx)
111{
112 //
113 // Another strange behaviour: TAxis::Set deletes the axis title!
114 //
115 TAxis &x = *h->GetXaxis();
116
117#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
118 TString xtitle = x.GetTitle();
119#endif
120
121 //
122 // This is a necessary workaround if one wants to set
123 // non-equidistant bins after the initialization
124 // TH1D::fNcells must be set correctly.
125 //
126 h->SetBins(binsx->GetNumBins(), 0, 1);
127
128 //
129 // Set the binning of the current histogram to the binning
130 // in one of the two given histograms
131 //
132 x.Set(binsx->GetNumBins(), binsx->GetEdges());
133#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
134 x.SetTitle(xtitle);
135#endif
136}
137
138void MH::SetBinning(TH1 *h, const MBinning *binsx, const MBinning *binsy)
139{
140 TAxis &x = *h->GetXaxis();
141 TAxis &y = *h->GetYaxis();
142
143 //
144 // Another strange behaviour: TAxis::Set deletes the axis title!
145 //
146#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
147 TString xtitle = x.GetTitle();
148 TString ytitle = y.GetTitle();
149#endif
150
151 //
152 // This is a necessary workaround if one wants to set
153 // non-equidistant bins after the initialization
154 // TH1D::fNcells must be set correctly.
155 //
156 h->SetBins(binsx->GetNumBins(), 0, 1,
157 binsy->GetNumBins(), 0, 1);
158
159 //
160 // Set the binning of the current histogram to the binning
161 // in one of the two given histograms
162 //
163 x.Set(binsx->GetNumBins(), binsx->GetEdges());
164 y.Set(binsy->GetNumBins(), binsy->GetEdges());
165#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
166 x.SetTitle(xtitle);
167 y.SetTitle(ytitle);
168#endif
169}
170
171void MH::SetBinning(TH1 *h, const MBinning *binsx, const MBinning *binsy, const MBinning *binsz)
172{
173 //
174 // Another strange behaviour: TAxis::Set deletes the axis title!
175 //
176 TAxis &x = *h->GetXaxis();
177 TAxis &y = *h->GetYaxis();
178 TAxis &z = *h->GetZaxis();
179
180#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
181 TString xtitle = x.GetTitle();
182 TString ytitle = y.GetTitle();
183 TString ztitle = z.GetTitle();
184#endif
185
186 //
187 // This is a necessary workaround if one wants to set
188 // non-equidistant bins after the initialization
189 // TH1D::fNcells must be set correctly.
190 //
191 h->SetBins(binsx->GetNumBins(), 0, 1,
192 binsy->GetNumBins(), 0, 1,
193 binsz->GetNumBins(), 0, 1);
194
195 //
196 // Set the binning of the current histogram to the binning
197 // in one of the two given histograms
198 //
199 x.Set(binsx->GetNumBins(), binsx->GetEdges());
200 y.Set(binsy->GetNumBins(), binsy->GetEdges());
201 z.Set(binsz->GetNumBins(), binsz->GetEdges());
202#if ROOT_VERSION_CODE < ROOT_VERSION(3,03,03)
203 x.SetTitle(xtitle);
204 y.SetTitle(ytitle);
205 z.SetTitle(ztitle);
206#endif
207}
208
209void MH::SetBinning(TH1 *h, const TArrayD *binsx)
210{
211 MBinning bx;
212 bx.SetEdges(*binsx);
213 SetBinning(h, &bx);
214}
215
216void MH::SetBinning(TH1 *h, const TArrayD *binsx, const TArrayD *binsy)
217{
218 MBinning bx;
219 MBinning by;
220 bx.SetEdges(*binsx);
221 by.SetEdges(*binsy);
222 SetBinning(h, &bx, &by);
223}
224
225void MH::SetBinning(TH1 *h, const TArrayD *binsx, const TArrayD *binsy, const TArrayD *binsz)
226{
227 MBinning bx;
228 MBinning by;
229 MBinning bz;
230 bx.SetEdges(*binsx);
231 by.SetEdges(*binsy);
232 bz.SetEdges(*binsz);
233 SetBinning(h, &bx, &by, &bz);
234}
235
236void MH::SetBinning(TH1 *h, const TAxis *binsx)
237{
238 const Int_t nx = binsx->GetNbins();
239
240 TArrayD bx(nx+1);
241 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
242 bx[nx] = binsx->GetXmax();
243
244 SetBinning(h, &bx);
245}
246
247void MH::SetBinning(TH1 *h, const TAxis *binsx, const TAxis *binsy)
248{
249 const Int_t nx = binsx->GetNbins();
250 const Int_t ny = binsy->GetNbins();
251
252 TArrayD bx(nx+1);
253 TArrayD by(ny+1);
254 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
255 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
256 bx[nx] = binsx->GetXmax();
257 by[ny] = binsy->GetXmax();
258
259 SetBinning(h, &bx, &by);
260}
261
262void MH::SetBinning(TH1 *h, const TAxis *binsx, const TAxis *binsy, const TAxis *binsz)
263{
264 const Int_t nx = binsx->GetNbins();
265 const Int_t ny = binsy->GetNbins();
266 const Int_t nz = binsz->GetNbins();
267
268 TArrayD bx(nx+1);
269 TArrayD by(ny+1);
270 TArrayD bz(nz+1);
271 for (int i=0; i<nx; i++) bx[i] = binsx->GetBinLowEdge(i+1);
272 for (int i=0; i<ny; i++) by[i] = binsy->GetBinLowEdge(i+1);
273 for (int i=0; i<nz; i++) bz[i] = binsz->GetBinLowEdge(i+1);
274 bx[nx] = binsx->GetXmax();
275 by[ny] = binsy->GetXmax();
276 bz[nz] = binsz->GetXmax();
277
278 SetBinning(h, &bx, &by, &bz);
279}
280
281void MH::SetBinning(TH1 *h, TH1 *x)
282{
283 SetBinning(h, x->GetXaxis(), x->GetYaxis(), x->GetZaxis());
284}
Note: See TracBrowser for help on using the repository browser.