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

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