source: trunk/MagicSoft/Mars/mhistmc/MHMcCT1CollectionArea.cc@ 1997

Last change on this file since 1997 was 1992, checked in by tbretz, 22 years ago
*** empty log message ***
File size: 10.6 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): A. Moralejo 3/2003 <mailto:moralejo@pd.infn.it>
19!
20! Copyright: MAGIC Software Development, 2000-2003
21!
22!
23\* ======================================================================== */
24
25//////////////////////////////////////////////////////////////////////////////
26// //
27// MHMcCT1CollectionArea //
28// //
29//////////////////////////////////////////////////////////////////////////////
30
31#include "MHMcCT1CollectionArea.h"
32
33#include <TH2.h>
34#include <TCanvas.h>
35
36#include "MMcEvt.hxx"
37#include "MH.h"
38#include "MBinning.h"
39#include "MParList.h"
40#include "MLog.h"
41#include "MLogManip.h"
42
43
44ClassImp(MHMcCT1CollectionArea);
45
46// --------------------------------------------------------------------------
47//
48// Creates the three necessary histograms:
49// - selected showers (input)
50// - all showers (input)
51// - collection area (result)
52//
53MHMcCT1CollectionArea::MHMcCT1CollectionArea(const char *name, const char *title)
54{
55 //
56 // nbins, minEnergy, maxEnergy defaults:
57 // we set the energy range from 100 Gev to 30000 GeV (in log, 3.5 orders
58 // of magnitude) and for each order we take 10 subdivisions --> 35 xbins
59 // we set the theta range from 12.5 to 48 deg, with 6 bins (the latter
60 // choice has been done to make the bin centers as close as possible to
61 // the actual zenith angles in the CT1 MC sample).
62 //
63
64 fName = name ? name : "MHMcCT1CollectionArea";
65 fTitle = title ? title : "Collection Area vs. log10 Energy";
66
67 fHistAll = new TH2D;
68 fHistSel = new TH2D;
69 fHistCol = new TH2D;
70
71 fHistCol->SetName(fName);
72 fHistAll->SetName("AllEvents");
73 fHistSel->SetName("SelectedEvents");
74
75 fHistCol->SetTitle(fTitle);
76 fHistAll->SetTitle("All showers - Theta vs log10 Energy distribution");
77 fHistSel->SetTitle("Selected showers - Theta vs log10 Energy distribution");
78
79 fHistAll->SetDirectory(NULL);
80 fHistSel->SetDirectory(NULL);
81 fHistCol->SetDirectory(NULL);
82
83 fHistAll->SetXTitle("log10 E [GeV]");
84 fHistAll->SetYTitle("theta [deg]");
85 fHistAll->SetZTitle("N");
86
87 fHistSel->SetXTitle("log10 E [GeV]");
88 fHistSel->SetYTitle("theta [deg]");
89 fHistSel->SetZTitle("N");
90
91 fHistCol->SetXTitle("log10 E [GeV]");
92 fHistCol->SetYTitle("theta [deg]");
93 fHistCol->SetZTitle("A [m^{2}]");
94}
95
96// --------------------------------------------------------------------------
97//
98// Delete the three histograms
99//
100MHMcCT1CollectionArea::~MHMcCT1CollectionArea()
101{
102 delete fHistAll;
103 delete fHistSel;
104 delete fHistCol;
105}
106
107// --------------------------------------------------------------------------
108//
109// Set the binnings and prepare the filling of the histograms
110//
111Bool_t MHMcCT1CollectionArea::SetupFill(const MParList *plist)
112{
113 const MBinning* binsenergy = (MBinning*)plist->FindObject("BinningE");
114 const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
115
116 if (!binsenergy || !binstheta)
117 {
118 *fLog << err << dbginf << "At least one MBinning not found... aborting.";
119 *fLog << endl;
120 return kFALSE;
121 }
122
123 SetBinning(fHistAll, binsenergy, binstheta);
124 SetBinning(fHistSel, binsenergy, binstheta);
125 SetBinning(fHistCol, binsenergy, binstheta);
126
127 fHistAll->Sumw2();
128 fHistSel->Sumw2();
129 fHistCol->Sumw2();
130
131 return kTRUE;
132}
133
134
135// --------------------------------------------------------------------------
136//
137// Fill data into the histogram which contains the selected showers
138//
139Bool_t MHMcCT1CollectionArea::Fill(const MParContainer *par, Double_t w)
140{
141 MMcEvt &mcevt = *(MMcEvt*)par;
142
143 fHistSel->Fill(log10(mcevt.GetEnergy()), kRad2Deg*mcevt.GetTelescopeTheta());
144 return kTRUE;
145}
146
147// --------------------------------------------------------------------------
148//
149// Draw the histogram with all showers
150//
151void MHMcCT1CollectionArea::DrawAll(Option_t* option)
152{
153 if (!gPad)
154 MH::MakeDefCanvas(fHistAll);
155
156 fHistAll->Draw(option);
157
158 gPad->Modified();
159 gPad->Update();
160}
161
162// --------------------------------------------------------------------------
163//
164// Draw the histogram with the selected showers only.
165//
166void MHMcCT1CollectionArea::DrawSel(Option_t* option)
167{
168 if (!gPad)
169 MH::MakeDefCanvas(fHistSel);
170
171 fHistSel->Draw(option);
172
173 gPad->Modified();
174 gPad->Update();
175}
176
177// --------------------------------------------------------------------------
178//
179// Creates a new canvas and draws the histogram into it.
180// Be careful: The histogram belongs to this object and won't get deleted
181// together with the canvas.
182//
183TObject *MHMcCT1CollectionArea::DrawClone(Option_t* option) const
184{
185 TCanvas *c = MH::MakeDefCanvas(fHistCol);
186
187 //
188 // This is necessary to get the expected behaviour of DrawClone
189 //
190 gROOT->SetSelectedPad(NULL);
191
192 fHistCol->DrawCopy(option);
193
194 c->Modified();
195 c->Update();
196
197 return c;
198}
199
200void MHMcCT1CollectionArea::Draw(Option_t* option)
201{
202 if (!gPad)
203 MH::MakeDefCanvas(fHistCol);
204
205 fHistCol->Draw(option);
206
207 gPad->Modified();
208 gPad->Update();
209}
210
211//
212// Calculate the Efficiency (collection area) for the CT1 MC sample
213// and set the 'ReadyToSave' flag
214//
215void MHMcCT1CollectionArea::CalcEfficiency()
216{
217 //
218 // Here we estimate the total number of showers in each energy bin
219 // from the known the energy range and spectral index of the generated
220 // showers. This procedure is intended for the CT1 MC files. The total
221 // number of generated events, collection area, spectral index etc will be
222 // set here by hand, so make sure that the MC sample you are using is the
223 // right one (check all these quantities in your files and compare with
224 // what is written below. In some theta bins, there are two different
225 // productions, with different energy limits but with the same spectral
226 // slope. We account for this when calculating the original number of
227 // events in each energy bin.
228 //
229
230 for (Int_t thetabin = 1; thetabin <= fHistAll->GetNbinsY(); thetabin++)
231 {
232 // This theta is not exactly the one of the MC events, just about
233 // the same:
234 Float_t theta = fHistAll->GetYaxis()->GetBinCenter(thetabin);
235
236 Float_t emin1, emax1, emin2, emax2;
237 Float_t index, expo, k1, k2;
238 Float_t numevts1, numevts2;
239 Float_t r1, r2; // Impact parameter range (on ground).
240
241 emin1 = 0; emax1 = 0; emin2 = 0; emax2 = 0;
242 expo = 0.; k1 = 0.; k2 = 0.; r1 = 0.; r2 = 0.;
243 numevts1 = 0; numevts2 = 0;
244
245 if (theta > 14 && theta < 16) // 15 deg
246 {
247 r1 = 0.;
248 r2 = 250.; //meters
249 emin1 = 300.;
250 emax1 = 400.; // Energies in GeV.
251 emin2 = 400.;
252 emax2 = 30000.;
253 numevts1 = 4000.;
254 numevts2 = 25740.;
255 }
256 else if (theta > 20 && theta < 21) // 20.5 deg
257 {
258 r1 = 0.;
259 r2 = 263.; //meters
260 emin1 = 300.;
261 emax1 = 400.; // Energies in GeV.
262 emin2 = 400.;
263 emax2 = 30000.;
264 numevts1 = 6611.;
265 numevts2 = 24448.;
266 }
267 else if (theta > 26 && theta < 27) // 26.5 degrees
268 {
269 r1 = 0.;
270 r2 = 290.; //meters
271 emin1 = 300.;
272 emax1 = 400.; // Energies in GeV.
273 emax2 = emax1; emin2 = 400.;
274 emax2 = 30000.;
275 numevts1 = 4000.;
276 numevts2 = 26316.;
277 }
278 else if (theta > 32 && theta < 33) // 32.5 degrees
279 {
280 r1 = 0.;
281 r2 = 350.; //meters
282 emin1 = 300.;
283 emax1 = 30000.; // Energies in GeV.
284 emax2 = emax1;
285 numevts1 = 33646.;
286 }
287 else if (theta > 38 && theta < 39) // 38.75 degrees
288 {
289 r1 = 0.;
290 r2 = 380.; //meters
291 emin1 = 300.;
292 emax1 = 30000.; // Energies in GeV.
293 emax2 = emax1;
294 numevts1 = 38415.;
295 }
296 else if (theta > 44 && theta < 46) // 45 degrees
297 {
298 r1 = 0.;
299 r2 = 565.; //meters
300 emin1 = 300.;
301 emax1 = 50000.; // Energies in GeV.
302 emax2 = emax1;
303 numevts1 = 30197.;
304 }
305
306 index = 1.5; // Differential spectral Index.
307 expo = 1.-index;
308 k1 = numevts1 / (pow(emax1,expo) - pow(emin1,expo));
309 k2 = numevts2 / (pow(emax2,expo) - pow(emin2,expo));
310
311 for (Int_t i=1; i <= fHistAll->GetNbinsX(); i++)
312 {
313 const Float_t e1 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i));
314 const Float_t e2 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i+1));
315
316 if (e1 < emin1 || e2 > emax2)
317 continue;
318
319 Float_t events;
320
321 if (e2 <= emax1)
322 events = k1 * (pow(e2, expo) - pow(e1, expo));
323 else if (e1 >= emin2)
324 events = k2 * (pow(e2, expo) - pow(e1, expo));
325 else
326 events =
327 k1 * (pow(emax1, expo) - pow(e1, expo))+
328 k2 * (pow(e2, expo) - pow(emin2, expo));
329
330 fHistAll->SetBinContent(i, thetabin, events);
331 fHistAll->SetBinError(i, thetabin, sqrt(events));
332 }
333
334 // -----------------------------------------------------------
335
336 const Float_t dr = TMath::Pi() * (r2*r2 - r1*r1);
337
338 for (Int_t ix = 1; ix <= fHistAll->GetNbinsX(); ix++)
339 {
340 const Float_t Na = fHistAll->GetBinContent(ix,thetabin);
341
342 if (Na <= 0)
343 {
344 //
345 // If energy is large, this case means that no or very few events
346 // were generated at this energy bin. In this case we assign it
347 // the effective area of the bin below it in energy. If energy is
348 // below 1E4, it means that no events triggered -> eff area = 0
349 //
350
351 if (fHistSel->GetXaxis()->GetBinLowEdge(ix) > 4.)
352 {
353 fHistCol->SetBinContent(ix, thetabin, fHistCol->GetBinContent(ix-1, thetabin));
354 fHistCol->SetBinError(ix, thetabin, fHistCol->GetBinError(ix-1, thetabin));
355 }
356 continue;
357 }
358
359 const Float_t Ns = fHistSel->GetBinContent(ix,thetabin);
360
361 // Since Na is an estimate of the total number of showers generated
362 // in the energy bin, it may happen that Ns (triggered showers) is
363 // larger than Na. In that case, the bin is skipped:
364
365 if (Na < Ns)
366 continue;
367
368 const Double_t eff = Ns/Na;
369 const Double_t efferr = sqrt((1.-eff)*Ns)/Na;
370
371
372 const Float_t area = dr * cos(theta*TMath::Pi()/180.);
373
374 fHistCol->SetBinContent(ix, thetabin, eff*area);
375 fHistCol->SetBinError(ix, thetabin, efferr*area);
376
377 }
378 }
379
380 SetReadyToSave();
381}
Note: See TracBrowser for help on using the repository browser.