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

Last change on this file since 2919 was 2414, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 7.5 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#include "MHHillasExt.h"
33
34#include <math.h>
35
36#include <TH1.h>
37#include <TPad.h>
38#include <TLegend.h>
39#include <TCanvas.h>
40
41#include "MLog.h"
42#include "MLogManip.h"
43
44#include "MGeomCam.h"
45
46#include "MParList.h"
47
48#include "MBinning.h"
49#include "MHillasExt.h"
50#include "MHillasSrc.h"
51
52ClassImp(MHHillasExt);
53
54using namespace std;
55
56// --------------------------------------------------------------------------
57//
58// Setup four histograms for Width, Length
59//
60MHHillasExt::MHHillasExt(const char *name, const char *title)
61 : fMm2Deg(1), fUseMmScale(kTRUE), fHilName("MHillasExt")
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 //
70 // loop over all Pixels and create two histograms
71 // one for the Low and one for the High gain
72 // connect all the histogram with the container fHist
73 //
74 fHAsym.SetDirectory(NULL);
75 fHM3Long.SetDirectory(NULL);
76 fHM3Trans.SetDirectory(NULL);
77
78 fHAsym.UseCurrentStyle();
79 fHM3Long.UseCurrentStyle();
80 fHM3Trans.UseCurrentStyle();
81
82 fHAsym.SetName("Asymmetry");
83 fHM3Long.SetName("3rd Mom Long");
84 fHM3Trans.SetName("3rd Mom Trans");
85
86 fHAsym.SetTitle("Asymmetry");
87 fHM3Long.SetTitle("3^{rd} Moment Longitudinal");
88 fHM3Trans.SetTitle("3^{rd} Moment Transverse");
89
90 fHAsym.SetXTitle("Asym [mm]");
91 fHM3Long.SetXTitle("3^{rd} M_{l} [mm]");
92 fHM3Trans.SetXTitle("3^{rd} M_{t} [mm]");
93
94 fHAsym.SetYTitle("Counts");
95 fHM3Long.SetYTitle("Counts");
96 fHM3Trans.SetYTitle("Counts");
97
98 fHAsym.SetFillStyle(4000);
99 fHM3Long.SetFillStyle(4000);
100 fHM3Trans.SetFillStyle(4000);
101
102 fHM3Trans.SetLineColor(kBlue);
103
104 MBinning bins;
105
106 bins.SetEdges(101, -326, 326);
107 bins.Apply(fHM3Long);
108 bins.Apply(fHM3Trans);
109
110 bins.SetEdges(101, -593, 593);
111 bins.Apply(fHAsym);
112}
113
114// --------------------------------------------------------------------------
115//
116// Setup the Binning for the histograms automatically if the correct
117// instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
118// are found in the parameter list
119// Use this function if you want to set the conversion factor which
120// is used to convert the mm-scale in the camera plain into the deg-scale
121// used for histogram presentations. The conversion factor is part of
122// the camera geometry. Please create a corresponding MGeomCam container.
123//
124Bool_t MHHillasExt::SetupFill(const MParList *plist)
125{
126 fHillasExt = (MHillasExt*)plist->FindObject(fHilName, "MHillasExt");
127 if (!fHillasExt)
128 {
129 *fLog << err << fHilName << "[MHillasExt] not found in parameter list... aborting." << endl;
130 return kFALSE;
131 }
132
133 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
134 if (!geom)
135 *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl;
136 else
137 {
138 fMm2Deg = geom->GetConvMm2Deg();
139 SetMmScale(kFALSE);
140 }
141
142 ApplyBinning(*plist, "Asym", &fHAsym);
143 ApplyBinning(*plist, "M3Long", &fHM3Long);
144 ApplyBinning(*plist, "M3Trans", &fHM3Trans);
145
146 return kTRUE;
147}
148
149// --------------------------------------------------------------------------
150//
151// Fill the four histograms with data from a MHillas-Container.
152// Be careful: Only call this with an object of type MHillas
153//
154Bool_t MHHillasExt::Fill(const MParContainer *par, const Stat_t w)
155{
156 const MHillasSrc *src = (MHillasSrc*)par;
157
158 const Double_t scale = TMath::Sign(fUseMmScale?1:fMm2Deg, (src ? src->GetCosDeltaAlpha() : 1));
159
160 fHAsym.Fill(scale*fHillasExt->GetAsym(), w);
161 fHM3Long.Fill(scale*fHillasExt->GetM3Long(), w);
162 fHM3Trans.Fill(scale*fHillasExt->GetM3Trans(), w);
163 //fHAsymna.Fill(scale*ext.GetAsymna());
164 //fHAsym0.Fill(scale*ext.GetAsym0());
165
166 return kTRUE;
167}
168
169// --------------------------------------------------------------------------
170//
171// With this function you can convert the histogram ('on the fly') between
172// degrees and millimeters.
173//
174void MHHillasExt::SetMmScale(Bool_t mmscale)
175{
176 if (fUseMmScale == mmscale)
177 return;
178
179 if (fMm2Deg<0)
180 {
181 *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
182 return;
183 }
184
185 const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
186 MH::ScaleAxis(&fHAsym, scale);
187 MH::ScaleAxis(&fHM3Long, scale);
188 MH::ScaleAxis(&fHM3Trans, scale);
189
190 if (mmscale)
191 {
192 fHAsym.SetXTitle("Asym [mm]");
193 fHM3Long.SetXTitle("3^{rd} M_{l} [mm]");
194 fHM3Trans.SetXTitle("3^{rd} M_{t} [mm]");
195 }
196 else
197 {
198 fHAsym.SetXTitle("Asym [\\circ]");
199 fHM3Long.SetXTitle("3^{rd} M_{l} [\\circ]");
200 fHM3Trans.SetXTitle("3^{rd} M_{t} [\\circ]");
201 }
202
203 fUseMmScale = mmscale;
204}
205
206// --------------------------------------------------------------------------
207//
208// Use this function to setup your own conversion factor between degrees
209// and millimeters. The conversion factor should be the one calculated in
210// MGeomCam. Use this function with Caution: You could create wrong values
211// by setting up your own scale factor.
212//
213void MHHillasExt::SetMm2Deg(Float_t mmdeg)
214{
215 if (mmdeg<0)
216 {
217 *fLog << warn << dbginf << "Warning - Conversion factor < 0 - nonsense. Ignored." << endl;
218 return;
219 }
220
221 if (fMm2Deg>=0)
222 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
223
224 fMm2Deg = mmdeg;
225}
226
227// --------------------------------------------------------------------------
228//
229// Creates a new canvas and draws the four histograms into it.
230// Be careful: The histograms belongs to this object and won't get deleted
231// together with the canvas.
232//
233void MHHillasExt::Draw(Option_t *)
234{
235 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
236 pad->SetBorderMode(0);
237
238 AppendPad("");
239
240 pad->Divide(2, 1);
241
242 pad->cd(1);
243 gPad->SetBorderMode(0);
244 MH::DrawSame(fHM3Long, fHM3Trans, "3^{rd} Moments");
245
246 pad->cd(2);
247 gPad->SetBorderMode(0);
248 fHAsym.Draw();
249
250 pad->Modified();
251 pad->Update();
252}
253
254TH1 *MHHillasExt::GetHistByName(const TString name)
255{
256 if (name.Contains("Asym", TString::kIgnoreCase))
257 return &fHAsym;
258 if (name.Contains("M3Long", TString::kIgnoreCase))
259 return &fHM3Long;
260 if (name.Contains("M3Trans", TString::kIgnoreCase))
261 return &fHM3Trans;
262
263 return NULL;
264}
Note: See TracBrowser for help on using the repository browser.