FEDRA emulsion software from the OPERA Collaboration
EdbCell1.h
Go to the documentation of this file.
1#ifndef ROOT_EdbCell1
2#define ROOT_EdbCell1
5
11
12#include "TObject.h"
13#include "TH1.h"
14#include "TObjArray.h"
15
17class EdbH1 : public TObject {
18
19 protected:
20
21 Int_t eN;
22 Float_t eMin;
23 Float_t eMax;
24 Float_t eBin;
25
26 Int_t eNcell;
27 Int_t *eNC;
28
29 public:
30 EdbH1();
31 EdbH1(int n, float min, float max) { Set0(); InitH1(n,min,max); }
32 EdbH1( const EdbH1 &h );
33 ~EdbH1();
34
35 void Set0();
36 void Copy( const EdbH1 &h );
37
38 int InitH1( const EdbH1 &h ) {return InitH1( h.N(), h.Xmin(), h.Xmax() ); }
39 int InitH1( int n, float min, float max );
40 void CleanCells();
41 void PrintStat();
42 void Print();
43 void Delete();
44
45 int DiscardHighCells(int nmax);
46
47 int Ncell() const {return eNcell;}
48 int N() const {return eN;}
49
50 int IX(float x) const { return (int)((x-eMin)/eBin); }
51 int Jcell(int ix) const { if(ix>=0&&ix<eN) return ix; else return -1;}
52 int Jcell(float x) const { return Jcell( IX(x) ); }
53 float X(int i) const { return eMin+eBin*(i+0.5); }
54 float Xmin() const { return eMin; }
55 float Xmax() const { return eMax; }
56 float Xbin() const { return eBin; }
57 int Bin(int ix) const { if(Jcell(ix)>-1) return eNC[Jcell(ix)]; else return 0; }
58 int MaxBin();
59 float XminA(float level=0);
60 float XmaxA(float level=0);
61
62 void AddBin(int jcell, int n) { if(jcell>=0&&jcell<eNcell) eNC[jcell]+=n; }
63 int Fill(float x) { return Fill(x,1); }
64 int Fill(float x, int n);
65 void SetBin(int ix, int n) { if(Jcell(ix)>-1) eNC[Jcell(ix)] = n; }
66
67 Long_t Integral();
68 Long_t Integral(int iv, int ir);
69
70 Float_t Mean() { return 1.*Integral()/eNcell; }
71
72 TH1I *DrawSpectrum( const char *name="EdbH1spectrun" );
73 TH1F *DrawH1( const char *name="EdbH1plot" , const char *title="EdbH1plot1D");
74
75 ClassDef(EdbH1,1) // fast 2-dim histogram class (used as a basis for EdbCell1)
76};
77
78
79/*
80//__________________________________________________________________________
81class EdbPeak1 : public EdbH1 {
82
83 public:
84 Int_t eNpeaks; // number of found peaks
85 TArrayF ePeak; //
86 TArrayF eMean3; //
87 TArrayF eMean; //
88 Float_t eNorm; // the norm-factor in case of the smoothing applied
89
90 public:
91 EdbPeak1() { InitPeaks(10); }
92 EdbPeak1( const EdbH1 &h ) : EdbH1( h ) { InitPeaks(10); }
93 ~EdbPeak1() {}
94
95 void Print();
96 void InitPeaks(int npeaks);
97 int FindPeak(int iv[2]);
98 int FindPeak(float v[2]);
99 int FindPeak(float &x, float &y);
100 float FindGlobalPeak(float &x, float &y, float ratio=0.1);
101 float ProbPeak();
102 float ProbPeak(int iv[2], int ir[2]);
103 int WipePeak(int iv[2], int ir[2]);
104 int ProbPeaks(int npeak);
105 int ProbPeaks(int npeak, int ir[2]);
106 float EstimatePeakVolume(int ipeak);
107 float EstimatePeakVolumeSafe(int ipeak);
108 float Smooth(Option_t *option="k5a");
109 float Xmean();
110 float Ymean();
111
112 ClassDef(EdbPeak1,1) // peak analyser for EdbH1
113};
114*/
115
117class EdbCell1 : public EdbH1 {
118
119 private:
120
121 Int_t eCellLim;
122 TObject **epO;
123 TObject ***epC;
124
125 public:
126 EdbCell1();
127 ~EdbCell1();
128
129 int InitCell(int maxpercell, int n, float min, float max );
130 void Delete();
131 void Reset() {CleanCells(); Delete();}
132
133 bool AddObject( float x, TObject *obj );
134 bool AddObject( int j, TObject *obj );
135
136 int SelectObjects(TObjArray &arr);
137 int SelectObjects(int min, int max, TObjArray &arr);
138 int SelectObjects(float min, float max, TObjArray &arr);
139 int SelectObjectsC(int iv, int ir, TObjArray &arr);
140 int SelectObjectsC(float v, int ir, TObjArray &arr) { return SelectObjectsC( IX(v), ir, arr); }
141
142 TObject *GetObject(int j, int ientr) const {
143 if(j>=0&&j<eNcell&&ientr>=0&&ientr<eCellLim) return epC[j][ientr];
144 else return 0;
145 }
146 void PrintStat();
147
148 ClassDef(EdbCell1,1) // class to group 2-dim objects
149};
150
151#endif /* EdbCell1 */
float min(TClonesArray *t)
Definition: bitview.cxx:275
class to group 2-dim objects
Definition: EdbCell1.h:117
int SelectObjectsC(int iv, int ir, TObjArray &arr)
Definition: EdbCell1.cpp:538
EdbCell1()
Definition: EdbCell1.cpp:476
~EdbCell1()
Definition: EdbCell1.cpp:484
void PrintStat()
Definition: EdbCell1.cpp:490
int InitCell(int maxpercell, int n, float min, float max)
Definition: EdbCell1.cpp:504
int SelectObjects(TObjArray &arr)
Definition: EdbCell1.cpp:566
TObject ** epO
pointers to objects [eNcell*eCellLim]
Definition: EdbCell1.h:122
TObject * GetObject(int j, int ientr) const
Definition: EdbCell1.h:142
Int_t eCellLim
max number of entries into one cell (for memory allocation)
Definition: EdbCell1.h:121
void Reset()
Definition: EdbCell1.h:131
int SelectObjectsC(float v, int ir, TObjArray &arr)
Definition: EdbCell1.h:140
void Delete()
Definition: EdbCell1.cpp:497
bool AddObject(float x, TObject *obj)
Definition: EdbCell1.cpp:519
TObject *** epC
pointers to cells [eNcell]
Definition: EdbCell1.h:123
fast 2-dim histogram class (used as a basis for EdbCell1)
Definition: EdbCell1.h:17
float X(int i) const
Definition: EdbCell1.h:53
Int_t eNcell
eN
Definition: EdbCell1.h:26
float Xbin() const
Definition: EdbCell1.h:56
EdbH1(int n, float min, float max)
Definition: EdbCell1.h:31
int Bin(int ix) const
Definition: EdbCell1.h:57
int N() const
Definition: EdbCell1.h:48
float XminA(float level=0)
Definition: EdbCell1.cpp:109
void AddBin(int jcell, int n)
Definition: EdbCell1.h:62
Int_t * eNC
[eNcell] number of objects/cell
Definition: EdbCell1.h:27
void Set0()
Definition: EdbCell1.cpp:28
int Jcell(float x) const
Definition: EdbCell1.h:52
int Jcell(int ix) const
Definition: EdbCell1.h:51
void CleanCells()
Definition: EdbCell1.cpp:85
void SetBin(int ix, int n)
Definition: EdbCell1.h:65
int Fill(float x)
Definition: EdbCell1.h:63
int Ncell() const
Definition: EdbCell1.h:47
void Copy(const EdbH1 &h)
Definition: EdbCell1.cpp:43
Float_t eBin
bin size
Definition: EdbCell1.h:24
float XmaxA(float level=0)
Definition: EdbCell1.cpp:118
int MaxBin()
Definition: EdbCell1.cpp:187
float Xmax() const
Definition: EdbCell1.h:55
void Print()
Definition: EdbCell1.cpp:161
TH1F * DrawH1(const char *name="EdbH1plot", const char *title="EdbH1plot1D")
Definition: EdbCell1.cpp:136
int IX(float x) const
Definition: EdbCell1.h:50
void PrintStat()
Definition: EdbCell1.cpp:146
Float_t eMax
max
Definition: EdbCell1.h:23
TH1I * DrawSpectrum(const char *name="EdbH1spectrun")
Definition: EdbCell1.cpp:127
int InitH1(const EdbH1 &h)
Definition: EdbCell1.h:38
void Delete()
Definition: EdbCell1.cpp:64
EdbH1()
Definition: EdbCell1.cpp:23
Float_t eMin
min
Definition: EdbCell1.h:22
~EdbH1()
Definition: EdbCell1.cpp:58
Int_t eN
divisions
Definition: EdbCell1.h:21
float Xmin() const
Definition: EdbCell1.h:54
Float_t Mean()
Definition: EdbCell1.h:70
Long_t Integral()
Definition: EdbCell1.cpp:170
int DiscardHighCells(int nmax)
Definition: EdbCell1.cpp:100
int max
Definition: check_shower.C:41
const char * name
Definition: merge_Energy_SytematicSources_Electron.C:24