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

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