source: tags/Mars-V0.9.1/mimage/MHNewImagePar.cc

Last change on this file was 6907, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 10.6 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): Wolfgang Wittek, 03/2003 <mailto:wittek@mppmu.mpg.de>
19! Author(s): Thomas Bretz, 04/2003 <mailto:tbretz@astro.uni-wuerzburg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MHNewImagePar
29//
30////////////////////////////////////////////////////////////////////////////
31#include "MHNewImagePar.h"
32
33#include <math.h>
34
35#include <TH1.h>
36#include <TPad.h>
37#include <TCanvas.h>
38
39#include "MLog.h"
40#include "MLogManip.h"
41
42#include "MGeomCam.h"
43#include "MBinning.h"
44#include "MParList.h"
45
46#include "MHillas.h"
47#include "MNewImagePar.h"
48
49ClassImp(MHNewImagePar);
50
51using namespace std;
52
53// --------------------------------------------------------------------------
54//
55// Setup histograms
56//
57MHNewImagePar::MHNewImagePar(const char *name, const char *title)
58 : fMm2Deg(1), fUseMmScale(kTRUE)
59{
60 fName = name ? name : "MHNewImagePar";
61 fTitle = title ? title : "Histograms of new image parameters";
62
63 fHistLeakage1.SetName("Leakage1");
64 fHistLeakage1.SetTitle("Leakage_{1}");
65 fHistLeakage1.SetXTitle("Leakage");
66 fHistLeakage1.SetYTitle("Counts");
67 fHistLeakage1.SetDirectory(NULL);
68 fHistLeakage1.UseCurrentStyle();
69 fHistLeakage1.SetFillStyle(4000);
70
71 fHistLeakage2.SetName("Leakage2");
72 fHistLeakage2.SetTitle("Leakage_{2}");
73 fHistLeakage2.SetXTitle("Leakage");
74 fHistLeakage2.SetYTitle("Counts");
75 fHistLeakage2.SetDirectory(NULL);
76 fHistLeakage2.UseCurrentStyle();
77 fHistLeakage2.SetLineColor(kBlue);
78 fHistLeakage2.SetFillStyle(4000);
79
80 fHistUsedPix.SetName("UsedPix");
81 fHistUsedPix.SetTitle("Number of used pixels");
82 fHistUsedPix.SetXTitle("Number of Pixels");
83 fHistUsedPix.SetYTitle("Counts");
84 fHistUsedPix.SetDirectory(NULL);
85 fHistUsedPix.UseCurrentStyle();
86 fHistUsedPix.SetLineColor(kBlue);
87 fHistUsedPix.SetFillStyle(4000);
88
89 fHistCorePix.SetName("CorePix");
90 fHistCorePix.SetTitle("Number of core pixels");
91 fHistCorePix.SetXTitle("Number of Pixels");
92 fHistCorePix.SetYTitle("Counts");
93 fHistCorePix.SetDirectory(NULL);
94 fHistCorePix.UseCurrentStyle();
95 fHistCorePix.SetLineColor(kBlack);
96 fHistCorePix.SetFillStyle(4000);
97
98 fHistUsedArea.SetName("UsedArea");
99 fHistUsedArea.SetTitle("Area of used pixels");
100 fHistUsedArea.SetXTitle("Area [m^2]");
101 fHistUsedArea.SetYTitle("Counts");
102 fHistUsedArea.SetDirectory(NULL);
103 fHistUsedArea.UseCurrentStyle();
104 fHistUsedArea.SetLineColor(kBlue);
105 fHistUsedArea.SetFillStyle(4000);
106
107 fHistCoreArea.SetName("CoreArea");
108 fHistCoreArea.SetTitle("Area of core pixels");
109 fHistCoreArea.SetXTitle("Area [m^2]");
110 fHistCoreArea.SetYTitle("Counts");
111 fHistCoreArea.SetDirectory(NULL);
112 fHistCoreArea.UseCurrentStyle();
113 fHistCoreArea.SetLineColor(kBlack);
114 fHistCoreArea.SetFillStyle(4000);
115
116 fHistConc.SetDirectory(NULL);
117 fHistConc1.SetDirectory(NULL);
118 fHistConc.SetName("Conc2");
119 fHistConc1.SetName("Conc1");
120 fHistConc.SetTitle("Ratio: Conc");
121 fHistConc1.SetTitle("Ratio: Conc1");
122 fHistConc.SetXTitle("Ratio");
123 fHistConc1.SetXTitle("Ratio");
124 fHistConc.SetYTitle("Counts");
125 fHistConc1.SetYTitle("Counts");
126 fHistConc1.UseCurrentStyle();
127 fHistConc.UseCurrentStyle();
128 fHistConc.SetFillStyle(4000);
129 fHistConc1.SetFillStyle(4000);
130 fHistConc1.SetLineColor(kBlue);
131
132
133 MBinning bins;
134
135 bins.SetEdges(100, 0, 1);
136 bins.Apply(fHistLeakage1);
137 bins.Apply(fHistLeakage2);
138 bins.Apply(fHistConc);
139 bins.Apply(fHistConc1);
140
141 bins.SetEdges(75, 0.5, 150.5);
142 bins.Apply(fHistUsedPix);
143 bins.Apply(fHistCorePix);
144
145 //bins.SetEdges(75, 0, 0.249);
146 //bins.Apply(fHistUsedArea);
147 //bins.Apply(fHistCoreArea);
148}
149
150// --------------------------------------------------------------------------
151//
152// Setup the Binning for the histograms automatically if the correct
153// instances of MBinning
154//
155Bool_t MHNewImagePar::SetupFill(const MParList *plist)
156{
157 MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
158 if (!geom)
159 *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl;
160 else
161 {
162 fMm2Deg = geom->GetConvMm2Deg();
163 SetMmScale(kFALSE);
164 }
165
166 const MBinning *bins = (MBinning*)plist->FindObject("BinningArea");
167 if (!bins)
168 {
169 float r = geom ? 1.5 : 0.133;
170
171 MBinning b;
172 b.SetEdges(50, 0, r);
173 b.Apply(fHistUsedArea);
174 b.Apply(fHistCoreArea);
175 }
176 else
177 {
178 bins->Apply(fHistUsedArea);
179 bins->Apply(fHistCoreArea);
180 }
181
182 ApplyBinning(*plist, "Leakage", &fHistLeakage1);
183 ApplyBinning(*plist, "Leakage", &fHistLeakage2);
184
185 ApplyBinning(*plist, "Pixels", &fHistUsedPix);
186 ApplyBinning(*plist, "Pixels", &fHistCorePix);
187
188 //ApplyBinning(*plist, "Area", &fHistUsedArea);
189 //ApplyBinning(*plist, "Area", &fHistCoreArea);
190
191 ApplyBinning(*plist, "Conc", &fHistConc);
192 ApplyBinning(*plist, "Conc1", &fHistConc1);
193
194 return kTRUE;
195}
196
197
198// --------------------------------------------------------------------------
199//
200// Fill the histograms with data from a MNewImagePar container.
201//
202Bool_t MHNewImagePar::Fill(const MParContainer *par, const Stat_t w)
203{
204 if (!par)
205 {
206 *fLog << err << "MHNewImagePar::Fill: Pointer (!=NULL) expected." << endl;
207 return kFALSE;
208 }
209
210 const Double_t scale = fUseMmScale ? 1e-6 : fMm2Deg*fMm2Deg;
211
212 const MNewImagePar &h = *(MNewImagePar*)par;
213
214 fHistLeakage1.Fill(h.GetLeakage1(), w);
215 fHistLeakage2.Fill(h.GetLeakage2(), w);
216
217 fHistUsedPix.Fill(h.GetNumUsedPixels(), w);
218 fHistCorePix.Fill(h.GetNumCorePixels(), w);
219
220 fHistUsedArea.Fill(h.GetUsedArea()*scale, w);
221 fHistCoreArea.Fill(h.GetCoreArea()*scale, w);
222
223 fHistConc.Fill(h.GetConc(), w);
224 fHistConc1.Fill(h.GetConc1(), w);
225
226 return kTRUE;
227}
228
229// --------------------------------------------------------------------------
230//
231// With this function you can convert the histogram ('on the fly') between
232// degrees and millimeters.
233//
234void MHNewImagePar::SetMmScale(Bool_t mmscale)
235{
236 if (fUseMmScale == mmscale)
237 return;
238
239 if (fMm2Deg<0)
240 {
241 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
242 return;
243 }
244
245 const Double_t scale = mmscale ? 1./(fMm2Deg*fMm2Deg) : (fMm2Deg*fMm2Deg);
246 MH::ScaleAxis(&fHistUsedArea, scale);
247 MH::ScaleAxis(&fHistCoreArea, scale);
248
249 if (mmscale)
250 {
251 fHistUsedArea.SetXTitle("A [m^{2}]");
252 fHistCoreArea.SetXTitle("A [m^{2}]");
253 }
254 else
255 {
256 fHistUsedArea.SetXTitle("A [deg^{2}]");
257 fHistCoreArea.SetXTitle("A [deg^{2}]");
258 }
259
260 fUseMmScale = mmscale;
261}
262
263// --------------------------------------------------------------------------
264//
265// Use this function to setup your own conversion factor between degrees
266// and millimeters. The conversion factor should be the one calculated in
267// MGeomCam. Use this function with Caution: You could create wrong values
268// by setting up your own scale factor.
269//
270void MHNewImagePar::SetMm2Deg(Float_t mmdeg)
271{
272 if (mmdeg<0)
273 {
274 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
275 return;
276 }
277
278 if (fMm2Deg>=0)
279 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
280
281 fMm2Deg = mmdeg;
282}
283
284void MHNewImagePar::Paint(Option_t *o)
285{
286 if (fHistLeakage1.GetMaximum()>0 && gPad->GetPad(1))
287 gPad->GetPad(1)->SetLogy();
288}
289
290// --------------------------------------------------------------------------
291//
292// Creates a new canvas and draws the two histograms into it.
293// Be careful: The histograms belongs to this object and won't get deleted
294// together with the canvas.
295//
296void MHNewImagePar::Draw(Option_t *o)
297{
298 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
299 pad->SetBorderMode(0);
300
301 AppendPad("");
302
303 // FIXME: If same-option given make two independant y-axis!
304 const TString opt(o);
305 const Bool_t same = opt.Contains("same");
306
307 if (!same)
308 pad->Divide(2,2);
309 else
310 {
311 fHistLeakage1.SetLineColor(kMagenta);
312 fHistLeakage1.SetLineColor(kCyan);
313
314 fHistCorePix.SetLineColor(kMagenta);
315 fHistUsedPix.SetLineColor(kCyan);
316
317 fHistConc1.SetLineColor(kMagenta);
318 fHistConc.SetLineColor(kCyan);
319
320 fHistCoreArea.SetLineColor(kMagenta);
321 fHistUsedArea.SetLineColor(kCyan);
322 }
323
324 pad->cd(1);
325 gPad->SetBorderMode(0);
326 TAxis &x = *fHistLeakage1.GetXaxis();
327 x.SetRangeUser(0.0, x.GetXmax());
328 MH::DrawSame(fHistLeakage1, fHistLeakage2, "Leakage1 and Leakage2", same);
329 fHistLeakage1.SetMinimum();
330 fHistLeakage2.SetMinimum();
331 fHistLeakage2.SetMaximum(0.1); // dummy value to allow log-scale
332
333 pad->cd(2);
334 gPad->SetBorderMode(0);
335 MH::DrawSame(fHistCorePix, fHistUsedPix, "Number of core/used Pixels", same);
336
337 pad->cd(3);
338 gPad->SetBorderMode(0);
339 MH::DrawSame(fHistConc1, fHistConc, "Concentrations", same);
340
341 pad->cd(4);
342 gPad->SetBorderMode(0);
343 MH::DrawSame(fHistCoreArea, fHistUsedArea, "Area of core/used Pixels", same);
344}
345
346TH1 *MHNewImagePar::GetHistByName(const TString name)
347{
348 if (name.Contains("Leakage1", TString::kIgnoreCase))
349 return &fHistLeakage1;
350 if (name.Contains("Leakage2", TString::kIgnoreCase))
351 return &fHistLeakage2;
352 if (name.Contains("Conc", TString::kIgnoreCase))
353 return &fHistConc;
354 if (name.Contains("Conc1", TString::kIgnoreCase))
355 return &fHistConc1;
356 if (name.Contains("UsedPix", TString::kIgnoreCase))
357 return &fHistUsedPix;
358 if (name.Contains("CorePix", TString::kIgnoreCase))
359 return &fHistCorePix;
360 if (name.Contains("UsedArea", TString::kIgnoreCase))
361 return &fHistUsedArea;
362 if (name.Contains("CoreArea", TString::kIgnoreCase))
363 return &fHistCoreArea;
364
365 return NULL;
366}
Note: See TracBrowser for help on using the repository browser.