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

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