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

Last change on this file since 915 was 867, checked in by tbretz, 23 years ago
*** empty log message ***
  • Property svn:executable set to *
File size: 4.7 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 <TH2.h>
29#include <TCanvas.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("All Events", "All showers - Radius vs log(E) distribution",
52 50, 0., 5.,
53 50, 0., 500.);
54
55 fHistAll->SetXTitle("log(E/GeV)");
56 fHistAll->SetYTitle("r/m");
57 fHistAll->SetZTitle("N");
58
59 fHistSel = new TH2D("Selected Events", "Selected showers - Radius vs log(E) distribution",
60 50, 0., 5.,
61 50, 0., 500.);
62
63 fHistSel->SetXTitle("log(E/GeV)");
64 fHistSel->SetYTitle("r/m");
65 fHistSel->SetYTitle("N");
66
67 fHistCol = new TH1D("Area", "Collection Area vs. log(E)",
68 50, 0., 5.);
69
70 fHistCol->SetXTitle("log(E/GeV)");
71 fHistCol->SetYTitle("A/m^2");
72}
73
74MHMcCollectionArea::~MHMcCollectionArea()
75{
76 delete fHistAll;
77 delete fHistSel;
78 delete fHistCol;
79}
80
81void MHMcCollectionArea::FillAll(Float_t log10E, Float_t radius)
82{
83 fHistAll->Fill(log10E, radius);
84}
85
86void MHMcCollectionArea::FillSel(Float_t log10E, Float_t radius)
87{
88 fHistSel->Fill(log10E, radius);
89}
90
91void MHMcCollectionArea::DrawAll(Option_t* option)
92{
93 TCanvas *c=new TCanvas(fHistAll->GetName(), fHistAll->GetTitle());
94
95 fHistSel->Draw(option);
96
97 c->Modified();
98 c->Update();
99}
100
101void MHMcCollectionArea::DrawSel(Option_t* option)
102{
103 TCanvas *c=new TCanvas(fHistSel->GetName(), fHistSel->GetTitle());
104
105 fHistSel->Draw(option);
106
107 c->Modified();
108 c->Update();
109}
110
111void MHMcCollectionArea::Draw(Option_t* option)
112{
113 TCanvas *c=new TCanvas(fHistCol->GetName(), fHistCol->GetTitle());
114
115 fHistCol->Draw(option);
116
117 c->Modified();
118 c->Update();
119}
120
121void MHMcCollectionArea::CalcEfficiency()
122{
123 // Description!
124
125 //
126 // first of all calculate the efficency
127 // do it here by hand to get the right error of efficency
128 //
129 const Int_t nbinx = fHistSel->GetXaxis()->GetNbins();
130 const Int_t nbiny = fHistSel->GetYaxis()->GetNbins();
131
132 for (Int_t ix=1; ix<=nbinx; ix++)
133 {
134 for (Int_t iy=1; iy<=nbiny; iy++)
135 {
136 const Float_t Ns = fHistSel->GetCellContent(ix, iy);
137 const Float_t Na = fHistAll->GetCellContent(ix, iy);
138
139 if (Na <= 0)
140 continue;
141
142 const Double_t eff = Ns/Na;
143 const Double_t err = sqrt(Na + Na*Ns - Ns*Ns - Ns) / (Na*Na);
144
145 fHistSel->SetCellContent(ix, iy, eff);
146 fHistSel->SetCellError(ix, iy, err);
147 }
148 }
149
150 //
151 // now calculate the Collection area for different
152 // energies
153 //
154 for (Int_t ix=1; ix<=nbinx; ix++)
155 {
156 Double_t errA = 0;
157 Double_t colA = 0;
158
159 for (Int_t iy=1; iy<=nbiny; iy++)
160 {
161 const Double_t r1 = fHistSel->GetYaxis()->GetBinLowEdge(iy);
162 const Double_t r2 = fHistSel->GetYaxis()->GetBinLowEdge(iy+1);
163
164 const Double_t A = kPI * (r2*r2 - r1*r1);
165
166 const Double_t eff = fHistSel->GetCellContent(ix, iy);
167 const Double_t err = fHistSel->GetCellError(ix, iy);
168
169 colA += eff*A;
170 errA += A*A * err*err;
171 }
172
173 fHistCol->SetBinContent(ix, colA);
174 fHistCol->SetBinError(ix, sqrt(errA));
175 }
176}
Note: See TracBrowser for help on using the repository browser.