source: trunk/MagicSoft/Mars/mimage/MHVsSize.cc@ 6942

Last change on this file since 6942 was 6890, checked in by tbretz, 20 years ago
*** empty log message ***
File size: 8.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): Thomas Bretz 2001 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Wolfgang Wittek 2002 <mailto:wittek@mppmu.mpg.de>
20!
21! Copyright: MAGIC Software Development, 2000-2004
22!
23!
24\* ======================================================================== */
25
26/////////////////////////////////////////////////////////////////////////////
27//
28// MHVsSize
29//
30// This class contains histograms for the source independent image parameters
31//
32/////////////////////////////////////////////////////////////////////////////
33#include "MHVsSize.h"
34
35#include <TH2.h>
36#include <TPad.h>
37#include <TStyle.h>
38#include <TCanvas.h>
39
40#include "MLog.h"
41#include "MLogManip.h"
42
43#include "MParList.h"
44
45#include "MHillas.h"
46#include "MHillasSrc.h"
47#include "MNewImagePar.h"
48#include "MGeomCam.h"
49#include "MBinning.h"
50
51ClassImp(MHVsSize);
52
53using namespace std;
54
55// --------------------------------------------------------------------------
56//
57// Setup four histograms for Width, Length
58//
59MHVsSize::MHVsSize(const char *name, const char *title)
60 : fHillas(NULL), fNewImagePar(NULL), fMm2Deg(1), fUseMmScale(kTRUE)
61{
62 //
63 // set the name and title of this object
64 //
65 fName = name ? name : "MHVsSize";
66 fTitle = title ? title : "Source independent image parameters";
67
68 fLength.SetNameTitle("Length", "Length vs. Size");
69 fWidth.SetNameTitle( "Width", "Width vs. Size");
70 fDist.SetNameTitle( "Dist", "Dist vs. Size");
71 fConc1.SetNameTitle( "Conc1", "Conc1 vs. Size");
72
73 fLength.SetDirectory(NULL);
74 fWidth.SetDirectory(NULL);
75 fDist.SetDirectory(NULL);
76 fConc1.SetDirectory(NULL);
77
78 fLength.SetXTitle("Size [phe]");
79 fWidth.SetXTitle("Size [phe]");
80 fDist.SetXTitle("Size [phe]");
81 fConc1.SetXTitle("Size [phe]");
82
83 fLength.SetYTitle("Length [mm]");
84 fWidth.SetYTitle("Width [mm]");
85 fDist.SetYTitle("Distance [mm]");
86 fConc1.SetYTitle("Conc1 [ratio]");
87
88 MBinning binse, binsl, binsd, binsc;
89 binse.SetEdgesLog(50, 10, 1e5);
90 binsl.SetEdges(100, 0, 296.7/2);
91 binsd.SetEdges(100, 0, 445);
92 binsc.SetEdgesLog(100, 1e-5, 5e-3);
93
94 MH::SetBinning(&fLength, &binse, &binsl);
95 MH::SetBinning(&fWidth, &binse, &binsl);
96 MH::SetBinning(&fDist, &binse, &binsd);
97 MH::SetBinning(&fConc1, &binse, &binsc);
98
99 fLength.UseCurrentStyle();
100 fWidth.UseCurrentStyle();
101 fDist.UseCurrentStyle();
102 fConc1.UseCurrentStyle();
103}
104
105// --------------------------------------------------------------------------
106//
107// Setup the Binning for the histograms automatically if the correct
108// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
109// are found in the parameter list
110// Use this function if you want to set the conversion factor which
111// is used to convert the mm-scale in the camera plain into the deg-scale
112// used for histogram presentations. The conversion factor is part of
113// the camera geometry. Please create a corresponding MGeomCam container.
114//
115Bool_t MHVsSize::SetupFill(const MParList *plist)
116{
117 MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
118 if (!geom)
119 *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl;
120 else
121 {
122 fMm2Deg = geom->GetConvMm2Deg();
123 SetMmScale(kFALSE);
124 }
125
126 fHillas = (MHillas*)plist->FindObject("MHillas");
127 if (!fHillas)
128 {
129 *fLog << err << "MHillas not found... abort." << endl;
130 return kFALSE;
131 }
132
133 fNewImagePar = (MNewImagePar*)plist->FindObject("MNewImagePar");
134 if (!fNewImagePar)
135 {
136 *fLog << err << "MNewImagePar not found... abort." << endl;
137 return kFALSE;
138 }
139
140 /*
141 ApplyBinning(*plist, "Width", fWidth);
142 ApplyBinning(*plist, "Length", fLength);
143 ApplyBinning(*plist, "Dist", fDistC);
144 ApplyBinning(*plist, "Delta", fDelta);
145 ApplyBinning(*plist, "Size", fSize);
146
147 const MBinning *bins = (MBinning*)plist->FindObject("BinningCamera");
148 if (!bins)
149 {
150 float r = fGeomCam ? fGeomCam->GetMaxRadius() : 600;
151 r *= 0.9;
152 if (!fUseMmScale)
153 r *= fMm2Deg;
154
155 MBinning b;
156 b.SetEdges(61, -r, r);
157 SetBinning(fCenter, &b, &b);
158 }
159 else
160 SetBinning(fCenter, bins, bins);
161 */
162
163 return kTRUE;
164}
165
166// --------------------------------------------------------------------------
167//
168// Use this function to setup your own conversion factor between degrees
169// and millimeters. The conversion factor should be the one calculated in
170// MGeomCam. Use this function with Caution: You could create wrong values
171// by setting up your own scale factor.
172//
173void MHVsSize::SetMm2Deg(Float_t mmdeg)
174{
175 if (mmdeg<0)
176 {
177 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
178 return;
179 }
180
181 if (fMm2Deg>=0)
182 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
183
184 fMm2Deg = mmdeg;
185}
186
187// --------------------------------------------------------------------------
188//
189// With this function you can convert the histogram ('on the fly') between
190// degrees and millimeters.
191//
192void MHVsSize::SetMmScale(Bool_t mmscale)
193{
194 if (fUseMmScale == mmscale)
195 return;
196
197 if (fMm2Deg<0)
198 {
199 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
200 return;
201 }
202
203 const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
204 MH::ScaleAxis(&fLength, 1, scale);
205 MH::ScaleAxis(&fWidth, 1, scale);
206 MH::ScaleAxis(&fDist, 1, scale);
207
208 if (mmscale)
209 {
210 fLength.SetYTitle("Length [mm]");
211 fWidth.SetYTitle("Width [mm]");
212 fDist.SetYTitle("Distance [mm]");
213 }
214 else
215 {
216 fLength.SetYTitle("Length [\\circ]");
217 fWidth.SetYTitle("Width [\\circ]");
218 fDist.SetYTitle("Distance [\\circ]");
219 }
220
221 fUseMmScale = mmscale;
222}
223
224// --------------------------------------------------------------------------
225//
226// Fill the histograms with data from a MHillas-Container.
227// Be careful: Only call this with an object of type MHillas
228//
229Bool_t MHVsSize::Fill(const MParContainer *par, const Stat_t w)
230{
231 const MHillasSrc *src = dynamic_cast<const MHillasSrc*>(par);
232 if (!src)
233 {
234 *fLog << err << "MHVsSize::Fill: Pointer (!=NULL) expected." << endl;
235 return kFALSE;
236 }
237
238 const Double_t scale = fUseMmScale ? 1 : fMm2Deg;
239
240 fLength.Fill(fHillas->GetSize(), scale*fHillas->GetLength(), w);
241 fWidth.Fill( fHillas->GetSize(), scale*fHillas->GetWidth(), w);
242 fDist.Fill( fHillas->GetSize(), scale*src->GetDist(), w);
243 fConc1.Fill( fHillas->GetSize(), scale*fNewImagePar->GetConc1(), w);
244
245 return kTRUE;
246}
247
248// --------------------------------------------------------------------------
249//
250// Creates a new canvas and draws the four histograms into it.
251// Be careful: The histograms belongs to this object and won't get deleted
252// together with the canvas.
253//
254void MHVsSize::Draw(Option_t *o)
255{
256 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
257 pad->SetBorderMode(0);
258
259 AppendPad("");
260
261 TString opt(o);
262 opt.ToLower();
263
264 // FIXME: If same-option given make two independant y-axis!
265 const Bool_t same = opt.Contains("same");
266
267 if (!same)
268 pad->Divide(2,2);
269 else
270 {
271 fDist.SetMarkerColor(kBlue);
272 fConc1.SetMarkerColor(kBlue);
273 fWidth.SetMarkerColor(kBlue);
274 fLength.SetMarkerColor(kBlue);
275 }
276
277 pad->cd(1);
278 gPad->SetBorderMode(0);
279 gPad->SetLogx();
280 fLength.Draw(same?"same":"");
281
282 pad->cd(2);
283 gPad->SetBorderMode(0);
284 gPad->SetLogx();
285 fWidth.Draw(same?"same":"");
286
287 pad->cd(3);
288 gPad->SetBorderMode(0);
289 gPad->SetLogx();
290 fDist.Draw(same?"same":"");
291
292 pad->cd(4);
293 gPad->SetBorderMode(0);
294 gPad->SetLogx();
295 gPad->SetLogy();
296 fConc1.Draw(same?"same":"");
297}
298
299void MHVsSize::Paint(Option_t *opt)
300{
301 //SetColors();
302 //MH::Paint();
303}
Note: See TracBrowser for help on using the repository browser.