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