source: trunk/MagicSoft/Mars/mimage/MHHillasSrc.cc@ 2445

Last change on this file since 2445 was 2438, checked in by tbretz, 21 years ago
*** empty log message ***
File size: 6.8 KB
Line 
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
50ClassImp(MHHillasSrc);
51
52using namespace std;
53
54// --------------------------------------------------------------------------
55//
56// Setup four histograms for Alpha, and Dist
57//
58MHHillasSrc::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", 181, -90, 90);
73 fDist = new TH1F("Dist", "Dist of Ellipse", 100, 0, 445);
74 fCosDA = new TH1F("CosDA", "cos(Delta,Alpha) of Ellipse", 101, -1, 1);
75
76 fAlpha->SetDirectory(NULL);
77 fDist->SetDirectory(NULL);
78 fCosDA->SetDirectory(NULL);
79
80 fAlpha->SetXTitle("\\alpha [\\circ]");
81 fDist->SetXTitle("Dist [mm]");
82 fCosDA->SetXTitle("cos(\\delta,\\alpha)");
83
84 fAlpha->SetYTitle("Counts");
85 fDist->SetYTitle("Counts");
86 fCosDA->SetYTitle("Counts");
87}
88
89// --------------------------------------------------------------------------
90//
91// Delete the four histograms
92//
93MHHillasSrc::~MHHillasSrc()
94{
95 delete fAlpha;
96 delete fDist;
97 delete fCosDA;
98}
99
100// --------------------------------------------------------------------------
101//
102// Setup the Binning for the histograms automatically if the correct
103// instances of MBinning (with the names 'BinningAlpha' and 'BinningDist')
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//
110Bool_t MHHillasSrc::SetupFill(const MParList *plist)
111{
112 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
113 if (!geom)
114 *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
115 else
116 {
117 fMm2Deg = geom->GetConvMm2Deg();
118 SetMmScale(kFALSE);
119 }
120
121 ApplyBinning(*plist, "Alpha", fAlpha);
122 ApplyBinning(*plist, "Dist", fDist);
123
124 return kTRUE;
125}
126
127// --------------------------------------------------------------------------
128//
129// Fill the four histograms with data from a MHillas-Container.
130// Be careful: Only call this with an object of type MHillas
131//
132Bool_t MHHillasSrc::Fill(const MParContainer *par, const Stat_t w)
133{
134 if (!par)
135 {
136 *fLog << err << "MHHillasSrc::Fill: Pointer (!=NULL) expected." << endl;
137 return kFALSE;
138 }
139
140 const MHillasSrc &h = *(MHillasSrc*)par;
141
142 fAlpha->Fill(h.GetAlpha(), w);
143 fDist ->Fill(fUseMmScale ? h.GetDist() : fMm2Deg*h.GetDist(), w);
144 fCosDA->Fill(h.GetCosDeltaAlpha(), w);
145
146 return kTRUE;
147}
148
149// --------------------------------------------------------------------------
150//
151// Use this function to setup your own conversion factor between degrees
152// and millimeters. The conversion factor should be the one calculated in
153// MGeomCam. Use this function with Caution: You could create wrong values
154// by setting up your own scale factor.
155//
156void MHHillasSrc::SetMm2Deg(Float_t mmdeg)
157{
158 if (mmdeg<=0)
159 {
160 *fLog << warn << dbginf << "Warning - Conversion factor <= 0 - nonsense. Ignored." << endl;
161 return;
162 }
163
164 if (fMm2Deg>0)
165 *fLog << warn << dbginf << "Warning - Conversion factor already set. Overwriting" << endl;
166
167 fMm2Deg = mmdeg;
168}
169
170// --------------------------------------------------------------------------
171//
172// With this function you can convert the histogram ('on the fly') between
173// degrees and millimeters.
174//
175void MHHillasSrc::SetMmScale(Bool_t mmscale)
176{
177 if (fUseMmScale == mmscale)
178 return;
179
180 if (fMm2Deg<0)
181 {
182 *fLog << warn << GetDescriptor() << ": Warning - Sorry, no conversion factor for conversion available." << endl;
183 return;
184 }
185
186 const Double_t scale = mmscale ? 1./fMm2Deg : fMm2Deg;
187 MH::ScaleAxis(fDist, scale);
188
189 fDist->SetXTitle(mmscale ? "Dist [mm]" : "Dist [\\circ]");
190
191 fUseMmScale = mmscale;
192}
193
194// --------------------------------------------------------------------------
195//
196// Creates a new canvas and draws the two histograms into it.
197// Be careful: The histograms belongs to this object and won't get deleted
198// together with the canvas.
199//
200void MHHillasSrc::Draw(Option_t *)
201{
202 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
203 pad->SetBorderMode(0);
204
205 AppendPad("");
206
207 // FIXME: Display Source position
208
209 pad->Divide(2, 2);
210
211 pad->cd(1);
212 gPad->SetBorderMode(0);
213 fAlpha->Draw();
214
215 pad->cd(2);
216 gPad->SetBorderMode(0);
217 fDist->Draw();
218
219 delete pad->GetPad(3);
220
221 pad->cd(4);
222 gPad->SetBorderMode(0);
223 //gPad->SetLogy();
224 fCosDA->Draw();
225
226 pad->Modified();
227 pad->Update();
228}
229
230void MHHillasSrc::Paint(Option_t *opt)
231{
232 if (fCosDA->GetEntries()==0)
233 return;
234
235 TVirtualPad *savepad = gPad;
236 savepad->cd(4);
237 gPad->SetLogy();
238 gPad = savepad;
239}
240
241TH1 *MHHillasSrc::GetHistByName(const TString name)
242{
243 if (name.Contains("Alpha", TString::kIgnoreCase))
244 return fAlpha;
245 if (name.Contains("Dist", TString::kIgnoreCase))
246 return fDist;
247 if (name.Contains("CosDA", TString::kIgnoreCase))
248 return fCosDA;
249
250 return NULL;
251}
Note: See TracBrowser for help on using the repository browser.