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

Last change on this file since 9439 was 9153, checked in by tbretz, 16 years ago
*** empty log message ***
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 MBinning b;
156 b.SetEdges(100, -r, r);
157 SetBinning(fStarMap, &b, &b);
158 }
159 else
160 SetBinning(fStarMap, bins, bins);
161
162 if (!geom)
163 {
164 *fLog << warn << "No Camera Geometry available. Using mm-scale for histograms." << endl;
165 return kTRUE;
166 }
167
168 *fLog << warn << "WARNING - Using MHStarMap doesn't take care of the Source Position!" << endl;
169
170 return kTRUE;
171}
172
173// --------------------------------------------------------------------------
174//
175// Get the observatory location "MObservatory" from parameter list
176//
177Bool_t MHStarMap::ReInit(MParList *pList)
178{
179 fObservatory = (MObservatory*)pList->FindObject(AddSerialNumber("MObservatory"));
180 if (!fObservatory)
181 *fLog << err << "MObservatory not found... no derotation." << endl;
182
183 return kTRUE;
184}
185
186// --------------------------------------------------------------------------
187//
188// Fill the four histograms with data from a MHillas-Container.
189// Be careful: Only call this with an object of type MHillas
190//
191Int_t MHStarMap::Fill(const MParContainer *par, const Stat_t w)
192{
193 const MHillas &h = *(MHillas*)par;
194
195 const Float_t delta = h.GetDelta();
196
197 const Float_t m = TMath::Tan(delta);
198
199 const Float_t cosd = 1.0/TMath::Sqrt(m*m+1);
200 const Float_t sind = TMath::Sqrt(1-cosd*cosd);
201
202 Float_t t = h.GetMeanY() - m*h.GetMeanX();
203
204 TVector2 src;
205 if (fSrcPos)
206 src = fSrcPos->GetXY();
207
208 if (!fUseMmScale)
209 {
210 t *= fMm2Deg;
211 src *= fMm2Deg;
212 }
213
214 // get step size ds along the main axis of the ellipse
215 const TAxis &axex = *fStarMap->GetXaxis();
216 const Float_t xmin = axex.GetXmin();
217 const Float_t xmax = axex.GetXmax();
218
219 // FIXME: Fixed number?
220 const Float_t ds = (xmax-xmin) / 200.0;
221
222 // Fill points along the main axis of the shower into the histogram;
223 // before filling
224 // - perform a translation by (xSource, ySource)
225 // - and perform a rotation to compensate the rotation of the
226 // sky image in the camera
227 const Double_t rho = fPointPos && fObservatory ?
228 fPointPos->RotationAngle(*fObservatory) : 0;
229
230 if (m>-1 && m<1)
231 {
232 const Float_t dx = ds * cosd;
233
234 for (float x=xmin+dx/2; x<(xmax-xmin)+dx; x+=dx)
235 {
236 TVector2 v(x, m*x+t);
237 v -= src;
238 v=v.Rotate(-rho);
239
240 fStarMap->Fill(v.X(), v.Y(), w);
241 }
242 }
243 else
244 {
245 const TAxis &axey = *fStarMap->GetYaxis();
246 const Float_t ymin = axey.GetXmin();
247 const Float_t ymax = axey.GetXmax();
248
249 const float dy = ds * sind;
250
251 for (float y=ymin+dy/2; y<(ymax-ymin)+dy; y+=dy)
252 {
253 TVector2 v((y-t)/m, y);
254 v -= src;
255 v=v.Rotate(-rho);
256
257 fStarMap->Fill(v.X(), v.Y(), w);
258 }
259 }
260
261 return kTRUE;
262}
263
264// --------------------------------------------------------------------------
265//
266// Use this function to setup your own conversion factor between degrees
267// and millimeters. The conversion factor should be the one calculated in
268// MGeomCam. Use this function with Caution: You could create wrong values
269// by setting up your own scale factor.
270//
271void MHStarMap::SetMm2Deg(Float_t mmdeg)
272{
273 if (mmdeg<0)
274 {
275 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
276 return;
277 }
278
279 if (fMm2Deg>=0)
280 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
281
282 fMm2Deg = mmdeg;
283}
284
285// --------------------------------------------------------------------------
286//
287// With this function you can convert the histogram ('on the fly') between
288// degrees and millimeters.
289//
290void MHStarMap::SetMmScale(Bool_t mmscale)
291{
292 if (fUseMmScale == mmscale)
293 return;
294
295 if (fMm2Deg<0)
296 {
297 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
298 return;
299 }
300
301 if (fUseMmScale)
302 {
303 fStarMap->SetXTitle("x [mm]");
304 fStarMap->SetYTitle("y [mm]");
305
306 fStarMap->Scale(1./fMm2Deg);
307 }
308 else
309 {
310 fStarMap->SetXTitle("x [\\circ]");
311 fStarMap->SetYTitle("y [\\circ]");
312
313 fStarMap->Scale(1./fMm2Deg);
314 }
315
316 fUseMmScale = mmscale;
317}
318
319// --------------------------------------------------------------------------
320//
321// Set the palette you wanna use:
322// - you could set the root "Pretty Palette Violet->Red" by
323// gStyle->SetPalette(1, 0), but in some cases this may look
324// confusing
325// - The maximum colors root allowes us to set by ourself
326// is 50 (idx: 51-100). This colors are set to a grayscaled
327// palette
328// - the number of contours must be two less than the number
329// of palette entries
330//
331void MHStarMap::PrepareDrawing() const
332{
333 const Int_t numg = 32; // number of gray scaled colors
334 const Int_t numw = 32; // number of white
335
336 Int_t palette[numg+numw];
337
338 //
339 // The first half of the colors are white.
340 // This is some kind of optical background supression
341 //
342 gROOT->GetColor(51)->SetRGB(1, 1, 1);
343
344 Int_t i;
345 for (i=0; i<numw; i++)
346 palette[i] = 51;
347
348 //
349 // now the (gray) scaled part is coming
350 //
351 for (;i<numw+numg; i++)
352 {
353 const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
354
355 gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
356 palette[i] = 52+i;
357 }
358
359 //
360 // Set the palette and the number of contour levels
361 //
362 gStyle->SetPalette(numg+numw, palette);
363 fStarMap->SetContour(numg+numw-2);
364}
365
366// --------------------------------------------------------------------------
367//
368// Draw clones of the histograms, so that the object can be deleted
369// and the histogram is still visible in the canvas.
370// The cloned object is deleted together with the canvas if the canvas is
371// destroyed. If you want to handle destroying the canvas you can get a
372// pointer to it from this function
373//
374TObject *MHStarMap::DrawClone(Option_t *opt) const
375{
376 return MH::DrawClone(opt, 500, 500);
377}
378
379// --------------------------------------------------------------------------
380//
381// Creates a new canvas and draw the histogram into it.
382// Be careful: The histogram belongs to this object and won't get deleted
383// together with the canvas.
384//
385void MHStarMap::Draw(Option_t *)
386{
387 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this, 500, 500);
388 pad->SetBorderMode(0);
389
390 pad->Divide(1,1);
391
392 pad->cd(1);
393 gPad->SetBorderMode(0);
394
395 AppendPad("");
396}
397
398void MHStarMap::Paint(Option_t *opt)
399{
400 //
401 // Maintain aspect ratio
402 //
403 const float w = gPad->GetWw();
404 const float h = gPad->GetWh();
405
406 if (h<w)
407 gPad->SetPad((1.-h/w)/2, 0, (h/w+1)/2, 1);
408 else
409 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1)/2);
410
411 //
412 // Maintain colors
413 //
414 PrepareDrawing();
415
416 //
417 // Paint Histogram
418 //
419 fStarMap->Paint("colz");
420}
Note: See TracBrowser for help on using the repository browser.