source: trunk/MagicSoft/Mars/mhist/MHHillasExt.cc@ 1765

Last change on this file since 1765 was 1762, checked in by wittek, 22 years ago
*** empty log message ***
File size: 9.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@uni-sw.gwdg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2002
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHHillasExt
28//
29// This class contains histograms for every Hillas parameter
30//
31/////////////////////////////////////////////////////////////////////////////
32
33#include "MHHillasExt.h"
34
35#include <math.h>
36
37#include <TH1.h>
38#include <TPad.h>
39#include <TLegend.h>
40#include <TCanvas.h>
41
42#include "MLog.h"
43#include "MLogManip.h"
44
45#include "MGeomCam.h"
46
47#include "MParList.h"
48
49#include "MBinning.h"
50#include "MHillasExt.h"
51#include "MHillasSrc.h"
52
53ClassImp(MHHillasExt);
54
55// --------------------------------------------------------------------------
56//
57// Setup four histograms for Width, Length
58//
59MHHillasExt::MHHillasExt(const char *name, const char *title,
60 const char *hil)
61 : fMm2Deg(1), fUseMmScale(kTRUE)
62{
63 //
64 // set the name and title of this object
65 //
66 fName = name ? name : "MHHillasExt";
67 fTitle = title ? title : "Container for extended Hillas histograms";
68
69 fHilName = hil;
70 //*fLog << "MHHillasExt : fHilName = " << fHilName << 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 fHConc1.SetLineColor(kBlue);
78 fHConc.SetFillStyle(0);
79
80 fHConc.SetDirectory(NULL);
81 fHConc1.SetDirectory(NULL);
82 fHAsym.SetDirectory(NULL);
83 fHM3Long.SetDirectory(NULL);
84 fHM3Trans.SetDirectory(NULL);
85
86 fHConc.SetName("Conc2");
87 fHConc1.SetName("Conc1");
88 fHAsym.SetName("Asymmetry");
89 fHM3Long.SetName("3rd Mom Long");
90 fHM3Trans.SetName("3rd Mom Trans");
91
92 fHConc.SetTitle("Ratio: Conc");
93 fHConc1.SetTitle("Ratio: Conc1");
94 fHAsym.SetTitle("Asymmetry");
95 fHM3Long.SetTitle("3^{rd} Moment Longitudinal");
96 fHM3Trans.SetTitle("3^{rd} Moment Transverse");
97
98 fHConc.SetXTitle("Ratio");
99 fHConc1.SetXTitle("Ratio");
100 fHAsym.SetXTitle("Asym [mm]");
101 fHM3Long.SetXTitle("3^{rd} M_{l} [mm]");
102 fHM3Trans.SetXTitle("3^{rd} M_{t} [mm]");
103
104 fHConc.SetYTitle("Counts");
105 fHConc1.SetYTitle("Counts");
106 fHAsym.SetYTitle("Counts");
107 fHM3Long.SetYTitle("Counts");
108 fHM3Trans.SetYTitle("Counts");
109
110
111 MBinning bins;
112
113 bins.SetEdges(100, 0, 1);
114 bins.Apply(fHConc);
115 bins.Apply(fHConc1);
116
117 bins.SetEdges(101, -326, 326);
118 bins.Apply(fHM3Long);
119 bins.Apply(fHM3Trans);
120
121 bins.SetEdges(101, -593, 593);
122 bins.Apply(fHAsym);
123}
124
125// --------------------------------------------------------------------------
126//
127// Delete the four histograms
128//
129MHHillasExt::~MHHillasExt()
130{
131}
132
133// --------------------------------------------------------------------------
134//
135// Setup the Binning for the histograms automatically if the correct
136// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
137// are found in the parameter list
138// Use this function if you want to set the conversion factor which
139// is used to convert the mm-scale in the camera plain into the deg-scale
140// used for histogram presentations. The conversion factor is part of
141// the camera geometry. Please create a corresponding MGeomCam container.
142//
143Bool_t MHHillasExt::SetupFill(const MParList *plist)
144{
145 TObject *obj = plist->FindObject(fHilName, "MHillas");
146 if (!obj)
147 {
148 *fLog << err << dbginf << "Sorry '" << fHilName
149 << "' not found in parameter list... aborting." << endl;
150 return kFALSE;
151 }
152 if (!obj->InheritsFrom(MHillasExt::Class()))
153 {
154 *fLog << err << dbginf << "Sorry 'MHillas' doesn't inherit from MHillasExt... aborting." << endl;
155 return kFALSE;
156 }
157 fHillasExt = (MHillasExt*)obj;
158
159 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
160 if (!geom)
161 *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl;
162 else
163 {
164 fMm2Deg = geom->GetConvMm2Deg();
165 SetMmScale(kFALSE);
166 }
167
168 ApplyBinning(*plist, "Conc", &fHConc);
169 ApplyBinning(*plist, "Conc1", &fHConc1);
170 ApplyBinning(*plist, "Asym", &fHAsym);
171 ApplyBinning(*plist, "M3Long", &fHM3Long);
172 ApplyBinning(*plist, "M3Trans", &fHM3Trans);
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//
182Bool_t MHHillasExt::Fill(const MParContainer *par)
183{
184 const MHillasSrc *src = (MHillasSrc*)par;
185
186 const Double_t scale = TMath::Sign(fUseMmScale?1:fMm2Deg, src ? src->GetCosDeltaAlpha() : 1);
187
188 fHConc.Fill(fHillasExt->GetConc());
189 fHConc1.Fill(fHillasExt->GetConc1());
190
191 fHAsym.Fill(scale*fHillasExt->GetAsym());
192 fHM3Long.Fill(scale*fHillasExt->GetM3Long());
193 fHM3Trans.Fill(scale*fHillasExt->GetM3Trans());
194 //fHAsymna.Fill(scale*ext.GetAsymna());
195 //fHAsym0.Fill(scale*ext.GetAsym0());
196
197 return kTRUE;
198}
199
200// --------------------------------------------------------------------------
201//
202// With this function you can convert the histogram ('on the fly') between
203// degrees and millimeters.
204//
205void MHHillasExt::SetMmScale(Bool_t mmscale)
206{
207 if (fUseMmScale == mmscale)
208 return;
209
210 if (fMm2Deg<0)
211 {
212 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
213 return;
214 }
215
216 const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
217 MH::ScaleAxis(&fHAsym, scale);
218 MH::ScaleAxis(&fHM3Long, scale);
219 MH::ScaleAxis(&fHM3Trans, scale);
220
221 if (mmscale)
222 {
223 fHAsym.SetXTitle("Asym [mm]");
224 fHM3Long.SetXTitle("3^{rd} M_{l} [mm]");
225 fHM3Trans.SetXTitle("3^{rd} M_{t} [mm]");
226 }
227 else
228 {
229 fHAsym.SetXTitle("Asym [\\circ]");
230 fHM3Long.SetXTitle("3^{rd} M_{l} [\\circ]");
231 fHM3Trans.SetXTitle("3^{rd} M_{t} [\\circ]");
232 }
233
234 fUseMmScale = mmscale;
235}
236
237// --------------------------------------------------------------------------
238//
239// Use this function to setup your own conversion factor between degrees
240// and millimeters. The conversion factor should be the one calculated in
241// MGeomCam. Use this function with Caution: You could create wrong values
242// by setting up your own scale factor.
243//
244void MHHillasExt::SetMm2Deg(Float_t mmdeg)
245{
246 if (mmdeg<0)
247 {
248 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
249 return;
250 }
251
252 if (fMm2Deg>=0)
253 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
254
255 fMm2Deg = mmdeg;
256}
257
258// --------------------------------------------------------------------------
259//
260// Draw clones of all four histograms. So that the object can be deleted
261// and the histograms are still visible in the canvas.
262// The cloned object are deleted together with the canvas if the canvas is
263// destroyed. If you want to handle dostroying the canvas you can get a
264// pointer to it from this function
265//
266TObject *MHHillasExt::DrawClone(Option_t *opt) const
267{
268 TCanvas &c = *MakeDefCanvas(this, 720, 540);
269 c.Divide(2, 2);
270
271 gROOT->SetSelectedPad(NULL);
272
273 c.cd(1);
274 DrawCopy(fHConc1, fHConc, "Concentrations");
275
276 c.cd(2);
277 ((TH1&)fHAsym).DrawCopy();
278
279 c.cd(3);
280 ((TH1&)fHM3Long).DrawCopy();
281
282 c.cd(4);
283 ((TH1&)fHM3Trans).DrawCopy();
284
285 c.Modified();
286 c.Update();
287
288 return &c;
289}
290
291// --------------------------------------------------------------------------
292//
293// Creates a new canvas and draws the four histograms into it.
294// Be careful: The histograms belongs to this object and won't get deleted
295// together with the canvas.
296//
297void MHHillasExt::Draw(Option_t *)
298{
299 if (!gPad)
300 MakeDefCanvas(this, 720, 540);
301
302 gPad->Divide(2, 2);
303
304 gPad->cd(1);
305 MH::Draw(fHConc1, fHConc, "Concentrations");
306
307 gPad->cd(2);
308 fHAsym.Draw();
309
310 gPad->cd(3);
311 fHM3Long.Draw();
312
313 gPad->cd(4);
314 fHM3Trans.Draw();
315
316 gPad->Modified();
317 gPad->Update();
318}
319
320TH1 *MHHillasExt::GetHistByName(const TString name)
321{
322 if (name.Contains("Conc", TString::kIgnoreCase))
323 return &fHConc;
324 if (name.Contains("Conc1", TString::kIgnoreCase))
325 return &fHConc1;
326 if (name.Contains("Asym", TString::kIgnoreCase))
327 return &fHAsym;
328 if (name.Contains("M3Long", TString::kIgnoreCase))
329 return &fHM3Long;
330 if (name.Contains("M3Trans", TString::kIgnoreCase))
331 return &fHM3Trans;
332
333 return NULL;
334}
Note: See TracBrowser for help on using the repository browser.