source: trunk/MagicSoft/Mars/mhistmc/MHMcCollectionArea.cc@ 3378

Last change on this file since 3378 was 2606, checked in by moralejo, 21 years ago
*** empty log message ***
File size: 9.1 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 12/2000 <mailto:tbretz@astro.uni-wuerzburg.de>
19! Author(s): Harald Kornmayer 1/2001
20!
21! Copyright: MAGIC Software Development, 2000-2002
22!
23!
24\* ======================================================================== */
25
26//////////////////////////////////////////////////////////////////////////////
27// //
28// MHMcCollectionArea //
29// //
30//////////////////////////////////////////////////////////////////////////////
31
32#include "MHMcCollectionArea.h"
33
34#include <TH2.h>
35#include <TCanvas.h>
36
37#include "MH.h"
38#include "MBinning.h"
39#include "MHMcEfficiency.h"
40#include "MHMcEnergyImpact.h"
41
42#include "MLog.h"
43#include "MLogManip.h"
44
45ClassImp(MHMcCollectionArea);
46
47using namespace std;
48
49// --------------------------------------------------------------------------
50//
51// Creates the three necessary histograms:
52// - selected showers (input)
53// - all showers (input)
54// - collection area (result)
55//
56MHMcCollectionArea::MHMcCollectionArea(const char *name, const char *title)
57{
58 // initialize the histogram for the distribution r vs E
59 //
60 // we set the energy range from 2 Gev to 20000 GeV (in log 4 orders
61 // of magnitude) and for each order we take 25 subdivision --> 100 xbins
62 //
63 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
64 //
65 fName = name ? name : "MHMcCollectionArea";
66 fTitle = title ? title : "Collection Area vs. Energy";
67
68 MBinning binsx;
69 MBinning binsy;
70 binsx.SetEdgesLog(100, 2., 20000);
71 binsy.SetEdges (50, 0, 500);
72
73 fHistAll = new TH2D;
74 fHistSel = new TH2D;
75 fHistCol = new TH1D;
76
77 MH::SetBinning(fHistAll, &binsx, &binsy);
78 MH::SetBinning(fHistSel, &binsx, &binsy);
79
80 fHistCol->SetName(fName);
81 fHistAll->SetName("AllEvents");
82 fHistSel->SetName("SelectedEvents");
83
84 fHistCol->SetTitle(fTitle);
85 fHistAll->SetTitle("All showers - Radius vs Energy distribution");
86 fHistSel->SetTitle("Selected showers - Radius vs Energy distribution");
87
88 fHistAll->SetDirectory(NULL);
89 fHistSel->SetDirectory(NULL);
90 fHistCol->SetDirectory(NULL);
91
92 fHistAll->UseCurrentStyle();
93 fHistSel->UseCurrentStyle();
94 fHistCol->UseCurrentStyle();
95
96 fHistAll->SetXTitle("E [GeV]");
97 fHistAll->SetYTitle("r [m]");
98 fHistAll->SetZTitle("N");
99
100 fHistSel->SetXTitle("E [GeV]");
101 fHistSel->SetYTitle("r [m]");
102 fHistSel->SetYTitle("N");
103
104 fHistCol->SetXTitle("E [GeV]");
105 fHistCol->SetYTitle("A [m^{2}]");
106}
107
108// --------------------------------------------------------------------------
109//
110// Delete the three histograms
111//
112MHMcCollectionArea::~MHMcCollectionArea()
113{
114 delete fHistAll;
115 delete fHistSel;
116 delete fHistCol;
117}
118
119// --------------------------------------------------------------------------
120//
121// Fill data into the histogram which contains all showers
122//
123void MHMcCollectionArea::FillAll(Double_t energy, Double_t radius)
124{
125 fHistAll->Fill(energy, radius);
126}
127
128// --------------------------------------------------------------------------
129//
130// Fill data into the histogram which contains the selected showers
131//
132void MHMcCollectionArea::FillSel(Double_t energy, Double_t radius)
133{
134 fHistSel->Fill(energy, radius);
135}
136
137// --------------------------------------------------------------------------
138//
139// Draw the histogram with all showers
140//
141void MHMcCollectionArea::DrawAll(Option_t* option)
142{
143 if (!gPad)
144 MH::MakeDefCanvas(fHistAll);
145
146 fHistAll->Draw(option);
147
148 gPad->SetLogx();
149
150 gPad->Modified();
151 gPad->Update();
152}
153
154// --------------------------------------------------------------------------
155//
156// Draw the histogram with the selected showers only.
157//
158void MHMcCollectionArea::DrawSel(Option_t* option)
159{
160 if (!gPad)
161 MH::MakeDefCanvas(fHistSel);
162
163 fHistSel->Draw(option);
164
165 gPad->SetLogx();
166
167 gPad->Modified();
168 gPad->Update();
169}
170
171void MHMcCollectionArea::Draw(Option_t* option)
172{
173 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
174 pad->SetBorderMode(0);
175
176 fHistCol->Draw(option);
177
178 pad->SetLogx();
179
180 pad->Modified();
181 pad->Update();
182}
183
184// --------------------------------------------------------------------------
185//
186// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
187// flag
188//
189void MHMcCollectionArea::CalcEfficiency()
190{
191 Calc(*fHistSel, *fHistAll);
192}
193
194// --------------------------------------------------------------------------
195//
196// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
197// flag
198//
199void MHMcCollectionArea::Calc(const MHMcEnergyImpact &mcsel, const MHMcEnergyImpact &mcall)
200{
201 Calc((TH2D&)*mcsel.GetHist(), (TH2D&)*mcall.GetHist());
202}
203
204// --------------------------------------------------------------------------
205//
206// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
207// flag
208//
209void MHMcCollectionArea::Calc(TH2D &hsel, TH2D &hall)
210{
211 MH::SetBinning(fHistCol, hsel.GetXaxis());
212
213 fHistCol->Sumw2();
214
215 TH1D &psel = *hsel.ProjectionX();
216 TH1D &pall = *hall.ProjectionX();
217
218 TH1D &pally = *hall.ProjectionY();
219
220 Double_t max = 0.;
221 for (Int_t i = hall.GetNbinsY(); i > 0; i--)
222 if (pally.GetBinContent(i) > 0)
223 {
224 max = pally.GetBinLowEdge(i+1);
225 break;
226 }
227
228 fHistCol->Divide(&psel, &pall, TMath::Pi()*max*max, 1);
229 fHistCol->SetEntries(hsel.GetEntries());
230
231 delete &psel;
232 delete &pall;
233 delete &pally;
234
235 SetReadyToSave();
236}
237
238// --------------------------------------------------------------------------
239//
240// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
241// flag
242//
243void MHMcCollectionArea::CalcEfficiency2()
244{
245 TH1D &histsel = *fHistSel->ProjectionX();
246 TH1D &histall = *fHistAll->ProjectionX();
247
248 TAxis &xaxis = *histsel.GetXaxis();
249 MH::SetBinning(fHistCol, &xaxis);
250 const Int_t nbinx = xaxis.GetNbins();
251
252 // -----------------------------------------------------------
253 //
254 // Impact parameter range: TO BE FIXED! Impact parameter shoud be
255 // read from run header, but it is not yet in!!
256 //
257 const Float_t r1 = 0;
258 const Float_t r2 = 300;
259
260 *fLog << warn << endl << dbginf << "WARNING! I will assume a maximum impact parameter of 300 meters for the MC events. Check that this is the true one!" <<endl<<endl;
261 const Float_t total_area = TMath::Pi() * (r2*r2 - r1*r1);
262
263 for (Int_t ix=1; ix<=nbinx; ix++)
264 {
265 const Float_t Na = histall.GetBinContent(ix);
266
267 if (Na <= 0)
268 continue;
269
270 const Float_t Ns = histsel.GetBinContent(ix);
271
272 // Since Na is an estimate of the total number of showers generated
273 // in the energy bin, it may happen that Ns (triggered showers) is
274 // larger than Na. In that case, the bin is skipped:
275
276 if (Na < Ns)
277 continue;
278
279 const Double_t eff = Ns/Na;
280
281 const Double_t efferr = sqrt((1.-eff)*Ns)/Na;
282
283 const Float_t col_area = eff * total_area;
284 const Float_t col_area_error = efferr * total_area;
285
286 fHistCol->SetBinContent(ix, col_area);
287 fHistCol->SetBinError(ix, col_area_error);
288 }
289
290 delete &histsel;
291
292 SetReadyToSave();
293}
294
295// --------------------------------------------------------------------------
296//
297// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
298// flag
299//
300void MHMcCollectionArea::Calc(const MHMcEfficiency &heff)
301{
302 //
303 // now calculate the Collection area for different
304 // energies
305 //
306 TH2D &h = (TH2D&)*heff.GetHist();
307
308 MH::SetBinning(fHistCol, h.GetXaxis());
309
310 const Int_t nbinx = h.GetXaxis()->GetNbins();
311 const Int_t nbiny = h.GetYaxis()->GetNbins();
312
313 for (Int_t ix=1; ix<=nbinx; ix++)
314 {
315 Double_t errA = 0;
316 Double_t colA = 0;
317
318 for (Int_t iy=1; iy<=nbiny; iy++)
319 {
320 TAxis *yaxis = h.GetYaxis();
321
322 const Double_t r1 = yaxis->GetBinLowEdge(iy);
323 const Double_t r2 = yaxis->GetBinLowEdge(iy+1);
324
325 const Double_t A = TMath::Pi() * (r2*r2 - r1*r1);
326
327 const Double_t eff = h.GetCellContent(ix, iy);
328 const Double_t efferr = h.GetCellError(ix, iy);
329
330 colA += eff*A;
331 errA += A*A * efferr*efferr;
332 }
333
334 fHistCol->SetBinContent(ix, colA);
335 fHistCol->SetBinError(ix, sqrt(errA));
336 }
337
338 SetReadyToSave();
339}
Note: See TracBrowser for help on using the repository browser.