source: trunk/MagicSoft/Mars/mreflector/MHReflector.cc@ 9484

Last change on this file since 9484 was 9320, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 8.6 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, 2/2007 <mailto:tbretz@astro.uni-wuerzburg.de>
19!
20! Copyright: MAGIC Software Development, 2000-2007
21!
22!
23\* ======================================================================== */
24
25/////////////////////////////////////////////////////////////////////////////
26//
27// MHReflector
28//
29/////////////////////////////////////////////////////////////////////////////
30#include "MHReflector.h"
31
32#include <TMath.h>
33#include <TLegend.h>
34#include <TCanvas.h>
35
36#include "MLog.h"
37#include "MLogManip.h"
38
39#include "MParList.h"
40#include "MBinning.h"
41#include "MString.h"
42
43#include "MGeomCam.h"
44
45#include "MRflEvtData.h"
46#include "MRflEvtHeader.h"
47
48
49ClassImp(MHReflector);
50
51using namespace std;
52
53// --------------------------------------------------------------------------
54//
55// Default constructor. Creates and initializes the histograms
56//
57MHReflector::MHReflector(const char *name, const char *title)
58{
59 fName = name ? name : "MHReflector";
60 fTitle = title ? title : "Histogram for the reflected photons";
61
62 fHistXY.SetName("ReflXY");
63 fHistXY.SetTitle("Histogram vs X/Y and Energy");
64 fHistXY.SetXTitle("X [\\circ]");
65 fHistXY.SetYTitle("Y [\\circ]");
66 fHistXY.SetZTitle("E [GeV]");
67 fHistXY.SetDirectory(NULL);
68 fHistXY.Sumw2();
69
70 fHistRad.SetName("ReflRad");
71 fHistRad.SetTitle("Histogram vs Radius and Energy");
72 fHistRad.SetXTitle("E [GeV]");
73 fHistRad.SetYTitle("R [\\circ]");
74 fHistRad.SetZTitle("Cnts/deg^{2}");
75 fHistRad.SetDirectory(NULL);
76 fHistRad.Sumw2();
77
78 fHistSize.SetName("ReflSize");
79 fHistSize.SetTitle("Histogram vs Size and Energy");
80 fHistSize.SetXTitle("E [GeV]");
81 fHistSize.SetYTitle("N\\gamma");
82 fHistSize.SetZTitle("Cnts");
83 fHistSize.SetDirectory(NULL);
84 fHistSize.Sumw2();
85
86 MBinning binse, binsd, binsr;
87 binse.SetEdgesLog(2*5, 10, 1000000);
88 binsd.SetEdges(50, -2.5, 2.5);
89 binsr.SetEdges(15, 0, 2.5);
90
91 SetBinning(&fHistXY, &binsd, &binsd, &binse);
92 SetBinning(&fHistRad, &binse, &binsr);
93 SetBinning(&fHistSize, &binse, &binse);
94}
95
96// --------------------------------------------------------------------------
97//
98// Search for MRflEvtData, MRflEvtHeader and MGeomCam
99//
100Bool_t MHReflector::SetupFill(const MParList *pList)
101{
102 fHeader = (MRflEvtHeader*)pList->FindObject("MRflEvtHeader");
103 if (!fHeader)
104 {
105 *fLog << err << "MRflEvtHeader not found... abort." << endl;
106 return kFALSE;
107 }
108 fData = (MRflEvtData*)pList->FindObject("MRflEvtData");
109 if (!fData)
110 {
111 *fLog << err << "MRflEvtData not found... abort." << endl;
112 return kFALSE;
113 }
114
115 MGeomCam *geom = (MGeomCam*)pList->FindObject("MGeomCam");
116 if (!geom)
117 {
118 *fLog << err << "MGeomCam not found... abort." << endl;
119 return kFALSE;
120 }
121
122 fMm2Deg = geom->GetConvMm2Deg();
123
124 return kTRUE;
125}
126
127// --------------------------------------------------------------------------
128//
129// Fill position and radius versus energy
130//
131#include "MRflSinglePhoton.h"
132Int_t MHReflector::Fill(const MParContainer *par, const Stat_t weight)
133{
134 const Double_t energy = fHeader->GetEnergy();
135
136 const Int_t n = fData->GetNumPhotons();
137
138 fHistSize.Fill(energy, n);
139
140 for (int i=0; i<n; i++)
141 {
142 const MRflSinglePhoton &ph = fData->GetPhoton(i);
143
144 const Double_t x = ph.GetX()*fMm2Deg;
145 const Double_t y = ph.GetY()*fMm2Deg;
146 const Double_t r = TMath::Hypot(x, y);
147 const Double_t U = r*TMath::TwoPi();
148
149 if (U==0)
150 continue;
151
152 fHistRad.Fill(energy, r, weight/U);
153 //fHistXY.Fill(x, y, energy);
154 }
155
156 //fData->FillRad(fHistRad, energy/*x*/, fMm2Deg);
157 //fData->Fill(fHistXY, energy/*z*/, fMm2Deg);
158
159 return kTRUE;
160}
161
162// --------------------------------------------------------------------------
163//
164// Scale radial histogram with area
165//
166Bool_t MHReflector::Finalize()
167{
168 /*
169 const TAxis &axey = *fHistRad.GetYaxis();
170 for (int y=1; y<=fHistRad.GetNbinsY(); y++)
171 {
172 const Double_t lo = axey.GetBinLowEdge(y);
173 const Double_t hi = axey.GetBinLowEdge(y+1);
174
175 const Double_t A = (hi*hi-lo*lo)*TMath::Pi()*TMath::Pi();
176
177 for (int x=1; x<=fHistRad.GetNbinsX(); x++)
178 fHistRad.SetBinContent(x, y, fHistRad.GetBinContent(x, y)/A);
179 }
180 */
181 return kTRUE;
182}
183
184#include "../mmc/MMcEvt.hxx"
185void MHReflector::Draw(Option_t *o)
186{
187 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
188 pad->SetBorderMode(0);
189
190 AppendPad();
191
192 pad->Divide(3,2);
193
194 pad->cd(1);
195 gPad->SetBorderMode(0);
196 gPad->SetLogx();
197 fHistRad.Draw("colz");
198
199 pad->cd(2);
200 gPad->SetBorderMode(0);
201 TH1 *h = fHistRad.ProjectionY("ProfRad", -1, -1, "e");
202 h->SetTitle("RadialProfile");
203 h->SetDirectory(NULL);
204 h->SetBit(kCanDelete);
205 h->SetLineWidth(2);
206 h->SetLineColor(kBlue);
207 h->SetStats(kFALSE);
208 h->Draw("");
209
210 pad->cd(3);
211 gPad->SetBorderMode(0);
212 gPad->SetLogx();
213 gPad->SetLogy();
214 h = fHistSize.ProjectionY("ProfSize", -1, -1, "e");
215 h->SetTitle("Size distribution");
216 h->SetDirectory(NULL);
217 h->SetBit(kCanDelete);
218 h->SetLineColor(kBlue);
219 h->SetLineWidth(2);
220 h->SetStats(kFALSE);
221 h->Draw("");
222
223 pad->cd(5);
224 gPad->SetBorderMode(0);
225
226 TLegend *leg = new TLegend(0.01, 0.01, 0.99, 0.99, "Energy Range");
227 for (int i=1; i<=fHistRad.GetNbinsX(); i++)
228 {
229 h = fHistRad.ProjectionY(MString::Format("ProjRad%d", i), i, i, "e");
230 h->SetDirectory(NULL);
231 h->SetBit(kCanDelete);
232 h->SetLineWidth(2);
233 h->SetStats(kFALSE);
234 h->Draw(i==1?"":"same");
235 if (i==1)
236 h->SetTitle("Radial Profile");
237
238 const Float_t xlo = fHistRad.GetXaxis()->GetBinLowEdge(i);
239 const Float_t xhi = fHistRad.GetXaxis()->GetBinUpEdge(i);
240
241 const TString lo = MMcEvt::GetEnergyStr(xlo);
242 const TString hi = MMcEvt::GetEnergyStr(xhi);
243
244 leg->AddEntry(h, MString::Format("%s - %s", lo.Data(), hi.Data()));
245 }
246
247 pad->cd(4);
248 gPad->SetBorderMode(0);
249 leg->SetBit(kCanDelete);
250 leg->Draw();
251
252 pad->cd(6);
253 gPad->SetBorderMode(0);
254
255 gPad->SetLogx();
256 gPad->SetLogy();
257
258 for (int i=1; i<=fHistSize.GetNbinsX(); i++)
259 {
260 h = fHistSize.ProjectionY(MString::Format("ProjSize%d", i), i, i, "e");
261 h->SetDirectory(NULL);
262 h->SetBit(kCanDelete);
263 h->SetLineWidth(2);
264 h->SetStats(kFALSE);
265 h->Draw(i==1?"":"same");
266 if (i==1)
267 h->SetTitle("Size distribution");
268 }
269}
270
271void MHReflector::Paint(Option_t *opt)
272{
273 TVirtualPad *pad = gPad;
274
275 pad->cd(1);
276 SetPalette("pretty");
277
278 TH1 *h=0;
279
280 pad->cd(2);
281 gPad->SetBorderMode(0);
282 if (gPad->FindObject("ProfRad"))
283 {
284 h = fHistRad.ProjectionY("ProfRad", -1, -1, "e");
285 h->Scale(1./h->Integral());
286 }
287
288 pad->cd(3);
289 gPad->SetBorderMode(0);
290 if (gPad->FindObject("ProfSize"))
291 {
292 h = fHistSize.ProjectionY("ProfSize", -1, -1, "e");
293 h->Scale(1./h->Integral());
294 }
295
296 pad->cd(5);
297
298 for (int i=1; i<=fHistRad.GetNbinsX(); i++)
299 {
300 const TString name = MString::Format("ProjRad%d", i);
301 if (!gPad->FindObject(name))
302 continue;
303
304 h = fHistRad.ProjectionY(name, i, i, "e");
305 h->SetLineColor(i);
306 h->Scale(1./fHistRad.Integral());
307 if (i==1)
308 h->SetMaximum(fHistRad.GetMaximum()/fHistRad.Integral()*1.05);
309 }
310
311 pad->cd(6);
312
313 for (int i=1; i<=fHistRad.GetNbinsX(); i++)
314 {
315 const TString name = MString::Format("ProjSize%d", i);
316 if (!gPad->FindObject(name))
317 continue;
318
319 h = fHistSize.ProjectionY(name, i, i, "e");
320 h->SetLineColor(i);
321 h->Scale(1./fHistSize.Integral());
322 if (i==1)
323 {
324 h->SetMaximum(TMath::Power(fHistSize.GetMaximum()/fHistSize.Integral(), 1.05));
325 h->SetMinimum(TMath::Power(fHistSize.GetMinimum(0)/fHistSize.Integral(), 0.95));
326 }
327 }
328
329}
Note: See TracBrowser for help on using the repository browser.