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