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