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

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