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 | // --------------------------------------------------------------------------
|
---|
53 | //
|
---|
54 | // Setup four histograms for Alpha, and Dist
|
---|
55 | //
|
---|
56 | MHHillasSrc::MHHillasSrc(const char *name, const char *title)
|
---|
57 | {
|
---|
58 | //
|
---|
59 | // set the name and title of this object
|
---|
60 | //
|
---|
61 | fName = name ? name : "MHHillasSrc";
|
---|
62 | fTitle = title ? title : "Container for Hillas histograms";
|
---|
63 |
|
---|
64 | //
|
---|
65 | // loop over all Pixels and create two histograms
|
---|
66 | // one for the Low and one for the High gain
|
---|
67 | // connect all the histogram with the container fHist
|
---|
68 | //
|
---|
69 | fAlpha = new TH1F("Alpha", "Alpha of Ellipse", 90, 0, 90);
|
---|
70 | fDist = new TH1F("Dist", "Dist of Ellipse", 100, 0, 600);
|
---|
71 |
|
---|
72 | fAlpha->SetDirectory(NULL);
|
---|
73 | fDist->SetDirectory(NULL);
|
---|
74 |
|
---|
75 | fAlpha->GetXaxis()->SetTitle("\\alpha [\\circ]");
|
---|
76 | fDist->GetXaxis()->SetTitle("Dist [mm]");
|
---|
77 |
|
---|
78 | fAlpha->GetYaxis()->SetTitle("Counts");
|
---|
79 | fDist->GetYaxis()->SetTitle("Counts");
|
---|
80 | }
|
---|
81 |
|
---|
82 | // --------------------------------------------------------------------------
|
---|
83 | //
|
---|
84 | // Delete the four histograms
|
---|
85 | //
|
---|
86 | MHHillasSrc::~MHHillasSrc()
|
---|
87 | {
|
---|
88 | delete fAlpha;
|
---|
89 | delete fDist;
|
---|
90 | }
|
---|
91 |
|
---|
92 | // --------------------------------------------------------------------------
|
---|
93 | //
|
---|
94 | // Setup the Binning for the histograms automatically if the correct
|
---|
95 | // instances of MBinning (with the names 'BinningAlpha' and 'BinningDist')
|
---|
96 | // are found in the parameter list
|
---|
97 | // Use this function if you want to set the conversion factor which
|
---|
98 | // is used to convert the mm-scale in the camera plain into the deg-scale
|
---|
99 | // used for histogram presentations. The conversion factor is part of
|
---|
100 | // the camera geometry. Please create a corresponding MGeomCam container.
|
---|
101 | //
|
---|
102 | Bool_t MHHillasSrc::SetupFill(const MParList *plist)
|
---|
103 | {
|
---|
104 | const MBinning* binsa = (MBinning*)plist->FindObject("BinningAlpha");
|
---|
105 | const MBinning* binsd = (MBinning*)plist->FindObject("BinningDist");
|
---|
106 | if (!binsa || !binsd)
|
---|
107 | {
|
---|
108 | *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl;
|
---|
109 | return kFALSE;
|
---|
110 | }
|
---|
111 |
|
---|
112 | SetBinning(fAlpha, binsa);
|
---|
113 | SetBinning(fDist, binsd);
|
---|
114 |
|
---|
115 | const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
|
---|
116 | if (!geom)
|
---|
117 | {
|
---|
118 | *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
|
---|
119 | return kTRUE;
|
---|
120 | }
|
---|
121 |
|
---|
122 | fDist->GetXaxis()->SetTitle("Dist [\\circ]");
|
---|
123 |
|
---|
124 | fMm2Deg = geom->GetConvMm2Deg();
|
---|
125 |
|
---|
126 | return kTRUE;
|
---|
127 | }
|
---|
128 |
|
---|
129 | // --------------------------------------------------------------------------
|
---|
130 | //
|
---|
131 | // Fill the four histograms with data from a MHillas-Container.
|
---|
132 | // Be careful: Only call this with an object of type MHillas
|
---|
133 | //
|
---|
134 | Bool_t MHHillasSrc::Fill(const MParContainer *par)
|
---|
135 | {
|
---|
136 | const MHillasSrc &h = *(MHillasSrc*)par;
|
---|
137 |
|
---|
138 | fAlpha->Fill(fabs(h.GetAlpha()));
|
---|
139 | fDist ->Fill(fUseMmScale ? h.GetDist() : fMm2Deg*h.GetDist());
|
---|
140 |
|
---|
141 | return kTRUE;
|
---|
142 | }
|
---|
143 |
|
---|
144 | // --------------------------------------------------------------------------
|
---|
145 | //
|
---|
146 | // Draw clones of all two histograms. So that the object can be deleted
|
---|
147 | // and the histograms are still visible in the canvas.
|
---|
148 | // The cloned object are deleted together with the canvas if the canvas is
|
---|
149 | // destroyed. If you want to handle dostroying the canvas you can get a
|
---|
150 | // pointer to it from this function
|
---|
151 | //
|
---|
152 | TObject *MHHillasSrc::DrawClone(Option_t *opt) const
|
---|
153 | {
|
---|
154 | TCanvas *c = MakeDefCanvas("Hillas", "Histograms of Source dependant Parameters",
|
---|
155 | 350, 500);
|
---|
156 | c->Divide(1, 2);
|
---|
157 |
|
---|
158 | // FIXME: Display Source position
|
---|
159 |
|
---|
160 | gROOT->SetSelectedPad(NULL);
|
---|
161 |
|
---|
162 | //
|
---|
163 | // This is necessary to get the expected bahviour of DrawClone
|
---|
164 | //
|
---|
165 | c->cd(1);
|
---|
166 | fAlpha->DrawCopy();
|
---|
167 |
|
---|
168 | c->cd(2);
|
---|
169 | fDist->DrawCopy();
|
---|
170 |
|
---|
171 | c->Modified();
|
---|
172 | c->Update();
|
---|
173 |
|
---|
174 | return c;
|
---|
175 | }
|
---|
176 |
|
---|
177 | // --------------------------------------------------------------------------
|
---|
178 | //
|
---|
179 | // Creates a new canvas and draws the two histograms into it.
|
---|
180 | // Be careful: The histograms belongs to this object and won't get deleted
|
---|
181 | // together with the canvas.
|
---|
182 | //
|
---|
183 | void MHHillasSrc::Draw(Option_t *)
|
---|
184 | {
|
---|
185 | if (!gPad)
|
---|
186 | MakeDefCanvas("Hillas", "Histograms of Src dependant Parameters", 350, 500);
|
---|
187 |
|
---|
188 | // FIXME: Display Source position
|
---|
189 |
|
---|
190 | gPad->Divide(1, 2);
|
---|
191 |
|
---|
192 | gPad->cd(1);
|
---|
193 | fAlpha->Draw();
|
---|
194 |
|
---|
195 | gPad->cd(2);
|
---|
196 | fDist->Draw();
|
---|
197 |
|
---|
198 | gPad->Modified();
|
---|
199 | gPad->Update();
|
---|
200 | }
|
---|