source: trunk/MagicSoft/Mars/mhist/MHHillasSrc.cc@ 1442

Last change on this file since 1442 was 1442, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 6.3 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", 90, 0, 90);
71 fDist = new TH1F("Dist", "Dist of Ellipse", 100, 0, 600);
72 fHeadTail = new TH1F("HeadTail", "HeadTail of Ellipse", 45, 0, 90);
73 fCosDA = new TH1F("CosDA", "cos(Delta,Alpha) of Ellipse", 100, -1, 1);
74
75 fAlpha->SetDirectory(NULL);
76 fDist->SetDirectory(NULL);
77 fHeadTail->SetDirectory(NULL);
78 fCosDA->SetDirectory(NULL);
79
80 fAlpha->SetXTitle("\\alpha [\\circ]");
81 fDist->SetXTitle("Dist [mm]");
82 fHeadTail->SetXTitle("Head-Tail [\\circ]");
83 fCosDA->SetXTitle("cos(\\delta,\\alpha) [mm]");
84
85 fAlpha->SetYTitle("Counts");
86 fDist->SetYTitle("Counts");
87 fHeadTail->SetYTitle("Counts");
88 fCosDA->SetYTitle("Counts");
89}
90
91// --------------------------------------------------------------------------
92//
93// Delete the four histograms
94//
95MHHillasSrc::~MHHillasSrc()
96{
97 delete fAlpha;
98 delete fDist;
99 delete fHeadTail;
100 delete fCosDA;
101}
102
103// --------------------------------------------------------------------------
104//
105// Setup the Binning for the histograms automatically if the correct
106// instances of MBinning (with the names 'BinningAlpha' and 'BinningDist')
107// are found in the parameter list
108// Use this function if you want to set the conversion factor which
109// is used to convert the mm-scale in the camera plain into the deg-scale
110// used for histogram presentations. The conversion factor is part of
111// the camera geometry. Please create a corresponding MGeomCam container.
112//
113Bool_t MHHillasSrc::SetupFill(const MParList *plist)
114{
115 const MBinning* binsa = (MBinning*)plist->FindObject("BinningAlpha");
116 const MBinning* binsd = (MBinning*)plist->FindObject("BinningDist");
117 if (!binsa || !binsd)
118 {
119 *fLog << err << dbginf << "At least one MBinning not found... aborting." << endl;
120 return kFALSE;
121 }
122
123 SetBinning(fAlpha, binsa);
124 SetBinning(fDist, binsd);
125
126 const MGeomCam *geom = (MGeomCam*)plist->FindObject("MGeomCam");
127 if (!geom)
128 {
129 *fLog << warn << dbginf << "No Camera Geometry available. Using mm-scale for histograms." << endl;
130 return kTRUE;
131 }
132
133 fDist->SetXTitle("Dist [\\circ]");
134
135 fMm2Deg = geom->GetConvMm2Deg();
136 fUseMmScale = kFALSE;
137
138 return kTRUE;
139}
140
141// --------------------------------------------------------------------------
142//
143// Fill the four histograms with data from a MHillas-Container.
144// Be careful: Only call this with an object of type MHillas
145//
146Bool_t MHHillasSrc::Fill(const MParContainer *par)
147{
148 const MHillasSrc &h = *(MHillasSrc*)par;
149
150 fAlpha ->Fill(fabs(h.GetAlpha()));
151 fDist ->Fill(fUseMmScale ? h.GetDist() : fMm2Deg*h.GetDist());
152 fHeadTail->Fill(h.GetHeadTail());
153 fCosDA ->Fill(h.GetCosDeltaAlpha());
154
155 return kTRUE;
156}
157
158// --------------------------------------------------------------------------
159//
160// Draw clones of all two histograms. So that the object can be deleted
161// and the histograms are still visible in the canvas.
162// The cloned object are deleted together with the canvas if the canvas is
163// destroyed. If you want to handle dostroying the canvas you can get a
164// pointer to it from this function
165//
166TObject *MHHillasSrc::DrawClone(Option_t *opt) const
167{
168 TCanvas *c = MakeDefCanvas("HillasSrc", "Histograms of Source dependant Parameters",
169 700, 500);
170 c->Divide(2, 2);
171
172 // FIXME: Display Source position
173
174 gROOT->SetSelectedPad(NULL);
175
176 //
177 // This is necessary to get the expected bahviour of DrawClone
178 //
179 c->cd(1);
180 fAlpha->DrawCopy();
181
182 c->cd(2);
183 fDist->DrawCopy();
184
185 c->cd(3);
186 fHeadTail->DrawCopy();
187
188 c->cd(4);
189 fCosDA->DrawCopy();
190
191 c->Modified();
192 c->Update();
193
194 return c;
195}
196
197// --------------------------------------------------------------------------
198//
199// Creates a new canvas and draws the two histograms into it.
200// Be careful: The histograms belongs to this object and won't get deleted
201// together with the canvas.
202//
203void MHHillasSrc::Draw(Option_t *)
204{
205 if (!gPad)
206 MakeDefCanvas("HillasSrc", "Histograms of Src dependant Parameters",
207 700, 500);
208
209 // FIXME: Display Source position
210
211 gPad->Divide(2, 2);
212
213 gPad->cd(1);
214 fAlpha->Draw();
215
216 gPad->cd(2);
217 fDist->Draw();
218
219 gPad->cd(1);
220 fHeadTail->Draw();
221
222 gPad->cd(2);
223 fCosDA->Draw();
224
225 gPad->Modified();
226 gPad->Update();
227}
Note: See TracBrowser for help on using the repository browser.