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

Last change on this file since 2219 was 2173, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 9.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 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->SetXTitle("E [GeV]");
93 fHistAll->SetYTitle("r [m]");
94 fHistAll->SetZTitle("N");
95
96 fHistSel->SetXTitle("E [GeV]");
97 fHistSel->SetYTitle("r [m]");
98 fHistSel->SetYTitle("N");
99
100 fHistCol->SetXTitle("E [GeV]");
101 fHistCol->SetYTitle("A [m^{2}]");
102}
103
104// --------------------------------------------------------------------------
105//
106// Delete the three histograms
107//
108MHMcCollectionArea::~MHMcCollectionArea()
109{
110 delete fHistAll;
111 delete fHistSel;
112 delete fHistCol;
113}
114
115// --------------------------------------------------------------------------
116//
117// Fill data into the histogram which contains all showers
118//
119void MHMcCollectionArea::FillAll(Double_t energy, Double_t radius)
120{
121 fHistAll->Fill(energy, radius);
122}
123
124// --------------------------------------------------------------------------
125//
126// Fill data into the histogram which contains the selected showers
127//
128void MHMcCollectionArea::FillSel(Double_t energy, Double_t radius)
129{
130 fHistSel->Fill(energy, radius);
131}
132
133// --------------------------------------------------------------------------
134//
135// Draw the histogram with all showers
136//
137void MHMcCollectionArea::DrawAll(Option_t* option)
138{
139 if (!gPad)
140 MH::MakeDefCanvas(fHistAll);
141
142 fHistAll->Draw(option);
143
144 gPad->SetLogx();
145
146 gPad->Modified();
147 gPad->Update();
148}
149
150// --------------------------------------------------------------------------
151//
152// Draw the histogram with the selected showers only.
153//
154void MHMcCollectionArea::DrawSel(Option_t* option)
155{
156 if (!gPad)
157 MH::MakeDefCanvas(fHistSel);
158
159 fHistSel->Draw(option);
160
161 gPad->SetLogx();
162
163 gPad->Modified();
164 gPad->Update();
165}
166
167void MHMcCollectionArea::Draw(Option_t* option)
168{
169 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
170 pad->SetBorderMode(0);
171
172 fHistCol->Draw(option);
173
174 pad->SetLogx();
175
176 pad->Modified();
177 pad->Update();
178}
179
180// --------------------------------------------------------------------------
181//
182// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
183// flag
184//
185void MHMcCollectionArea::CalcEfficiency()
186{
187 Calc(*fHistSel, *fHistAll);
188}
189
190// --------------------------------------------------------------------------
191//
192// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
193// flag
194//
195void MHMcCollectionArea::Calc(const MHMcEnergyImpact &mcsel, const MHMcEnergyImpact &mcall)
196{
197 Calc((TH2D&)*mcsel.GetHist(), (TH2D&)*mcall.GetHist());
198}
199
200// --------------------------------------------------------------------------
201//
202// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
203// flag
204//
205void MHMcCollectionArea::Calc(TH2D &hsel, TH2D &hall)
206{
207 MH::SetBinning(fHistCol, hsel.GetXaxis());
208
209 fHistCol->Sumw2();
210
211 TH1D &psel = *hsel.ProjectionX();
212 TH1D &pall = *hall.ProjectionX();
213
214 TH1D &pally = *hall.ProjectionY();
215
216 Double_t max = 0.;
217 for (Int_t i = hall.GetNbinsY(); i > 0; i--)
218 if (pally.GetBinContent(i) > 0)
219 {
220 max = pally.GetBinLowEdge(i+1);
221 break;
222 }
223
224 fHistCol->Divide(&psel, &pall, TMath::Pi()*max*max, 1);
225 fHistCol->SetEntries(hsel.GetEntries());
226
227 delete &psel;
228 delete &pall;
229 delete &pally;
230
231 SetReadyToSave();
232}
233
234// --------------------------------------------------------------------------
235//
236// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
237// flag
238//
239void MHMcCollectionArea::CalcEfficiency(UInt_t numevts, Float_t emin, Float_t emax, Float_t index)
240{
241 // Here we estimate the total number of showers in each energy bin
242 // known the energy range and spectral index of the generated showers
243 // (set now by hand since the info is not available in the header!)
244
245 TH1D &histsel = *fHistSel->ProjectionX();
246
247 TH1D histall;
248
249 TAxis &xaxis = *histsel.GetXaxis();
250
251 MH::SetBinning(&histall, &xaxis);
252 MH::SetBinning(fHistCol, &xaxis);
253
254 const Float_t expo = 1+index;
255
256 const Float_t k = (Float_t)numevts / (pow(emax,expo) - pow(emin,expo));
257
258 const Int_t nbinx = xaxis.GetNbins();
259
260 for (Int_t i=1; i<=nbinx; i++)
261 {
262 const Float_t e1 = histall.GetBinLowEdge(i);
263 const Float_t e2 = histall.GetBinLowEdge(i+1);
264
265 if (e1 < emin || e2 > emax)
266 continue;
267
268 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
269
270 histall.SetBinContent(i, events);
271 histall.SetBinError(i, sqrt(events));
272
273 }
274
275 // -----------------------------------------------------------
276 //
277 // Impact parameter range: TO BE FIXED! Impact parameter shoud be
278 // read from run header, but it is not yet in!!
279 //
280 const Float_t r1 = 0;
281 const Float_t r2 = 300;
282
283 *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;
284 const Float_t area = TMath::Pi() * (r2*r2 - r1*r1);
285
286 for (Int_t ix=1; ix<=nbinx; ix++)
287 {
288 const Float_t Na = histall.GetBinContent(ix);
289
290 if (Na <= 0)
291 continue;
292
293 const Float_t Ns = histsel.GetBinContent(ix);
294
295 // Since Na is an estimate of the total number of showers generated
296 // in the energy bin, it may happen that Ns (triggered showers) is
297 // larger than Na. In that case, the bin is skipped:
298
299 if (Na < Ns)
300 continue;
301
302 const Double_t eff = Ns/Na;
303
304 const Double_t efferr = sqrt((1.-eff)*Ns)/Na;
305
306 fHistCol->SetBinContent(ix, eff*area);
307 fHistCol->SetBinError(ix, efferr*area);
308 }
309
310 delete &histsel;
311
312 SetReadyToSave();
313}
314
315// --------------------------------------------------------------------------
316//
317// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
318// flag
319//
320void MHMcCollectionArea::Calc(const MHMcEfficiency &heff)
321{
322 //
323 // now calculate the Collection area for different
324 // energies
325 //
326 TH2D &h = (TH2D&)*heff.GetHist();
327
328 MH::SetBinning(fHistCol, h.GetXaxis());
329
330 const Int_t nbinx = h.GetXaxis()->GetNbins();
331 const Int_t nbiny = h.GetYaxis()->GetNbins();
332
333 for (Int_t ix=1; ix<=nbinx; ix++)
334 {
335 Double_t errA = 0;
336 Double_t colA = 0;
337
338 for (Int_t iy=1; iy<=nbiny; iy++)
339 {
340 TAxis *yaxis = h.GetYaxis();
341
342 const Double_t r1 = yaxis->GetBinLowEdge(iy);
343 const Double_t r2 = yaxis->GetBinLowEdge(iy+1);
344
345 const Double_t A = TMath::Pi() * (r2*r2 - r1*r1);
346
347 const Double_t eff = h.GetCellContent(ix, iy);
348 const Double_t efferr = h.GetCellError(ix, iy);
349
350 colA += eff*A;
351 errA += A*A * efferr*efferr;
352 }
353
354 fHistCol->SetBinContent(ix, colA);
355 fHistCol->SetBinError(ix, sqrt(errA));
356 }
357
358 SetReadyToSave();
359}
Note: See TracBrowser for help on using the repository browser.