source: trunk/MagicSoft/Mars/mhflux/MHCollectionArea.cc@ 9375

Last change on this file since 9375 was 9375, checked in by tbretz, 16 years ago
*** empty log message ***
File size: 13.2 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-2009
22!
23!
24\* ======================================================================== */
25
26//////////////////////////////////////////////////////////////////////////////
27//
28// MHCollectionArea
29//
30// Class Version 2:
31// ----------------
32// + added //! to fMcEvt which was missing before
33// - removed obsolete data member fEnergy
34//
35//////////////////////////////////////////////////////////////////////////////
36#include "MHCollectionArea.h"
37
38#include <TMath.h>
39
40#include <TLatex.h>
41#include <TCanvas.h>
42#include <TPaveStats.h>
43
44#include "MLog.h"
45#include "MLogManip.h"
46
47#include "MString.h"
48#include "MBinning.h"
49
50#include "MMcEvt.hxx"
51#include "MMcRunHeader.hxx"
52#include "MMcCorsikaRunHeader.h"
53
54#include "MParList.h"
55#include "MParameters.h"
56
57ClassImp(MHCollectionArea);
58
59using namespace std;
60
61// --------------------------------------------------------------------------
62//
63// Creates the three necessary histograms:
64// - selected showers (input)
65// - all showers (input)
66// - collection area (result)
67//
68MHCollectionArea::MHCollectionArea(const char *name, const char *title)
69 : fMcEvt(0), fHeader(0), fMcAreaRadius(-1), fIsExtern(kFALSE)
70{
71 // initialize the histogram for the distribution r vs E
72 //
73 // we set the energy range from 2 Gev to 20000 GeV (in log 4 orders
74 // of magnitude) and for each order we take 25 subdivision --> 100 xbins
75 //
76 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
77 //
78 fName = name ? name : "MHCollectionArea";
79 fTitle = title ? title : "Collection Area vs. Energy/Theta";
80
81 fHistSel.SetName("SelEvts");
82 fHistSel.SetTitle("Number of Events after cuts");
83 fHistSel.SetXTitle("\\Theta [deg]");
84 fHistSel.SetYTitle("E [GeV]");
85 fHistSel.SetDirectory(NULL);
86 fHistSel.UseCurrentStyle();
87 fHistSel.SetLineColor(kBlue);
88
89 fHistAll.SetName("AllEvts");
90 fHistAll.SetTitle("Number of events produced");
91 fHistAll.SetXTitle("\\Theta [deg]");
92 fHistAll.SetYTitle("E_{mc} [GeV]");
93 fHistAll.SetDirectory(NULL);
94 fHistAll.UseCurrentStyle();
95
96 fHEnergy.SetName("CollEnergy");
97 fHEnergy.SetTitle("Collection Area");
98 fHEnergy.SetXTitle("E [GeV]");
99 fHEnergy.SetYTitle("A_{eff} [m^{2}]");
100 fHEnergy.SetDirectory(NULL);
101 fHEnergy.UseCurrentStyle();
102
103 MBinning binsa, binse, binst;
104 binse.SetEdgesLog(21, 6.3, 100000);
105 binst.SetEdgesASin(67, -0.005, 0.665);
106
107 binse.Apply(fHEnergy);
108
109 MH::SetBinning(&fHistSel, &binst, &binse);
110 MH::SetBinning(&fHistAll, &binst, &binse);
111
112 // For some unknown reasons this must be called after
113 // the binning has been initialized at least once
114 fHistSel.Sumw2();
115 fHistAll.Sumw2();
116 fHEnergy.Sumw2();
117}
118
119// --------------------------------------------------------------------------
120//
121// Return the Area defined by fMcAreaRadius
122//
123Double_t MHCollectionArea::GetCollectionAreaAbs() const
124{
125 return TMath::Pi()*fMcAreaRadius*fMcAreaRadius;
126}
127
128// --------------------------------------------------------------------------
129//
130// Calculate the Efficiency (collection area) and set the 'ReadyToSave'
131// flag
132//
133void MHCollectionArea::CalcEfficiency()
134{
135 TH1D *hsel = fHistSel.ProjectionY("Spy", -1, -1, "E");;
136 TH1D *hall = fHistAll.ProjectionY("Apy", -1, -1, "E");
137
138 //
139 // Impact parameter range.
140 //
141 const Float_t totalarea = GetCollectionAreaAbs();//TMath::Pi() * (r2*r2 - r1*r1);
142
143 // "b" option: calculate binomial errors
144 // Do not use totalarea inside the binomial error calculation:
145 // it is not a weight.
146 fHEnergy.Divide(hsel, hall, 1, 1, "b");
147#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
148 MH::SetBinomialErrors(fHEnergy, *hsel, *hall);
149#endif
150 if (fMcAreaRadius>0)
151 fHEnergy.Scale(totalarea);
152
153 delete hsel;
154 delete hall;
155}
156
157
158Bool_t MHCollectionArea::SetupFill(const MParList *pl)
159{
160 fHistSel.Reset();
161 if (!fIsExtern)
162 fHistAll.Reset();
163
164 fHeader = (MMcRunHeader*)pl->FindObject("MMcRunHeader");
165 if (!fHeader)
166 {
167 *fLog << err << "MMcRunHeader not found... abort." << endl;
168 return kFALSE;
169 }
170
171 fMcEvt = (MMcEvt*)pl->FindObject("MMcEvt");
172 if (!fMcEvt)
173 {
174 *fLog << err << "MMcEvt not found... abort." << endl;
175 return kFALSE;
176 }
177 /*
178 fEnergy = (MParameterD*)pl->FindObject("MEnergyEst", "MParameterD");
179 if (!fEnergy)
180 {
181 *fLog << err << "MEnergyEst [MParameterD] not found... abort." << endl;
182 return kFALSE;
183 }
184 */
185 MBinning binst, binse;
186 binst.SetEdges(fHistAll, 'x');
187 binse.SetEdges(fHistAll, 'y');
188
189 //if (!fIsExtern)
190 {
191 MBinning *bins = (MBinning*)pl->FindObject("BinningTheta", "MBinning");
192 if (bins)
193 binst.SetEdges(*bins);
194
195 bins = (MBinning*)pl->FindObject("BinningEnergyEst", "MBinning");
196 if (bins)
197 binse.SetEdges(*bins);
198 }
199
200 binse.Apply(fHEnergy);
201
202 MH::SetBinning(&fHistSel, &binst, &binse);
203 MH::SetBinning(&fHistAll, &binst, &binse);
204
205 fMcAreaRadius = -1;
206 fCorsikaVersion = 0;
207
208 return kTRUE;
209}
210
211void MHCollectionArea::GetImpactMax()
212{
213 if (fHeader->GetImpactMax()<=fMcAreaRadius*100 || fHeader->GetImpactMax()<0)
214 return;
215
216 fMcAreaRadius = 0.01*fHeader->GetImpactMax(); // cm->m
217 *fLog << inf << "Maximum simulated impact: " << fMcAreaRadius << "m" << endl;
218}
219
220Bool_t MHCollectionArea::ReInit(MParList *plist)
221{
222 GetImpactMax();
223
224 if (fCorsikaVersion!=0 && fCorsikaVersion!=fHeader->GetCorsikaVersion())
225 {
226 *fLog << warn;
227 *fLog << "Warning - Read files have different Corsika versions..." << endl;
228 *fLog << " Last file=" << fCorsikaVersion << " New file=" << fHeader->GetCorsikaVersion() << endl;
229 }
230 fCorsikaVersion = fHeader->GetCorsikaVersion();
231
232 if (fIsExtern)
233 return kTRUE;
234
235 fTotalNumSimulatedShowers += fHeader->GetNumSimulatedShowers();
236 *fLog << inf << "Total Number of Simulated showers: " << fTotalNumSimulatedShowers << endl;
237
238 fAllEvtsTriggered |= fHeader->GetAllEvtsTriggered();
239 *fLog << inf << "Only triggered events avail: " << (fAllEvtsTriggered?"yes":"no") << endl;
240
241 MMcCorsikaRunHeader *crh = (MMcCorsikaRunHeader*)plist->FindObject("MMcCorsikaRunHeader");
242 if (!crh)
243 {
244 *fLog << err << "MMcCorsikaRunHeader not found... abort." << endl;
245 return kFALSE;
246 }
247
248 //
249 // Calculate approximately the original number of events in each
250 // energy bin:
251 //
252 const Float_t emin = crh->GetELowLim();
253 const Float_t emax = crh->GetEUppLim();
254 const Float_t expo = 1 + crh->GetSlopeSpec();
255 const Float_t k = fHeader->GetNumSimulatedShowers() /
256 (pow(emax,expo) - pow(emin,expo));
257
258 const Int_t nbiny = fHistAll.GetNbinsY();
259
260 TAxis &axe = *fHistAll.GetYaxis();
261 for (Int_t i = 1; i <= nbiny; i++)
262 {
263 const Float_t e1 = axe.GetBinLowEdge(i);
264 const Float_t e2 = axe.GetBinLowEdge(i+1);
265
266 if (e1 < emin || e2 > emax)
267 continue;
268
269 const Float_t events = k * (pow(e2, expo) - pow(e1, expo));
270 //
271 // We fill the i-th energy bin, with the total number of events
272 // Second argument of Fill would be impact parameter of each
273 // event, but we don't really need it for the collection area,
274 // so we just put a dummy value (1.)
275 //
276
277 const Float_t energy = (e1+e2)/2.;
278 for (int j=0; j<TMath::Nint(events); j++)
279 fHistAll.Fill(0., energy);
280 // you have MMcRunHeader.fShowerThetaMin and MMcRunHeader.fShowerThetaMax
281 }
282
283 return kTRUE;
284}
285
286void MHCollectionArea::Paint(Option_t *option)
287{
288 // This is a workaround to support also older files
289 if (TString(option)=="paint3")
290 return;
291
292 if (TString(option)=="paint4")
293 {
294 if (fMcAreaRadius<=0)
295 return;
296
297 const TString txt = MString::Format("r_{max}=%.0fm --> A_{max}=%.0fm^{2}",
298 fMcAreaRadius, GetCollectionAreaAbs());
299
300 TLatex text(0.31, 0.95, txt);
301 text.SetBit(TLatex::kTextNDC);
302 text.SetTextSize(0.04);
303 text.Paint();
304 return;
305 }
306
307 TVirtualPad *pad;
308
309 TPaveStats *st=0;
310 for (int x=0; x<4; x++)
311 {
312 pad=gPad->GetPad(x+1);
313 if (!pad || !(st = (TPaveStats*)pad->GetPrimitive("stats")))
314 continue;
315
316 if (st->GetOptStat()==11)
317 continue;
318
319 const Double_t y1 = st->GetY1NDC();
320 const Double_t y2 = st->GetY2NDC();
321 const Double_t x1 = st->GetX1NDC();
322 const Double_t x2 = st->GetX2NDC();
323
324 st->SetY1NDC((y2-y1)/3+y1);
325 st->SetX1NDC((x2-x1)/3+x1);
326 st->SetOptStat(11);
327 }
328
329 pad = gPad;
330
331 TH1 *h1=0, *h2=0;
332
333 pad->cd(1);
334 if (gPad->FindObject("ProjSelX"))
335 fHistSel.ProjectionX("ProjSelX", -1, -1, "E");
336
337 pad->cd(2);
338 if (gPad->FindObject("ProjAllY"))
339 h1=fHistAll.ProjectionY("ProjAllY", -1, -1, "E");
340 if (gPad->FindObject("ProjSelY"))
341 h2=fHistSel.ProjectionY("ProjSelY", -1, -1, "E");
342
343 if (h1 && h1->GetMaximum()>0)
344 {
345 gPad->SetLogx();
346 gPad->SetLogy();
347 }
348
349 pad->cd(3);
350 TH1 *h=dynamic_cast<TH1*>(gPad->FindObject("Efficiency"));
351 if (h1 && h2 && h)
352 {
353 h->Divide(h2, h1, 1, 1, "b");
354#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
355 MH::SetBinomialErrors(*h, *h2, *h1);
356#endif
357 h->SetMinimum(0);
358 }
359
360 pad->cd(4);
361 CalcEfficiency();
362 if (fHEnergy.GetMaximum()>0)
363 {
364 gPad->SetLogx();
365 gPad->SetLogy();
366 }
367}
368
369void MHCollectionArea::Draw(Option_t *option)
370{
371 TVirtualPad *pad = gPad ? gPad : MakeDefCanvas(this);
372
373 // Do the projection before painting the histograms into
374 // the individual pads
375 AppendPad();
376
377 pad->SetBorderMode(0);
378 pad->Divide(2,2);
379
380 TH1 *h=0, *h1=0, *h2=0;
381
382 if (fHistSel.GetNbinsX()>1)
383 {
384 pad->cd(1);
385 gPad->SetBorderMode(0);
386 gPad->SetGridx();
387 gPad->SetGridy();
388 /*
389 h = fHistAll.ProjectionX("ProjAllX", -1, -1, "E");
390 h->SetXTitle("\\Theta [\\circ]");
391 h->SetDirectory(NULL);
392 h->SetLineColor(kGreen);
393 h->SetBit(kCanDelete);
394 h->Draw();
395 */
396 h = fHistSel.ProjectionX("ProjSelX", -1, -1, "E");
397 h->SetXTitle("\\Theta [\\circ]");
398 h->SetDirectory(NULL);
399 h->SetLineColor(kRed);
400 h->SetBit(kCanDelete);
401 h->Draw("hist"/*"same"*/);
402 }
403 else
404 delete pad->GetPad(1);
405
406 if (fHistSel.GetNbinsY()>1)
407 {
408 pad->cd(2);
409 gPad->SetBorderMode(0);
410 gPad->SetGridx();
411 gPad->SetGridy();
412
413 h1 = fHistAll.ProjectionY("ProjAllY", -1, -1, "E");
414 h1->SetDirectory(NULL);
415 h1->SetLineColor(kGreen);
416 h1->SetXTitle("E [GeV]");
417 h1->SetBit(kCanDelete);
418 h1->Draw();
419
420 h2 = fHistSel.ProjectionY("ProjSelY", -1, -1, "E");
421 h2->SetDirectory(NULL);
422 h2->SetLineColor(kRed);
423 h2->SetBit(kCanDelete);
424 h2->Draw("same");
425 }
426 else
427 delete pad->GetPad(2);
428
429 if (h1 && h2)
430 {
431 pad->cd(3);
432 gPad->SetBorderMode(0);
433 gPad->SetGridx();
434 gPad->SetGridy();
435 gPad->SetLogx();
436 h = h2->DrawCopy();
437 h->Divide(h2, h1, 1, 1, "b");
438#if ROOT_VERSION_CODE < ROOT_VERSION(5,13,04)
439 MH::SetBinomialErrors(*h, *h2, *h1);
440#endif
441 h->SetNameTitle("Efficiency", "Efficiency");
442 h->SetDirectory(NULL);
443 //AppendPad("paint3");
444 }
445 else
446 delete pad->GetPad(4);
447
448 if (fHEnergy.GetNbinsX()>1)
449 {
450 pad->cd(4);
451 gPad->SetBorderMode(0);
452 gPad->SetGridx();
453 gPad->SetGridy();
454 fHEnergy.Draw();
455 AppendPad("paint4");
456 }
457 else
458 delete pad->GetPad(4);
459}
460
461Int_t MHCollectionArea::Fill(const MParContainer *par, const Stat_t weight)
462{
463 // This is not perfect because it selects the maximum impact only
464 // from the selected events. Hoever, it will get overwritten
465 // in finalize anyway.
466 const Double_t impact = fMcEvt->GetImpact()*0.01; // cm->m
467 if (impact>fMcAreaRadius)
468 fMcAreaRadius = impact;
469
470 const Double_t energy = fMcEvt->GetEnergy();
471 const Double_t theta = fMcEvt->GetTelescopeTheta()*TMath::RadToDeg();
472
473 fHistSel.Fill(theta, energy, weight);
474
475 return kTRUE;
476}
477
478Bool_t MHCollectionArea::Finalize()
479{
480 GetImpactMax();
481
482 *fLog << all << "Maximum simulated impact found: " << fMcAreaRadius << "m" << endl;
483
484 CalcEfficiency();
485
486 return kTRUE;
487}
Note: See TracBrowser for help on using the repository browser.