]> git.uio.no Git - u/mrichter/AliRoot.git/blob - HLT/hough/AliL3HoughMaxFinder.h
Merged Bergen, mergen Cvetan TransformerRow and
[u/mrichter/AliRoot.git] / HLT / hough / AliL3HoughMaxFinder.h
1 // @(#) $Id$
2
3 #ifndef ALIL3_HOUGH_MaxFinder
4 #define ALIL3_HOUGH_MaxFinder
5
6 #include "AliL3RootTypes.h"
7 #include "AliL3StandardIncludes.h"
8
9 class AliL3Histogram;
10 class AliL3TrackArray;
11 class AliL3HoughTrack;
12 class TNtuple;
13
14 struct AxisWindow
15 {
16   Int_t ymin;
17   Int_t ymax;
18   Int_t xbin;
19   Int_t weight;
20 };
21
22 class AliL3HoughMaxFinder {
23   
24  private:
25
26   Int_t fThreshold;
27   Int_t fCurrentEtaSlice;
28   AliL3Histogram *fCurrentHisto;  //!
29   
30   Float_t fGradX;
31   Float_t fGradY;
32   Float_t *fXPeaks; //!
33   Float_t *fYPeaks; //!
34   Int_t *fSTARTXPeaks; //!
35   Int_t *fSTARTYPeaks; //!
36   Int_t *fENDXPeaks; //!
37   Int_t *fENDYPeaks; //!
38   Int_t *fSTARTETAPeaks; //!
39   Int_t *fENDETAPeaks; //!
40   Int_t *fWeight;   //!
41   Int_t fN1PeaksPrevEtaSlice;
42   Int_t fN2PeaksPrevEtaSlice;
43   Int_t fNPeaks;
44   Int_t fNMax;
45   
46   Char_t fHistoType;
47
48 #ifndef no_root
49   TNtuple *fNtuppel; //!
50 #endif
51
52  public:
53   AliL3HoughMaxFinder(); 
54   AliL3HoughMaxFinder(Char_t *histotype,Int_t nmax,AliL3Histogram *hist=0);
55   virtual ~AliL3HoughMaxFinder();
56   void Reset();
57
58   void CreateNtuppel();
59   void WriteNtuppel(Char_t *filename);
60
61   //Simple maxima finders:
62   void FindAbsMaxima();
63   void FindBigMaxima();
64   void FindMaxima(Int_t threshold=0);
65   void FindAdaptedPeaks(Int_t nkappawindow,Float_t cut_ratio);
66   //Peak finder for HoughTransformerRow
67   void FindAdaptedRowPeaks(Int_t kappawindow,Int_t xsize,Int_t ysize);
68   //More sophisticated peak finders:
69   void FindPeak(Int_t t1,Double_t t2,Int_t t3);
70   void FindPeak1(Int_t y_window=2,Int_t x_bin_sides=1);
71   void SortPeaks(struct AxisWindow **a,Int_t first,Int_t last);
72   Int_t PeakCompare(struct AxisWindow *a,struct AxisWindow *b);
73   
74   //Setters:
75   void SetGradient(Float_t x,Float_t y) {fGradX=x; fGradY=y;}
76   void SetThreshold(Int_t f) {fThreshold = f;}
77   void SetHistogram(AliL3Histogram *hist) {fCurrentHisto = hist;}
78   void SetEtaSlice(Int_t etaslice) {fCurrentEtaSlice = etaslice;}
79   
80   //Getters:
81   Float_t GetXPeak(Int_t i);
82   Float_t GetYPeak(Int_t i);
83   Float_t GetXPeakSize(Int_t i);
84   Float_t GetYPeakSize(Int_t i);
85   Int_t GetWeight(Int_t i);
86   Int_t GetStartEta(Int_t i);
87   Int_t GetEndEta(Int_t i);
88   Int_t GetEntries() {return fNPeaks;}
89
90   ClassDef(AliL3HoughMaxFinder,1) //Maximum finder class
91
92 };
93
94 inline Float_t AliL3HoughMaxFinder::GetXPeak(Int_t i)
95 {
96   if(i<0 || i>fNMax)
97     {
98       STDCERR<<"AliL3HoughMaxFinder::GetXPeak : Invalid index "<<i<<STDENDL;
99       return 0;
100     }
101   return fXPeaks[i];
102 }
103
104 inline Float_t AliL3HoughMaxFinder::GetYPeak(Int_t i)
105 {
106   if(i<0 || i>fNMax)
107     {
108       STDCERR<<"AliL3HoughMaxFinder::GetYPeak : Invalid index "<<i<<STDENDL;
109       return 0;
110     }
111   return fYPeaks[i];
112
113 }
114
115 inline Int_t AliL3HoughMaxFinder::GetWeight(Int_t i)
116 {
117   if(i<0 || i>fNMax)
118     {
119       STDCERR<<"AliL3HoughMaxFinder::GetWeight : Invalid index "<<i<<STDENDL;
120       return 0;
121     }
122   return fWeight[i];
123 }
124
125 inline Int_t AliL3HoughMaxFinder::GetStartEta(Int_t i)
126 {
127   if(i<0 || i>fNMax)
128     {
129       STDCERR<<"AliL3HoughMaxFinder::GetStartEta : Invalid index "<<i<<STDENDL;
130       return 0;
131     }
132   return fSTARTETAPeaks[i];
133 }
134
135 inline Int_t AliL3HoughMaxFinder::GetEndEta(Int_t i)
136 {
137   if(i<0 || i>fNMax)
138     {
139       STDCERR<<"AliL3HoughMaxFinder::GetStartEta : Invalid index "<<i<<STDENDL;
140       return 0;
141     }
142   return fENDETAPeaks[i];
143 }
144
145 #endif
146