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

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