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

Last change on this file since 3367 was 3366, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 12.2 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
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
64ClassImp(MHStarMap);
65
66using namespace std;
67
68// --------------------------------------------------------------------------
69//
70// Create the star map histogram
71//
72MHStarMap::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//
106MHStarMap::~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//
121Bool_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//
178Bool_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 fSinLat = TMath::Sin(obs->GetPhi());
188 fCosLat = TMath::Cos(obs->GetPhi());
189
190 return kTRUE;
191}
192
193// --------------------------------------------------------------------------
194//
195// Calls MObservatory::RotationAngle
196//
197void MHStarMap::GetRotationAngle(Double_t &sin, Double_t &cos)
198{
199 fObservatory->RotationAngle(fMcEvt->GetTelescopeTheta(),
200 fMcEvt->GetTelescopePhi(), sin, cos);
201}
202
203// --------------------------------------------------------------------------
204//
205// Fill the four histograms with data from a MHillas-Container.
206// Be careful: Only call this with an object of type MHillas
207//
208Bool_t MHStarMap::Fill(const MParContainer *par, const Stat_t w)
209{
210 const MHillas &h = *(MHillas*)par;
211
212 const Float_t delta = h.GetDelta();
213
214 const Float_t m = TMath::Tan(delta);
215
216 const Float_t cosd = 1.0/TMath::Sqrt(m*m+1);
217 const Float_t sind = TMath::Sqrt(1-cosd*cosd);
218
219 Float_t t = h.GetMeanY() - m*h.GetMeanX();
220
221 Float_t xSource = fSrcPos->GetX();
222 Float_t ySource = fSrcPos->GetY();
223
224 if (!fUseMmScale)
225 {
226 t *= fMm2Deg;
227 xSource *= fMm2Deg;
228 ySource *= fMm2Deg;
229 }
230
231 // get step size ds along the main axis of the ellipse
232 const TAxis &axex = *fStarMap->GetXaxis();
233 const Float_t xmin = axex.GetXmin();
234 const Float_t xmax = axex.GetXmax();
235
236 // FIXME: Fixed number?
237 const Float_t ds = (xmax-xmin) / 200.0;
238
239 // Fill points along the main axis of the shower into the histogram;
240 // before filling
241 // - perform a translation by (xSource, ySource)
242 // - and perform a rotation to compensate the rotation of the
243 // sky image in the camera
244
245 const Double_t theta = fMcEvt->GetTelescopeTheta();
246 const Double_t phi = fMcEvt->GetTelescopePhi();
247
248 Double_t cosal;
249 Double_t sinal;
250 GetRotationAngle(sinal, cosal);
251
252 if (m>-1 && m<1)
253 {
254 const Float_t dx = ds * cosd;
255
256 for (float x=xmin+dx/2; x<(xmax-xmin)+dx; x+=dx)
257 {
258 const Float_t xorig = x - xSource;
259 const Float_t yorig = m*x+t - ySource;
260
261 const Float_t xfill = cosal*xorig + sinal*yorig;
262 const Float_t yfill = -sinal*xorig + cosal*yorig;
263 fStarMap->Fill(xfill, yfill, w);
264 }
265 }
266 else
267 {
268 const TAxis &axey = *fStarMap->GetYaxis();
269 const Float_t ymin = axey.GetXmin();
270 const Float_t ymax = axey.GetXmax();
271
272 const float dy = ds * sind;
273
274 for (float y=ymin+dy/2; y<(ymax-ymin)+dy; y+=dy)
275 {
276 const Float_t xorig = (y-t)/m - xSource;
277 const Float_t yorig = y - ySource;
278
279 const Float_t xfill = cosal*xorig + sinal*yorig;
280 const Float_t yfill = -sinal*xorig + cosal*yorig;
281 fStarMap->Fill(xfill, yfill, w);
282 }
283 }
284
285 return kTRUE;
286}
287
288// --------------------------------------------------------------------------
289//
290// Use this function to setup your own conversion factor between degrees
291// and millimeters. The conversion factor should be the one calculated in
292// MGeomCam. Use this function with Caution: You could create wrong values
293// by setting up your own scale factor.
294//
295void MHStarMap::SetMm2Deg(Float_t mmdeg)
296{
297 if (mmdeg<0)
298 {
299 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
300 return;
301 }
302
303 if (fMm2Deg>=0)
304 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
305
306 fMm2Deg = mmdeg;
307}
308
309// --------------------------------------------------------------------------
310//
311// With this function you can convert the histogram ('on the fly') between
312// degrees and millimeters.
313//
314void MHStarMap::SetMmScale(Bool_t mmscale)
315{
316 if (fUseMmScale == mmscale)
317 return;
318
319 if (fMm2Deg<0)
320 {
321 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
322 return;
323 }
324
325 if (fUseMmScale)
326 {
327 fStarMap->SetXTitle("x [mm]");
328 fStarMap->SetYTitle("y [mm]");
329
330 fStarMap->Scale(1./fMm2Deg);
331 }
332 else
333 {
334 fStarMap->SetXTitle("x [\\circ]");
335 fStarMap->SetYTitle("y [\\circ]");
336
337 fStarMap->Scale(1./fMm2Deg);
338 }
339
340 fUseMmScale = mmscale;
341}
342
343// --------------------------------------------------------------------------
344//
345// Set the palette you wanna use:
346// - you could set the root "Pretty Palette Violet->Red" by
347// gStyle->SetPalette(1, 0), but in some cases this may look
348// confusing
349// - The maximum colors root allowes us to set by ourself
350// is 50 (idx: 51-100). This colors are set to a grayscaled
351// palette
352// - the number of contours must be two less than the number
353// of palette entries
354//
355void MHStarMap::PrepareDrawing() const
356{
357 const Int_t numg = 32; // number of gray scaled colors
358 const Int_t numw = 32; // number of white
359
360 Int_t palette[numg+numw];
361
362 //
363 // The first half of the colors are white.
364 // This is some kind of optical background supression
365 //
366 gROOT->GetColor(51)->SetRGB(1, 1, 1);
367
368 Int_t i;
369 for (i=0; i<numw; i++)
370 palette[i] = 51;
371
372 //
373 // now the (gray) scaled part is coming
374 //
375 for (;i<numw+numg; i++)
376 {
377 const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
378
379 gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
380 palette[i] = 52+i;
381 }
382
383 //
384 // Set the palette and the number of contour levels
385 //
386 gStyle->SetPalette(numg+numw, palette);
387 fStarMap->SetContour(numg+numw-2);
388}
389
390// --------------------------------------------------------------------------
391//
392// Draw clones of the histograms, so that the object can be deleted
393// and the histogram is still visible in the canvas.
394// The cloned object is deleted together with the canvas if the canvas is
395// destroyed. If you want to handle destroying the canvas you can get a
396// pointer to it from this function
397//
398TObject *MHStarMap::DrawClone(Option_t *opt) const
399{
400 return MH::DrawClone(opt, 500, 500);
401}
402
403// --------------------------------------------------------------------------
404//
405// Creates a new canvas and draw the histogram into it.
406// Be careful: The histogram belongs to this object and won't get deleted
407// together with the canvas.
408//
409void MHStarMap::Draw(Option_t *)
410{
411 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this, 500, 500);
412 pad->SetBorderMode(0);
413
414 pad->Divide(1,1);
415
416 pad->cd(1);
417 gPad->SetBorderMode(0);
418
419 AppendPad("");
420}
421
422void MHStarMap::Paint(Option_t *opt)
423{
424 //
425 // Maintain aspect ratio
426 //
427 const float w = gPad->GetWw();
428 const float h = gPad->GetWh();
429
430 if (h<w)
431 gPad->SetPad((1.-h/w)/2, 0, (h/w+1)/2, 1);
432 else
433 gPad->SetPad(0, (1.-w/h)/2, 1, (w/h+1)/2);
434
435 //
436 // Maintain colors
437 //
438 PrepareDrawing();
439
440 //
441 // Paint Histogram
442 //
443 fStarMap->Paint("colz");
444}
Note: See TracBrowser for help on using the repository browser.