]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/hough/AliL3Histogram.cxx
ac5ee804e40029f0ff3d87a29363d76837c79fca
[u/mrichter/AliRoot.git] / HLT / hough / AliL3Histogram.cxx
1 //Author:        Anders Strand Vestbo
2 //Last Modified: 28.6.01
3
4 #include "AliL3Logging.h"
5 #include "AliL3Histogram.h"
6
7 //2D histogram class.
8
9 ClassImp(AliL3Histogram)
10
11 AliL3Histogram::AliL3Histogram()
12 {
13   fNxbins = 0;
14   fNybins = 0;
15   fNcells = 0;
16   fXmin = 0;
17   fYmin = 0;
18   fXmax = 0;
19   fYmax = 0;
20   fFirstXbin = 0;
21   fLastXbin = 0;
22   fFirstYbin = 0;
23   fLastYbin = 0;
24   fEntries = 0;
25   fContent = 0;
26 }
27
28   
29 AliL3Histogram::AliL3Histogram(Char_t *name,Char_t *id,Int_t nxbin,Double_t xmin,Double_t xmax,Int_t nybin,Double_t ymin,Double_t ymax) 
30 {
31   
32   strcpy(fName,name);
33   fNxbins = nxbin;
34   fNybins = nybin;
35   fNcells = (nxbin+2)*(nybin+2);
36   
37   fXmin = xmin;
38   fYmin = ymin;
39   fXmax = xmax;
40   fYmax = ymax;
41   fEntries = 0;
42   fFirstXbin = 1;
43   fFirstYbin = 1;
44   fLastXbin = nxbin;
45   fLastYbin = nybin;
46   fRootHisto = 0;
47
48   fContent = new Double_t[fNcells];
49   Reset();
50 }
51
52 AliL3Histogram::~AliL3Histogram()
53 {
54   //Destructor
55   if(fContent)
56     delete [] fContent;
57   if(fRootHisto)
58     delete fRootHisto;
59 }
60
61
62 void AliL3Histogram::Reset()
63 {
64   
65   for(Int_t i=0; i<fNcells; i++)
66     fContent[i] = 0;
67   fEntries=0;
68 }
69
70 void AliL3Histogram::Fill(Double_t x,Double_t y,Int_t weight)
71 {
72   Int_t bin = FindBin(x,y);
73   AddBinContent(bin,weight);
74
75 }
76
77 Int_t AliL3Histogram::FindBin(Double_t x,Double_t y)
78 {
79   
80   Int_t xbin = FindXbin(x);
81   Int_t ybin = FindYbin(y);
82   
83   return GetBin(xbin,ybin);
84 }
85
86 Int_t AliL3Histogram::FindXbin(Double_t x)
87 {
88   if(x < fXmin || x > fXmax)
89     return 0;
90   
91   return 1 + (Int_t)(fNxbins*(x-fXmin)/(fXmax-fXmin));
92
93 }
94
95 Int_t AliL3Histogram::FindYbin(Double_t y)
96 {
97   if(y < fYmin || y > fYmax)
98     return 0;
99   
100   return 1 + (Int_t)(fNybins*(y-fYmin)/(fYmax-fYmin));
101
102 }
103
104 Int_t AliL3Histogram::GetBin(Int_t xbin,Int_t ybin)
105 {
106   if(xbin < 0 || xbin > GetLastXbin())
107     {
108       LOG(AliL3Log::kError,"AliL3Histogram::GetBin","array")<<AliL3Log::kDec<<
109         "xbin out of range "<<xbin<<ENDLOG;
110       return 0;
111     }
112   if(ybin < 0 || ybin > GetLastYbin())
113     {
114       LOG(AliL3Log::kError,"AliL3Histogram::FindYbin","array")<<AliL3Log::kDec<<
115         "ybin out of range "<<xbin<<ENDLOG;
116       return 0;
117     }
118     
119   return xbin + ybin*(fNxbins+2);
120 }
121
122 Double_t AliL3Histogram::GetBinContent(Int_t bin)
123 {
124   if(bin >= fNcells)
125     {
126       LOG(AliL3Log::kError,"AliL3Histogram::GetBinContent","array")<<AliL3Log::kDec<<
127         "bin out of range "<<bin<<ENDLOG;
128       return 0;
129     }
130   
131   return fContent[bin];
132 }
133
134 void AliL3Histogram::SetBinContent(Int_t xbin,Int_t ybin,Int_t value)
135 {
136   Int_t bin = GetBin(xbin,ybin);
137   if(bin == 0) 
138     return;
139   SetBinContent(bin,value);
140 }
141
142 void AliL3Histogram::SetBinContent(Int_t bin,Int_t value)
143 {
144
145   if(bin >= fNcells)
146     {
147       LOG(AliL3Log::kError,"AliL3Histogram::SetBinContent","array")<<AliL3Log::kDec<<
148         "bin out of range "<<bin<<ENDLOG;
149       return;
150     }
151   if(bin == 0)
152     return;
153   fContent[bin]=value;
154   
155 }
156
157 void AliL3Histogram::AddBinContent(Int_t xbin,Int_t ybin,Int_t weight)
158 {
159   Int_t bin = GetBin(xbin,ybin);
160   if(bin == 0)
161     return;
162   AddBinContent(bin,weight);
163
164 }
165
166 void AliL3Histogram::AddBinContent(Int_t bin,Int_t weight)
167 {
168   if(bin < 0 || bin > fNcells)
169     {
170       LOG(AliL3Log::kError,"AliL3Histogram::AddBinContent","array")<<AliL3Log::kDec<<
171         "bin-value out of range "<<bin<<ENDLOG;
172       return;
173     }
174   if(bin == 0)
175     return;
176   fEntries++;
177   fContent[bin] += weight;
178 }
179
180 Double_t AliL3Histogram::GetBinCenterX(Int_t xbin)
181 {
182   
183   Double_t binwidth = (fXmax - fXmin) / fNxbins;
184   return fXmin + (xbin-1) * binwidth + 0.5*binwidth;
185   
186 }
187
188 Double_t AliL3Histogram::GetBinCenterY(Int_t ybin)
189 {
190   
191   Double_t binwidth = (fYmax - fYmin) / fNybins;
192   return fYmin + (ybin-1) * binwidth + 0.5*binwidth;
193   
194 }
195
196
197 void AliL3Histogram::Draw(Char_t *option)
198 {
199   fRootHisto = new TH2F(fName,"",fNxbins,fXmin,fXmax,fNybins,fYmin,fYmax);
200   for(Int_t bin=0; bin<fNcells; bin++)
201     {
202       fRootHisto->AddBinContent(bin,fContent[bin]);
203     }
204   //printf("ncells %d %d\n",(hist->GetNbinsX()+2)*(hist->GetNbinsY()+2),fNcells);
205   
206   fRootHisto->Draw(option);
207   
208 }