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

Last change on this file since 1956 was 1668, checked in by tbretz, 22 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 10.0 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* option) const
169{
170 TCanvas *c = MH::MakeDefCanvas(fHistCol);
171
172 //
173 // This is necessary to get the expected bahviour of DrawClone
174 //
175 gROOT->SetSelectedPad(NULL);
176
177 fHistCol->DrawCopy(option);
178
179 gPad->SetLogx();
180
181 c->Modified();
182 c->Update();
183
184 return c;
185}
186
187void MHMcCollectionArea::Draw(Option_t* option)
188{
189 if (!gPad)
190 MH::MakeDefCanvas(fHistCol);
191
192 fHistCol->Draw(option);
193
194 gPad->SetLogx();
195
196 gPad->Modified();
197 gPad->Update();
198}
199
200// --------------------------------------------------------------------------
201//
202// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
203// flag
204//
205void MHMcCollectionArea::CalcEfficiency()
206{
207 Calc(*fHistSel, *fHistAll);
208}
209
210// --------------------------------------------------------------------------
211//
212// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
213// flag
214//
215void MHMcCollectionArea::Calc(const MHMcEnergyImpact &mcsel, const MHMcEnergyImpact &mcall)
216{
217 Calc((TH2D&)*mcsel.GetHist(), (TH2D&)*mcall.GetHist());
218}
219
220// --------------------------------------------------------------------------
221//
222// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
223// flag
224//
225void MHMcCollectionArea::Calc(TH2D &hsel, TH2D &hall)
226{
227 MH::SetBinning(fHistCol, hsel.GetXaxis());
228
229 fHistCol->Sumw2();
230
231 TH1D &psel = *hsel.ProjectionX();
232 TH1D &pall = *hall.ProjectionX();
233
234 const Double_t max = psel.GetYaxis()->GetXmax();
235
236 fHistCol->Divide(&psel, &pall, TMath::Pi()*max*max, 1);
237 fHistCol->SetEntries(hsel.GetEntries());
238
239 delete &psel;
240 delete &pall;
241
242 SetReadyToSave();
243}
244
245// --------------------------------------------------------------------------
246//
247// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
248// flag
249//
250void MHMcCollectionArea::CalcEfficiency(UInt_t numevts, Float_t angle)
251{
252 // Here we estimate the total number of showers in each energy bin
253 // known the energy range and spectral index of the generated showers
254 // (set now by hand since the info is not available in the header!)
255
256 TH1D &histsel = *fHistSel->ProjectionX();
257
258 TH1D histall;
259
260 TAxis &xaxis = *histsel.GetXaxis();
261
262 MH::SetBinning(&histall, &xaxis);
263 MH::SetBinning(fHistCol, &xaxis);
264
265 // TO BE FIXED! In forthcoming camera versions emin, emax and index
266 // will be available in a run header, and should be read from it.
267
268 const Float_t emin = 10.;
269 const Float_t emax = 30000.; // Energies in GeV.
270 const Float_t index = 2.6; // Differential spectral Index.
271
272 const Float_t expo = 1.-index;
273
274 const Float_t k = (Float_t)numevts / (pow(emax,expo) - pow(emin,expo));
275
276 const Int_t nbinx = xaxis.GetNbins();
277
278 for (Int_t i=1; i<=nbinx; i++)
279 {
280 const Float_t e1 = histall.GetBinLowEdge(i);
281 const Float_t e2 = histall.GetBinLowEdge(i+1);
282
283 if (e1 < emin || e2 > emax)
284 continue;
285
286 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
287
288 histall.SetBinContent(i, events);
289 histall.SetBinError(i, sqrt(events));
290
291 }
292
293 // -----------------------------------------------------------
294
295 // Impact parameter range:
296 const Float_t r1 = 0;
297 const Float_t r2 = 400;
298
299 const Float_t dr = TMath::Pi() * (r2*r2 - r1*r1);
300
301 angle *= TMath::Pi()/180;
302
303 for (Int_t ix=1; ix<=nbinx; ix++)
304 {
305 const Float_t Na = histall.GetBinContent(ix);
306
307 if (Na <= 0)
308 continue;
309
310 const Float_t Ns = histsel.GetBinContent(ix);
311
312 // Since Na is an estimate of the total number of showers generated
313 // in the energy bin, it may happen that Ns (triggered showers) is
314 // larger than Na. In that case, the bin is skipped:
315
316 if (Na < Ns)
317 continue;
318
319 const Double_t eff = Ns/Na;
320
321 const Double_t err = sqrt((1.-eff)*Ns)/Na;
322
323 const Float_t area = dr * cos(angle);
324
325 fHistCol->SetBinContent(ix, eff*area);
326 fHistCol->SetBinError(ix, err*area);
327 }
328
329 delete &histsel;
330
331 SetReadyToSave();
332}
333
334// --------------------------------------------------------------------------
335//
336// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
337// flag
338//
339void MHMcCollectionArea::Calc(const MHMcEfficiency &heff)
340{
341 //
342 // now calculate the Collection area for different
343 // energies
344 //
345 TH2D &h = (TH2D&)*heff.GetHist();
346
347 MH::SetBinning(fHistCol, h.GetXaxis());
348
349 const Int_t nbinx = h.GetXaxis()->GetNbins();
350 const Int_t nbiny = h.GetYaxis()->GetNbins();
351
352 for (Int_t ix=1; ix<=nbinx; ix++)
353 {
354 Double_t errA = 0;
355 Double_t colA = 0;
356
357 for (Int_t iy=1; iy<=nbiny; iy++)
358 {
359 TAxis *yaxis = h.GetYaxis();
360
361 const Double_t r1 = yaxis->GetBinLowEdge(iy);
362 const Double_t r2 = yaxis->GetBinLowEdge(iy+1);
363
364 const Double_t A = TMath::Pi() * (r2*r2 - r1*r1);
365
366 const Double_t eff = h.GetCellContent(ix, iy);
367 const Double_t err = h.GetCellError(ix, iy);
368
369 colA += eff*A;
370 errA += A*A * err*err;
371 }
372
373 fHistCol->SetBinContent(ix, colA);
374 fHistCol->SetBinError(ix, sqrt(errA));
375 }
376
377 SetReadyToSave();
378}
Note: See TracBrowser for help on using the repository browser.