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

Last change on this file since 1376 was 1330, checked in by tbretz, 23 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 //
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 'BinningWidth' and 'BinningLength')
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 *= 2./3;
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 return kTRUE;
139}
140
141// --------------------------------------------------------------------------
142//
143// Fill the four histograms with data from a MHillas-Container.
144// Be careful: Only call this with an object of type MHillas
145//
146Bool_t MHStarMap::Fill(const MParContainer *par)
147{
148 const MHillas &h = *(MHillas*)par;
149
150 const float delta = h.GetDelta();
151
152 const float m = tan(delta);
153
154 float t = h.GetMeanY() - m*h.GetMeanX();
155
156 if (!fUseMmScale)
157 t *= fMm2Deg;
158
159 if (m>-1 && m<1)
160 {
161 TAxis &axe = *fStarMap->GetXaxis();
162
163 const int N = axe.GetXbins()->GetSize();
164 for (int i=0; i<N; i++)
165 {
166 const float x = axe.GetBinCenter(i);
167 const float y = m*x+t;
168
169 fStarMap->Fill(x, y);
170 }
171 }
172 else
173 {
174 TAxis &axe = *fStarMap->GetYaxis();
175
176 const int N = axe.GetXbins()->GetSize();
177 for (int i=0; i<N; i++)
178 {
179 const float y = axe.GetBinCenter(i);
180 const float x = (y-t)/m;
181
182 fStarMap->Fill(x, y);
183 }
184 }
185
186 return kTRUE;
187}
188
189// --------------------------------------------------------------------------
190//
191// Use this function to setup your own conversion factor between degrees
192// and millimeters. The conversion factor should be the one calculated in
193// MGeomCam. Use this function with Caution: You could create wrong values
194// by setting up your own scale factor.
195//
196void MHStarMap::SetMm2Deg(Float_t mmdeg)
197{
198 if (mmdeg<0)
199 {
200 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
201 return;
202 }
203
204 if (fMm2Deg>=0)
205 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
206
207 fMm2Deg = mmdeg;
208}
209
210// --------------------------------------------------------------------------
211//
212// With this function you can convert the histogram ('on the fly') between
213// degrees and millimeters.
214//
215void MHStarMap::SetMmScale(Bool_t mmscale)
216{
217 if (fUseMmScale == mmscale)
218 return;
219
220 if (fMm2Deg<0)
221 {
222 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
223 return;
224 }
225
226 if (fUseMmScale)
227 {
228 fStarMap->SetXTitle("x [mm]");
229 fStarMap->SetYTitle("y [mm]");
230
231 fStarMap->Scale(1./fMm2Deg);
232 }
233 else
234 {
235 fStarMap->SetXTitle("x [\\circ]");
236 fStarMap->SetYTitle("y [\\circ]");
237
238 fStarMap->Scale(1./fMm2Deg);
239 }
240
241 fUseMmScale = mmscale;
242}
243
244// --------------------------------------------------------------------------
245//
246// Set the palette you wanna use:
247// - you could set the root "Pretty Palette Violet->Red" by
248// gStyle->SetPalette(1, 0), but in some cases this may look
249// confusing
250// - The maximum colors root allowes us to set by ourself
251// is 50 (idx: 51-100). This colors are set to a grayscaled
252// palette
253// - the number of contours must be two less than the number
254// of palette entries
255//
256void MHStarMap::PrepareDrawing() const
257{
258 const Int_t numg = 32; // number of gray scaled colors
259 const Int_t numw = 32; // number of white
260
261 Int_t palette[numg+numw];
262
263 //
264 // The first half of the colors are white.
265 // This is some kind of optical background supression
266 //
267 gROOT->GetColor(51)->SetRGB(1, 1, 1);
268
269 Int_t i;
270 for (i=0; i<numw; i++)
271 palette[i] = 51;
272
273 //
274 // now the (gray) scaled part is coming
275 //
276 for (;i<numw+numg; i++)
277 {
278 const Float_t gray = 1.0-(float)(i-numw)/(numg-1.0);
279
280 gROOT->GetColor(52+i)->SetRGB(gray, gray, gray);
281 palette[i] = 52+i;
282 }
283
284 //
285 // Set the palette and the number of contour levels
286 //
287 gStyle->SetPalette(numg+numw, palette);
288 fStarMap->SetContour(numg+numw-2);
289}
290
291
292// --------------------------------------------------------------------------
293//
294// Draw clones of the histograms, so that the object can be deleted
295// and the histogram is still visible in the canvas.
296// The cloned object is deleted together with the canvas if the canvas is
297// destroyed. If you want to handle destroying the canvas you can get a
298// pointer to it from this function
299//
300TObject *MHStarMap::DrawClone(Option_t *opt) const
301{
302 TCanvas *c=MakeDefCanvas(fStarMap, 500, 500);
303
304 //
305 // This is necessary to get the expected bahviour of DrawClone
306 //
307 gROOT->SetSelectedPad(NULL);
308
309 PrepareDrawing();
310
311 fStarMap->DrawCopy("colz");
312
313 c->Modified();
314 c->Update();
315
316 return c;
317}
318
319// --------------------------------------------------------------------------
320//
321// Creates a new canvas and draw the histogram into it.
322// Be careful: The histogram belongs to this object and won't get deleted
323// together with the canvas.
324//
325void MHStarMap::Draw(Option_t *)
326{
327 if (!gPad)
328 MakeDefCanvas(fStarMap, 500, 500);
329
330 PrepareDrawing();
331
332 fStarMap->Draw("colz");
333
334 gPad->Modified();
335 gPad->Update();
336}
Note: See TracBrowser for help on using the repository browser.