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

Last change on this file since 1463 was 1336, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 9.8 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
61#include "MH3.h"
62
63#include <TH2.h>
64#include <TH3.h>
65#include <TProfile.h>
66#include <TProfile2D.h>
67#include <TPad.h>
68#include <TCanvas.h>
69
70#include "MLog.h"
71#include "MLogManip.h"
72
73#include "MParList.h"
74
75#include "MDataChain.h"
76
77ClassImp(MH3);
78
79// --------------------------------------------------------------------------
80//
81// Creates an TH1F. memberx is filled into the X-bins. For a more detailed
82// description see the class description above.
83//
84MH3::MH3(const char *memberx)
85 : fDimension(1)
86{
87 fHist = new TH1F;
88
89 fData[0] = new MDataChain(memberx);
90 fData[1] = NULL;
91 fData[2] = NULL;
92
93 fName = "MH3";
94 fTitle = "Container for a 1D Mars Histogram";
95
96 fHist->SetDirectory(NULL);
97 fHist->SetYTitle("Counts");
98
99 fScale[0] = 1;
100 fScale[1] = 1;
101 fScale[2] = 1;
102}
103
104// --------------------------------------------------------------------------
105//
106// Creates an TH2F. memberx is filled into the X-bins. membery is filled
107// into the Y-bins. For a more detailed description see the class
108// description above.
109//
110MH3::MH3(const char *memberx, const char *membery)
111 : fDimension(2)
112{
113 fHist = new TH2F;
114
115 fData[0] = new MDataChain(memberx);
116 fData[1] = new MDataChain(membery);
117 fData[2] = NULL;
118
119 fName = "MH3";
120 fTitle = "Container for a 2D Mars Histogram";
121
122 fHist->SetDirectory(NULL);
123 fHist->SetZTitle("Counts");
124
125 fScale[0] = 1;
126 fScale[1] = 1;
127 fScale[2] = 1;
128}
129
130// --------------------------------------------------------------------------
131//
132// Creates an TH3F. memberx is filled into the X-bins. membery is filled
133// into the Y-bins. membery is filled into the Z-bins. For a more detailed
134// description see the class description above.
135//
136MH3::MH3(const char *memberx, const char *membery, const char *memberz)
137 : fDimension(3)
138{
139 fHist = new TH3F;
140
141 fData[0] = new MDataChain(memberx);
142 fData[1] = new MDataChain(membery);
143 fData[2] = new MDataChain(memberz);
144
145 fName = "MH3";
146 fTitle = "Container for a 3D Mars Histogram";
147
148 fHist->SetDirectory(NULL);
149
150 fScale[0] = 1;
151 fScale[1] = 1;
152 fScale[2] = 1;
153}
154
155// --------------------------------------------------------------------------
156//
157// Deletes the histogram
158//
159MH3::~MH3()
160{
161 delete fHist;
162
163 for (int i=0; i<3; i++)
164 if (fData[i])
165 delete fData[i];
166}
167
168// --------------------------------------------------------------------------
169//
170// Setup the Binning for the histograms automatically if the correct
171// instances of MBinning are found in the parameter list
172// For a more detailed description see class description above.
173//
174Bool_t MH3::SetupFill(const MParList *plist)
175{
176 TString bname("Binning");
177 bname += fName;
178
179 MBinning *binsx = NULL;
180 MBinning *binsy = NULL;
181 MBinning *binsz = NULL;
182 switch (fDimension)
183 {
184 case 3:
185 binsz = (MBinning*)plist->FindObject(bname+"Z");
186 if (!binsz)
187 {
188 *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
189 return kFALSE;
190 }
191 fHist->SetZTitle(fData[2]->GetTitle());
192 if (!fData[2]->PreProcess(plist))
193 return kFALSE;
194 case 2:
195 binsy = (MBinning*)plist->FindObject(bname+"Y");
196 if (!binsy)
197 {
198 *fLog << err << dbginf << "MBinning '" << bname << "Y' not found... aborting." << endl;
199 return kFALSE;
200 }
201 fHist->SetYTitle(fData[1]->GetTitle());
202 if (!fData[1]->PreProcess(plist))
203 return kFALSE;
204 case 1:
205 binsx = (MBinning*)plist->FindObject(bname+"X");
206 if (!binsx)
207 {
208 *fLog << err << dbginf << "MBinning '" << bname << "X' not found... aborting." << endl;
209 return kFALSE;
210 }
211 fHist->SetXTitle(fData[0]->GetTitle());
212 if (!fData[0]->PreProcess(plist))
213 return kFALSE;
214 }
215
216 fHist->SetName(fName);
217
218 TString title("Histogram for ");
219 title += fName;
220
221 switch (fDimension)
222 {
223 case 1:
224 fHist->SetTitle(title+" (1D)");
225 SetBinning(fHist, binsx);
226 return kTRUE;
227 case 2:
228 fHist->SetTitle(title+" (2D)");
229 SetBinning((TH2*)fHist, binsx, binsy);
230 return kTRUE;
231 case 3:
232 fHist->SetTitle(title+" (3D)");
233 SetBinning((TH3*)fHist, binsx, binsy, binsz);
234 return kTRUE;
235 }
236
237 return kTRUE;
238}
239
240// --------------------------------------------------------------------------
241//
242// Set the name of the histogram ant the MH3 container
243//
244void MH3::SetName(const char *name)
245{
246 fHist->SetName(name);
247 MParContainer::SetName(name);
248}
249
250// --------------------------------------------------------------------------
251//
252// Set the title of the histogram ant the MH3 container
253//
254void MH3::SetTitle(const char *title)
255{
256 fHist->SetTitle(title);
257 MParContainer::SetTitle(title);
258}
259
260// --------------------------------------------------------------------------
261//
262// Fills the one, two or three data members into our histogram
263//
264Bool_t MH3::Fill(const MParContainer *par)
265{
266 Double_t x=0;
267 Double_t y=0;
268 Double_t z=0;
269
270 switch (fDimension)
271 {
272 case 3:
273 z = fData[2]->GetValue()*fScale[2];
274 case 2:
275 y = fData[1]->GetValue()*fScale[1];
276 case 1:
277 x = fData[0]->GetValue()*fScale[0];
278 }
279
280 switch (fDimension)
281 {
282 case 3:
283 ((TH3*)fHist)->Fill(x, y, z);
284 return kTRUE;
285 case 2:
286 ((TH2*)fHist)->Fill(x, y);
287 return kTRUE;
288 case 1:
289 fHist->Fill(x);
290 return kTRUE;
291 }
292
293 return kFALSE;
294}
295
296// --------------------------------------------------------------------------
297//
298// Draw clone of histogram. So that the object can be deleted
299// and the histogram is still visible in the canvas.
300// The cloned object are deleted together with the canvas if the canvas is
301// destroyed. If you want to handle destroying the canvas you can get a
302// pointer to it from this function
303//
304TObject *MH3::DrawClone(Option_t *opt) const
305{
306 TCanvas &c = *MH::MakeDefCanvas(fHist);
307
308 //
309 // This is necessary to get the expected bahviour of DrawClone
310 //
311 gROOT->SetSelectedPad(NULL);
312
313 fHist->DrawCopy(opt);
314
315 TString str(opt);
316 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
317 {
318 TProfile *p = ((TH2*)fHist)->ProfileX();
319 p->Draw("same");
320 p->SetBit(kCanDelete);
321 p->SetDirectory(NULL);
322 }
323 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
324 {
325 TProfile *p = ((TH2*)fHist)->ProfileY();
326 p->Draw("same");
327 p->SetBit(kCanDelete);
328 p->SetDirectory(NULL);
329 }
330
331 c.Modified();
332 c.Update();
333
334 return &c;
335}
336
337// --------------------------------------------------------------------------
338//
339// Creates a new canvas and draws the histogram into it.
340// Be careful: The histogram belongs to this object and won't get deleted
341// together with the canvas.
342//
343void MH3::Draw(Option_t *opt)
344{
345 if (!gPad)
346 MH::MakeDefCanvas(fHist);
347
348 fHist->Draw(opt);
349
350 TString str(opt);
351 if (str.Contains("PROFX", TString::kIgnoreCase) && fDimension==2)
352 {
353 TProfile *p = ((TH2*)fHist)->ProfileX();
354 p->Draw("same");
355 p->SetBit(kCanDelete);
356 p->SetDirectory(NULL);
357 }
358 if (str.Contains("PROFY", TString::kIgnoreCase) && fDimension==2)
359 {
360 TProfile *p = ((TH2*)fHist)->ProfileY();
361 p->Draw("same");
362 p->SetBit(kCanDelete);
363 p->SetDirectory(NULL);
364 }
365
366 gPad->Modified();
367 gPad->Update();
368}
Note: See TracBrowser for help on using the repository browser.