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

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