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 <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 |
|
---|
52 | ClassImp(MHHillasExt);
|
---|
53 |
|
---|
54 | // --------------------------------------------------------------------------
|
---|
55 | //
|
---|
56 | // Setup four histograms for Width, Length
|
---|
57 | //
|
---|
58 | MHHillasExt::MHHillasExt(const char *name, const char *title)
|
---|
59 | {
|
---|
60 | //
|
---|
61 | // set the name and title of this object
|
---|
62 | //
|
---|
63 | fName = name ? name : "MHHillasExt";
|
---|
64 | fTitle = title ? title : "Container for Hillas (ext) histograms";
|
---|
65 |
|
---|
66 | //
|
---|
67 | // loop over all Pixels and create two histograms
|
---|
68 | // one for the Low and one for the High gain
|
---|
69 | // connect all the histogram with the container fHist
|
---|
70 | //
|
---|
71 | fHConc.SetDirectory(NULL);
|
---|
72 | fHConc1.SetDirectory(NULL);
|
---|
73 | fHAsym.SetDirectory(NULL);
|
---|
74 | fHM3Long.SetDirectory(NULL);
|
---|
75 | fHM3Trans.SetDirectory(NULL);
|
---|
76 |
|
---|
77 | fHConc.SetTitle("Ratio: Conc");
|
---|
78 | fHConc1.SetTitle("Ratio: Conc1");
|
---|
79 | fHAsym.SetTitle("Asymmetry");
|
---|
80 | fHM3Long.SetTitle("3^{rd} Moment Longitudinal");
|
---|
81 | fHM3Trans.SetTitle("3^{rd} Moment Transverse");
|
---|
82 |
|
---|
83 | fHConc.SetXTitle("Ratio");
|
---|
84 | fHConc1.SetXTitle("Ratio");
|
---|
85 | fHAsym.SetXTitle("Asym [mm]");
|
---|
86 | fHM3Long.SetXTitle("3^{rd} M_{l} [mm]");
|
---|
87 | fHM3Trans.SetXTitle("3^{rd} M_{t} [mm]");
|
---|
88 |
|
---|
89 | fHConc.SetYTitle("Counts");
|
---|
90 | fHConc1.SetYTitle("Counts");
|
---|
91 | fHAsym.SetYTitle("Counts");
|
---|
92 | fHM3Long.SetYTitle("Counts");
|
---|
93 | fHM3Trans.SetYTitle("Counts");
|
---|
94 |
|
---|
95 |
|
---|
96 | MBinning bins;
|
---|
97 |
|
---|
98 | bins.SetEdges(100, 0, 1);
|
---|
99 | bins.Apply(fHConc);
|
---|
100 | bins.Apply(fHConc1);
|
---|
101 |
|
---|
102 | bins.SetEdges(101, -326, 326);
|
---|
103 | bins.Apply(fHM3Long);
|
---|
104 | bins.Apply(fHM3Trans);
|
---|
105 |
|
---|
106 | bins.SetEdges(89, -297, 593);
|
---|
107 | bins.Apply(fHAsym);
|
---|
108 | }
|
---|
109 |
|
---|
110 | // --------------------------------------------------------------------------
|
---|
111 | //
|
---|
112 | // Delete the four histograms
|
---|
113 | //
|
---|
114 | MHHillasExt::~MHHillasExt()
|
---|
115 | {
|
---|
116 | }
|
---|
117 |
|
---|
118 | // --------------------------------------------------------------------------
|
---|
119 | //
|
---|
120 | // Setup the Binning for the histograms automatically if the correct
|
---|
121 | // instances of MBinning (with the names 'BinningWidth' and 'BinningLength')
|
---|
122 | // are found in the parameter list
|
---|
123 | // Use this function if you want to set the conversion factor which
|
---|
124 | // is used to convert the mm-scale in the camera plain into the deg-scale
|
---|
125 | // used for histogram presentations. The conversion factor is part of
|
---|
126 | // the camera geometry. Please create a corresponding MGeomCam container.
|
---|
127 | //
|
---|
128 | Bool_t MHHillasExt::SetupFill(const MParList *plist)
|
---|
129 | {
|
---|
130 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
|
---|
131 | if (!geom)
|
---|
132 | *fLog << warn << GetDescriptor() << ": No Camera Geometry available. Using mm-scale for histograms." << endl;
|
---|
133 | else
|
---|
134 | {
|
---|
135 | fMm2Deg = geom->GetConvMm2Deg();
|
---|
136 | SetMmScale(kFALSE);
|
---|
137 | }
|
---|
138 |
|
---|
139 | ApplyBinning(*plist, "Conc", &fHConc);
|
---|
140 | ApplyBinning(*plist, "Conc1", &fHConc1);
|
---|
141 | ApplyBinning(*plist, "Asym", &fHAsym);
|
---|
142 | ApplyBinning(*plist, "M3Long", &fHM3Long);
|
---|
143 | ApplyBinning(*plist, "M3Trans", &fHM3Trans);
|
---|
144 |
|
---|
145 | fHillasSrc = (MHillasSrc*)plist->FindObject("HillasSource");
|
---|
146 | if (!fHillasSrc)
|
---|
147 | *fLog << warn << dbginf << "Sorry 'HillasSource' [MHillasSrc] not found. Histograms will be source independant." << endl;
|
---|
148 |
|
---|
149 | return MHHillas::SetupFill(plist);
|
---|
150 | }
|
---|
151 |
|
---|
152 | // --------------------------------------------------------------------------
|
---|
153 | //
|
---|
154 | // Fill the four histograms with data from a MHillas-Container.
|
---|
155 | // Be careful: Only call this with an object of type MHillas
|
---|
156 | //
|
---|
157 | Bool_t MHHillasExt::Fill(const MParContainer *par)
|
---|
158 | {
|
---|
159 | const MHillasExt &ext = *(MHillasExt*)par;
|
---|
160 |
|
---|
161 | const Double_t scale = fHillasSrc ? TMath::Sign(fUseMmScale?1:fMm2Deg, fHillasSrc->GetCosDeltaAlpha()) : 1;
|
---|
162 |
|
---|
163 | fHConc.Fill(ext.GetConc());
|
---|
164 | fHConc1.Fill(ext.GetConc1());
|
---|
165 |
|
---|
166 | fHAsym.Fill(scale*ext.GetAsym());
|
---|
167 | fHM3Long.Fill(scale*ext.GetM3Long());
|
---|
168 | fHM3Trans.Fill(scale*ext.GetM3Trans());
|
---|
169 | //fHAsymna.Fill(scale*ext.GetAsymna());
|
---|
170 | //fHAsym0.Fill(scale*ext.GetAsym0());
|
---|
171 |
|
---|
172 | return MHHillas::Fill(par);
|
---|
173 | }
|
---|
174 |
|
---|
175 | // --------------------------------------------------------------------------
|
---|
176 | //
|
---|
177 | // With this function you can convert the histogram ('on the fly') between
|
---|
178 | // degrees and millimeters.
|
---|
179 | //
|
---|
180 | void MHHillasExt::SetMmScale(Bool_t mmscale)
|
---|
181 | {
|
---|
182 | if (fUseMmScale == mmscale)
|
---|
183 | return;
|
---|
184 |
|
---|
185 | if (fMm2Deg<0)
|
---|
186 | {
|
---|
187 | *fLog << warn << dbginf << "Warning - Sorry, no conversion factor for conversion available." << endl;
|
---|
188 | return;
|
---|
189 | }
|
---|
190 |
|
---|
191 | const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
|
---|
192 | MH::ScaleAxis(&fHAsym, scale);
|
---|
193 | MH::ScaleAxis(&fHM3Long, scale);
|
---|
194 | MH::ScaleAxis(&fHM3Trans, scale);
|
---|
195 |
|
---|
196 | if (mmscale)
|
---|
197 | {
|
---|
198 | fHAsym.SetXTitle("Asym [mm]");
|
---|
199 | fHM3Long.SetXTitle("3^{rd} M_{l}[mm]");
|
---|
200 | fHM3Trans.SetXTitle("3^{rd} M_{t} [mm]");
|
---|
201 | }
|
---|
202 | else
|
---|
203 | {
|
---|
204 | fHAsym.SetXTitle("Asym [\\circ]");
|
---|
205 | fHM3Long.SetXTitle("3^{rd} M_{l} [\\circ]");
|
---|
206 | fHM3Trans.SetXTitle("3^{rd} M_{t} [\\circ]");
|
---|
207 | }
|
---|
208 |
|
---|
209 | MHHillas::SetMmScale(mmscale);
|
---|
210 | }
|
---|
211 |
|
---|
212 | // --------------------------------------------------------------------------
|
---|
213 | //
|
---|
214 | // Draw clones of all four histograms. So that the object can be deleted
|
---|
215 | // and the histograms are still visible in the canvas.
|
---|
216 | // The cloned object are deleted together with the canvas if the canvas is
|
---|
217 | // destroyed. If you want to handle dostroying the canvas you can get a
|
---|
218 | // pointer to it from this function
|
---|
219 | //
|
---|
220 | TObject *MHHillasExt::DrawClone(Option_t *opt) const
|
---|
221 | {
|
---|
222 | TCanvas &c = *MakeDefCanvas("Hillas", "Histograms of Hillas Parameters",
|
---|
223 | 720, 810);
|
---|
224 | c.Divide(2, 3);
|
---|
225 |
|
---|
226 | gROOT->SetSelectedPad(NULL);
|
---|
227 |
|
---|
228 | //
|
---|
229 | // This is necessary to get the expected bahviour of DrawClone
|
---|
230 | //
|
---|
231 | c.cd(1);
|
---|
232 | ((TH1F&)fHConc).DrawCopy();
|
---|
233 |
|
---|
234 | c.cd(2);
|
---|
235 | ((TH1F&)fHConc1).DrawCopy();
|
---|
236 |
|
---|
237 | c.cd(5);
|
---|
238 | ((TH1F&)fHAsym).DrawCopy();
|
---|
239 |
|
---|
240 | c.cd(3);
|
---|
241 | ((TH1F&)fHM3Long).DrawCopy();
|
---|
242 |
|
---|
243 | c.cd(4);
|
---|
244 | ((TH1F&)fHM3Trans).DrawCopy();
|
---|
245 |
|
---|
246 | c.Modified();
|
---|
247 | c.Update();
|
---|
248 |
|
---|
249 | MHHillas::DrawClone();
|
---|
250 |
|
---|
251 | return &c;
|
---|
252 | }
|
---|
253 |
|
---|
254 | // --------------------------------------------------------------------------
|
---|
255 | //
|
---|
256 | // Creates a new canvas and draws the four histograms into it.
|
---|
257 | // Be careful: The histograms belongs to this object and won't get deleted
|
---|
258 | // together with the canvas.
|
---|
259 | //
|
---|
260 | void MHHillasExt::Draw(Option_t *)
|
---|
261 | {
|
---|
262 | if (!gPad)
|
---|
263 | MakeDefCanvas("Hillas", "Histograms of Hillas Parameters",
|
---|
264 | 720, 810);
|
---|
265 |
|
---|
266 | gPad->Divide(2, 3);
|
---|
267 |
|
---|
268 | gPad->cd(1);
|
---|
269 | fHConc.DrawCopy();
|
---|
270 |
|
---|
271 | gPad->cd(2);
|
---|
272 | fHConc1.DrawCopy();
|
---|
273 |
|
---|
274 | gPad->cd(5);
|
---|
275 | fHAsym.DrawCopy();
|
---|
276 |
|
---|
277 | gPad->cd(3);
|
---|
278 | fHM3Long.DrawCopy();
|
---|
279 |
|
---|
280 | gPad->cd(4);
|
---|
281 | fHM3Trans.DrawCopy();
|
---|
282 |
|
---|
283 | gPad->Modified();
|
---|
284 | gPad->Update();
|
---|
285 |
|
---|
286 | MHHillas::DrawClone();
|
---|
287 | }
|
---|