source: trunk/Mars/mimage/MHHillas.cc@ 19876

Last change on this file since 19876 was 19679, checked in by tbretz, 5 years ago
Fixed a type in the units of the COG
File size: 9.1 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@astro.uni-wuerzburg.de>
19! Author(s): Wolfgang Wittek 2002 <mailto:wittek@mppmu.mpg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2009
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MHHillas
29//
30// This class contains histograms for the source independent image parameters
31//
32// ClassVersion 2:
33// ---------------
34// - fMm2Deg
35// - fUseMmScale
36//
37/////////////////////////////////////////////////////////////////////////////
38#include "MHHillas.h"
39
40#include <math.h>
41
42#include <TH1.h>
43#include <TH2.h>
44#include <TPad.h>
45#include <TStyle.h>
46#include <TCanvas.h>
47
48#include "MLog.h"
49#include "MLogManip.h"
50
51#include "MParList.h"
52
53#include "MHillas.h"
54#include "MGeomCam.h"
55#include "MBinning.h"
56
57#include "MHCamera.h"
58
59ClassImp(MHHillas);
60
61using namespace std;
62
63// --------------------------------------------------------------------------
64//
65// Setup four histograms for Width, Length
66//
67MHHillas::MHHillas(const char *name, const char *title)
68 : fGeomCam(0)
69{
70 //
71 // set the name and title of this object
72 //
73 fName = name ? name : "MHHillas";
74 fTitle = title ? title : "Source independent image parameters";
75
76 fLength = new TH1F("Length", "Length of Ellipse", 100, 0, 1.0);
77 fWidth = new TH1F("Width", "Width of Ellipse", 100, 0, 1.0);
78 fDistC = new TH1F("DistC", "Distance from center of camera", 100, 0, 1.5);
79 fDelta = new TH1F("Delta", "Angle (Main axis - x-axis)", 101, -90, 90);
80
81 fLength->SetLineColor(kBlue);
82
83 fLength->SetDirectory(NULL);
84 fWidth->SetDirectory(NULL);
85 fDistC->SetDirectory(NULL);
86 fDelta->SetDirectory(NULL);
87
88 fLength->SetXTitle("Length [#circ]");
89 fWidth->SetXTitle("Width [#circ]");
90 fDistC->SetXTitle("Distance [#circ]");
91 fDelta->SetXTitle("Delta [#circ]");
92
93 fLength->SetYTitle("Counts");
94 fWidth->SetYTitle("Counts");
95 fDistC->SetYTitle("Counts");
96 fDelta->SetYTitle("Counts");
97
98 fDelta->SetMinimum(0);
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 fSize->UseCurrentStyle();
113
114 bins.Apply(*fSize);
115
116 fCenter = new TH2F("Center", "Center of gravity", 51, -1.5, 1.5, 51, -1.5, 1.5);
117 fCenter->SetDirectory(NULL);
118 fCenter->SetXTitle("x [#circ]");
119 fCenter->SetYTitle("y [#circ]");
120 fCenter->SetZTitle("Counts");
121}
122
123// --------------------------------------------------------------------------
124//
125// Delete the histograms
126//
127MHHillas::~MHHillas()
128{
129 delete fLength;
130 delete fWidth;
131
132 delete fDistC;
133 delete fDelta;
134
135 delete fSize;
136 delete fCenter;
137}
138
139// --------------------------------------------------------------------------
140//
141// Setup the Binning for the histograms automatically if the correct
142// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
143// are found in the parameter list
144// Use this function if you want to set the conversion factor which
145// is used to convert the mm-scale in the camera plain into the deg-scale
146// used for histogram presentations. The conversion factor is part of
147// the camera geometry. Please create a corresponding MGeomCam container.
148//
149Bool_t MHHillas::SetupFill(const MParList *plist)
150{
151 fGeomCam = (MGeomCam*)plist->FindObject("MGeomCam");
152 if (!fGeomCam)
153 {
154 *fLog << err << "MGeomCam not found... abort." << endl;
155 return kFALSE;
156 }
157
158
159 ApplyBinning(*plist, "Width", *fWidth);
160 ApplyBinning(*plist, "Length", *fLength);
161 ApplyBinning(*plist, "Dist", *fDistC);
162 ApplyBinning(*plist, "Delta", *fDelta);
163 ApplyBinning(*plist, "Size", *fSize);
164
165 const MBinning *bins = (MBinning*)plist->FindObject("BinningCamera");
166 if (!bins)
167 {
168 const Float_t r = fGeomCam->GetMaxRadius()*fGeomCam->GetConvMm2Deg();
169
170 const MBinning b(61, -r, r);
171 SetBinning(*fCenter, b, b);
172 }
173 else
174 SetBinning(*fCenter, *bins, *bins);
175
176 return kTRUE;
177}
178
179// --------------------------------------------------------------------------
180//
181// Fill the histograms with data from a MHillas-Container.
182// Be careful: Only call this with an object of type MHillas
183//
184Int_t MHHillas::Fill(const MParContainer *par, const Stat_t w)
185{
186 if (!par)
187 {
188 *fLog << err << "MHHillas::Fill: Pointer (!=NULL) expected." << endl;
189 return kERROR;
190 }
191
192 const MHillas &h = *(MHillas*)par;
193
194 const Double_t scale = fGeomCam->GetConvMm2Deg();
195
196 fLength->Fill(scale*h.GetLength(), w);
197 fWidth ->Fill(scale*h.GetWidth(), w);
198 fDistC ->Fill(scale*h.GetDist0(), w);
199 fCenter->Fill(scale*h.GetMeanX(), scale*h.GetMeanY(), w);
200 fDelta ->Fill(kRad2Deg*h.GetDelta(), w);
201 fSize ->Fill(h.GetSize(), w);
202
203 return kTRUE;
204}
205
206// --------------------------------------------------------------------------
207//
208// Creates a new canvas and draws the four histograms into it.
209// Be careful: The histograms belongs to this object and won't get deleted
210// together with the canvas.
211//
212void MHHillas::Draw(Option_t *o)
213{
214 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
215 pad->SetBorderMode(0);
216
217 AppendPad("");
218
219 TString opt(o);
220 opt.ToLower();
221
222 // FIXME: If same-option given make two independant y-axis!
223 const Bool_t same = opt.Contains("same");
224
225 if (!same)
226 pad->Divide(2,3);
227 else
228 {
229 fLength->SetName("LengthSame");
230 fWidth->SetName("WidthSame");
231 fDistC->SetName("DistCSame");
232 fDelta->SetName("DeltaSame");
233 fSize->SetName("SizeSame");
234 fCenter->SetName("CenterSame");
235
236 fLength->SetDirectory(0);
237 fWidth->SetDirectory(0);
238 fDistC->SetDirectory(0);
239 fDelta->SetDirectory(0);
240 fSize->SetDirectory(0);
241 fCenter->SetDirectory(0);
242
243 fDistC->SetLineColor(kBlue);
244 fSize->SetLineColor(kBlue);
245 fDelta->SetLineColor(kBlue);
246 fWidth->SetLineColor(kMagenta);
247 fLength->SetLineColor(kCyan);
248 }
249
250 pad->cd(1);
251 gPad->SetBorderMode(0);
252 gPad->SetGridx();
253 gPad->SetGridy();
254 RemoveFromPad("WidthSame");
255 RemoveFromPad("LengthSame");
256 MH::DrawSame(*fWidth, *fLength, "Width'n'Length", same);
257
258 pad->cd(2);
259 gPad->SetBorderMode(0);
260 gPad->SetGridx();
261 gPad->SetGridy();
262 RemoveFromPad("DistCSame");
263 fDistC->Draw(same?"same":"");
264
265 pad->cd(3);
266 gPad->SetGridx();
267 gPad->SetGridy();
268 gPad->SetBorderMode(0);
269 gPad->SetLogx();
270 gPad->SetLogy();
271 RemoveFromPad("SizeSame");
272 fSize->Draw(same?"same":"");
273
274 pad->cd(4);
275 gPad->SetBorderMode(0);
276 gPad->SetGridx();
277 gPad->SetGridy();
278 gPad->SetPad(0.51, 0.01, 0.99, 0.65);
279 if (same)
280 {
281 TH2 *h=dynamic_cast<TH2*>(gPad->FindObject("Center"));
282 if (h)
283 {
284 // This causes crashes in THistPainter::PaintTable
285 // if the z-axis is not kept. No idea why...
286 h->SetDrawOption("z");
287 h->SetMarkerColor(kBlack);
288 }
289
290 RemoveFromPad("CenterSame");
291 fCenter->SetMarkerColor(kBlue);
292 fCenter->Draw("same");
293 }
294 else
295 fCenter->Draw("colz");
296
297 if (fGeomCam)
298 {
299// MHCamera *cam = new MHCamera(*fGeomCam);
300// cam->Draw("same");
301// cam->SetBit(kCanDelete);
302 }
303
304 pad->cd(5);
305 gPad->SetBorderMode(0);
306 gPad->SetGridx();
307 gPad->SetGridy();
308 RemoveFromPad("DeltaSame");
309 fDelta->Draw(same?"same":"");
310
311 pad->cd(6);
312 if (gPad && !same)
313 delete gPad;
314}
315
316TH1 *MHHillas::GetHistByName(const TString name) const
317{
318 if (name.Contains("Width", TString::kIgnoreCase))
319 return fWidth;
320 if (name.Contains("Length", TString::kIgnoreCase))
321 return fLength;
322 if (name.Contains("Size", TString::kIgnoreCase))
323 return fSize;
324 if (name.Contains("Delta", TString::kIgnoreCase))
325 return fDelta;
326 if (name.Contains("DistC", TString::kIgnoreCase))
327 return fDistC;
328 if (name.Contains("Center", TString::kIgnoreCase))
329 return fCenter;
330
331 return NULL;
332}
333
334void MHHillas::Paint(Option_t *opt)
335{
336 MH::SetPalette("pretty");
337 MH::Paint();
338}
Note: See TracBrowser for help on using the repository browser.