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