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

Last change on this file since 1441 was 1441, checked in by tbretz, 22 years ago
*** empty log message ***
  • Property svn:executable set to *
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 2001 <mailto:tbretz@uni-sw.gwdg.de>
19! Wolfgang Wittek 2002 <mailto:wittek@mppmu.mpg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2002
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MHHillas
29//
30// This class contains histograms for the source independent image parameters
31//
32/////////////////////////////////////////////////////////////////////////////
33
34#include "MHHillas.h"
35
36#include <math.h>
37
38#include <TH1.h>
39#include <TH2.h>
40#include <TPad.h>
41#include <TStyle.h>
42#include <TCanvas.h>
43
44#include "MLog.h"
45#include "MLogManip.h"
46
47#include "MParList.h"
48
49#include "MHillas.h"
50#include "MGeomCam.h"
51#include "MBinning.h"
52
53
54ClassImp(MHHillas);
55
56// --------------------------------------------------------------------------
57//
58// Setup four histograms for Width, Length
59//
60MHHillas::MHHillas(const char *name, const char *title)
61 : fMm2Deg(1), fUseMmScale(kFALSE)
62{
63 //
64 // set the name and title of this object
65 //
66 fName = name ? name : "MHHillas";
67 fTitle = title ? title : "Source independent image parameters";
68
69 fLength = new TH1F("Length", "Length of Ellipse", 100, 0, 300);
70 fLength->SetDirectory(NULL);
71 fLength->SetXTitle("Length [mm]");
72 fLength->SetYTitle("Counts");
73
74 fWidth = new TH1F("Width", "Width of Ellipse", 100, 0, 300);
75 fWidth->SetDirectory(NULL);
76 fWidth->SetXTitle("Width [mm]");
77 fWidth->SetYTitle("Counts");
78
79 fDistC = new TH1F("DistC", "Distance from center of camera", 100, 0, 600);
80 fDistC->SetDirectory(NULL);
81 fDistC->SetXTitle("Distance [mm]");
82 fDistC->SetYTitle("Counts");
83
84 fDelta = new TH1F("Delta", "Angle (Main axis - x-axis)", 100, -90, 90);
85 fDelta->SetDirectory(NULL);
86 fDelta->SetXTitle("Delta [\\circ]");
87 fDelta->SetYTitle("Counts");
88
89 MBinning bins;
90 bins.SetEdgesLog(50, 1, 1e7);
91
92 fSize = new TH1F;
93 fSize->SetName("Size");
94 fSize->SetTitle("Number of Photons");
95 fSize->SetDirectory(NULL);
96 fSize->SetXTitle("Size");
97 fSize->SetYTitle("Counts");
98 fSize->GetXaxis()->SetTitleOffset(1.2);
99 fSize->GetXaxis()->SetLabelOffset(-0.015);
100
101 MH::SetBinning(fSize, &bins);
102
103 fCenter = new TH2F("Center", "Center of Ellipse", 60, -600, 600, 60, -600, 600);
104 fCenter->SetDirectory(NULL);
105 fCenter->SetXTitle("x [mm]");
106 fCenter->SetYTitle("y [mm]");
107 fCenter->SetZTitle("Counts");
108}
109
110// --------------------------------------------------------------------------
111//
112// Delete the histograms
113//
114MHHillas::~MHHillas()
115{
116 delete fLength;
117 delete fWidth;
118
119 delete fDistC;
120 delete fDelta;
121
122 delete fSize;
123 delete fCenter;
124}
125
126// --------------------------------------------------------------------------
127//
128// Setup the Binning for the histograms automatically if the correct
129// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
130// are found in the parameter list
131// Use this function if you want to set the conversion factor which
132// is used to convert the mm-scale in the camera plain into the deg-scale
133// used for histogram presentations. The conversion factor is part of
134// the camera geometry. Please create a corresponding MGeomCam container.
135//
136Bool_t MHHillas::SetupFill(const MParList *plist)
137{
138 const MBinning *binsw = (MBinning*)plist->FindObject("BinningWidth");
139 const MBinning *binsl = (MBinning*)plist->FindObject("BinningLength");
140 const MBinning *binsd = (MBinning*)plist->FindObject("BinningDist");
141 const MBinning *binsc = (MBinning*)plist->FindObject("BinningCamera");
142
143 if (!binsw || !binsl || !binsd || !binsc)
144 {
145 *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl;
146 return kFALSE;
147 }
148
149 SetBinning(fWidth, binsw);
150 SetBinning(fLength, binsl);
151 SetBinning(fDistC, binsd);
152 SetBinning(fCenter, binsc, binsc);
153
154 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
155 if (!geom)
156 {
157 *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
158 return kTRUE;
159 }
160
161 fMm2Deg = geom->GetConvMm2Deg();
162 fUseMmScale = kFALSE;
163
164 fLength->SetXTitle("Length [\\circ]");
165 fWidth->SetXTitle("Width [\\circ]");
166 fDistC->SetXTitle("Distance [\\circ]");
167 fCenter->SetXTitle("x [\\circ]");
168 fCenter->SetYTitle("y [\\circ]");
169
170 return kTRUE;
171}
172
173// --------------------------------------------------------------------------
174//
175// Use this function to setup your own conversion factor between degrees
176// and millimeters. The conversion factor should be the one calculated in
177// MGeomCam. Use this function with Caution: You could create wrong values
178// by setting up your own scale factor.
179//
180void MHHillas::SetMm2Deg(Float_t mmdeg)
181{
182 if (mmdeg<0)
183 {
184 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
185 return;
186 }
187
188 if (fMm2Deg>=0)
189 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
190
191 fMm2Deg = mmdeg;
192}
193
194// --------------------------------------------------------------------------
195//
196// With this function you can convert the histogram ('on the fly') between
197// degrees and millimeters.
198//
199void MHHillas::SetMmScale(Bool_t mmscale)
200{
201 if (fUseMmScale == mmscale)
202 return;
203
204 if (fMm2Deg<0)
205 {
206 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
207 return;
208 }
209
210 if (fUseMmScale)
211 {
212 MH::ScaleAxis(fLength, 1./fMm2Deg);
213 MH::ScaleAxis(fWidth, 1./fMm2Deg);
214 MH::ScaleAxis(fDistC, 1./fMm2Deg);
215 MH::ScaleAxis(fCenter, 1./fMm2Deg, 1./fMm2Deg);
216
217 fLength->SetXTitle("Length [mm]");
218 fWidth->SetXTitle("Width [mm]");
219 fCenter->SetXTitle("x [mm]");
220 fCenter->SetYTitle("y [mm]");
221 }
222 else
223 {
224 MH::ScaleAxis(fLength, fMm2Deg);
225 MH::ScaleAxis(fWidth, fMm2Deg);
226 MH::ScaleAxis(fDistC, fMm2Deg);
227 MH::ScaleAxis(fCenter, fMm2Deg, fMm2Deg);
228
229 fLength->SetXTitle("Length [\\circ]");
230 fWidth->SetXTitle("Width [\\circ]");
231 fDistC->SetXTitle("Distance [\\circ]");
232 fCenter->SetXTitle("x [\\circ]");
233 fCenter->SetYTitle("y [\\circ]");
234 }
235
236 fUseMmScale = mmscale;
237}
238
239// --------------------------------------------------------------------------
240//
241// Fill the histograms with data from a MHillas-Container.
242// Be careful: Only call this with an object of type MHillas
243//
244Bool_t MHHillas::Fill(const MParContainer *par)
245{
246 const MHillas &h = *(MHillas*)par;
247
248 const Double_t d = sqrt(h.GetMeanX()*h.GetMeanX() + h.GetMeanY()*h.GetMeanY());
249
250 if (fUseMmScale)
251 {
252 fLength->Fill(h.GetLength());
253 fWidth ->Fill(h.GetWidth());
254 fDistC ->Fill(d);
255 fCenter->Fill(h.GetMeanX(), h.GetMeanY());
256 }
257 else
258 {
259 fLength->Fill(fMm2Deg*h.GetLength());
260 fWidth ->Fill(fMm2Deg*h.GetWidth());
261 fDistC ->Fill(fMm2Deg*d);
262 fCenter->Fill(fMm2Deg*h.GetMeanX(), fMm2Deg*h.GetMeanY());
263 }
264
265 fDelta->Fill(kRad2Deg*h.GetDelta());
266 fSize->Fill(h.GetSize());
267
268 return kTRUE;
269}
270
271// --------------------------------------------------------------------------
272//
273// Setup a inversed deep blue sea palette for the fCenter histogram.
274//
275void MHHillas::SetColors() const
276{
277 gStyle->SetPalette(51, NULL);
278 Int_t c[50];
279 for (int i=0; i<50; i++)
280 c[49-i] = gStyle->GetColorPalette(i);
281 gStyle->SetPalette(50, c);
282}
283
284// --------------------------------------------------------------------------
285//
286// Draw clones of four histograms. So that the object can be deleted
287// and the histograms are still visible in the canvas.
288// The cloned object are deleted together with the canvas if the canvas is
289// destroyed. If you want to handle dostroying the canvas you can get a
290// pointer to it from this function
291//
292TObject *MHHillas::DrawClone(Option_t *opt) const
293{
294 TCanvas *c = MakeDefCanvas("Hillas", fTitle, 700, 750);
295 c->Divide(2,3);
296
297 gROOT->SetSelectedPad(NULL);
298
299 c->cd(1);
300 fLength->DrawCopy();
301
302 c->cd(2);
303 fWidth->DrawCopy();
304
305 c->cd(3);
306 gPad->SetLogx();
307 fSize->DrawCopy();
308
309 c->cd(4);
310 fDelta->DrawCopy();
311
312 c->cd(5);
313 fDistC->DrawCopy();
314
315 c->cd(6);
316 SetColors();
317 fCenter->DrawCopy("colz");
318
319 c->Modified();
320 c->Update();
321
322 return c;
323}
324
325// --------------------------------------------------------------------------
326//
327// Creates a new canvas and draws the four histograms into it.
328// Be careful: The histograms belongs to this object and won't get deleted
329// together with the canvas.
330//
331void MHHillas::Draw(Option_t *)
332{
333 if (!gPad)
334 MakeDefCanvas("Hillas", fTitle, 700, 750);
335
336 gPad->Divide(2,3);
337
338 gPad->cd(1);
339 fLength->Draw();
340
341 gPad->cd(2);
342 fWidth->Draw();
343
344 gPad->cd(3);
345 gPad->SetLogx();
346 fSize->Draw();
347
348 gPad->cd(4);
349 fDelta->Draw();
350
351 gPad->cd(5);
352 fDistC->Draw();
353
354 gPad->cd(6);
355 SetColors();
356 fCenter->Draw("colz");
357
358 gPad->Modified();
359 gPad->Update();
360}
Note: See TracBrowser for help on using the repository browser.