FEDRA emulsion software from the OPERA Collaboration
EdbTestAl Class Reference

alignment test class More...

#include <EdbTestAl.h>

Inheritance diagram for EdbTestAl:
Collaboration diagram for EdbTestAl:

Public Member Functions

int CheckMaxBin ()
 
int CheckMaxBin (float dz, float phi, float &meanbin, float &xmax, float &ymax)
 
int DubletsFilterOut (EdbPattern &p, float xbin, float ybin, float dMin=5, float dtMin=0.05)
 
 EdbTestAl ()
 
int FillTree (float gdz=0)
 
void HDistance (EdbPattern &p1, EdbPattern &p2)
 
int MakeTrans (EdbAffine2D &aff, float dz=0, const char *ccut="abs(dy)<400&&abs(dx)<400")
 
void PositionPlot (EdbPattern &p1, EdbPattern &p2, float xbin, float ybin, TObjArray &arr1, TObjArray &arr2)
 
void PositionPlot (EdbPattern &p1, EdbPattern &p2, TH2F &hd, float xbin=100, float ybin=100, TTree *posnt=0)
 
virtual ~EdbTestAl ()
 

Public Attributes

float eBinSize
 microns More...
 
TNtuple * eBinTree
 put bins value for all attempts More...
 
Float_t eD0 [4]
 limits, and found value of the peak for dx,dy,dz,phi More...
 
Float_t eDmax [4]
 
Float_t eDmin [4]
 
TFile * eFile
 
Int_t eITMAX
 angular step (def=50) More...
 
Float_t eMaxBin
 the max value corresponding to eD0 More...
 
Int_t eN [4]
 nx,ny,nz,nphi - number of divisions More...
 
Int_t eOCMAX
 occupancy (def=100) More...
 
float eOffset
 
TObjArray * eS1
 
TObjArray * eS2
 pointers to segments selected by HDistance More...
 
TTree * eT
 
TH2F * HD
 
TH2F * HDF
 
TH2F * HDF2
 

Detailed Description

alignment test class

Constructor & Destructor Documentation

◆ EdbTestAl()

EdbTestAl::EdbTestAl ( )
22{
23 eITMAX=50; // angular normalization (def=50) (the step is ~1./eITMAX)
24 eOCMAX=100; // occupancy (def=100)
25
26 eOffset=5000;
27 eBinSize=50.;//microns
28
29 eS1=0;
30 eS2=0;
31
32 eBinTree=0;
33 eT = 0;
34 eFile= 0;
35 HD=0;
36 HDF=0;
37 HDF2=0;
38
39 eMaxBin=0;
40 for(int i=0; i<4; i++) {eN[i]=0; eDmin[i]=0; eDmax[i]=0; eD0[i]=0;}
41
42}
float eBinSize
microns
Definition: EdbTestAl.h:21
TObjArray * eS1
Definition: EdbTestAl.h:23
Float_t eD0[4]
limits, and found value of the peak for dx,dy,dz,phi
Definition: EdbTestAl.h:34
TFile * eFile
Definition: EdbTestAl.h:28
TObjArray * eS2
pointers to segments selected by HDistance
Definition: EdbTestAl.h:23
TH2F * HDF2
Definition: EdbTestAl.h:31
Int_t eITMAX
angular step (def=50)
Definition: EdbTestAl.h:17
Float_t eMaxBin
the max value corresponding to eD0
Definition: EdbTestAl.h:35
Int_t eOCMAX
occupancy (def=100)
Definition: EdbTestAl.h:18
Int_t eN[4]
nx,ny,nz,nphi - number of divisions
Definition: EdbTestAl.h:33
Float_t eDmin[4]
Definition: EdbTestAl.h:34
Float_t eDmax[4]
Definition: EdbTestAl.h:34
TTree * eT
Definition: EdbTestAl.h:27
float eOffset
Definition: EdbTestAl.h:20
TH2F * HD
Definition: EdbTestAl.h:29
TH2F * HDF
Definition: EdbTestAl.h:30
TNtuple * eBinTree
put bins value for all attempts
Definition: EdbTestAl.h:25

◆ ~EdbTestAl()

EdbTestAl::~EdbTestAl ( )
virtual
46{
47 if(eT) eT->AutoSave("SaveSelf");
48 if(eFile) eFile->Close();
49}

Member Function Documentation

◆ CheckMaxBin() [1/2]

int EdbTestAl::CheckMaxBin ( )

Output: vmax - x,y,z,phi of the highest bin

231{
233
234 printf("CheckMaxBin:(z:nphi): %d %f %f %d %f %f\n",eN[3],eDmin[3],eDmax[3],eN[2],eDmin[2],eDmax[2]);
235
236 TH2F *h2 = new TH2F("z_phi","z_phi",eN[3],eDmin[3],eDmax[3],eN[2],eDmin[2],eDmax[2]);
237
238 float binz = (eDmax[2]-eDmin[2])/eN[2];
239 float binp = (eDmax[3]-eDmin[3])/eN[3];
240 float dz=0, dphi=0, mean=0;
241 float xmax=0,ymax=0;
242 int mbin=0;
243 for(int iz=0; iz<eN[2]; iz++) {
244 dz = eDmin[2] + iz*binz +binz/2;
245 for(int ip=0; ip<eN[3]; ip++) {
246 dphi = eDmin[3] + ip*binp + binp/2;
247 mbin = CheckMaxBin(dz,dphi,mean,xmax,ymax);
248 h2->Fill(dphi,dz,mbin);
249 printf("dx,dy,dz,dphi = %f %f %f %f; maxbin/mean= %d/%f = %f\n",xmax,ymax,dz,dphi, mbin,mean, mbin/mean);
250 if(mbin>=eMaxBin) {
251 eD0[0] = xmax; eD0[1]=ymax; eD0[2]=dz; eD0[3]=dphi;
252 eMaxBin=mbin;
253 }
254 }
255 }
256 h2->Draw("colz");
257
258 int mx,my,mz;
259 h2->GetMaximumBin(mx,my,mz);
260 float phi = h2->GetXaxis()->GetBinCenter(mx);
261 float z = h2->GetYaxis()->GetBinCenter(my);
262 mbin = (int)(h2->GetMaximum());
263 mean = h2->GetSum()/h2->GetNbinsX()/h2->GetNbinsY();
264 printf("max(%f,%f) = %d ; mean = %f\n",phi,z,mbin,mean);
265
266 return mbin;
267}
brick dz
Definition: RecDispMC.C:107
int CheckMaxBin()
Definition: EdbTestAl.cxx:230
TH1F * h2
Definition: energy.C:19

◆ CheckMaxBin() [2/2]

int EdbTestAl::CheckMaxBin ( float  dz,
float  phi,
float &  meanbin,
float &  xmax,
float &  ymax 
)

Output: vmax - x,y,z,phi of the highest bin

271{
273
274 if(!eS1) return 0;
275 EdbSegP *s1=0, *s2=0;
276 float x2p,y2p,dx,dy;
277 TH2F h2("H2","H2",eN[0],eDmin[0],eDmax[0],eN[1],eDmin[1],eDmax[1]);
278 int n = eS1->GetEntries();
279 for(int i=0; i<n; i++) {
280 s1 = (EdbSegP *)(eS1->UncheckedAt(i));
281 s2 = (EdbSegP *)(eS2->UncheckedAt(i));
282
283 x2p = s2->X()*Cos(phi)-s2->Y()*Sin(phi);
284 y2p = s2->X()*Sin(phi)+s2->Y()*Cos(phi);
285 // todo - transform angles?
286
287 dx= (x2p-s2->TX()*gdz) - s1->X();
288 dy= (y2p-s2->TY()*gdz) - s1->Y();
289 h2.Fill(dx,dy);
290 }
291
292 if(eBinTree) { // disabled test
293 for(int i=0; i<eN[0]; i++)
294 for(int j=0; j<eN[1]; j++) {
295 float bin=h2.GetBinContent(h2.GetBin(i,j));
296 if(bin<1) continue;
297 dx = eDmin[0] + (i+0.5)*(eDmax[0]-eDmin[0])/eN[0];
298 dy = eDmin[1] + (j+0.5)*(eDmax[1]-eDmin[1])/eN[1];
299 eBinTree->Fill(gdz,phi,dx,dy,bin);
300 }
301 }
302 meanbin = h2.GetSum()/h2.GetNbinsX()/h2.GetNbinsY();
303 int mx,my,mz;
304 h2.GetMaximumBin(mx,my,mz);
305 xmax = h2.GetXaxis()->GetBinCenter(mx);
306 ymax = h2.GetYaxis()->GetBinCenter(my);
307 return (int)(h2.GetMaximum());
308}
Definition: EdbSegP.h:21
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
EdbSegP * s1
Definition: tlg2couples.C:29
EdbSegP * s2
Definition: tlg2couples.C:30

◆ DubletsFilterOut()

int EdbTestAl::DubletsFilterOut ( EdbPattern p,
float  xbin,
float  ybin,
float  dMin = 5,
float  dtMin = 0.05 
)

assign flag = -10 for the duplicated segments

53{
55 int nout=0;
56
57 float dxMin=dMin,dyMin=dMin, dtxMin=dtMin, dtyMin=dtMin;
58
59 TObjArray arr1(10000), arr2(10000);
60 PositionPlot(p,p, xbin,ybin, arr1,arr2);
61 int n1 = arr1.GetEntries();
62 EdbSegP *s1,*s2;
63 for(int i=0; i<n1; i++) {
64 s1 = (EdbSegP*)arr1.At(i);
65 s2 = (EdbSegP*)arr2.At(i);
66 if( Abs(s2->X()-s1->X()) > dxMin ) continue;
67 if( Abs(s2->Y()-s1->Y()) > dyMin ) continue;
68 if( Abs(s2->TX()-s1->TX()) > dtxMin ) continue;
69 if( Abs(s2->TY()-s1->TY()) > dtyMin ) continue;
70 if( s2->Aid(0)==s1->Aid(0) && s2->Aid(1)==s1->Aid(1) && s1->Side()==s2->Side() ) continue;
71 if( s2->W()>s1->W() ) s1->SetFlag(-10);
72 else s2->SetFlag(-10);
73 nout++;
74 }
75 Log(2,"DubletsFilterOut","%d segments discarded",nout);
76 return nout;
77}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
Int_t Side() const
mandatory virtual functions:
Definition: EdbSegP.h:170
Float_t W() const
Definition: EdbSegP.h:151
Int_t Aid(int i) const
Definition: EdbSegP.h:169
void SetFlag(int flag)
Definition: EdbSegP.h:130
void PositionPlot(EdbPattern &p1, EdbPattern &p2, float xbin, float ybin, TObjArray &arr1, TObjArray &arr2)
Definition: EdbTestAl.cxx:106
p
Definition: testBGReduction_AllMethods.C:8

◆ FillTree()

int EdbTestAl::FillTree ( float  gdz = 0)
312{
313 if(!eS1) return 0;
314
315 eFile = new TFile("alt.root","RECREATE");
316 eT = new TTree("T","Alignment tree");
317 float dx,dy;
318 static DTR dtr;
319 eT->Branch("dtr",&dtr,"it1/I:it2/I:x1/F:y1/F:z1/F:x2/F:y2/F:z2/F:tx1/F:ty1/F:tx2/F:ty2/F:chi1/F:chi2/F");
320 EdbSegP *s1,*s2;
321 int steps=int(2.*eOffset/eBinSize);
322 if(HD) HD->Reset();
323 else HD=new TH2F("HD","HD",steps,-eOffset,eOffset,steps,-eOffset,eOffset);
324 if(HDF) HDF->Reset();
325 else HDF=new TH2F("HDF","HDF",100,-eBinSize,eBinSize,100,-eBinSize,eBinSize);
326 HDF2=new TH2F("HDF2","HDF2",100,-3*eBinSize,3*eBinSize,100,-3*eBinSize,3*eBinSize);
327
328 int n = eS1->GetEntries();
329 for(int i=0; i<n; i++) {
330 s1 = (EdbSegP *)(eS1->UncheckedAt(i));
331 s2 = (EdbSegP *)(eS2->UncheckedAt(i));
332 dtr.x1=s1->X();
333 dtr.y1=s1->Y();
334 dtr.z1=s1->Z();
335 dtr.x2=s2->X();
336 dtr.y2=s2->Y();
337 dtr.z2=s2->Z();
338 dtr.chi1=s1->Chi2();
339 dtr.chi2=s2->Chi2();
340
341 dtr.tx1=s1->TX();
342 dtr.ty1=s1->TY();
343 dtr.tx2=s2->TX();
344 dtr.ty2=s2->TY();
345
346 dx= (dtr.x2+s2->TX()*gdz) - dtr.x1;
347 dy= (dtr.y2+s2->TY()*gdz) - dtr.y1;
348 eT->Fill();
349 HD->Fill(dx,dy);
350 HDF->Fill(dx,dy);
351 HDF2->Fill(dx,dy);
352 }
353 HD->Write();
354 HDF->Write();
355 HDF2->Write();
356 eT->Write();
357 return n;
358}
Float_t Chi2() const
Definition: EdbSegP.h:157
Float_t Z() const
Definition: EdbSegP.h:153
Definition: EdbTestAl.cxx:13
Float_t ty1
Definition: EdbTestAl.cxx:16
Float_t y1
Definition: EdbTestAl.cxx:15
Float_t x2
Definition: EdbTestAl.cxx:15
Float_t chi1
Definition: EdbTestAl.cxx:17
Float_t chi2
Definition: EdbTestAl.cxx:17
Float_t tx1
Definition: EdbTestAl.cxx:16
Float_t z2
Definition: EdbTestAl.cxx:15
Float_t y2
Definition: EdbTestAl.cxx:15
Float_t z1
Definition: EdbTestAl.cxx:15
Float_t ty2
Definition: EdbTestAl.cxx:16
Float_t x1
Definition: EdbTestAl.cxx:15
Float_t tx2
Definition: EdbTestAl.cxx:16

◆ HDistance()

void EdbTestAl::HDistance ( EdbPattern p1,
EdbPattern p2 
)
177{
178 int nx=200, ny=200, no=100;
179 Long_t *ind1 = new Long_t[nx*ny*no];
180 Long_t *ind2 = new Long_t[nx*ny*no];
181 Int_t *oc1 = new Int_t[nx*ny];
182 Int_t *oc2 = new Int_t[nx*ny];
183 memset(ind1,0,nx*ny*no*sizeof(Long_t));
184 memset(ind2,0,nx*ny*no*sizeof(Long_t));
185 memset( oc1,0,nx*ny*sizeof(Int_t));
186 memset( oc2,0,nx*ny*sizeof(Int_t));
187
188 int itx,ity;
189
190 for(int i=0;i<eITMAX;i++) for(int j=0;j<eITMAX;j++) { oc1[i*nx+j]=0; oc2[i*nx+j]=0;}
191 for(int i=0;i<p1.N();i++){
192 itx=int( (p1.GetSegment(i)->TX()+1.)*eITMAX/2. );
193 ity=int( (p1.GetSegment(i)->TY()+1.)*eITMAX/2. );
194 if(itx<eITMAX)
195 if(ity<eITMAX)
196 if(itx>=0)
197 if(ity>=0)
198 if(oc1[itx*nx+ity]<eOCMAX-1)
199 { ind1[itx*nx*ny+ity*ny+oc1[itx*nx+ity]]=i; (oc1[itx*nx+ity])++;}
200 }
201 for(int i=0;i<p2.N();i++){
202 itx=int( (p2.GetSegment(i)->TX()+1.)*eITMAX/2 );
203 ity=int( (p2.GetSegment(i)->TY()+1.)*eITMAX/2 );
204 if(itx<eITMAX)
205 if(ity<eITMAX)
206 if(itx>=0)
207 if(ity>=0)
208 if(oc2[itx*nx+ity]<eOCMAX-1)
209 { ind2[itx*nx*ny+ity*ny+oc2[itx*nx+ity]]=i; (oc2[itx*nx+ity])++;}
210 }
211
212 if(!eS1) eS1=new TObjArray(1000); else eS1->Clear();
213 if(!eS2) eS2=new TObjArray(1000); else eS2->Clear();
214
215 for(int itx=0;itx<eITMAX;itx++) for(int ity=0;ity<eITMAX;ity++)
216 for(int i=0;i<oc1[itx*nx+ity];i++)
217 for(int j=0;j<oc2[itx*nx+ity];j++)
218 {
219 eS1->Add( p1.GetSegment(ind1[itx*nx*ny+ity*ny+i]) );
220 eS2->Add( p2.GetSegment(ind2[itx*nx*ny+ity*ny+j]) );
221 }
222
223 delete [] ind1;
224 delete [] ind2;
225 delete [] oc1;
226 delete [] oc2;
227}
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66

◆ MakeTrans()

int EdbTestAl::MakeTrans ( EdbAffine2D aff,
float  dz = 0,
const char *  ccut = "abs(dy)<400&&abs(dx)<400" 
)
383{
384 static DTR dtr;
385 eT->SetBranchAddress("dtr",&dtr);
386
387 printf("MakeTrans: propagate the first pattern to %f before...\n",dz);
388 char str[64];
389 eT->SetAlias("dtx","tx2-tx1");
390 eT->SetAlias("dty","ty2-ty1");
391 sprintf(str,"x2-(%f*tx1+x1)",dz); eT->SetAlias("dx",str);
392 sprintf(str,"y2-(%f*ty1+y1)",dz); eT->SetAlias("dy",str);
393
394 TCanvas *c = new TCanvas("alcan","alignment check");
395 c->Divide(2,2);
396
397 TCut cut=ccut;
398
399 c->cd(1); eT->Draw("dy:dx",cut);
400 c->cd(2); eT->Draw("dy:x2",cut);
401 c->cd(4); eT->Draw("dx:y2",cut);
402
403 TEventList *lst =0;
404 eT->Draw(">>lst", cut );
405 lst = (TEventList*)gDirectory->GetList()->FindObject("lst");
406 int nlst =lst->GetN();
407 int entr=0;
408
409 TArrayF x1(nlst);
410 TArrayF x2(nlst);
411 TArrayF y1(nlst);
412 TArrayF y2(nlst);
413
414 for(int i=0; i<nlst; i++ ) {
415 entr = lst->GetEntry(i);
416 eT->GetEntry(entr);
417 x1[i]=dtr.x1;
418 y1[i]=dtr.y1;
419 x2[i]=dtr.x2;
420 y2[i]=dtr.y2;
421 }
422
423 aff.Calculate(nlst,x1.GetArray(),y1.GetArray(),x2.GetArray(),y2.GetArray());
424 return nlst;
425}
Int_t Calculate(EdbPointsBox2D *b1, EdbPointsBox2D *b2)
Definition: EdbAffine.cxx:231
TCut cut
Definition: check_shower.C:6
new TCanvas()

◆ PositionPlot() [1/2]

void EdbTestAl::PositionPlot ( EdbPattern p1,
EdbPattern p2,
float  xbin,
float  ybin,
TObjArray &  arr1,
TObjArray &  arr2 
)

speed-optimized position plot when xbin,ybin << (xmax-xmin),(ymax-ymin)

107{
109
110 EdbSegP *s = p1.GetSegment(0);
111 float xmin = p1.Xmin(), xmax = p1.Xmax();
112 float ymin = p1.Ymin(), ymax = p1.Ymax();
113 int nx = (int)((xmax-xmin)/xbin+2);
114 int ny = (int)((ymax-ymin)/ybin+2);
115 //printf("xbin = %f, ybin=%f nx = %d ny = %d\n", xbin, ybin, nx,ny);
116
117 int meancell1 = (int)Sqrt(p1.N()/nx/ny);
118 int meancell2 = (int)Sqrt(p2.N()/nx/ny);
119 int meancell = (int)Max(meancell1,meancell2)+10;
120
121 TObjArray as1(nx*ny), as2(nx*ny);
122 for(int i=0; i<nx*ny; i++) {
123 as1.Add(new TObjArray(meancell));
124 as2.Add(new TObjArray(meancell));
125 }
126
127 int kx,ky, k;
128 for(int i=0; i<p1.N(); i++) {
129 s = p1.GetSegment(i);
130 if(s->Flag()==-10) continue;
131 ky = (int)((s->Y()-ymin)/ybin);
132 kx = (int)((s->X()-xmin)/xbin);
133 k = nx*ky + kx;
134 ((TObjArray*)(as1.At(k)))->Add(s);
135 }
136 for(int i=0; i<p2.N(); i++) {
137 s = p2.GetSegment(i);
138 if(s->Flag()==-10) continue;
139 if( s->X() <= xmin-xbin/2) continue;
140 if( s->X() >= xmax+xbin/2) continue;
141 if( s->Y() <= ymin-ybin/2) continue;
142 if( s->Y() >= ymax+ybin/2) continue;
143
144 kx = (int)((s->X()-xmin - xbin/2)/xbin);
145 ky = (int)((s->Y()-ymin - ybin/2)/ybin);
146 if(kx>nx-1||kx<0||ky>ny-1||ky<0) continue;
147 ((TObjArray*)(as2.At(nx*ky+kx)))->Add(s);
148 ((TObjArray*)(as2.At(nx*ky+kx+1)))->Add(s);
149 ((TObjArray*)(as2.At(nx*(ky+1)+kx)))->Add(s);
150 ((TObjArray*)(as2.At(nx*(ky+1)+kx+1)))->Add(s);
151 }
152
153 TObjArray *a1,*a2;
154 EdbSegP *s1,*s2;
155 for(int i=0; i<nx*ny; i++) {
156 a1 = (TObjArray*)as1.At(i);
157 a2 = (TObjArray*)as2.At(i);
158 for(int i1=0; i1<a1->GetEntriesFast(); i1++)
159 {
160 s1 = (EdbSegP*)a1->At(i1);
161 for(int i2=0; i2<a2->GetEntriesFast(); i2++)
162 {
163 s2 = (EdbSegP*)a2->At(i2);
164 if(s1==s2) continue; // useful in case if p1==p2
165 if( Abs(s2->TX()-s1->TX()) > 0.25) continue;
166 if( Abs(s2->TY()-s1->TY()) > 0.25) continue;
167 arr1.Add(s1);
168 arr2.Add(s2);
169 }
170 }
171 }
172
173}
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
s
Definition: check_shower.C:55

◆ PositionPlot() [2/2]

void EdbTestAl::PositionPlot ( EdbPattern p1,
EdbPattern p2,
TH2F &  hd,
float  xbin = 100,
float  ybin = 100,
TTree *  posnt = 0 
)
81{
82 TObjArray arr1(10000), arr2(10000);
83 PositionPlot(p1,p2,xbin,ybin,arr1,arr2);
84 int n1 = arr1.GetEntries();
85 EdbSegP *s1,*s2;
86 for(int i=0; i<n1; i++) {
87 s1 = (EdbSegP*)arr1.At(i);
88 s2 = (EdbSegP*)arr2.At(i);
89
90 hd.Fill(s2->X()-s1->X(), s2->Y()-s1->Y());
91 if(posnt) {
92 posnt->SetBranchAddress("s1.",&s1);
93 posnt->SetBranchAddress("s2.",&s2);
94 posnt->Fill();
95 }
96
97// if(posnt) posnt->Fill(
98// s1->X(),s1->Y(),s1->Z(),s1->TX(),s1->TY(),
99// s2->X(),s2->Y(),s2->Z(),s2->TX(),s2->TY(),
100// s1->Vid(0),s1->Vid(1),s2->Vid(0),s2->Vid(1)
101// );
102 }
103}

Member Data Documentation

◆ eBinSize

float EdbTestAl::eBinSize

microns

◆ eBinTree

TNtuple* EdbTestAl::eBinTree

put bins value for all attempts

◆ eD0

Float_t EdbTestAl::eD0[4]

limits, and found value of the peak for dx,dy,dz,phi

◆ eDmax

Float_t EdbTestAl::eDmax[4]

◆ eDmin

Float_t EdbTestAl::eDmin[4]

◆ eFile

TFile* EdbTestAl::eFile

◆ eITMAX

Int_t EdbTestAl::eITMAX

angular step (def=50)

◆ eMaxBin

Float_t EdbTestAl::eMaxBin

the max value corresponding to eD0

◆ eN

Int_t EdbTestAl::eN[4]

nx,ny,nz,nphi - number of divisions

◆ eOCMAX

Int_t EdbTestAl::eOCMAX

occupancy (def=100)

◆ eOffset

float EdbTestAl::eOffset

◆ eS1

TObjArray* EdbTestAl::eS1

◆ eS2

TObjArray * EdbTestAl::eS2

pointers to segments selected by HDistance

◆ eT

TTree* EdbTestAl::eT

◆ HD

TH2F* EdbTestAl::HD

◆ HDF

TH2F* EdbTestAl::HDF

◆ HDF2

TH2F* EdbTestAl::HDF2

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