source: trunk/Mars/mhist/MHStarMap.cc

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