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

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