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-2001
|
---|
21 | !
|
---|
22 | !
|
---|
23 | \* ======================================================================== */
|
---|
24 |
|
---|
25 | ///////////////////////////////////////////////////////////////////////
|
---|
26 | //
|
---|
27 | // MHHillasSrc
|
---|
28 | //
|
---|
29 | // This class contains histograms for every Hillas parameter
|
---|
30 | //
|
---|
31 | ///////////////////////////////////////////////////////////////////////
|
---|
32 | #include "MHHillasSrc.h"
|
---|
33 |
|
---|
34 | #include <math.h>
|
---|
35 |
|
---|
36 | #include <TH1.h>
|
---|
37 | #include <TPad.h>
|
---|
38 | #include <TCanvas.h>
|
---|
39 |
|
---|
40 | #include "MLog.h"
|
---|
41 | #include "MLogManip.h"
|
---|
42 |
|
---|
43 | #include "MGeomCam.h"
|
---|
44 |
|
---|
45 | #include "MParList.h"
|
---|
46 |
|
---|
47 | #include "MHillas.h"
|
---|
48 | #include "MHillasSrc.h"
|
---|
49 |
|
---|
50 | ClassImp(MHHillasSrc);
|
---|
51 |
|
---|
52 | using namespace std;
|
---|
53 |
|
---|
54 | // --------------------------------------------------------------------------
|
---|
55 | //
|
---|
56 | // Setup four histograms for Alpha, and Dist
|
---|
57 | //
|
---|
58 | MHHillasSrc::MHHillasSrc(const char *name, const char *title)
|
---|
59 | : fUseMmScale(kTRUE)
|
---|
60 | {
|
---|
61 | //
|
---|
62 | // set the name and title of this object
|
---|
63 | //
|
---|
64 | fName = name ? name : "MHHillasSrc";
|
---|
65 | fTitle = title ? title : "Container for Hillas histograms";
|
---|
66 |
|
---|
67 | //
|
---|
68 | // loop over all Pixels and create two histograms
|
---|
69 | // one for the Low and one for the High gain
|
---|
70 | // connect all the histogram with the container fHist
|
---|
71 | //
|
---|
72 | fAlpha = new TH1F("Alpha", "Alpha of Ellipse", 90, -90, 90);
|
---|
73 | fDist = new TH1F("Dist", "Dist of Ellipse", 70, 0, 623);
|
---|
74 | fCosDA = new TH1F("CosDA", "cos(Delta,Alpha) of Ellipse", 101, -1, 1);
|
---|
75 | fDCA = new TH1F("DCA", "Distance of closest aproach", 101, -500, 500);
|
---|
76 | fDCADelta = new TH1F("DCADelta", "Angle between shower and x-axis", 80, 0, 360);
|
---|
77 |
|
---|
78 | fAlpha->SetDirectory(NULL);
|
---|
79 | fDist->SetDirectory(NULL);
|
---|
80 | fCosDA->SetDirectory(NULL);
|
---|
81 | fDCA->SetDirectory(NULL);
|
---|
82 | fDCADelta->SetDirectory(NULL);
|
---|
83 |
|
---|
84 | fAlpha->SetXTitle("\\alpha [\\circ]");
|
---|
85 | fDist->SetXTitle("Dist [mm]");
|
---|
86 | fCosDA->SetXTitle("cos(\\delta,\\alpha)");
|
---|
87 | fDCA->SetXTitle("DCA [\\circ]");
|
---|
88 | fDCADelta->SetXTitle("DCADelta [0, 2\\pi]");
|
---|
89 |
|
---|
90 | fAlpha->SetYTitle("Counts");
|
---|
91 | fDist->SetYTitle("Counts");
|
---|
92 | fCosDA->SetYTitle("Counts");
|
---|
93 | fDCA->SetYTitle("Counts");
|
---|
94 | fDCADelta->SetYTitle("Counts");
|
---|
95 |
|
---|
96 | fAlpha->SetMinimum(0);
|
---|
97 | fCosDA->SetMinimum(0);
|
---|
98 | fDCADelta->SetMinimum(0);
|
---|
99 | }
|
---|
100 |
|
---|
101 | // --------------------------------------------------------------------------
|
---|
102 | //
|
---|
103 | // Delete the four histograms
|
---|
104 | //
|
---|
105 | MHHillasSrc::~MHHillasSrc()
|
---|
106 | {
|
---|
107 | delete fAlpha;
|
---|
108 | delete fDist;
|
---|
109 | delete fCosDA;
|
---|
110 | delete fDCA;
|
---|
111 | delete fDCADelta;
|
---|
112 | }
|
---|
113 |
|
---|
114 | // --------------------------------------------------------------------------
|
---|
115 | //
|
---|
116 | // Setup the Binning for the histograms automatically if the correct
|
---|
117 | // instances of MBinning (with the names 'BinningAlpha' and 'BinningDist')
|
---|
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 | //
|
---|
124 | Bool_t MHHillasSrc::SetupFill(const MParList *plist)
|
---|
125 | {
|
---|
126 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
|
---|
127 | if (!geom)
|
---|
128 | *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
|
---|
129 | else
|
---|
130 | {
|
---|
131 | fMm2Deg = geom->GetConvMm2Deg();
|
---|
132 | SetMmScale(kFALSE);
|
---|
133 | }
|
---|
134 |
|
---|
135 | ApplyBinning(*plist, "Alpha", fAlpha);
|
---|
136 | ApplyBinning(*plist, "Dist", fDist);
|
---|
137 | ApplyBinning(*plist, "DCA", fDCA);
|
---|
138 | ApplyBinning(*plist, "DCADelta", fDCADelta);
|
---|
139 |
|
---|
140 | return kTRUE;
|
---|
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 MHHillasSrc::Fill(const MParContainer *par, const Stat_t w)
|
---|
149 | {
|
---|
150 | if (!par)
|
---|
151 | {
|
---|
152 | *fLog << err << "MHHillasSrc::Fill: Pointer (!=NULL) expected." << endl;
|
---|
153 | return kFALSE;
|
---|
154 | }
|
---|
155 |
|
---|
156 | const MHillasSrc &h = *(MHillasSrc*)par;
|
---|
157 |
|
---|
158 | fAlpha->Fill(h.GetAlpha(), w);
|
---|
159 | fDist ->Fill(fUseMmScale ? h.GetDist() : fMm2Deg*h.GetDist(), w);
|
---|
160 | fCosDA->Fill(h.GetCosDeltaAlpha(), w);
|
---|
161 | fDCA ->Fill(fUseMmScale ? h.GetDCA() : fMm2Deg*h.GetDCA(), w);
|
---|
162 | fDCADelta->Fill(h.GetDCADelta(), w);
|
---|
163 |
|
---|
164 | return kTRUE;
|
---|
165 | }
|
---|
166 |
|
---|
167 | // --------------------------------------------------------------------------
|
---|
168 | //
|
---|
169 | // Use this function to setup your own conversion factor between degrees
|
---|
170 | // and millimeters. The conversion factor should be the one calculated in
|
---|
171 | // MGeomCam. Use this function with Caution: You could create wrong values
|
---|
172 | // by setting up your own scale factor.
|
---|
173 | //
|
---|
174 | void MHHillasSrc::SetMm2Deg(Float_t mmdeg)
|
---|
175 | {
|
---|
176 | if (mmdeg<=0)
|
---|
177 | {
|
---|
178 | *fLog << warn << dbginf << "Warning - Conversion factor <= 0 - nonsense. Ignored." << endl;
|
---|
179 | return;
|
---|
180 | }
|
---|
181 |
|
---|
182 | if (fMm2Deg>0)
|
---|
183 | *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
|
---|
184 |
|
---|
185 | fMm2Deg = mmdeg;
|
---|
186 | }
|
---|
187 |
|
---|
188 | // --------------------------------------------------------------------------
|
---|
189 | //
|
---|
190 | // With this function you can convert the histogram ('on the fly') between
|
---|
191 | // degrees and millimeters.
|
---|
192 | //
|
---|
193 | void MHHillasSrc::SetMmScale(Bool_t mmscale)
|
---|
194 | {
|
---|
195 | if (fUseMmScale == mmscale)
|
---|
196 | return;
|
---|
197 |
|
---|
198 | if (fMm2Deg<0)
|
---|
199 | {
|
---|
200 | *fLog << warn << GetDescriptor() << ": Warning - Sorry, no conversion factor for conversion available." << endl;
|
---|
201 | return;
|
---|
202 | }
|
---|
203 |
|
---|
204 | const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
|
---|
205 | MH::ScaleAxis(fDist, scale);
|
---|
206 | MH::ScaleAxis(fDCA, scale);
|
---|
207 |
|
---|
208 | fDist->SetXTitle(mmscale ? "Dist [mm]" : "Dist [\\circ]");
|
---|
209 | fDCA ->SetXTitle(mmscale ? "DCA [mm]" : "DCA [\\circ]");
|
---|
210 |
|
---|
211 | fUseMmScale = mmscale;
|
---|
212 | }
|
---|
213 |
|
---|
214 | // --------------------------------------------------------------------------
|
---|
215 | //
|
---|
216 | // Creates a new canvas and draws the two histograms into it.
|
---|
217 | // Be careful: The histograms belongs to this object and won't get deleted
|
---|
218 | // together with the canvas.
|
---|
219 | //
|
---|
220 | void MHHillasSrc::Draw(Option_t *o)
|
---|
221 | {
|
---|
222 | TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
|
---|
223 | pad->SetBorderMode(0);
|
---|
224 |
|
---|
225 | AppendPad("");
|
---|
226 |
|
---|
227 | // FIXME: Display Source position
|
---|
228 |
|
---|
229 | // FIXME: If same-option given make two independant y-axis!
|
---|
230 | const TString opt(o);
|
---|
231 | const Bool_t same = opt.Contains("same");
|
---|
232 |
|
---|
233 | if (!same)
|
---|
234 | pad->Divide(2,2);
|
---|
235 | else
|
---|
236 | {
|
---|
237 | fAlpha->SetName("AlphaSame");
|
---|
238 | fDist ->SetName("DistSame");
|
---|
239 | fCosDA->SetName("CosDASame");
|
---|
240 | fDCA ->SetName("DCASame");
|
---|
241 | fDCADelta->SetName("DCADeltaSame");
|
---|
242 |
|
---|
243 | fAlpha->SetDirectory(0);
|
---|
244 | fDist ->SetDirectory(0);
|
---|
245 | fCosDA->SetDirectory(0);
|
---|
246 | fDCA ->SetDirectory(0);
|
---|
247 | fDCADelta->SetDirectory(0);
|
---|
248 |
|
---|
249 | fAlpha->SetLineColor(kGreen);
|
---|
250 | fDist->SetLineColor(kGreen);
|
---|
251 | fDCA->SetLineColor(kGreen);
|
---|
252 | fCosDA->SetLineColor(kGreen);
|
---|
253 | fDCADelta->SetLineColor(kGreen);
|
---|
254 | }
|
---|
255 |
|
---|
256 | pad->cd(1);
|
---|
257 | gPad->SetBorderMode(0);
|
---|
258 | RemoveFromPad("AlphaSame");
|
---|
259 | fAlpha->Draw(same?"same":"");
|
---|
260 |
|
---|
261 | pad->cd(2);
|
---|
262 | gPad->SetBorderMode(0);
|
---|
263 | RemoveFromPad("DistSame");
|
---|
264 | fDist->Draw(same?"same":"");
|
---|
265 |
|
---|
266 | pad->cd(3);
|
---|
267 | gPad->SetBorderMode(0);
|
---|
268 | RemoveFromPad("DCASame");
|
---|
269 | fDCA->Draw(same?"same":"");
|
---|
270 |
|
---|
271 | pad->cd(4);
|
---|
272 | gPad->SetBorderMode(0);
|
---|
273 |
|
---|
274 | TVirtualPad *p = gPad;
|
---|
275 | if (!same)
|
---|
276 | p->Divide(1,2);
|
---|
277 | p->cd(1);
|
---|
278 | gPad->SetBorderMode(0);
|
---|
279 | RemoveFromPad("CosDASame");
|
---|
280 | fCosDA->Draw(same?"same":"");
|
---|
281 |
|
---|
282 | p->cd(2);
|
---|
283 | gPad->SetBorderMode(0);
|
---|
284 | RemoveFromPad("DCADeltaSame");
|
---|
285 | fDCADelta->Draw(same?"same":"");
|
---|
286 | }
|
---|
287 |
|
---|
288 | void MHHillasSrc::Paint(Option_t *opt)
|
---|
289 | {
|
---|
290 | if (fCosDA->GetEntries()==0)
|
---|
291 | return;
|
---|
292 |
|
---|
293 | TVirtualPad *savepad = gPad;
|
---|
294 | savepad->cd(4);
|
---|
295 | gPad->SetLogy();
|
---|
296 | gPad = savepad;
|
---|
297 | }
|
---|
298 |
|
---|
299 | TH1 *MHHillasSrc::GetHistByName(const TString name) const
|
---|
300 | {
|
---|
301 | if (name.Contains("Alpha", TString::kIgnoreCase))
|
---|
302 | return fAlpha;
|
---|
303 | if (name.Contains("Dist", TString::kIgnoreCase))
|
---|
304 | return fDist;
|
---|
305 | if (name.Contains("CosDA", TString::kIgnoreCase))
|
---|
306 | return fCosDA;
|
---|
307 |
|
---|
308 | return NULL;
|
---|
309 | }
|
---|