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

Last change on this file since 9488 was 9350, checked in by tbretz, 16 years ago
*** empty log message ***
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 [mm]");
119 fCenter->SetYTitle("y [mm]");
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 MBinning b;
171 b.SetEdges(61, -r, r);
172 SetBinning(fCenter, &b, &b);
173 }
174 else
175 SetBinning(fCenter, bins, bins);
176
177
178 return kTRUE;
179}
180
181// --------------------------------------------------------------------------
182//
183// Fill the histograms with data from a MHillas-Container.
184// Be careful: Only call this with an object of type MHillas
185//
186Int_t MHHillas::Fill(const MParContainer *par, const Stat_t w)
187{
188 if (!par)
189 {
190 *fLog << err << "MHHillas::Fill: Pointer (!=NULL) expected." << endl;
191 return kERROR;
192 }
193
194 const MHillas &h = *(MHillas*)par;
195
196 const Double_t scale = fGeomCam->GetConvMm2Deg();
197
198 fLength->Fill(scale*h.GetLength(), w);
199 fWidth ->Fill(scale*h.GetWidth(), w);
200 fDistC ->Fill(scale*h.GetDist0(), w);
201 fCenter->Fill(scale*h.GetMeanX(), scale*h.GetMeanY(), w);
202 fDelta ->Fill(kRad2Deg*h.GetDelta(), w);
203 fSize ->Fill(h.GetSize(), w);
204
205 return kTRUE;
206}
207
208// --------------------------------------------------------------------------
209//
210// Creates a new canvas and draws the four histograms into it.
211// Be careful: The histograms belongs to this object and won't get deleted
212// together with the canvas.
213//
214void MHHillas::Draw(Option_t *o)
215{
216 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
217 pad->SetBorderMode(0);
218
219 AppendPad("");
220
221 TString opt(o);
222 opt.ToLower();
223
224 // FIXME: If same-option given make two independant y-axis!
225 const Bool_t same = opt.Contains("same");
226
227 if (!same)
228 pad->Divide(2,3);
229 else
230 {
231 fLength->SetName("LengthSame");
232 fWidth->SetName("WidthSame");
233 fDistC->SetName("DistCSame");
234 fDelta->SetName("DeltaSame");
235 fSize->SetName("SizeSame");
236 fCenter->SetName("CenterSame");
237
238 fLength->SetDirectory(0);
239 fWidth->SetDirectory(0);
240 fDistC->SetDirectory(0);
241 fDelta->SetDirectory(0);
242 fSize->SetDirectory(0);
243 fCenter->SetDirectory(0);
244
245 fDistC->SetLineColor(kBlue);
246 fSize->SetLineColor(kBlue);
247 fDelta->SetLineColor(kBlue);
248 fWidth->SetLineColor(kMagenta);
249 fLength->SetLineColor(kCyan);
250 }
251
252 pad->cd(1);
253 gPad->SetBorderMode(0);
254 gPad->SetGridx();
255 gPad->SetGridy();
256 RemoveFromPad("WidthSame");
257 RemoveFromPad("LengthSame");
258 MH::DrawSame(*fWidth, *fLength, "Width'n'Length", same);
259
260 pad->cd(2);
261 gPad->SetBorderMode(0);
262 gPad->SetGridx();
263 gPad->SetGridy();
264 RemoveFromPad("DistCSame");
265 fDistC->Draw(same?"same":"");
266
267 pad->cd(3);
268 gPad->SetGridx();
269 gPad->SetGridy();
270 gPad->SetBorderMode(0);
271 gPad->SetLogx();
272 gPad->SetLogy();
273 RemoveFromPad("SizeSame");
274 fSize->Draw(same?"same":"");
275
276 pad->cd(4);
277 gPad->SetBorderMode(0);
278 gPad->SetGridx();
279 gPad->SetGridy();
280 gPad->SetPad(0.51, 0.01, 0.99, 0.65);
281 if (same)
282 {
283 TH2 *h=dynamic_cast<TH2*>(gPad->FindObject("Center"));
284 if (h)
285 {
286 // This causes crashes in THistPainter::PaintTable
287 // if the z-axis is not kept. No idea why...
288 h->SetDrawOption("z");
289 h->SetMarkerColor(kBlack);
290 }
291
292 RemoveFromPad("CenterSame");
293 fCenter->SetMarkerColor(kBlue);
294 fCenter->Draw("same");
295 }
296 else
297 fCenter->Draw("colz");
298
299 if (fGeomCam)
300 {
301 MHCamera *cam = new MHCamera(*fGeomCam);
302 cam->Draw("same");
303 cam->SetBit(kCanDelete);
304 }
305
306 pad->cd(5);
307 gPad->SetBorderMode(0);
308 gPad->SetGridx();
309 gPad->SetGridy();
310 RemoveFromPad("DeltaSame");
311 fDelta->Draw(same?"same":"");
312
313 pad->cd(6);
314 if (gPad && !same)
315 delete gPad;
316}
317
318TH1 *MHHillas::GetHistByName(const TString name) const
319{
320 if (name.Contains("Width", TString::kIgnoreCase))
321 return fWidth;
322 if (name.Contains("Length", TString::kIgnoreCase))
323 return fLength;
324 if (name.Contains("Size", TString::kIgnoreCase))
325 return fSize;
326 if (name.Contains("Delta", TString::kIgnoreCase))
327 return fDelta;
328 if (name.Contains("DistC", TString::kIgnoreCase))
329 return fDistC;
330 if (name.Contains("Center", TString::kIgnoreCase))
331 return fCenter;
332
333 return NULL;
334}
335
336void MHHillas::Paint(Option_t *opt)
337{
338 MH::SetPalette("pretty");
339 MH::Paint();
340}
Note: See TracBrowser for help on using the repository browser.