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 | ClassImp(MHMcCT1CollectionArea);
|
---|
44 |
|
---|
45 | using namespace std;
|
---|
46 |
|
---|
47 | // --------------------------------------------------------------------------
|
---|
48 | //
|
---|
49 | // Creates the three necessary histograms:
|
---|
50 | // - selected showers (input)
|
---|
51 | // - all showers (input)
|
---|
52 | // - collection area (result)
|
---|
53 | //
|
---|
54 | MHMcCT1CollectionArea::MHMcCT1CollectionArea(const char *name, const char *title)
|
---|
55 | {
|
---|
56 | //
|
---|
57 | // nbins, minEnergy, maxEnergy defaults:
|
---|
58 | // we set the energy range from 100 Gev to 30000 GeV (in log, 3.5 orders
|
---|
59 | // of magnitude) and for each order we take 10 subdivisions --> 35 xbins
|
---|
60 | // we set the theta range from 12.5 to 48 deg, with 6 bins (the latter
|
---|
61 | // choice has been done to make the bin centers as close as possible to
|
---|
62 | // the actual zenith angles in the CT1 MC sample).
|
---|
63 | //
|
---|
64 |
|
---|
65 | fName = name ? name : "MHMcCT1CollectionArea";
|
---|
66 | fTitle = title ? title : "Collection Area vs. log10 Energy";
|
---|
67 |
|
---|
68 | fEaxis = kLog10Energy;
|
---|
69 |
|
---|
70 | fHistAll = new TH2D;
|
---|
71 | fHistSel = new TH2D;
|
---|
72 | fHistCol = new TH2D;
|
---|
73 |
|
---|
74 | fHistCol->SetName(fName);
|
---|
75 | fHistAll->SetName("AllEvents");
|
---|
76 | fHistSel->SetName("SelectedEvents");
|
---|
77 |
|
---|
78 | fHistCol->SetTitle(fTitle);
|
---|
79 | fHistAll->SetTitle("All showers - Theta vs log10 Energy distribution");
|
---|
80 | fHistSel->SetTitle("Selected showers - Theta vs log10 Energy distribution");
|
---|
81 |
|
---|
82 | fHistAll->SetDirectory(NULL);
|
---|
83 | fHistSel->SetDirectory(NULL);
|
---|
84 | fHistCol->SetDirectory(NULL);
|
---|
85 |
|
---|
86 | fHistAll->SetXTitle("log10 E [GeV]");
|
---|
87 | fHistAll->SetYTitle("theta [deg]");
|
---|
88 | fHistAll->SetZTitle("N");
|
---|
89 |
|
---|
90 | fHistSel->SetXTitle("log10 E [GeV]");
|
---|
91 | fHistSel->SetYTitle("theta [deg]");
|
---|
92 | fHistSel->SetZTitle("N");
|
---|
93 |
|
---|
94 | fHistCol->SetXTitle("log10 E [GeV]");
|
---|
95 | fHistCol->SetYTitle("theta [deg]");
|
---|
96 | fHistCol->SetZTitle("A [m^{2}]");
|
---|
97 |
|
---|
98 | }
|
---|
99 |
|
---|
100 | // --------------------------------------------------------------------------
|
---|
101 | //
|
---|
102 | // Delete the three histograms
|
---|
103 | //
|
---|
104 | MHMcCT1CollectionArea::~MHMcCT1CollectionArea()
|
---|
105 | {
|
---|
106 | delete fHistAll;
|
---|
107 | delete fHistSel;
|
---|
108 | delete fHistCol;
|
---|
109 | }
|
---|
110 |
|
---|
111 | // --------------------------------------------------------------------------
|
---|
112 | //
|
---|
113 | // Set the binnings and prepare the filling of the histograms
|
---|
114 | //
|
---|
115 | Bool_t MHMcCT1CollectionArea::SetupFill(const MParList *plist)
|
---|
116 | {
|
---|
117 | const MBinning* binsenergy = (MBinning*)plist->FindObject("BinningE");
|
---|
118 | const MBinning* binstheta = (MBinning*)plist->FindObject("BinningTheta");
|
---|
119 |
|
---|
120 | if (!binsenergy || !binstheta)
|
---|
121 | {
|
---|
122 | *fLog << err << dbginf << "At least one MBinning not found... aborting.";
|
---|
123 | *fLog << endl;
|
---|
124 | return kFALSE;
|
---|
125 | }
|
---|
126 |
|
---|
127 | SetBinning(fHistAll, binsenergy, binstheta);
|
---|
128 | SetBinning(fHistSel, binsenergy, binstheta);
|
---|
129 | SetBinning(fHistCol, binsenergy, binstheta);
|
---|
130 |
|
---|
131 | fHistAll->Sumw2();
|
---|
132 | fHistSel->Sumw2();
|
---|
133 | fHistCol->Sumw2();
|
---|
134 |
|
---|
135 | if (fEaxis == kEnergy)
|
---|
136 | {
|
---|
137 | fTitle = "Collection Area vs. Energy";
|
---|
138 | fHistCol->SetTitle(fTitle);
|
---|
139 | fHistAll->SetTitle("All showers - Theta vs Energy distribution");
|
---|
140 | fHistSel->SetTitle("Selected showers - Theta vs Energy distribution");
|
---|
141 | fHistCol->SetXTitle("E [GeV]");
|
---|
142 | fHistAll->SetXTitle("E [GeV]");
|
---|
143 | fHistSel->SetXTitle("E [GeV]");
|
---|
144 | }
|
---|
145 |
|
---|
146 | return kTRUE;
|
---|
147 | }
|
---|
148 |
|
---|
149 |
|
---|
150 | // --------------------------------------------------------------------------
|
---|
151 | //
|
---|
152 | // Fill data into the histogram which contains the selected showers
|
---|
153 | //
|
---|
154 | Bool_t MHMcCT1CollectionArea::Fill(const MParContainer *par, const Stat_t w)
|
---|
155 | {
|
---|
156 | MMcEvt &mcevt = *(MMcEvt*)par;
|
---|
157 |
|
---|
158 | if (fEaxis == kLog10Energy)
|
---|
159 | fHistSel->Fill(log10(mcevt.GetEnergy()), kRad2Deg*mcevt.GetTelescopeTheta(), w);
|
---|
160 | else
|
---|
161 | fHistSel->Fill(mcevt.GetEnergy(), kRad2Deg*mcevt.GetTelescopeTheta(), w);
|
---|
162 |
|
---|
163 | return kTRUE;
|
---|
164 | }
|
---|
165 |
|
---|
166 | // --------------------------------------------------------------------------
|
---|
167 | //
|
---|
168 | // Draw the histogram with all showers
|
---|
169 | //
|
---|
170 | void MHMcCT1CollectionArea::DrawAll(Option_t* option)
|
---|
171 | {
|
---|
172 | if (!gPad)
|
---|
173 | MH::MakeDefCanvas(fHistAll);
|
---|
174 |
|
---|
175 | fHistAll->Draw(option);
|
---|
176 |
|
---|
177 | gPad->Modified();
|
---|
178 | gPad->Update();
|
---|
179 | }
|
---|
180 |
|
---|
181 | // --------------------------------------------------------------------------
|
---|
182 | //
|
---|
183 | // Draw the histogram with the selected showers only.
|
---|
184 | //
|
---|
185 | void MHMcCT1CollectionArea::DrawSel(Option_t* option)
|
---|
186 | {
|
---|
187 | if (!gPad)
|
---|
188 | MH::MakeDefCanvas(fHistSel);
|
---|
189 |
|
---|
190 | fHistSel->Draw(option);
|
---|
191 |
|
---|
192 | gPad->Modified();
|
---|
193 | gPad->Update();
|
---|
194 | }
|
---|
195 |
|
---|
196 | // --------------------------------------------------------------------------
|
---|
197 | //
|
---|
198 | // Creates a new canvas and draws the histogram into it.
|
---|
199 | // Be careful: The histogram belongs to this object and won't get deleted
|
---|
200 | // together with the canvas.
|
---|
201 | //
|
---|
202 | TObject *MHMcCT1CollectionArea::DrawClone(Option_t* option) const
|
---|
203 | {
|
---|
204 | TCanvas &c = *MakeDefCanvas("CollArea", "Collection area plots", 600, 600);
|
---|
205 | c.Divide(2,2);
|
---|
206 |
|
---|
207 | //
|
---|
208 | // This is necessary to get the expected behaviour of DrawClone
|
---|
209 | //
|
---|
210 | gROOT->SetSelectedPad(NULL);
|
---|
211 |
|
---|
212 | c.cd(1);
|
---|
213 | if (fEaxis == kEnergy)
|
---|
214 | gPad->SetLogx();
|
---|
215 | fHistCol->SetDirectory(NULL);
|
---|
216 | fHistCol->DrawCopy(option);
|
---|
217 |
|
---|
218 | c.cd(2);
|
---|
219 | if (fEaxis == kEnergy)
|
---|
220 | gPad->SetLogx();
|
---|
221 | fHistSel->SetDirectory(NULL);
|
---|
222 | fHistSel->DrawCopy(option);
|
---|
223 |
|
---|
224 | c.cd(3);
|
---|
225 | if (fEaxis == kEnergy)
|
---|
226 | gPad->SetLogx();
|
---|
227 | fHistAll->SetDirectory(NULL);
|
---|
228 | fHistAll->DrawCopy(option);
|
---|
229 |
|
---|
230 |
|
---|
231 | c.Modified();
|
---|
232 | c.Update();
|
---|
233 |
|
---|
234 | return &c;
|
---|
235 | }
|
---|
236 |
|
---|
237 | void MHMcCT1CollectionArea::Draw(Option_t* option)
|
---|
238 | {
|
---|
239 | if (!gPad)
|
---|
240 | MH::MakeDefCanvas(fHistCol);
|
---|
241 |
|
---|
242 | fHistCol->Draw(option);
|
---|
243 |
|
---|
244 | gPad->Modified();
|
---|
245 | gPad->Update();
|
---|
246 | }
|
---|
247 |
|
---|
248 | //
|
---|
249 | // Calculate the Efficiency (collection area) for the CT1 MC sample
|
---|
250 | // and set the 'ReadyToSave' flag
|
---|
251 | //
|
---|
252 | void MHMcCT1CollectionArea::CalcEfficiency()
|
---|
253 | {
|
---|
254 | //
|
---|
255 | // Here we estimate the total number of showers in each energy bin
|
---|
256 | // from the known the energy range and spectral index of the generated
|
---|
257 | // showers. This procedure is intended for the CT1 MC files. The total
|
---|
258 | // number of generated events, collection area, spectral index etc will be
|
---|
259 | // set here by hand, so make sure that the MC sample you are using is the
|
---|
260 | // right one (check all these quantities in your files and compare with
|
---|
261 | // what is written below. In some theta bins, there are two different
|
---|
262 | // productions, with different energy limits but with the same spectral
|
---|
263 | // slope. We account for this when calculating the original number of
|
---|
264 | // events in each energy bin.
|
---|
265 | //
|
---|
266 | // The theta angle with which the MC data (from D. Kranich) were produced
|
---|
267 | // is not exactly the center of the theta bins we are using (the bin limits
|
---|
268 | // should be 0.0, 17.5, 23.5, 29.5, 35.5, 42., 50.). The theta variable in
|
---|
269 | // the MC root file has nevertheless been changed (W.Wittek) to correspond
|
---|
270 | // to the centers of these bins. Only in the first bin is the difference big:
|
---|
271 | // the data were produced at theta = 15 degrees, whreas the bin center is at
|
---|
272 | // 8.75 degrees. Howeverm at such low z.a. the shower characteristics change
|
---|
273 | // very slowly with theta.
|
---|
274 | //
|
---|
275 | //
|
---|
276 | //
|
---|
277 | // Only for the binning taken from D. Kranich :
|
---|
278 | //
|
---|
279 |
|
---|
280 | for (Int_t thetabin = 1; thetabin <= fHistAll->GetNbinsY(); thetabin++)
|
---|
281 | {
|
---|
282 | // This theta is not exactly the one of the MC events, just about
|
---|
283 | // the same (bins have been selected so):
|
---|
284 |
|
---|
285 | Float_t theta = fHistAll->GetYaxis()->GetBinCenter(thetabin);
|
---|
286 |
|
---|
287 | Float_t emin[4]; // Minimum energy in MC sample
|
---|
288 | Float_t emax[4]; // Maximum energy in MC sample
|
---|
289 | Float_t index[4]; // Spectral index
|
---|
290 | Float_t numevts[4]; // Number of events
|
---|
291 | Float_t multfactor[4]; // Factor by which the original number of events in an MC
|
---|
292 | // sample has been multiplied to account for the differences
|
---|
293 | // in the generation areas of the various samples.
|
---|
294 | Float_t rmax; // Maximum impact parameter range (on ground up to 45 degrees,
|
---|
295 | // on a plane perpendicular to Shower axis for 55 and 65 deg).
|
---|
296 |
|
---|
297 | memset(emin, 0, 4*sizeof(Float_t));
|
---|
298 | memset(emax, 0, 4*sizeof(Float_t));
|
---|
299 | memset(index, 0, 4*sizeof(Float_t));
|
---|
300 | memset(numevts, 0, 4*sizeof(Float_t));
|
---|
301 | rmax = 0.;
|
---|
302 |
|
---|
303 | multfactor[0] = 1.;
|
---|
304 | multfactor[1] = 1.;
|
---|
305 | multfactor[2] = 1.;
|
---|
306 | multfactor[3] = 1.;
|
---|
307 |
|
---|
308 | //
|
---|
309 | // rmin and rmax are the minimum and maximum values of the impact
|
---|
310 | // parameter of the shower on the ground (horizontal plane).
|
---|
311 | //
|
---|
312 |
|
---|
313 | Int_t num_MC_samples = 0;
|
---|
314 |
|
---|
315 | if (theta > 8 && theta < 9) // 8.75 deg
|
---|
316 | {
|
---|
317 | emin[0] = 300.;
|
---|
318 | emax[0] = 400.; // Energies in GeV.
|
---|
319 | index[0] = 1.5;
|
---|
320 | numevts[0] = 4000.;
|
---|
321 |
|
---|
322 | emin[1] = 400.;
|
---|
323 | emax[1] = 30000.;
|
---|
324 | index[1] = 1.5;
|
---|
325 | numevts[1] = 25740.;
|
---|
326 |
|
---|
327 | rmax = 250.; //meters
|
---|
328 | num_MC_samples = 2;
|
---|
329 | }
|
---|
330 | else if (theta > 20 && theta < 21) // 20.5 deg
|
---|
331 | {
|
---|
332 | emin[0] = 300.;
|
---|
333 | emax[0] = 400.; // Energies in GeV.
|
---|
334 | index[0] = 1.5;
|
---|
335 | numevts[0] = 6611.;
|
---|
336 |
|
---|
337 | emin[1] = 400.;
|
---|
338 | emax[1] = 30000.;
|
---|
339 | index[1] = 1.5;
|
---|
340 | numevts[1] = 24448.;
|
---|
341 |
|
---|
342 | rmax = 263.;
|
---|
343 | num_MC_samples = 2;
|
---|
344 | }
|
---|
345 | else if (theta > 26 && theta < 27) // 26.5 degrees
|
---|
346 | {
|
---|
347 | emin[0] = 300.;
|
---|
348 | emax[0] = 400.; // Energies in GeV.
|
---|
349 | index[0] = 1.5;
|
---|
350 | numevts[0] = 4000.;
|
---|
351 |
|
---|
352 | emin[1] = 400.;
|
---|
353 | emax[1] = 30000.;
|
---|
354 | index[1] = 1.5;
|
---|
355 | numevts[1] = 26316.;
|
---|
356 |
|
---|
357 | rmax = 290.; //meters
|
---|
358 | num_MC_samples = 2;
|
---|
359 | }
|
---|
360 | else if (theta > 32 && theta < 33) // 32.5 degrees
|
---|
361 | {
|
---|
362 | emin[0] = 300.;
|
---|
363 | emax[0] = 30000.; // Energies in GeV.
|
---|
364 | index[0] = 1.5;
|
---|
365 | numevts[0] = 33646.;
|
---|
366 |
|
---|
367 | rmax = 350.; //meters
|
---|
368 | num_MC_samples = 1;
|
---|
369 | }
|
---|
370 | else if (theta > 38 && theta < 39) // 38.75 degrees
|
---|
371 | {
|
---|
372 | emin[0] = 300.;
|
---|
373 | emax[0] = 30000.; // Energies in GeV.
|
---|
374 | index[0] = 1.5;
|
---|
375 | numevts[0] = 38415.;
|
---|
376 |
|
---|
377 | rmax = 380.; //meters
|
---|
378 | num_MC_samples = 1;
|
---|
379 | }
|
---|
380 | else if (theta > 45 && theta < 47) // 46 degrees
|
---|
381 | {
|
---|
382 | emin[0] = 300.;
|
---|
383 | emax[0] = 50000.; // Energies in GeV.
|
---|
384 | index[0] = 1.5;
|
---|
385 | numevts[0] = 30197.;
|
---|
386 |
|
---|
387 | rmax = 565.; //meters
|
---|
388 | num_MC_samples = 1;
|
---|
389 | }
|
---|
390 | else if (theta > 54 && theta < 56) // 55 degrees
|
---|
391 | {
|
---|
392 | //
|
---|
393 | // The value of numevts in the first sample (below) has been
|
---|
394 | // changed to simplify calculations. We have multiplied it
|
---|
395 | // times 1.2808997 to convert it to the number it would be if
|
---|
396 | // the generation area was equal to that of the other samples
|
---|
397 | // at 55 degrees (pi*600**2 m2). This has to be taken into account
|
---|
398 | // in the error in the number of events.
|
---|
399 | //
|
---|
400 |
|
---|
401 | emin[0] = 500.;
|
---|
402 | emax[0] = 50000.; // Energies in GeV.
|
---|
403 | index[0] = 1.5;
|
---|
404 | numevts[0] = 3298.;
|
---|
405 | multfactor[0] = 1.2808997;
|
---|
406 |
|
---|
407 | emin[1] = 1500.;
|
---|
408 | emax[1] = 50000.; // Energies in GeV.
|
---|
409 | index[1] = 1.5;
|
---|
410 | numevts[1] = 22229.;
|
---|
411 |
|
---|
412 | emin[2] = 1500.;
|
---|
413 | emax[2] = 50000.; // Energies in GeV.
|
---|
414 | index[2] = 1.7;
|
---|
415 | numevts[2] = 7553.;
|
---|
416 |
|
---|
417 | rmax = 600; //meters
|
---|
418 | num_MC_samples = 3;
|
---|
419 | }
|
---|
420 |
|
---|
421 | else if (theta > 64 && theta < 66) // 65 degrees
|
---|
422 | {
|
---|
423 | emin[0] = 2000.;
|
---|
424 | emax[0] = 50000.; // Energies in GeV.
|
---|
425 | index[0] = 1.5;
|
---|
426 | numevts[0] = 16310.;
|
---|
427 |
|
---|
428 | emin[1] = 2000.;
|
---|
429 | emax[1] = 50000.; // Energies in GeV.
|
---|
430 | index[1] = 1.7;
|
---|
431 | numevts[1] = 3000.;
|
---|
432 |
|
---|
433 | //
|
---|
434 | // The value of numevts in the next two samples (below) has been
|
---|
435 | // changed to simplify calculations. We have converted them to the
|
---|
436 | // number it would be if the generation area was equal to that of
|
---|
437 | // the first two samples at 65 degrees (pi*800**2 m2) (four times
|
---|
438 | // as many, since the original maximum impact parameter was 400
|
---|
439 | // instead of 800. This is taken into account in the error too.
|
---|
440 | //
|
---|
441 |
|
---|
442 | emin[2] = 5000.;
|
---|
443 | emax[2] = 50000.; // Energies in GeV.
|
---|
444 | index[2] = 1.5;
|
---|
445 | numevts[2] = 56584.;
|
---|
446 | multfactor[2] = 4;
|
---|
447 |
|
---|
448 | emin[3] = 5000.;
|
---|
449 | emax[3] = 50000.; // Energies in GeV.
|
---|
450 | index[3] = 1.7;
|
---|
451 | numevts[3] = 11464;
|
---|
452 | multfactor[3] = 4;
|
---|
453 |
|
---|
454 | rmax = 800; // meters
|
---|
455 | num_MC_samples = 4;
|
---|
456 | }
|
---|
457 |
|
---|
458 |
|
---|
459 | for (Int_t i=1; i <= fHistAll->GetNbinsX(); i++)
|
---|
460 | {
|
---|
461 | Float_t e1;
|
---|
462 | Float_t e2;
|
---|
463 |
|
---|
464 | if (fEaxis == kLog10Energy)
|
---|
465 | {
|
---|
466 | e1 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i));
|
---|
467 | e2 = pow(10.,fHistAll->GetXaxis()->GetBinLowEdge(i+1));
|
---|
468 | }
|
---|
469 | else
|
---|
470 | {
|
---|
471 | e1 = fHistAll->GetXaxis()->GetBinLowEdge(i);
|
---|
472 | e2 = fHistAll->GetXaxis()->GetBinLowEdge(i+1);
|
---|
473 | }
|
---|
474 |
|
---|
475 | Float_t events = 0.;
|
---|
476 | Float_t errevents = 0.;
|
---|
477 |
|
---|
478 | for (Int_t sample = 0; sample < num_MC_samples; sample++)
|
---|
479 | {
|
---|
480 | Float_t expo = 1.-index[sample];
|
---|
481 | Float_t k = numevts[sample] / (pow(emax[sample],expo) - pow(emin[sample],expo));
|
---|
482 |
|
---|
483 | if (e2 < emin[sample] || e1 > emax[sample])
|
---|
484 | continue;
|
---|
485 |
|
---|
486 | if (emin[sample] > e1)
|
---|
487 | e1 = emin[sample];
|
---|
488 |
|
---|
489 | if (emax[sample] < e2)
|
---|
490 | e2 = emax[sample];
|
---|
491 |
|
---|
492 | events += k * (pow(e2, expo) - pow(e1, expo));
|
---|
493 | errevents += multfactor[sample] * events;
|
---|
494 | }
|
---|
495 |
|
---|
496 | errevents= sqrt(errevents);
|
---|
497 |
|
---|
498 | fHistAll->SetBinContent(i, thetabin, events);
|
---|
499 | fHistAll->SetBinError(i, thetabin, errevents);
|
---|
500 | }
|
---|
501 |
|
---|
502 | // -----------------------------------------------------------
|
---|
503 |
|
---|
504 | const Float_t dr = TMath::Pi() * rmax * rmax;
|
---|
505 |
|
---|
506 | for (Int_t ix = 1; ix <= fHistAll->GetNbinsX(); ix++)
|
---|
507 | {
|
---|
508 | const Float_t Na = fHistAll->GetBinContent(ix,thetabin);
|
---|
509 |
|
---|
510 | if (Na <= 0)
|
---|
511 | {
|
---|
512 | //
|
---|
513 | // If energy is large, this case means that no or very few events
|
---|
514 | // were generated at this energy bin. In this case we assign it
|
---|
515 | // the effective area of the bin below it in energy. If energy is
|
---|
516 | // below 1E4, it means that no events triggered -> eff area = 0
|
---|
517 | //
|
---|
518 | // NOW DISABLED: because collection area after analysis does not
|
---|
519 | // saturate at high E!
|
---|
520 | //
|
---|
521 |
|
---|
522 | /*
|
---|
523 | if (fHistSel->GetXaxis()->GetBinLowEdge(ix) > 4.)
|
---|
524 | {
|
---|
525 | fHistCol->SetBinContent(ix, thetabin, fHistCol->GetBinContent(ix-1, thetabin));
|
---|
526 | fHistCol->SetBinError(ix, thetabin, fHistCol->GetBinError(ix-1, thetabin));
|
---|
527 | }
|
---|
528 | */
|
---|
529 | continue;
|
---|
530 | }
|
---|
531 |
|
---|
532 | const Float_t Ns = fHistSel->GetBinContent(ix,thetabin);
|
---|
533 |
|
---|
534 | // Since Na is an estimate of the total number of showers generated
|
---|
535 | // in the energy bin, it may happen that Ns (triggered showers) is
|
---|
536 | // larger than Na. In that case, the bin is skipped:
|
---|
537 |
|
---|
538 | if (Na < Ns)
|
---|
539 | continue;
|
---|
540 |
|
---|
541 | const Double_t eff = Ns/Na;
|
---|
542 | const Double_t efferr = sqrt((1.-eff)*Ns)/Na;
|
---|
543 |
|
---|
544 | //
|
---|
545 | // Now we get the total area, perpendicular to the observation direction
|
---|
546 | // in which the events were generated (correct for cos theta):
|
---|
547 | //
|
---|
548 |
|
---|
549 | Float_t area;
|
---|
550 |
|
---|
551 | if (theta < 50)
|
---|
552 | area = dr * cos(theta*TMath::Pi()/180.);
|
---|
553 | else
|
---|
554 | area = dr;
|
---|
555 |
|
---|
556 | // Above 50 degrees MC was generated with Corsika 6.xx, and the cores
|
---|
557 | // were distributed on a circle perpendicular to the observation direction,
|
---|
558 | // and not on ground, hence the correction for cos(theta) is not necessary.
|
---|
559 | //
|
---|
560 |
|
---|
561 |
|
---|
562 | fHistCol->SetBinContent(ix, thetabin, eff*area);
|
---|
563 | fHistCol->SetBinError(ix, thetabin, efferr*area);
|
---|
564 |
|
---|
565 | }
|
---|
566 | }
|
---|
567 |
|
---|
568 | SetReadyToSave();
|
---|
569 | }
|
---|