]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/hough/AliL3Histogram1D.cxx
Compilation on Sun
[u/mrichter/AliRoot.git] / HLT / hough / AliL3Histogram1D.cxx
1 // @(#) $Id$
2
3 // Author: Anders Vestbo <mailto:vestbo@fi.uib.no>
4 //*-- Copyright &copy ALICE HLT Group
5
6 #include <strings.h>
7 #include "AliL3StandardIncludes.h"
8
9 #include "AliL3Logging.h"
10 #include "AliL3Histogram1D.h"
11
12 #if __GNUC__ == 3
13 using namespace std;
14 #endif
15
16 //_____________________________________________________________
17 // AliL3Histogram1D
18 //
19 // 1D histogram class.
20
21 ClassImp(AliL3Histogram1D)
22
23 AliL3Histogram1D::AliL3Histogram1D()
24 {
25   fNbins = 0;
26   fNcells = 0;
27   fEntries = 0;
28   fXmin = 0;
29   fXmax = 0;
30 #ifdef use_root
31   fRootHisto = 0;
32 #endif
33   fThreshold = 0;
34   fContent = 0;
35   
36 }
37   
38 AliL3Histogram1D::AliL3Histogram1D(Char_t *name,Char_t *id,Int_t nxbin,Double_t xmin,Double_t xmax)
39
40 {
41   
42   strcpy(fName,name);
43   fNbins = nxbin;
44   fNcells = fNbins + 2;
45   fEntries = 0;
46   fXmin = xmin;
47   fXmax = xmax;
48 #ifdef use_root
49   fRootHisto = 0;
50 #endif
51   fThreshold = 0;
52   
53   fContent = new Double_t[fNcells];
54   Reset();
55 }
56
57 AliL3Histogram1D::~AliL3Histogram1D()
58 {
59   //Destructor
60   if(fContent)
61     delete [] fContent;
62 #ifdef use_root
63   if(fRootHisto)
64     delete fRootHisto;
65 #endif
66 }
67
68
69 void AliL3Histogram1D::Reset()
70 {
71   bzero(fContent,fNcells*sizeof(Double_t));
72   fEntries=0;
73 }
74
75 void AliL3Histogram1D::Fill(Double_t x,Int_t weight)
76 {
77   Int_t bin = FindBin(x);
78   AddBinContent(bin,weight);
79 }
80
81
82 Int_t AliL3Histogram1D::FindBin(Double_t x)
83 {
84   if(x < fXmin || x > fXmax)
85     return 0;
86   
87   return 1 + (Int_t)(fNbins*(x-fXmin)/(fXmax-fXmin));
88
89 }
90
91 Int_t AliL3Histogram1D::GetMaximumBin()
92 {
93   Double_t max_value=0;
94   Int_t max_bin=0;
95   for(Int_t i=0; i<fNcells; i++)
96     {
97       if(fContent[i] > max_value)
98         {
99           max_value=fContent[i];
100           max_bin = i;
101         }
102     }
103   return max_bin;
104 }
105
106 Double_t AliL3Histogram1D::GetBinContent(Int_t bin)
107 {
108   if(bin >= fNcells)
109     {
110       LOG(AliL3Log::kError,"AliL3Histogram::GetBinContent","array")<<AliL3Log::kDec<<
111         "bin out of range "<<bin<<ENDLOG;
112       return 0;
113     }
114   
115   if(fContent[bin] < fThreshold)
116     return 0;
117   return fContent[bin];
118 }
119
120
121 void AliL3Histogram1D::SetBinContent(Int_t bin,Int_t value)
122 {
123
124   if(bin >= fNcells)
125     {
126       LOG(AliL3Log::kError,"AliL3Histogram::SetBinContent","array")<<AliL3Log::kDec<<
127         "bin out of range "<<bin<<ENDLOG;
128       return;
129     }
130   if(bin == 0)
131     return;
132   fContent[bin]=value;
133   
134 }
135
136 void AliL3Histogram1D::AddBinContent(Int_t bin,Int_t weight)
137 {
138   if(bin < 0 || bin > fNcells)
139     {
140       LOG(AliL3Log::kError,"AliL3Histogram::AddBinContent","array")<<AliL3Log::kDec<<
141         "bin-value out of range "<<bin<<ENDLOG;
142       return;
143     }
144   if(bin == 0)
145     return;
146   fEntries++;
147   fContent[bin] += weight;
148 }
149
150 Double_t AliL3Histogram1D::GetBinCenter(Int_t bin)
151 {
152   
153   Double_t binwidth = (fXmax - fXmin) / fNbins;
154   return fXmin + (bin-1) * binwidth + 0.5*binwidth;
155   
156 }
157
158 #ifdef use_root
159 void AliL3Histogram1D::Draw(Char_t *option)
160 {
161   fRootHisto = new TH1F(fName,"",fNbins,fXmin,fXmax);
162   for(Int_t bin=0; bin<fNcells; bin++)
163     fRootHisto->AddBinContent(bin,GetBinContent(bin));
164   
165   fRootHisto->Draw(option);
166   
167 }
168 #endif