source: trunk/MagicSoft/Mars/mmontecarlo/MCollArea.cc@ 744

Last change on this file since 744 was 710, checked in by tbretz, 24 years ago
*** empty log message ***
File size: 3.3 KB
Line 
1#include "MCollArea.h"
2
3#include <MLog.h>
4#include <TH2.h>
5
6ClassImp(MCollArea)
7
8
9
10MCollArea::MCollArea(const char *name, const char *title)
11{
12 //
13 // default constructor
14 //
15
16 // initialize the histogram for the distribution r vs E
17 //
18 // we set the energy range from 1 Gev to 10000 GeV (in log 5 orders
19 // of magnitude) and for each order we take 10 subdivision --> 50 xbins
20 //
21 // we set the radius range from 0 m to 500 m with 10 m bin --> 50 ybins
22
23
24 *fName = name ? name : "MCollArea";
25 *fTitle = title ? title : "Data to Calculate Coll-Area";
26
27
28 fHistAll = new TH2D("collAreaAll", "all showers - Radius vs log(E) distribution",
29 50, 0., 5.,
30 50, 0., 500. ) ;
31
32 fHistSel = new TH2D("collAreaSel", "selected showers - Radius vs log(E) distribution",
33 50, 0., 5.,
34 50, 0., 500. ) ;
35
36 fHistCol = new TH1D("collArea", "Collection Area",
37 50, 0., 5.) ;
38
39}
40
41MCollArea::~MCollArea()
42{
43 delete fHistAll ;
44 delete fHistSel ;
45 delete fHistCol ;
46}
47
48void MCollArea::FillAll(Float_t log10E, Float_t radius)
49{
50 fHistAll->Fill(log10E, radius ) ;
51}
52
53void MCollArea::FillSel(Float_t log10E, Float_t radius)
54{
55 fHistSel->Fill(log10E, radius ) ;
56}
57
58void MCollArea::DrawAll()
59{
60 fHistAll->Draw() ;
61}
62
63void MCollArea::DrawSel()
64{
65 fHistSel->Draw() ;
66}
67
68void MCollArea::Draw(Option_t* option)
69{
70 fHistCol->Draw(option) ;
71}
72
73void MCollArea::CalcEfficiency()
74{
75 // Description!
76
77 //
78 // first of all calculate the efficency
79 // do it here by hand to get the right error of efficency
80 //
81 const Int_t iBinx = ((TAxis *)fHistSel->GetXaxis())->GetNbins();
82 const Int_t iBiny = ((TAxis *)fHistSel->GetYaxis())->GetNbins();
83
84 for (Int_t ix=1; ix<=iBiny; ix++ )
85 {
86 for (Int_t iy=1; iy<=iBiny; iy++ )
87 {
88 const Float_t N = fHistSel->GetCellContent(ix, iy);
89 const Float_t Nall = fHistAll->GetCellContent(ix, iy);
90
91 if ( Nall <= 0 ) {
92 // cout << ix << " " << iy << endl ;
93 continue;
94 }
95
96 const Double_t eff = N / Nall ;
97 const Double_t err = sqrt(Nall + Nall*N - N*N - N) / (Nall*Nall);
98 /*
99 cout << ix << " " << iy
100 << " N " << N
101 << " Nall " << Nall
102 << " effi " << eff
103 << " error " << err
104 << endl ;
105 */
106 fHistSel->SetCellContent(ix, iy, eff);
107 fHistSel->SetCellError(ix, iy, err);
108 }
109 }
110
111 //
112 // now calculate the Collection area for different
113 // energies
114 //
115 for (Int_t ix=1; ix<=iBiny; ix++ )
116 {
117 Double_t errA = 0;
118 Double_t colA = 0;
119
120 for (Int_t iy=1; iy<=iBiny; iy++ )
121 {
122 const Double_t r1 = ((TAxis *)fHistSel->GetYaxis())->GetBinLowEdge(iy);
123 const Double_t r2 = ((TAxis *)fHistSel->GetYaxis())->GetBinLowEdge(iy+1);
124
125 const Double_t A = kPI * (r2*r2 - r1*r1);
126
127 const Double_t eff = fHistSel->GetCellContent(ix, iy);
128 const Double_t err = fHistSel->GetCellError(ix, iy);
129
130 colA += eff*A;
131 errA += (A*A) * err * err;
132 }
133
134 fHistCol->SetBinContent(ix, colA);
135 fHistCol->SetBinError(ix, sqrt(errA));
136 }
137}
Note: See TracBrowser for help on using the repository browser.