source: trunk/MagicSoft/Mars/mhist/MHMcCT1CollectionArea.cc@ 1824

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