source: trunk/MagicSoft/Mars/mhist/MHMcCollectionArea.cc@ 1502

Last change on this file since 1502 was 1330, checked in by tbretz, 23 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 9.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 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 "MHMcEfficiency.h"
39#include "MHMcEnergyImpact.h"
40
41ClassImp(MHMcCollectionArea);
42
43// --------------------------------------------------------------------------
44//
45// Creates the three necessary histograms:
46// - selected showers (input)
47// - all showers (input)
48// - collection area (result)
49//
50MHMcCollectionArea::MHMcCollectionArea(const char *name, const char *title)
51{
52 // initialize the histogram for the distribution r vs E
53 //
54 // we set the energy range from 1 Gev to 10000 GeV (in log 5 orders
55 // of magnitude) and for each order we take 10 subdivision --> 50 xbins
56 //
57 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
58
59
60 fName = name ? name : "MHMcCollectionArea";
61 fTitle = title ? title : "Collection Area vs. Energy";
62
63 const Int_t nbins = 50;
64 const Float_t maxx = 5;
65 const Float_t maxy = 500;
66
67 Float_t *binsx = new Float_t[nbins+1];
68 for (int i=0; i<nbins+1; i++)
69 binsx[i] = pow(10, maxx*i/nbins);
70
71 Float_t *binsy = new Float_t[nbins+1];
72 for (int i=0; i<nbins+1; i++)
73 binsy[i] = maxy*i/nbins;
74
75 fHistAll = new TH2D("AllEvents", "All showers - Radius vs Energy distribution",
76 nbins, binsx, nbins, binsy);
77 fHistSel = new TH2D("SelectedEvents", "Selected showers - Radius vs Energy distribution",
78 nbins, binsx, nbins, binsy);
79
80 delete binsx;
81 delete binsy;
82
83 fHistCol = new TH1D;
84 fHistCol->SetName(fName);
85 fHistCol->SetTitle(fTitle);
86
87 fHistAll->SetDirectory(NULL);
88 fHistSel->SetDirectory(NULL);
89 fHistCol->SetDirectory(NULL);
90
91 fHistAll->SetXTitle("E [GeV]");
92 fHistAll->SetYTitle("r [m]");
93 fHistAll->SetZTitle("N");
94
95 fHistSel->SetXTitle("E [GeV]");
96 fHistSel->SetYTitle("r [m]");
97 fHistSel->SetYTitle("N");
98
99 fHistCol->SetXTitle("E [GeV]");
100 fHistCol->SetYTitle("A [m^{2}]");
101}
102
103// --------------------------------------------------------------------------
104//
105// Delete the three histograms
106//
107MHMcCollectionArea::~MHMcCollectionArea()
108{
109 delete fHistAll;
110 delete fHistSel;
111 delete fHistCol;
112}
113
114// --------------------------------------------------------------------------
115//
116// Fill data into the histogram which contains all showers
117//
118void MHMcCollectionArea::FillAll(Float_t energy, Float_t radius)
119{
120 fHistAll->Fill(energy, radius);
121}
122
123// --------------------------------------------------------------------------
124//
125// Fill data into the histogram which contains the selected showers
126//
127void MHMcCollectionArea::FillSel(Float_t energy, Float_t radius)
128{
129 fHistSel->Fill(energy, radius);
130}
131
132// --------------------------------------------------------------------------
133//
134// Draw the histogram with all showers
135//
136void MHMcCollectionArea::DrawAll(Option_t* option)
137{
138 if (!gPad)
139 MH::MakeDefCanvas(fHistAll);
140
141 fHistAll->Draw(option);
142
143 gPad->SetLogx();
144
145 gPad->Modified();
146 gPad->Update();
147}
148
149// --------------------------------------------------------------------------
150//
151// Draw the histogram with the selected showers only.
152//
153void MHMcCollectionArea::DrawSel(Option_t* option)
154{
155 if (!gPad)
156 MH::MakeDefCanvas(fHistSel);
157
158 fHistSel->Draw(option);
159
160 gPad->SetLogx();
161
162 gPad->Modified();
163 gPad->Update();
164}
165
166// --------------------------------------------------------------------------
167//
168// Creates a new canvas and draws the histogram into it.
169// Be careful: The histogram belongs to this object and won't get deleted
170// together with the canvas.
171//
172TObject *MHMcCollectionArea::DrawClone(Option_t* option) const
173{
174 TCanvas *c = MH::MakeDefCanvas(fHistCol);
175
176 //
177 // This is necessary to get the expected bahviour of DrawClone
178 //
179 gROOT->SetSelectedPad(NULL);
180
181 fHistCol->DrawCopy(option);
182
183 gPad->SetLogx();
184
185 c->Modified();
186 c->Update();
187
188 return c;
189}
190
191void MHMcCollectionArea::Draw(Option_t* option)
192{
193 if (!gPad)
194 MH::MakeDefCanvas(fHistCol);
195
196 fHistCol->Draw(option);
197
198 gPad->SetLogx();
199
200 gPad->Modified();
201 gPad->Update();
202}
203
204// --------------------------------------------------------------------------
205//
206// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
207// flag
208//
209void MHMcCollectionArea::CalcEfficiency()
210{
211 MHMcEfficiency heff;
212 heff.Calc(*fHistSel, *fHistAll);
213
214 Calc(heff);
215}
216
217void MHMcCollectionArea::CalcEfficiency(UInt_t numevts, Float_t angle)
218{
219 // Here we estimate the total number of showers in each energy bin
220 // known the energy range and spectral index of the generated showers
221 // (set now by hand since the info is not available in the header!)
222
223 TH1D &histsel = *fHistSel->ProjectionX();
224
225 TH1D histall;
226
227 TAxis &xaxis = *histsel.GetXaxis();
228
229 MH::SetBinning(&histall, &xaxis);
230 MH::SetBinning(fHistCol, &xaxis);
231
232 const Float_t emin = 10.;
233 const Float_t emax = 30000.; // Energies in GeV.
234 const Float_t index = 2.6; // Differential spectral Index.
235
236 const Float_t expo = 1.-index;
237
238 const Float_t k = (Float_t)numevts / (pow(emax,expo) - pow(emin,expo));
239
240 const Int_t nbinx = xaxis.GetNbins();
241
242 for (Int_t i=1; i<=nbinx; i++)
243 {
244 const Float_t e1 = histall.GetBinLowEdge(i);
245 const Float_t e2 = histall.GetBinLowEdge(i+1);
246
247 if (e1 < emin || e2 > emax)
248 continue;
249
250 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
251
252 histall.SetBinContent(i, events);
253 histall.SetBinError(i, sqrt(events));
254
255 }
256
257 // -----------------------------------------------------------
258
259 // Impact parameter range:
260 const Float_t r1 = 0;
261 const Float_t r2 = 400;
262
263 const Float_t dr = TMath::Pi() * (r2*r2 - r1*r1);
264
265 angle *= TMath::Pi()/180;
266
267 for (Int_t ix=1; ix<=nbinx; ix++)
268 {
269 const Float_t Na = histall.GetBinContent(ix);
270
271 if (Na <= 0)
272 continue;
273
274 const Float_t Ns = histsel.GetBinContent(ix);
275
276 // Since Na is an estimate of the total number of showers generated
277 // in the energy bin, it may happen that Ns (triggered showers) is
278 // larger than Na. In that case, the bin is skipped:
279
280 if (Na < Ns)
281 continue;
282
283 const Double_t eff = Ns/Na;
284
285 const Double_t err = sqrt((1.-eff)*Ns)/Na;
286
287 const Float_t area = dr * cos(angle);
288
289 fHistCol->SetBinContent(ix, eff*area);
290 fHistCol->SetBinError(ix, err*area);
291 }
292
293 delete &histsel;
294
295 SetReadyToSave();
296}
297
298void MHMcCollectionArea::Calc(const MHMcEnergyImpact &mcsel, const MHMcEnergyImpact &mcall)
299{
300 MHMcEfficiency heff;
301 heff.Calc(*mcsel.GetHist(), *mcall.GetHist());
302
303 Calc(heff);
304}
305
306void MHMcCollectionArea::Calc(const MHMcEfficiency &heff)
307{
308 //
309 // now calculate the Collection area for different
310 // energies
311 //
312 TH2D &h = (TH2D&)*heff.GetHist();
313
314 MH::SetBinning(fHistCol, h.GetXaxis());
315
316 const Int_t nbinx = h.GetXaxis()->GetNbins();
317 const Int_t nbiny = h.GetYaxis()->GetNbins();
318
319 for (Int_t ix=1; ix<=nbinx; ix++)
320 {
321 Double_t errA = 0;
322 Double_t colA = 0;
323
324 for (Int_t iy=1; iy<=nbiny; iy++)
325 {
326 TAxis *yaxis = h.GetYaxis();
327
328 const Double_t r1 = yaxis->GetBinLowEdge(iy);
329 const Double_t r2 = yaxis->GetBinLowEdge(iy+1);
330
331 const Double_t A = TMath::Pi() * (r2*r2 - r1*r1);
332
333 const Double_t eff = h.GetCellContent(ix, iy);
334 const Double_t err = h.GetCellError(ix, iy);
335
336 colA += eff*A;
337 errA += A*A * err*err;
338 }
339
340 fHistCol->SetBinContent(ix, colA);
341 fHistCol->SetBinError(ix, sqrt(errA));
342 }
343
344 SetReadyToSave();
345}
346
347// --------------------------------------------------------------------------
348//
349// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
350// flag
351//
352/*
353void MHMcCollectionArea::Calc(MHMcEnergyImpact &mcsel, UInt_t numevents, Float_t angle)
354{
355}
356*/
Note: See TracBrowser for help on using the repository browser.