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 12/2000 <mailto:tbretz@astro.uni-wuerzburg.de>
|
---|
19 | ! Author(s): Wolfgang Wittek 02/2004 <mailto:wittek@mppmu.mpg.de>
|
---|
20 | !
|
---|
21 | ! Copyright: MAGIC Software Development, 2000-2004
|
---|
22 | !
|
---|
23 | !
|
---|
24 | \* ======================================================================== */
|
---|
25 |
|
---|
26 | /////////////////////////////////////////////////////////////////////////////
|
---|
27 | //
|
---|
28 | // MHStarMap
|
---|
29 | //
|
---|
30 | // This class contains a 2-dimensional histogram. It should show some
|
---|
31 | // kind of star map. The algorith which calculates the star map
|
---|
32 | // from the Hillas parameters (Fill) can be enhanced.
|
---|
33 | //
|
---|
34 | // For a given a shower, a series of points along its main axis are filled
|
---|
35 | // into the 2-dim histogram (x, y) of the camera plane.
|
---|
36 | //
|
---|
37 | // Before filling a point (x, y) into the histogram it is
|
---|
38 | // - shifted by (xSrc, ySrc) (the expected source position)
|
---|
39 | // - and rotated in order to compensate the rotation of the
|
---|
40 | // sky image in the camera
|
---|
41 | //
|
---|
42 | //
|
---|
43 | /////////////////////////////////////////////////////////////////////////////
|
---|
44 | #include "MHStarMap.h"
|
---|
45 |
|
---|
46 | #include <TH2.h>
|
---|
47 | #include <TStyle.h> // gStyle
|
---|
48 | #include <TColor.h> // SetRGB
|
---|
49 | #include <TCanvas.h> // TCanvas
|
---|
50 |
|
---|
51 | #include "MLog.h"
|
---|
52 | #include "MLogManip.h"
|
---|
53 |
|
---|
54 | #include "MParList.h"
|
---|
55 |
|
---|
56 | #include "MAstro.h"
|
---|
57 | #include "MGeomCam.h"
|
---|
58 | #include "MHillas.h"
|
---|
59 | #include "MBinning.h"
|
---|
60 | #include "MMcEvt.hxx"
|
---|
61 | #include "MSrcPosCam.h"
|
---|
62 | #include "MObservatory.h"
|
---|
63 |
|
---|
64 | ClassImp(MHStarMap);
|
---|
65 |
|
---|
66 | using namespace std;
|
---|
67 |
|
---|
68 | // --------------------------------------------------------------------------
|
---|
69 | //
|
---|
70 | // Create the star map histogram
|
---|
71 | //
|
---|
72 | MHStarMap::MHStarMap(const char *name, const char *title)
|
---|
73 | : fMm2Deg(-1), fUseMmScale(kTRUE)
|
---|
74 | {
|
---|
75 | //
|
---|
76 | // default constructor
|
---|
77 | // creates an a list of histograms for all pixels and both gain channels
|
---|
78 | //
|
---|
79 |
|
---|
80 | //
|
---|
81 | // set the name and title of this object
|
---|
82 | //
|
---|
83 | fName = name ? name : "MHStarMap" ;
|
---|
84 | fTitle = title ? title : "Container for a Star Map" ;
|
---|
85 |
|
---|
86 | //
|
---|
87 | // loop over all Pixels and create two histograms
|
---|
88 | // one for the Low and one for the High gain
|
---|
89 | // connect all the histogram with the container fHist
|
---|
90 | //
|
---|
91 | fStarMap = new TH2F("StarMap", " 2D Hillas Star Map ",
|
---|
92 | 150, -300, 300,
|
---|
93 | 150, -300, 300);
|
---|
94 |
|
---|
95 | fStarMap->SetDirectory(NULL);
|
---|
96 |
|
---|
97 | fStarMap->SetXTitle("x [mm]");
|
---|
98 | fStarMap->SetYTitle("y [mm]");
|
---|
99 | fStarMap->SetZTitle("Counts");
|
---|
100 | }
|
---|
101 |
|
---|
102 | // --------------------------------------------------------------------------
|
---|
103 | //
|
---|
104 | // delete the histogram instance
|
---|
105 | //
|
---|
106 | MHStarMap::~MHStarMap()
|
---|
107 | {
|
---|
108 | delete fStarMap;
|
---|
109 | }
|
---|
110 |
|
---|
111 | // --------------------------------------------------------------------------
|
---|
112 | //
|
---|
113 | // Setup the Binning for the histograms automatically if the correct
|
---|
114 | // instances of MBinning (with the names 'BinningCamera')
|
---|
115 | // are found in the parameter list
|
---|
116 | // Use this function if you want to set the conversion factor which
|
---|
117 | // is used to convert the mm-scale in the camera plain into the deg-scale
|
---|
118 | // used for histogram presentations. The conversion factor is part of
|
---|
119 | // the camera geometry. Please create a corresponding MGeomCam container.
|
---|
120 | //
|
---|
121 | Bool_t MHStarMap::SetupFill(const MParList *plist)
|
---|
122 | {
|
---|
123 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
|
---|
124 | if (geom)
|
---|
125 | {
|
---|
126 | fMm2Deg = geom->GetConvMm2Deg();
|
---|
127 | fUseMmScale = kFALSE;
|
---|
128 |
|
---|
129 | fStarMap->SetXTitle("x [\\circ]");
|
---|
130 | fStarMap->SetYTitle("y [\\circ]");
|
---|
131 | }
|
---|
132 |
|
---|
133 | fMcEvt = (MMcEvt*)plist->FindObject("MMcEvt");
|
---|
134 | if (!fMcEvt)
|
---|
135 | {
|
---|
136 | *fLog << err << "MMcEvt not found... aborting." << endl;
|
---|
137 | return kFALSE;
|
---|
138 | }
|
---|
139 |
|
---|
140 | fSrcPos = (MSrcPosCam*)plist->FindObject(AddSerialNumber("MSrcPosCam"));
|
---|
141 | if (!fSrcPos)
|
---|
142 | {
|
---|
143 | *fLog << err << "MSrcPosCam not found... aborting" << endl;
|
---|
144 | return kFALSE;
|
---|
145 | }
|
---|
146 |
|
---|
147 |
|
---|
148 | const MBinning *bins = (MBinning*)plist->FindObject("BinningStarMap");
|
---|
149 | if (!bins)
|
---|
150 | {
|
---|
151 | float r = geom ? geom->GetMaxRadius() : 600;
|
---|
152 | r *= 5./6;
|
---|
153 | if (!fUseMmScale)
|
---|
154 | r *= fMm2Deg;
|
---|
155 |
|
---|
156 | MBinning b;
|
---|
157 | b.SetEdges(100, -r, r);
|
---|
158 | SetBinning(fStarMap, &b, &b);
|
---|
159 | }
|
---|
160 | else
|
---|
161 | SetBinning(fStarMap, bins, bins);
|
---|
162 |
|
---|
163 | if (!geom)
|
---|
164 | {
|
---|
165 | *fLog << warn << "No Camera Geometry available. Using mm-scale for histograms." << endl;
|
---|
166 | return kTRUE;
|
---|
167 | }
|
---|
168 |
|
---|
169 | *fLog << warn << "WARNING - Using MHStarMap doesn't take care of the Source Position!" << endl;
|
---|
170 |
|
---|
171 | return kTRUE;
|
---|
172 | }
|
---|
173 |
|
---|
174 | // --------------------------------------------------------------------------
|
---|
175 | //
|
---|
176 | // Get the observatory location "MObservatory" from parameter list
|
---|
177 | //
|
---|
178 | Bool_t MHStarMap::ReInit(MParList *pList)
|
---|
179 | {
|
---|
180 | fObservatory = (MObservatory*)pList->FindObject(AddSerialNumber("MObservatory"));
|
---|
181 | if (!fObservatory)
|
---|
182 | {
|
---|
183 | *fLog << err << "MObservatory not found... aborting" << endl;
|
---|
184 | return kFALSE;
|
---|
185 | }
|
---|
186 |
|
---|
187 | return kTRUE;
|
---|
188 | }
|
---|
189 |
|
---|
190 | // --------------------------------------------------------------------------
|
---|
191 | //
|
---|
192 | // Calls MObservatory::RotationAngle
|
---|
193 | //
|
---|
194 | void MHStarMap::GetRotationAngle(Double_t &sin, Double_t &cos)
|
---|
195 | {
|
---|
196 | fObservatory->RotationAngle(fMcEvt->GetTelescopeTheta(),
|
---|
197 | fMcEvt->GetTelescopePhi(), sin, cos);
|
---|
198 | }
|
---|
199 |
|
---|
200 | // --------------------------------------------------------------------------
|
---|
201 | //
|
---|
202 | // Fill the four histograms with data from a MHillas-Container.
|
---|
203 | // Be careful: Only call this with an object of type MHillas
|
---|
204 | //
|
---|
205 | Bool_t MHStarMap::Fill(const MParContainer *par, const Stat_t w)
|
---|
206 | {
|
---|
207 | const MHillas &h = *(MHillas*)par;
|
---|
208 |
|
---|
209 | const Float_t delta = h.GetDelta();
|
---|
210 |
|
---|
211 | const Float_t m = TMath::Tan(delta);
|
---|
212 |
|
---|
213 | const Float_t cosd = 1.0/TMath::Sqrt(m*m+1);
|
---|
214 | const Float_t sind = TMath::Sqrt(1-cosd*cosd);
|
---|
215 |
|
---|
216 | Float_t t = h.GetMeanY() - m*h.GetMeanX();
|
---|
217 |
|
---|
218 | Float_t xSource = fSrcPos->GetX();
|
---|
219 | Float_t ySource = fSrcPos->GetY();
|
---|
220 |
|
---|
221 | if (!fUseMmScale)
|
---|
222 | {
|
---|
223 | t *= fMm2Deg;
|
---|
224 | xSource *= fMm2Deg;
|
---|
225 | ySource *= fMm2Deg;
|
---|
226 | }
|
---|
227 |
|
---|
228 | // get step size ds along the main axis of the ellipse
|
---|
229 | const TAxis &axex = *fStarMap->GetXaxis();
|
---|
230 | const Float_t xmin = axex.GetXmin();
|
---|
231 | const Float_t xmax = axex.GetXmax();
|
---|
232 |
|
---|
233 | // FIXME: Fixed number?
|
---|
234 | const Float_t ds = (xmax-xmin) / 200.0;
|
---|
235 |
|
---|
236 | // Fill points along the main axis of the shower into the histogram;
|
---|
237 | // before filling
|
---|
238 | // - perform a translation by (xSource, ySource)
|
---|
239 | // - and perform a rotation to compensate the rotation of the
|
---|
240 | // sky image in the camera
|
---|
241 | Double_t cosal;
|
---|
242 | Double_t sinal;
|
---|
243 | GetRotationAngle(sinal, cosal);
|
---|
244 |
|
---|
245 | if (m>-1 && m<1)
|
---|
246 | {
|
---|
247 | const Float_t dx = ds * cosd;
|
---|
248 |
|
---|
249 | for (float x=xmin+dx/2; x<(xmax-xmin)+dx; x+=dx)
|
---|
250 | {
|
---|
251 | const Float_t xorig = x - xSource;
|
---|
252 | const Float_t yorig = m*x+t - ySource;
|
---|
253 |
|
---|
254 | const Float_t xfill = cosal*xorig + sinal*yorig;
|
---|
255 | const Float_t yfill = -sinal*xorig + cosal*yorig;
|
---|
256 | fStarMap->Fill(xfill, yfill, w);
|
---|
257 | }
|
---|
258 | }
|
---|
259 | else
|
---|
260 | {
|
---|
261 | const TAxis &axey = *fStarMap->GetYaxis();
|
---|
262 | const Float_t ymin = axey.GetXmin();
|
---|
263 | const Float_t ymax = axey.GetXmax();
|
---|
264 |
|
---|
265 | const float dy = ds * sind;
|
---|
266 |
|
---|
267 | for (float y=ymin+dy/2; y<(ymax-ymin)+dy; y+=dy)
|
---|
268 | {
|
---|
269 | const Float_t xorig = (y-t)/m - xSource;
|
---|
270 | const Float_t yorig = y - ySource;
|
---|
271 |
|
---|
272 | const Float_t xfill = cosal*xorig + sinal*yorig;
|
---|
273 | const Float_t yfill = -sinal*xorig + cosal*yorig;
|
---|
274 | fStarMap->Fill(xfill, yfill, w);
|
---|
275 | }
|
---|
276 | }
|
---|
277 |
|
---|
278 | return kTRUE;
|
---|
279 | }
|
---|
280 |
|
---|
281 | // --------------------------------------------------------------------------
|
---|
282 | //
|
---|
283 | // Use this function to setup your own conversion factor between degrees
|
---|
284 | // and millimeters. The conversion factor should be the one calculated in
|
---|
285 | // MGeomCam. Use this function with Caution: You could create wrong values
|
---|
286 | // by setting up your own scale factor.
|
---|
287 | //
|
---|
288 | void MHStarMap::SetMm2Deg(Float_t mmdeg)
|
---|
289 | {
|
---|
290 | if (mmdeg<0)
|
---|
291 | {
|
---|
292 | *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
|
---|
293 | return;
|
---|
294 | }
|
---|
295 |
|
---|
296 | if (fMm2Deg>=0)
|
---|
297 | *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
|
---|
298 |
|
---|
299 | fMm2Deg = mmdeg;
|
---|
300 | }
|
---|
301 |
|
---|
302 | // --------------------------------------------------------------------------
|
---|
303 | //
|
---|
304 | // With this function you can convert the histogram ('on the fly') between
|
---|
305 | // degrees and millimeters.
|
---|
306 | //
|
---|
307 | void MHStarMap::SetMmScale(Bool_t mmscale)
|
---|
308 | {
|
---|
309 | if (fUseMmScale == mmscale)
|
---|
310 | return;
|
---|
311 |
|
---|
312 | if (fMm2Deg<0)
|
---|
313 | {
|
---|
314 | *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
|
---|
315 | return;
|
---|
316 | }
|
---|
317 |
|
---|
318 | if (fUseMmScale)
|
---|
319 | {
|
---|
320 | fStarMap->SetXTitle("x [mm]");
|
---|
321 | fStarMap->SetYTitle("y [mm]");
|
---|
322 |
|
---|
323 | fStarMap->Scale(1./fMm2Deg);
|
---|
324 | }
|
---|
325 | else
|
---|
326 | {
|
---|
327 | fStarMap->SetXTitle("x [\\circ]");
|
---|
328 | fStarMap->SetYTitle("y [\\circ]");
|
---|
329 |
|
---|
330 | fStarMap->Scale(1./fMm2Deg);
|
---|
331 | }
|
---|
332 |
|
---|
333 | fUseMmScale = mmscale;
|
---|
334 | }
|
---|
335 |
|
---|
336 | // --------------------------------------------------------------------------
|
---|
337 | //
|
---|
338 | // Set the palette you wanna use:
|
---|
339 | // - you could set the root "Pretty Palette Violet->Red" by
|
---|
340 | // gStyle->SetPalette(1, 0), but in some cases this may look
|
---|
341 | // confusing
|
---|
342 | // - The maximum colors root allowes us to set by ourself
|
---|
343 | // is 50 (idx: 51-100). This colors are set to a grayscaled
|
---|
344 | // palette
|
---|
345 | // - the number of contours must be two less than the number
|
---|
346 | // of palette entries
|
---|
347 | //
|
---|
348 | void MHStarMap::PrepareDrawing() const
|
---|
349 | {
|
---|
350 | const Int_t numg = 32; // number of gray scaled colors
|
---|
351 | const Int_t numw = 32; // number of white
|
---|
352 |
|
---|
353 | Int_t palette[numg+numw];
|
---|
354 |
|
---|
355 | //
|
---|
356 | // The first half of the colors are white.
|
---|
357 | // This is some kind of optical background supression
|
---|
358 | //
|
---|
359 | gROOT->GetColor(51)->SetRGB(1, 1, 1);
|
---|
360 |
|
---|
361 | Int_t i;
|
---|
362 | for (i=0; i<numw; i++)
|
---|
363 | palette[i] = 51;
|
---|
364 |
|
---|
365 | //
|
---|
366 | // now the (gray) scaled part is coming
|
---|
367 | //
|
---|
368 | for (;i<numw+numg; i++)
|
---|
369 | {
|
---|
370 | const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
|
---|
371 |
|
---|
372 | gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
|
---|
373 | palette[i] = 52+i;
|
---|
374 | }
|
---|
375 |
|
---|
376 | //
|
---|
377 | // Set the palette and the number of contour levels
|
---|
378 | //
|
---|
379 | gStyle->SetPalette(numg+numw, palette);
|
---|
380 | fStarMap->SetContour(numg+numw-2);
|
---|
381 | }
|
---|
382 |
|
---|
383 | // --------------------------------------------------------------------------
|
---|
384 | //
|
---|
385 | // Draw clones of the histograms, so that the object can be deleted
|
---|
386 | // and the histogram is still visible in the canvas.
|
---|
387 | // The cloned object is deleted together with the canvas if the canvas is
|
---|
388 | // destroyed. If you want to handle destroying the canvas you can get a
|
---|
389 | // pointer to it from this function
|
---|
390 | //
|
---|
391 | TObject *MHStarMap::DrawClone(Option_t *opt) const
|
---|
392 | {
|
---|
393 | return MH::DrawClone(opt, 500, 500);
|
---|
394 | }
|
---|
395 |
|
---|
396 | // --------------------------------------------------------------------------
|
---|
397 | //
|
---|
398 | // Creates a new canvas and draw the histogram into it.
|
---|
399 | // Be careful: The histogram belongs to this object and won't get deleted
|
---|
400 | // together with the canvas.
|
---|
401 | //
|
---|
402 | void MHStarMap::Draw(Option_t *)
|
---|
403 | {
|
---|
404 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this, 500, 500);
|
---|
405 | pad->SetBorderMode(0);
|
---|
406 |
|
---|
407 | pad->Divide(1,1);
|
---|
408 |
|
---|
409 | pad->cd(1);
|
---|
410 | gPad->SetBorderMode(0);
|
---|
411 |
|
---|
412 | AppendPad("");
|
---|
413 | }
|
---|
414 |
|
---|
415 | void MHStarMap::Paint(Option_t *opt)
|
---|
416 | {
|
---|
417 | //
|
---|
418 | // Maintain aspect ratio
|
---|
419 | //
|
---|
420 | const float w = gPad->GetWw();
|
---|
421 | const float h = gPad->GetWh();
|
---|
422 |
|
---|
423 | if (h<w)
|
---|
424 | gPad->SetPad((1.-h/w)/2, 0, (h/w+1)/2, 1);
|
---|
425 | else
|
---|
426 | gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1)/2);
|
---|
427 |
|
---|
428 | //
|
---|
429 | // Maintain colors
|
---|
430 | //
|
---|
431 | PrepareDrawing();
|
---|
432 |
|
---|
433 | //
|
---|
434 | // Paint Histogram
|
---|
435 | //
|
---|
436 | fStarMap->Paint("colz");
|
---|
437 | }
|
---|