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

Last change on this file since 1463 was 1463, checked in by tbretz, 22 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 9.0 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(kTRUE)
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, 296.7);
70 fWidth = new TH1F("Width", "Width of Ellipse", 100, 0, 296.7);
71 fDistC = new TH1F("DistC", "Distance from center of camera", 100, 0, 445);
72 fDelta = new TH1F("Delta", "Angle (Main axis - x-axis)", 101, -90, 90);
73
74 fLength->SetDirectory(NULL);
75 fWidth->SetDirectory(NULL);
76 fDistC->SetDirectory(NULL);
77 fDelta->SetDirectory(NULL);
78
79 fLength->SetXTitle("Length [mm]");
80 fWidth->SetXTitle("Width [mm]");
81 fDistC->SetXTitle("Distance [mm]");
82 fDelta->SetXTitle("Delta [\\circ]");
83
84 fLength->SetYTitle("Counts");
85 fWidth->SetYTitle("Counts");
86 fDistC->SetYTitle("Counts");
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", 51, -445, 445, 51, -445, 445);
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 MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
139 if (!geom)
140 *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
141 else
142 {
143 fMm2Deg = geom->GetConvMm2Deg();
144 SetMmScale(kFALSE);
145 }
146
147 ApplyBinning(*plist, "Width", fWidth);
148 ApplyBinning(*plist, "Length", fLength);
149 ApplyBinning(*plist, "Dist", fDistC);
150 ApplyBinning(*plist, "Camera", fCenter);
151 ApplyBinning(*plist, "Delta", fDelta);
152 ApplyBinning(*plist, "Size", fSize);
153
154 return kTRUE;
155}
156
157// --------------------------------------------------------------------------
158//
159// Use this function to setup your own conversion factor between degrees
160// and millimeters. The conversion factor should be the one calculated in
161// MGeomCam. Use this function with Caution: You could create wrong values
162// by setting up your own scale factor.
163//
164void MHHillas::SetMm2Deg(Float_t mmdeg)
165{
166 if (mmdeg<0)
167 {
168 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
169 return;
170 }
171
172 if (fMm2Deg>=0)
173 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
174
175 fMm2Deg = mmdeg;
176}
177
178// --------------------------------------------------------------------------
179//
180// With this function you can convert the histogram ('on the fly') between
181// degrees and millimeters.
182//
183void MHHillas::SetMmScale(Bool_t mmscale)
184{
185 if (fUseMmScale == mmscale)
186 return;
187
188 if (fMm2Deg<0)
189 {
190 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
191 return;
192 }
193
194 const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
195 MH::ScaleAxis(fLength, scale);
196 MH::ScaleAxis(fWidth, scale);
197 MH::ScaleAxis(fDistC, scale);
198 MH::ScaleAxis(fCenter, scale, scale);
199
200 if (mmscale)
201 {
202 fLength->SetXTitle("Length [mm]");
203 fWidth->SetXTitle("Width [mm]");
204 fDistC->SetXTitle("Distance [mm]");
205 fCenter->SetXTitle("x [mm]");
206 fCenter->SetYTitle("y [mm]");
207 }
208 else
209 {
210 fLength->SetXTitle("Length [\\circ]");
211 fWidth->SetXTitle("Width [\\circ]");
212 fDistC->SetXTitle("Distance [\\circ]");
213 fCenter->SetXTitle("x [\\circ]");
214 fCenter->SetYTitle("y [\\circ]");
215 }
216
217 fUseMmScale = mmscale;
218}
219
220// --------------------------------------------------------------------------
221//
222// Fill the histograms with data from a MHillas-Container.
223// Be careful: Only call this with an object of type MHillas
224//
225Bool_t MHHillas::Fill(const MParContainer *par)
226{
227 const MHillas &h = *(MHillas*)par;
228
229 const Double_t d = sqrt(h.GetMeanX()*h.GetMeanX() + h.GetMeanY()*h.GetMeanY());
230 const Double_t scale = fUseMmScale ? 1 : fMm2Deg;
231
232 fLength->Fill(scale*h.GetLength());
233 fWidth ->Fill(scale*h.GetWidth());
234 fDistC ->Fill(scale*d);
235 fCenter->Fill(scale*h.GetMeanX(), scale*h.GetMeanY());
236 fDelta ->Fill(kRad2Deg*h.GetDelta());
237 fSize ->Fill(h.GetSize());
238
239 return kTRUE;
240}
241
242// --------------------------------------------------------------------------
243//
244// Setup a inversed deep blue sea palette for the fCenter histogram.
245//
246void MHHillas::SetColors() const
247{
248 // FIXME: This must be redone each time the canvas is repainted....
249 gStyle->SetPalette(51, NULL);
250 Int_t c[50];
251 for (int i=0; i<50; i++)
252 c[49-i] = gStyle->GetColorPalette(i);
253 gStyle->SetPalette(50, c);
254}
255
256// --------------------------------------------------------------------------
257//
258// Draw clones of four histograms. So that the object can be deleted
259// and the histograms are still visible in the canvas.
260// The cloned object are deleted together with the canvas if the canvas is
261// destroyed. If you want to handle dostroying the canvas you can get a
262// pointer to it from this function
263//
264TObject *MHHillas::DrawClone(Option_t *opt) const
265{
266 TCanvas *c = MakeDefCanvas("Hillas", fTitle, 720, 810);
267 c->Divide(2,3);
268
269 gROOT->SetSelectedPad(NULL);
270
271 c->cd(1);
272 fLength->DrawCopy();
273
274 c->cd(2);
275 fWidth->DrawCopy();
276
277 c->cd(3);
278 gPad->SetLogx();
279 fSize->DrawCopy();
280
281 c->cd(4);
282 fDelta->DrawCopy();
283
284 c->cd(5);
285 fDistC->DrawCopy();
286
287 c->cd(6);
288 SetColors();
289 fCenter->DrawCopy("colz");
290
291 c->Modified();
292 c->Update();
293
294 return c;
295}
296
297// --------------------------------------------------------------------------
298//
299// Creates a new canvas and draws the four histograms into it.
300// Be careful: The histograms belongs to this object and won't get deleted
301// together with the canvas.
302//
303void MHHillas::Draw(Option_t *)
304{
305 if (!gPad)
306 MakeDefCanvas("Hillas", fTitle, 720, 810);
307
308 gPad->Divide(2,3);
309
310 gPad->cd(1);
311 fLength->Draw();
312
313 gPad->cd(2);
314 fWidth->Draw();
315
316 gPad->cd(3);
317 gPad->SetLogx();
318 fSize->Draw();
319
320 gPad->cd(4);
321 fDelta->Draw();
322
323 gPad->cd(5);
324 fDistC->Draw();
325
326 gPad->cd(6);
327 SetColors();
328 fCenter->Draw("colz");
329
330 gPad->Modified();
331 gPad->Update();
332}
Note: See TracBrowser for help on using the repository browser.