FEDRA emulsion software from the OPERA Collaboration
EdbPatCell2 Class Reference

helper class to analyse the single pattern More...

#include <EdbPatCell2.h>

Inheritance diagram for EdbPatCell2:
Collaboration diagram for EdbPatCell2:

Public Member Functions

void AcceptCorrections (EdbLayer &layer)
 
void ApplyCorrections (EdbLayer &layer)
 
void ApplyCorrections (EdbPattern &pat)
 
 EdbPatCell2 ()
 
int FillCell (EdbPattern &pat)
 limits should be already defined More...
 
int FillCell (TObjArray &pat)
 limits should be already defined More...
 
int FillCombinations (EdbPatCell2 &cell, int ir2[2], TObjArray &arrC1, TObjArray &arrC2, bool doFill=1)
 
EdbH2FillSelectedH2 ()
 
void InitPat (EdbPattern &pat, int nx, int ny)
 define limits to fit the pattern More...
 
void InitPatBin (EdbPattern &pat, float binx, float biny)
 define limits to fit the pattern More...
 
void PrintCorr ()
 
void ResetCorr ()
 
float TXs (EdbSegP &s)
 
float TYs (EdbSegP &s)
 
float Xs (EdbSegP &s)
 
float XsPos (EdbSegP &s)
 
float XsRot (EdbSegP &s)
 
float Ys (EdbSegP &s)
 
float YsPos (EdbSegP &s)
 
float YsRot (EdbSegP &s)
 
virtual ~EdbPatCell2 ()
 
- Public Member Functions inherited from EdbCell2
bool AddObject (float v[2], TObject *obj)
 
bool AddObject (float x, float y, TObject *obj)
 
bool AddObject (int ix, int iy, TObject *obj)
 
bool AddObject (int j, TObject *obj)
 
int CellLim ()
 
void Copy (const EdbCell2 &cell)
 
void Delete ()
 
 EdbCell2 ()
 
 EdbCell2 (const EdbCell2 &cell)
 
TObject * GetObject (float x, float y, int ientr) const
 
TObject * GetObject (int ix, int iy, int ientr) const
 
TObject * GetObject (int j, int ientr) const
 
int InitCell (EdbCell2 &c)
 
int InitCell (int maxpercell, int n[2], float min[2], float max[2])
 
int InitCell (int nx, float minx, float maxx, int ny, float miny, float maxy, int maxpercell)
 
void InitEPC ()
 
void PrintStat ()
 
void Reset ()
 
int SelectObjects (float min[2], float max[2], TObjArray &arr)
 
int SelectObjects (int min[2], int max[2], TObjArray &arr)
 
int SelectObjects (TObjArray &arr)
 
int SelectObjectsC (float v[2], float r, TObjArray &arr)
 
int SelectObjectsC (float v[2], int ir[2], TObjArray &arr)
 
int SelectObjectsC (int iv[2], int ir[2], TObjArray &arr)
 
int SelectObjectsCJ (int j, int ir, TObjArray &arr)
 
void Set0 ()
 
 ~EdbCell2 ()
 
- Public Member Functions inherited from EdbH2
void AddBin (int jcell, int n)
 
int Bin (float x, float y) const
 
int Bin (int iv[2]) const
 
int Bin (int ix, int iy) const
 
int Bin (int j) const
 
void CleanCells ()
 
void Copy (const EdbH2 &h)
 
void Delete ()
 
int DiscardHighCells (int nmax)
 
TH2F * DrawH2 (const char *name="plot2d", const char *title="EdbH2plot2D")
 
TH1I * DrawSpectrum (const char *name="plot1d", const char *title="EdbH2 DrawSpectrun")
 
 EdbH2 ()
 
 EdbH2 (const EdbH2 &h)
 
 EdbH2 (int nx, float minx, float maxx, int ny, float miny, float maxy)
 
int Fill (float x, float y)
 
int Fill (float x, float y, int n)
 
int InitH2 (const EdbH2 &h)
 
int InitH2 (int n[2], float min[2], float max[2])
 
int InitH2 (int nx, float minx, float maxx, int ny, float miny, float maxy)
 
Long_t Integral ()
 
Long_t Integral (int iv[2], int ir[2])
 
int IX (float x) const
 
int IX (int jcell) const
 
int IY (float y) const
 
int IY (int jcell) const
 
int Jcell (float v[2]) const
 
int Jcell (float x, float y) const
 
int Jcell (int ix, int iy) const
 
int MaxBin ()
 
Float_t Mean ()
 
int Ncell () const
 
int NX () const
 
int NY () const
 
void PrintStat ()
 
EdbH1ProjectionX ()
 
EdbH1ProjectionY ()
 
void Set0 ()
 
void SetBin (int ix, int iy, int n)
 
void SetBin (int j, int n)
 
float X (int ix) const
 
float Xbin () const
 
float Xj (int j) const
 
float Xmax () const
 
float XmaxA (float level=0)
 
float Xmin () const
 
float XminA (float level=0)
 
float Y (int iy) const
 
float Ybin () const
 
float Yj (int j) const
 
float Ymax () const
 
float YmaxA (float level=0)
 
float Ymin () const
 
float YminA (float level=0)
 
 ~EdbH2 ()
 

Public Attributes

bool eApplyCorr
 
TH1I * eDoubletsRate
 to be filled in FillCombinations() More...
 
float eDTX
 
float eDTXlim
 
float eDTY
 corrections to be applied if eApplyCorr==true More...
 
float eDTYlim
 
float eDX
 
float eDXlim
 
float eDY
 
float eDYlim
 acceptance limits for the combinations selection More...
 
float eDZ
 corrections to be applied if eApplyCorr==true More...
 
float ePhi
 rotation angle (using dx,dy one can set-up the center of rotation) More...
 
float eShr
 corrections to be applied if eApplyCorr==true More...
 
float eXmarg
 
float eYmarg
 margins for the cell definition More...
 

Additional Inherited Members

- Protected Attributes inherited from EdbH2
Float_t eBin [2]
 bin size More...
 
Float_t eMax [2]
 max More...
 
Float_t eMin [2]
 min More...
 
Int_t eN [2]
 divisions More...
 
Int_t * eNC
 [eNcell] number of objects/cell More...
 
Int_t eNcell
 eNx*eNy More...
 

Detailed Description

helper class to analyse the single pattern

Constructor & Destructor Documentation

◆ EdbPatCell2()

EdbPatCell2::EdbPatCell2 ( )
13{
14 eDXlim=3; eDYlim=3; eDTXlim=0.003; eDTYlim=0.003;
15 eXmarg = eYmarg = 0.0001; // margins for the cell definition
17 ResetCorr();
18}
float eDYlim
acceptance limits for the combinations selection
Definition: EdbPatCell2.h:12
float eDXlim
Definition: EdbPatCell2.h:12
TH1I * eDoubletsRate
to be filled in FillCombinations()
Definition: EdbPatCell2.h:23
void ResetCorr()
Definition: EdbPatCell2.cxx:28
float eXmarg
Definition: EdbPatCell2.h:15
float eDTYlim
Definition: EdbPatCell2.h:13
float eDTXlim
Definition: EdbPatCell2.h:13
float eYmarg
margins for the cell definition
Definition: EdbPatCell2.h:15

◆ ~EdbPatCell2()

virtual EdbPatCell2::~EdbPatCell2 ( )
inlinevirtual
27{ResetCorr();}

Member Function Documentation

◆ AcceptCorrections()

void EdbPatCell2::AcceptCorrections ( EdbLayer layer)

accept corrections from the layer

38{
40 eDX = layer.GetAffineXY()->B1();
41 eDY = layer.GetAffineXY()->B1();
42 eDZ = layer.Zcorr();
43 eDTX = layer.GetAffineTXTY()->B1();
44 eDTY = layer.GetAffineTXTY()->B2();
45 eShr = layer.Shr();
46}
Float_t B2() const
Definition: EdbAffine.h:48
Float_t B1() const
Definition: EdbAffine.h:47
float Zcorr() const
Definition: EdbLayer.h:90
float Shr() const
Definition: EdbLayer.h:89
EdbAffine2D * GetAffineTXTY()
Definition: EdbLayer.h:120
EdbAffine2D * GetAffineXY()
Definition: EdbLayer.h:119
float eShr
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:19
float eDZ
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:18
float eDX
Definition: EdbPatCell2.h:18
float eDTX
Definition: EdbPatCell2.h:20
float eDY
Definition: EdbPatCell2.h:18
float eDTY
corrections to be applied if eApplyCorr==true
Definition: EdbPatCell2.h:20

◆ ApplyCorrections() [1/2]

void EdbPatCell2::ApplyCorrections ( EdbLayer layer)

propagate corrections to layer

50{
52 EdbAffine2D a2( 1,0,0,1, eDX, eDY);
53 EdbAffine2D at2(1,0,0,1, eDTX, eDTY);
54 layer.ApplyCorrections( eShr, eDZ, a2, at2 );
55 Log(3,"EdbPatCell2::ApplyCorrections","");
57
58}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
Definition: EdbAffine.h:17
void ApplyCorrections(EdbLayer &la)
Definition: EdbLayer.cxx:81
void PrintCorr()
Definition: EdbPatCell2.cxx:21
gEDBDEBUGLEVEL
Definition: energy.C:7

◆ ApplyCorrections() [2/2]

void EdbPatCell2::ApplyCorrections ( EdbPattern pat)

apply corrections to the segments of a given pat

62{
64 pat.TransformShr(eShr);
65 EdbAffine2D aff(1,0,0,1, eDX, eDY);
66 pat.Transform(&aff);
67 EdbAffine2D affa(1,0,0,1, eDTX, eDTY);
68 pat.TransformA(&affa);
69}
virtual void Transform(const EdbAffine2D *a)
Definition: EdbVirtual.cxx:155
void TransformA(const EdbAffine2D *affA)
Definition: EdbPattern.cxx:324
void TransformShr(const float shr)
Definition: EdbPattern.cxx:341

◆ FillCell() [1/2]

int EdbPatCell2::FillCell ( EdbPattern pat)

limits should be already defined

assume that the cell is already initialized

97{
99 EdbSegP *s=0;
100 int n = pat.N();
101 for(int i=0; i<n; i++) {
102 s = pat.GetSegment(i);
103 AddObject( Xs(*s), Ys(*s), (TObject*)s );
104 }
105 return eNcell;
106}
bool AddObject(float v[2], TObject *obj)
Definition: EdbCell2.h:173
Int_t eNcell
eNx*eNy
Definition: EdbCell2.h:28
float Xs(EdbSegP &s)
Definition: EdbPatCell2.h:40
float Ys(EdbSegP &s)
Definition: EdbPatCell2.h:41
Definition: EdbSegP.h:21
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
s
Definition: check_shower.C:55

◆ FillCell() [2/2]

int EdbPatCell2::FillCell ( TObjArray &  pat)

limits should be already defined

assume that the cell is already initialized

110{
112 EdbSegP *s=0;
113 int n = pat.GetEntriesFast();
114 for(int i=0; i<n; i++) {
115 s = (EdbSegP*)pat.UncheckedAt(i);
116 AddObject( Xs(*s), Ys(*s), (TObject*)s );
117 }
118 return eNcell;
119}

◆ FillCombinations()

int EdbPatCell2::FillCombinations ( EdbPatCell2 cell,
int  ir2[2],
TObjArray &  arrC1,
TObjArray &  arrC2,
bool  doFill = 1 
)
123{
124 int nout=0;
125
126 TObjArray arr1(10000);
127 TObjArray arr2(100);
128 float v[2];
129 EdbSegP *s1,*s2;
130
131 int ncomb; // douplications rate counter
132 int n1 = SelectObjects(arr1);
133 for(int i=0; i<n1; i++) {
134 ncomb=0;
135 s1 = (EdbSegP*)arr1.UncheckedAt(i);
136 if(s1->Flag()==-10) continue;
137 arr2.Clear();
138 v[0] = Xs(*s1);
139 v[1] = Ys(*s1);
140 int n2 = cell.SelectObjectsC(v,ir2,arr2);
141 if(n2<1) continue;
142 float tx1 = TXs(*s1);
143 float ty1 = TYs(*s1);
144 for(int i2=0; i2<n2; i2++) {
145 s2 = (EdbSegP*)arr2.UncheckedAt(i2);
146 if(s2->Flag()==-10) continue;
147 if( s2==s1 ) continue;
148 if( Abs(cell.Xs(*s2) - v[0]) > eDXlim ) continue;
149 if( Abs(cell.Ys(*s2) - v[1]) > eDYlim ) continue;
150 if( Abs(cell.TXs(*s2)- tx1) > eDTXlim ) continue;
151 if( Abs(cell.TYs(*s2)- ty1) > eDTYlim ) continue;
152
153 if(doFill) {
154 arrC1.Add(s1);
155 arrC2.Add(s2);
156 }
157 ncomb++;
158 }
159 nout+=ncomb;
160 if(eDoubletsRate) eDoubletsRate->Fill(ncomb);
161 }
162
163 Log(3,"EdbPatCell2::FillCombinations","%7d couples found with the acceptance:(dx,dy): %6.2f %6.2f (dtx,dty): %6.4f %6.4f",
166 if(gEDBDEBUGLEVEL>2&&cell.eApplyCorr) cell.PrintCorr();
167 return nout;
168}
int SelectObjects(TObjArray &arr)
Definition: EdbCell2.cpp:738
int SelectObjectsC(int iv[2], int ir[2], TObjArray &arr)
Definition: EdbCell2.cpp:707
float TYs(EdbSegP &s)
Definition: EdbPatCell2.h:50
float TXs(EdbSegP &s)
Definition: EdbPatCell2.h:49
bool eApplyCorr
Definition: EdbPatCell2.h:17
void Clear()
Definition: EdbSegP.h:86
Int_t Flag() const
Definition: EdbSegP.h:149
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ FillSelectedH2()

EdbH2 * EdbPatCell2::FillSelectedH2 ( )
172{
173 EdbH2 *h = new EdbH2();
174 h->InitH2(eN, eMin, eMax);
175 EdbSegP *s;
176 for(int i=0; i<eNcell; i++)
177 for(int j=0; j<eNC[i]; j++) {
178 s = (EdbSegP*)(GetObject(i,j));
179 if(s->Flag()==-10) continue;
180 h->AddBin(i,1);
181 }
182 return h;
183}
TObject * GetObject(float x, float y, int ientr) const
Definition: EdbCell2.h:195
fast 2-dim histogram class (used as a basis for EdbCell2)
Definition: EdbCell2.h:19
Int_t * eNC
[eNcell] number of objects/cell
Definition: EdbCell2.h:29
Float_t eMin[2]
min
Definition: EdbCell2.h:24
void AddBin(int jcell, int n)
Definition: EdbCell2.h:87
int InitH2(const EdbH2 &h)
Definition: EdbCell2.cpp:78
Int_t eN[2]
divisions
Definition: EdbCell2.h:23
Float_t eMax[2]
max
Definition: EdbCell2.h:25
EdbH2()
Definition: EdbCell2.cpp:24

◆ InitPat()

void EdbPatCell2::InitPat ( EdbPattern pat,
int  nx,
int  ny 
)

define limits to fit the pattern

73{
74 float min[2] = {pat.Xmin()-eXmarg, pat.Ymin()-eYmarg };
75 float max[2] = {pat.Xmax()+eXmarg, pat.Ymax()+eYmarg };
76 int n[2] = {nx,ny};
77 int maxcell = pat.N()/n[0]/n[1]+10;
78 maxcell += (int)(5*Sqrt(maxcell));
79 Log(2, "EdbPatCell2::InitPat", "maxcell = %d\n",maxcell);
80 InitCell(maxcell,n,min,max);
81}
float min(TClonesArray *t)
Definition: bitview.cxx:275
int InitCell(EdbCell2 &c)
Definition: EdbCell2.h:167
virtual Float_t Xmax() const
Definition: EdbVirtual.cxx:196
virtual Float_t Ymin() const
Definition: EdbVirtual.cxx:206
virtual Float_t Xmin() const
Definition: EdbVirtual.cxx:186
virtual Float_t Ymax() const
Definition: EdbVirtual.cxx:216
int max
Definition: check_shower.C:41

◆ InitPatBin()

void EdbPatCell2::InitPatBin ( EdbPattern pat,
float  binx,
float  biny 
)

define limits to fit the pattern

85{
86 float min[2] = {pat.Xmin()-eXmarg, pat.Ymin()-eYmarg };
87 float max[2] = {pat.Xmax()+eXmarg, pat.Ymax()+eYmarg };
88 int n[2] = { (int)((max[0]-min[0])/binx+1), (int)((max[1]-min[1])/biny+1) };
89 int maxcell = pat.N()/n[0]/n[1]+10;
90 maxcell += (int)(5*Sqrt(maxcell));
91 Log(2, "EdbPatCell2::InitPatBin", "maxcell = %d\n",maxcell);
92 InitCell(maxcell,n,min,max);
93}

◆ PrintCorr()

void EdbPatCell2::PrintCorr ( )
22{
23 printf("EdbPatCell2 corrections: (dx,dy,dz): %7.2f %7.2f %7.2f (shr): %6.4f (dtx,dty): %6.4f %6.4f\n",
25}

◆ ResetCorr()

void EdbPatCell2::ResetCorr ( )
29{
30 eApplyCorr = false;
31 ePhi = eDX = eDY = eDZ = eDTX = eDTY = 0.;
32 eShr = 1;
33 SafeDelete(eDoubletsRate);
34}
float ePhi
rotation angle (using dx,dy one can set-up the center of rotation)
Definition: EdbPatCell2.h:21

◆ TXs()

float EdbPatCell2::TXs ( EdbSegP s)
inline
49{ return eApplyCorr? s.eTX/eShr+eDTX : s.eTX; }

◆ TYs()

float EdbPatCell2::TYs ( EdbSegP s)
inline
50{ return eApplyCorr? s.eTY/eShr+eDTY : s.eTY; }

◆ Xs()

float EdbPatCell2::Xs ( EdbSegP s)
inline
40{ if(!eApplyCorr) return s.eX; if(ePhi<0.00001) return XsPos(s); else return XsRot(s); }
float XsPos(EdbSegP &s)
Definition: EdbPatCell2.h:43
float XsRot(EdbSegP &s)
Definition: EdbPatCell2.h:46

◆ XsPos()

float EdbPatCell2::XsPos ( EdbSegP s)
inline
43{ return eApplyCorr? eDX + s.eX + TXs(s)*eDZ : s.eX; }

◆ XsRot()

float EdbPatCell2::XsRot ( EdbSegP s)
inline
46{ return XsPos(s)*TMath::Cos(ePhi)-YsPos(s)*TMath::Sin(ePhi); }
float YsPos(EdbSegP &s)
Definition: EdbPatCell2.h:44

◆ Ys()

float EdbPatCell2::Ys ( EdbSegP s)
inline
41{ if(!eApplyCorr) return s.eY; if(ePhi<0.00001) return YsPos(s); else return YsRot(s); }
float YsRot(EdbSegP &s)
Definition: EdbPatCell2.h:47

◆ YsPos()

float EdbPatCell2::YsPos ( EdbSegP s)
inline
44{ return eApplyCorr? eDY + s.eY + TYs(s)*eDZ : s.eY; }

◆ YsRot()

float EdbPatCell2::YsRot ( EdbSegP s)
inline
47{ return XsPos(s)*TMath::Sin(ePhi)+YsPos(s)*TMath::Cos(ePhi); }

Member Data Documentation

◆ eApplyCorr

bool EdbPatCell2::eApplyCorr

◆ eDoubletsRate

TH1I* EdbPatCell2::eDoubletsRate

to be filled in FillCombinations()

◆ eDTX

float EdbPatCell2::eDTX

◆ eDTXlim

float EdbPatCell2::eDTXlim

◆ eDTY

float EdbPatCell2::eDTY

corrections to be applied if eApplyCorr==true

◆ eDTYlim

float EdbPatCell2::eDTYlim

◆ eDX

float EdbPatCell2::eDX

◆ eDXlim

float EdbPatCell2::eDXlim

◆ eDY

float EdbPatCell2::eDY

◆ eDYlim

float EdbPatCell2::eDYlim

acceptance limits for the combinations selection

◆ eDZ

float EdbPatCell2::eDZ

corrections to be applied if eApplyCorr==true

◆ ePhi

float EdbPatCell2::ePhi

rotation angle (using dx,dy one can set-up the center of rotation)

◆ eShr

float EdbPatCell2::eShr

corrections to be applied if eApplyCorr==true

◆ eXmarg

float EdbPatCell2::eXmarg

◆ eYmarg

float EdbPatCell2::eYmarg

margins for the cell definition


The documentation for this class was generated from the following files: