source: trunk/MagicSoft/Mars/mhist/MHStarMap.cc@ 3622

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