source: trunk/MagicSoft/Mars/mimage/MHHillas.cc@ 1985

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