source: trunk/MagicSoft/Mars/mhist/MHMcCollectionArea.cc@ 859

Last change on this file since 859 was 859, checked in by tbretz, 23 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 4.1 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): Harald Kornmayer 1/2001 (harald@mppmu.mpg.de)
19! Thomas Bretz 12/2000 (tbretz@uni-sw.gwdg.de)
20!
21! Copyright: MAGIC Software Development, 2000-2001
22!
23!
24\* ======================================================================== */
25
26#include "MHMcCollectionArea.h"
27
28#include <MLog.h>
29#include <TH2.h>
30
31ClassImp(MHMcCollectionArea);
32
33MHMcCollectionArea::MHMcCollectionArea(const char *name, const char *title)
34{
35 //
36 // default constructor
37 //
38
39 // initialize the histogram for the distribution r vs E
40 //
41 // we set the energy range from 1 Gev to 10000 GeV (in log 5 orders
42 // of magnitude) and for each order we take 10 subdivision --> 50 xbins
43 //
44 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
45
46
47 *fName = name ? name : "MHMcCollectionArea";
48 *fTitle = title ? title : "Data to Calculate Coll-Area";
49
50
51 fHistAll = new TH2D("collAreaAll", "all showers - Radius vs log(E) distribution",
52 50, 0., 5.,
53 50, 0., 500.);
54
55 fHistSel = new TH2D("collAreaSel", "selected showers - Radius vs log(E) distribution",
56 50, 0., 5.,
57 50, 0., 500.);
58
59 fHistCol = new TH1D("collArea", "Collection Area",
60 50, 0., 5.);
61}
62
63MHMcCollectionArea::~MHMcCollectionArea()
64{
65 delete fHistAll;
66 delete fHistSel;
67 delete fHistCol;
68}
69
70void MHMcCollectionArea::FillAll(Float_t log10E, Float_t radius)
71{
72 fHistAll->Fill(log10E, radius);
73}
74
75void MHMcCollectionArea::FillSel(Float_t log10E, Float_t radius)
76{
77 fHistSel->Fill(log10E, radius);
78}
79
80void MHMcCollectionArea::DrawAll()
81{
82 fHistAll->Draw();
83}
84
85void MHMcCollectionArea::DrawSel()
86{
87 fHistSel->Draw();
88}
89
90void MHMcCollectionArea::Draw(Option_t* option)
91{
92 fHistCol->Draw(option);
93}
94
95void MHMcCollectionArea::CalcEfficiency()
96{
97 // Description!
98
99 //
100 // first of all calculate the efficency
101 // do it here by hand to get the right error of efficency
102 //
103 const Int_t nbinx = fHistSel->GetXaxis()->GetNbins();
104 const Int_t nbiny = fHistSel->GetYaxis()->GetNbins();
105
106 for (Int_t ix=1; ix<=nbinx; ix++)
107 {
108 for (Int_t iy=1; iy<=nbiny; iy++)
109 {
110 const Float_t Ns = fHistSel->GetCellContent(ix, iy);
111 const Float_t Na = fHistAll->GetCellContent(ix, iy);
112
113 if (Na <= 0)
114 continue;
115
116 const Double_t eff = Ns/Na;
117 const Double_t err = sqrt(Na + Na*Ns - Ns*Ns - Ns) / (Na*Na);
118
119 fHistSel->SetCellContent(ix, iy, eff);
120 fHistSel->SetCellError(ix, iy, err);
121 }
122 }
123
124 //
125 // now calculate the Collection area for different
126 // energies
127 //
128 for (Int_t ix=1; ix<=nbinx; ix++)
129 {
130 Double_t errA = 0;
131 Double_t colA = 0;
132
133 for (Int_t iy=1; iy<=nbiny; iy++)
134 {
135 const Double_t r1 = fHistSel->GetYaxis()->GetBinLowEdge(iy);
136 const Double_t r2 = fHistSel->GetYaxis()->GetBinLowEdge(iy+1);
137
138 const Double_t A = kPI * (r2*r2 - r1*r1);
139
140 const Double_t eff = fHistSel->GetCellContent(ix, iy);
141 const Double_t err = fHistSel->GetCellError(ix, iy);
142
143 colA += eff*A;
144 errA += A*A * err*err;
145 }
146
147 fHistCol->SetBinContent(ix, colA);
148 fHistCol->SetBinError(ix, sqrt(errA));
149 }
150}
Note: See TracBrowser for help on using the repository browser.