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