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

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