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

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