source: trunk/MagicSoft/Mars/mimage/MHHillasExt.cc@ 9343

Last change on this file since 9343 was 9343, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 8.2 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!
20! Copyright: MAGIC Software Development, 2000-2007
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHHillasExt
28//
29// This class contains histograms for every Hillas parameter
30//
31// Class Version 2:
32// ----------------
33// - fHMaxDist
34// + fHSlopeL
35//
36// ClassVersion 3:
37// ---------------
38// - fMm2Deg
39// - fUseMmScale
40//
41/////////////////////////////////////////////////////////////////////////////
42#include "MHHillasExt.h"
43
44#include <math.h>
45
46#include <TPad.h>
47#include <TLegend.h>
48#include <TCanvas.h>
49
50#include "MLog.h"
51#include "MLogManip.h"
52
53#include "MGeomCam.h"
54
55#include "MParList.h"
56
57#include "MBinning.h"
58
59#include "MHillas.h"
60#include "MHillasExt.h"
61#include "MHillasSrc.h"
62
63ClassImp(MHHillasExt);
64
65using namespace std;
66
67// --------------------------------------------------------------------------
68//
69// Setup four histograms for Width, Length
70//
71MHHillasExt::MHHillasExt(const char *name, const char *title)
72 : fGeom(0), fHillas(0), fHillasExt(0), fHilName("MHillasExt")
73{
74 //
75 // set the name and title of this object
76 //
77 fName = name ? name : "MHHillasExt";
78 fTitle = title ? title : "Container for extended Hillas histograms";
79
80 //
81 // loop over all Pixels and create two histograms
82 // one for the Low and one for the High gain
83 // connect all the histogram with the container fHist
84 //
85 fHAsym.UseCurrentStyle();
86 fHM3Long.UseCurrentStyle();
87 fHM3Trans.UseCurrentStyle();
88 fHSlopeL.UseCurrentStyle();
89
90 fHAsym.SetName("Asymmetry");
91 fHM3Long.SetName("M3l");
92 fHM3Trans.SetName("M3t");
93 fHSlopeL.SetName("SlopeL");
94
95 fHAsym.SetTitle("Asymmetry");
96 fHM3Long.SetTitle("3^{rd} Moment Longitudinal");
97 fHM3Trans.SetTitle("3^{rd} Moment Transverse");
98 fHSlopeL.SetTitle("Longitudinal time-slope vs. Dist");
99
100 fHAsym.SetXTitle("Asym [\\circ]");
101 fHM3Long.SetXTitle("3^{rd} M_{l} [\\circ]");
102 fHM3Trans.SetXTitle("3^{rd} M_{t} [\\circ]");
103 fHSlopeL.SetXTitle("D [\\circ]");
104
105 fHAsym.SetYTitle("Counts");
106 fHM3Long.SetYTitle("Counts");
107 fHM3Trans.SetYTitle("Counts");
108 fHSlopeL.SetYTitle("S_{l} [ns/\\circ]");
109
110 fHAsym.SetFillStyle(4000);
111 fHM3Long.SetFillStyle(4000);
112 fHM3Trans.SetFillStyle(4000);
113 //fHSlopeL.SetFillStyle(4000);
114
115 fHAsym.SetDirectory(NULL);
116 fHM3Long.SetDirectory(NULL);
117 fHM3Trans.SetDirectory(NULL);
118 fHSlopeL.SetDirectory(NULL);
119
120 fHM3Trans.SetLineColor(kBlue);
121
122 MBinning binsx, binsy;
123
124 binsx.SetEdges(51, -326, 326);
125 binsx.Apply(fHM3Long);
126 binsx.Apply(fHM3Trans);
127
128 binsx.SetEdges(51, -593, 593);
129 binsx.Apply(fHAsym);
130
131 binsx.SetEdges(100, 0, 445);
132 binsy.SetEdges(100, -0.04, 0.04);
133 MH::SetBinning(&fHSlopeL, &binsx, &binsy);
134}
135
136// --------------------------------------------------------------------------
137//
138// Setup the Binning for the histograms automatically if the correct
139// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
140// are found in the parameter list
141// Use this function if you want to set the conversion factor which
142// is used to convert the mm-scale in the camera plain into the deg-scale
143// used for histogram presentations. The conversion factor is part of
144// the camera geometry. Please create a corresponding MGeomCam container.
145//
146Bool_t MHHillasExt::SetupFill(const MParList *plist)
147{
148 fHillasExt = (MHillasExt*)plist->FindObject(fHilName, "MHillasExt");
149 if (!fHillasExt)
150 {
151 *fLog << err << fHilName << "[MHillasExt] not found in parameter list... aborting." << endl;
152 return kFALSE;
153 }
154
155 fHillas = (MHillas*)plist->FindObject("MHillas");
156 if (!fHillas)
157 {
158 *fLog << err << "MHillas not found in parameter list... aborting." << endl;
159 return kFALSE;
160 }
161
162 fGeom = (MGeomCam*)plist->FindObject("MGeomCam");
163 if (!fGeom)
164 {
165 *fLog << err << "MGeomCam not found... abort." << endl;
166 return kFALSE;
167 }
168
169 ApplyBinning(*plist, "Asym", &fHAsym);
170 ApplyBinning(*plist, "M3Long", &fHM3Long);
171 ApplyBinning(*plist, "M3Trans", &fHM3Trans);
172 ApplyBinning(*plist, "Dist", "Slope", &fHSlopeL);
173
174 return kTRUE;
175}
176
177// --------------------------------------------------------------------------
178//
179// Fill the four histograms with data from a MHillas-Container.
180// Be careful: Only call this with an object of type MHillas
181//
182Int_t MHHillasExt::Fill(const MParContainer *par, const Stat_t w)
183{
184 const MHillasSrc *src = (MHillasSrc*)par;
185
186 const Double_t scale = TMath::Sign(fGeom->GetConvMm2Deg(), (src ? src->GetCosDeltaAlpha() : 1));
187 const Double_t dist = src ? src->GetDist() : fHillas->GetDist0();
188
189 fHAsym.Fill(scale*fHillasExt->GetAsym(), w);
190 fHM3Long.Fill(scale*fHillasExt->GetM3Long(), w);
191 fHM3Trans.Fill(scale*fHillasExt->GetM3Trans(), w);
192 fHSlopeL.Fill(scale*dist, fHillasExt->GetSlopeLong()/scale, w);
193
194 return kTRUE;
195}
196
197// --------------------------------------------------------------------------
198//
199// Creates a new canvas and draws the four histograms into it.
200// Be careful: The histograms belongs to this object and won't get deleted
201// together with the canvas.
202//
203void MHHillasExt::Draw(Option_t *o)
204{
205 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
206 pad->SetBorderMode(0);
207
208 AppendPad("");
209
210 // FIXME: If same-option given make two independant y-axis!
211 const TString opt(o);
212 const Bool_t same = opt.Contains("same");
213
214 if (!same)
215 pad->Divide(2,2);
216 else
217 {
218 fHAsym.SetName("AsymmetrySame");
219 fHM3Long.SetName("M3lSame");
220 fHM3Trans.SetName("M3tSame");
221 fHSlopeL.SetName("SlopeLSame");
222
223 fHAsym.SetDirectory(0);
224 fHM3Long.SetDirectory(0);
225 fHM3Trans.SetDirectory(0);
226 fHSlopeL.SetDirectory(0);
227
228 fHM3Long.SetLineColor(kMagenta);
229 fHM3Trans.SetLineColor(kCyan);
230 fHAsym.SetLineColor(kBlue);
231 fHSlopeL.SetMarkerColor(kBlue);
232 }
233
234 pad->cd(1);
235 gPad->SetBorderMode(0);
236 gPad->SetGridx();
237 gPad->SetGridy();
238 RemoveFromPad("M3lSame");
239 RemoveFromPad("M3tSame");
240 MH::DrawSame(fHM3Long, fHM3Trans, "3^{rd} Moments", same);
241
242 pad->cd(3);
243 gPad->SetBorderMode(0);
244 gPad->SetGridx();
245 gPad->SetGridy();
246 RemoveFromPad("AsymmetrySame");
247 fHAsym.Draw(same?"same":"");
248
249 pad->cd(2);
250 gPad->SetBorderMode(0);
251 gPad->SetGridx();
252 gPad->SetGridy();
253 //RemoveFromPad("SlopeLSame");
254 //fHSlopeL.Draw(same?"same":"");
255 if (same)
256 {
257 TH2 *h=dynamic_cast<TH2*>(gPad->FindObject("SlopeL"));
258 if (h)
259 {
260 // This causes crashes in THistPainter::PaintTable
261 // if the z-axis is not kept. No idea why...
262 h->SetDrawOption("z");
263 h->SetMarkerColor(kBlack);
264 }
265
266 RemoveFromPad("SlopeLSame");
267 fHSlopeL.SetMarkerColor(kGreen);
268 fHSlopeL.Draw("same");
269 }
270 else
271 fHSlopeL.Draw("colz");
272
273 delete pad->GetPad(4);
274}
275
276TH1 *MHHillasExt::GetHistByName(const TString name) const
277{
278 if (name.Contains("Asym", TString::kIgnoreCase))
279 return const_cast<TH1F*>(&fHAsym);
280 if (name.Contains("M3L", TString::kIgnoreCase))
281 return const_cast<TH1F*>(&fHM3Long);
282 if (name.Contains("M3T", TString::kIgnoreCase))
283 return const_cast<TH1F*>(&fHM3Trans);
284 if (name.Contains("SlopeL", TString::kIgnoreCase))
285 return const_cast<TH2F*>(&fHSlopeL);
286
287 return NULL;
288}
Note: See TracBrowser for help on using the repository browser.