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

Last change on this file since 1788 was 1715, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 8.9 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 *fLog << warn << "WARNING - Using MHStarMap doesn't take care of the Source Position!" << endl;
71
72 //
73 // loop over all Pixels and create two histograms
74 // one for the Low and one for the High gain
75 // connect all the histogram with the container fHist
76 //
77 fStarMap = new TH2F("StarMap", " 2D Hillas Star Map ",
78 150, -300, 300,
79 150, -300, 300);
80
81 fStarMap->SetDirectory(NULL);
82
83 fStarMap->SetXTitle("x [mm]");
84 fStarMap->SetYTitle("y [mm]");
85 fStarMap->SetZTitle("Counts");
86}
87
88// --------------------------------------------------------------------------
89//
90// delete the histogram instance
91//
92MHStarMap::~MHStarMap()
93{
94 delete fStarMap;
95}
96
97// --------------------------------------------------------------------------
98//
99// Setup the Binning for the histograms automatically if the correct
100// instances of MBinning (with the names 'BinningCamera')
101// are found in the parameter list
102// Use this function if you want to set the conversion factor which
103// is used to convert the mm-scale in the camera plain into the deg-scale
104// used for histogram presentations. The conversion factor is part of
105// the camera geometry. Please create a corresponding MGeomCam container.
106//
107Bool_t MHStarMap::SetupFill(const MParList *plist)
108{
109 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
110 if (geom)
111 {
112 fMm2Deg = geom->GetConvMm2Deg();
113 fUseMmScale = kFALSE;
114
115 fStarMap->SetXTitle("x [\\circ]");
116 fStarMap->SetYTitle("y [\\circ]");
117 }
118
119 const MBinning *bins = (MBinning*)plist->FindObject("BinningStarMap");
120 if (!bins)
121 {
122 float r = geom ? geom->GetMaxRadius() : 600;
123 r *= 5./6;
124 if (!fUseMmScale)
125 r *= fMm2Deg;
126
127 MBinning b;
128 b.SetEdges(100, -r, r);
129 SetBinning(fStarMap, &b, &b);
130 }
131 else
132 SetBinning(fStarMap, bins, bins);
133
134 if (!geom)
135 {
136 *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
137 return kTRUE;
138 }
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)
149{
150 const MHillas &h = *(MHillas*)par;
151
152 const float delta = h.GetDelta();
153
154 const float m = tan(delta);
155
156 float t = h.GetMeanY() - m*h.GetMeanX();
157
158 if (!fUseMmScale)
159 t *= fMm2Deg;
160
161 if (m>-1 && m<1)
162 {
163 TAxis &axe = *fStarMap->GetXaxis();
164
165 const int N = axe.GetXbins()->GetSize();
166 for (int i=0; i<N; i++)
167 {
168 const float x = axe.GetBinCenter(i);
169 const float y = m*x+t;
170
171 fStarMap->Fill(x, y);
172 }
173 }
174 else
175 {
176 TAxis &axe = *fStarMap->GetYaxis();
177
178 const int N = axe.GetXbins()->GetSize();
179 for (int i=0; i<N; i++)
180 {
181 const float y = axe.GetBinCenter(i);
182 const float x = (y-t)/m;
183
184 fStarMap->Fill(x, y);
185 }
186 }
187
188 return kTRUE;
189}
190
191// --------------------------------------------------------------------------
192//
193// Use this function to setup your own conversion factor between degrees
194// and millimeters. The conversion factor should be the one calculated in
195// MGeomCam. Use this function with Caution: You could create wrong values
196// by setting up your own scale factor.
197//
198void MHStarMap::SetMm2Deg(Float_t mmdeg)
199{
200 if (mmdeg<0)
201 {
202 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
203 return;
204 }
205
206 if (fMm2Deg>=0)
207 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
208
209 fMm2Deg = mmdeg;
210}
211
212// --------------------------------------------------------------------------
213//
214// With this function you can convert the histogram ('on the fly') between
215// degrees and millimeters.
216//
217void MHStarMap::SetMmScale(Bool_t mmscale)
218{
219 if (fUseMmScale == mmscale)
220 return;
221
222 if (fMm2Deg<0)
223 {
224 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
225 return;
226 }
227
228 if (fUseMmScale)
229 {
230 fStarMap->SetXTitle("x [mm]");
231 fStarMap->SetYTitle("y [mm]");
232
233 fStarMap->Scale(1./fMm2Deg);
234 }
235 else
236 {
237 fStarMap->SetXTitle("x [\\circ]");
238 fStarMap->SetYTitle("y [\\circ]");
239
240 fStarMap->Scale(1./fMm2Deg);
241 }
242
243 fUseMmScale = mmscale;
244}
245
246// --------------------------------------------------------------------------
247//
248// Set the palette you wanna use:
249// - you could set the root "Pretty Palette Violet->Red" by
250// gStyle->SetPalette(1, 0), but in some cases this may look
251// confusing
252// - The maximum colors root allowes us to set by ourself
253// is 50 (idx: 51-100). This colors are set to a grayscaled
254// palette
255// - the number of contours must be two less than the number
256// of palette entries
257//
258void MHStarMap::PrepareDrawing() const
259{
260 const Int_t numg = 32; // number of gray scaled colors
261 const Int_t numw = 32; // number of white
262
263 Int_t palette[numg+numw];
264
265 //
266 // The first half of the colors are white.
267 // This is some kind of optical background supression
268 //
269 gROOT->GetColor(51)->SetRGB(1, 1, 1);
270
271 Int_t i;
272 for (i=0; i<numw; i++)
273 palette[i] = 51;
274
275 //
276 // now the (gray) scaled part is coming
277 //
278 for (;i<numw+numg; i++)
279 {
280 const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
281
282 gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
283 palette[i] = 52+i;
284 }
285
286 //
287 // Set the palette and the number of contour levels
288 //
289 gStyle->SetPalette(numg+numw, palette);
290 fStarMap->SetContour(numg+numw-2);
291}
292
293
294// --------------------------------------------------------------------------
295//
296// Draw clones of the histograms, so that the object can be deleted
297// and the histogram is still visible in the canvas.
298// The cloned object is deleted together with the canvas if the canvas is
299// destroyed. If you want to handle destroying the canvas you can get a
300// pointer to it from this function
301//
302TObject *MHStarMap::DrawClone(Option_t *opt) const
303{
304 TCanvas *c=MakeDefCanvas(fStarMap, 500, 500);
305
306 //
307 // This is necessary to get the expected bahviour of DrawClone
308 //
309 gROOT->SetSelectedPad(NULL);
310
311 PrepareDrawing();
312
313 fStarMap->DrawCopy("colz");
314
315 c->Modified();
316 c->Update();
317
318 return c;
319}
320
321// --------------------------------------------------------------------------
322//
323// Creates a new canvas and draw the histogram into it.
324// Be careful: The histogram belongs to this object and won't get deleted
325// together with the canvas.
326//
327void MHStarMap::Draw(Option_t *)
328{
329 if (!gPad)
330 MakeDefCanvas(fStarMap, 500, 500);
331
332 PrepareDrawing();
333
334 fStarMap->Draw("colz");
335
336 gPad->Modified();
337 gPad->Update();
338}
Note: See TracBrowser for help on using the repository browser.