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

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