FEDRA emulsion software from the OPERA Collaboration
EdbPattern Class Reference

#include <EdbPattern.h>

Inheritance diagram for EdbPattern:
Collaboration diagram for EdbPattern:

Public Member Functions

TIndexCellCell () const
 
 EdbPattern ()
 
 EdbPattern (EdbPattern &p)
 
 EdbPattern (float x0, float y0, float z0, int n=0)
 
EdbPatternExtractSubPattern (float min[5], float max[5], int MCevt=-1)
 
void FillCell (float stepx, float stepy, float steptx, float stepty)
 
int FindCompliments (EdbSegP &s, TObjArray &arr, float nsig, float nsigt)
 
EdbSegPFindSegment (int id)
 
int ID () const
 
Int_t NAff () const
 
int PID () const
 
Int_t Plate () const
 
void Reset ()
 
EdbID ScanID () const
 
void Set0 ()
 
void SetID (int id)
 
void SetNAff (int n)
 
void SetPID (int pid)
 
void SetScanID (EdbID id)
 
void SetSegmentsPID ()
 
void SetSegmentsScanID (EdbID id)
 
void SetSide (int side)
 
void SetSigma (float sx, float sy, float stx, float sty)
 
void SetStep (float stepx, float stepy, float steptx, float stepty)
 
Int_t Side () const
 
float StepTX () const
 
float StepTY () const
 
float StepX () const
 
float StepY () const
 
float SummaryPath ()
 
float Xmean ()
 
float Ymean ()
 
virtual ~EdbPattern ()
 
- Public Member Functions inherited from EdbSegmentsBox
EdbSegPAddSegment (EdbSegP &s)
 
EdbSegPAddSegment (EdbSegP &s1, EdbSegP &s2)
 
EdbSegPAddSegment (int i, EdbSegP &s)
 
EdbSegPAddSegment (int id, float x, float y, float tx, float ty, float w=0, int flag=0)
 
EdbSegPAddSegmentNoDuplicate (EdbSegP &s)
 
EdbPointAt (int i) const
 
Int_t CalculateAXAY (EdbSegmentsBox *p, EdbAffine2D *affA)
 
Int_t CalculateXY (EdbSegmentsBox *p, EdbAffine2D *aff)
 
Float_t Diff (EdbSegmentsBox &p)
 
Float_t DiffAff (EdbAffine2D *aff)
 
Float_t DZ () const
 
 EdbSegmentsBox (EdbSegmentsBox &box)
 
 EdbSegmentsBox (float x0, float y0, float z0, int nseg=0)
 
 EdbSegmentsBox (int nseg=0)
 
Int_t GetN () const
 
EdbSegPGetSegment (int i) const
 
EdbSegPGetSegmentLast () const
 
TClonesArray * GetSegments () const
 
voidGetSegmentsAddr ()
 
Float_t GetSize (Int_t XorY)
 
Float_t GetSizeX ()
 
Float_t GetSizeXY ()
 
Float_t GetSizeY ()
 
Float_t GetTrackDensity ()
 
Float_t GetTrackDensitymm2 ()
 
Int_t N () const
 
void Print (Option_t *opt="") const
 other functions More...
 
void ProjectTo (const float dz)
 
void Reset ()
 
void Set0 ()
 
void SetSegmentsDZ (float dz)
 
void SetSegmentsPlate (int plate)
 
void SetSegmentsZ ()
 
void SetX (float x)
 mandatory virtual functions: More...
 
void SetY (float y)
 
void SetZ (float z)
 
void TransformA (const EdbAffine2D *affA)
 
void TransformARot (const EdbAffine2D *affA)
 
void TransformShr (const float shr)
 
Float_t X () const
 
Float_t Y () const
 
Float_t Z () const
 
virtual ~EdbSegmentsBox ()
 
- Public Member Functions inherited from EdbPointsBox2D
virtual EdbPointAt (int i) const =0
 
virtual Float_t DeltaX ()
 
virtual Float_t DeltaY ()
 
virtual void DrawPoints (int style=23, int col=4, float size=1.)
 
 EdbPointsBox2D ()
 
 EdbPointsBox2D (const EdbPointsBox2D &pb)
 
virtual const EdbAffine2DGetKeep () const
 
virtual void GetKeep (EdbAffine2D &aff)
 
virtual Int_t N () const =0
 
virtual void Print (Option_t *opt="") const
 
virtual void Retransform ()
 
virtual void Rotate (float angle)
 
virtual void ScaleX (float scaleFactor)
 
virtual void ScaleY (float scaleFactor)
 
virtual void SetKeep (float a11, float a12, float a21, float a22, float b1, float b2)
 
virtual void SetX (float x)
 
virtual void SetY (float y)
 
virtual void SetZ (float z)
 
virtual void ShiftX (float offset)
 
virtual void ShiftY (float offset)
 
virtual void SmearXY (float sigmax, float sigmay)
 
virtual void Substruct (EdbPointsBox2D *b)
 
virtual void Transform (const EdbAffine2D *a)
 
virtual Float_t X () const
 
virtual TH1F * Xhist ()
 
virtual Float_t Xmax () const
 
virtual Float_t Xmin () const
 
virtual TH2F * XYhist ()
 
virtual Float_t Y () const
 
virtual TH1F * Yhist ()
 
virtual Float_t Ymax () const
 
virtual Float_t Ymin () const
 
virtual Float_t Z () const
 
virtual ~EdbPointsBox2D ()
 
- Public Member Functions inherited from EdbPoint3D
virtual void Print (Option_t *opt="") const
 
virtual void SetZ (float z)=0
 
virtual void Substruct (EdbPoint *p)
 
virtual void Test () const
 
virtual void TestPoint3D () const
 
virtual void Transform (const EdbAffine3D *a)
 
virtual Float_t Z () const =0
 
virtual ~EdbPoint3D ()
 
- Public Member Functions inherited from EdbPoint2D
virtual void Print (Option_t *opt="") const
 
virtual void SetX (float x)=0
 
virtual void SetY (float y)=0
 
virtual void SetZ (float z)
 
virtual void Substruct (EdbPoint *p)
 
virtual void Test () const
 
virtual void TestPoint2D () const
 
virtual void Transform (const EdbAffine2D *a)
 
virtual Float_t X () const =0
 
virtual Float_t Y () const =0
 
virtual Float_t Z () const
 
virtual ~EdbPoint2D ()
 
- Public Member Functions inherited from EdbPoint
virtual void SetX (float x)=0
 
virtual void SetY (float y)=0
 
virtual void SetZ (float z)=0
 
virtual void Substruct (EdbPoint *p)=0
 
virtual void Test () const
 
virtual void Transform (const EdbAffine2D *a)
 
virtual void Transform (const EdbAffine3D *a)
 
virtual Float_t X () const =0
 
virtual Float_t Y () const =0
 
virtual Float_t Z () const =0
 
virtual ~EdbPoint ()
 

Private Attributes

TIndexCelleCell
 associated with eSegments More...
 
Int_t eFlag
 pattern flag More...
 
Int_t eID
 pattern id in the volume More...
 
Int_t eNAff
 number of segments selected for affine calculation More...
 
Int_t ePID
 correspond to the piece ID More...
 
EdbID eScanID
 main scanning ID for this pattern, defined when possible More...
 
Int_t eSide
 emulsion side 0/1/2, defined when possible More...
 
Float_t eSigma [4]
 accuracy in comparison to neibour pattern More...
 
Float_t eStepTX
 
Float_t eStepTY
 bin size for the index cell More...
 
Float_t eStepX
 
Float_t eStepY
 bin size for the index cell More...
 

Constructor & Destructor Documentation

◆ EdbPattern() [1/3]

EdbPattern::EdbPattern ( )
1255 : EdbSegmentsBox()
1256{
1257 eCell = new TIndexCell();
1258 Set0();
1259}
TIndexCell * eCell
associated with eSegments
Definition: EdbPattern.h:279
void Set0()
Definition: EdbPattern.cxx:1290
EdbSegmentsBox(int nseg=0)
Definition: EdbPattern.cxx:31
sort collection with attributes
Definition: TIndexCell.h:19

◆ EdbPattern() [2/3]

EdbPattern::EdbPattern ( EdbPattern p)
1270{
1271 eID = p.eID;
1272 ePID = p.ePID;
1273 eSide = p.eSide;
1274 eScanID = p.eScanID;
1275 eFlag = p.eFlag;
1276 eNAff = p.eNAff;
1277 eStepX = p.eStepX; eStepY = p.eStepY;
1278 eStepTX = p.eStepTX; eStepTY = p.eStepTY;
1279 for(int i=0; i<4; i++) eSigma[i] = p.eSigma[i];
1280 if(p.eCell) eCell = new TIndexCell(*(p.eCell));
1281}
Int_t eID
pattern id in the volume
Definition: EdbPattern.h:277
Float_t eStepY
bin size for the index cell
Definition: EdbPattern.h:280
Int_t ePID
correspond to the piece ID
Definition: EdbPattern.h:278
Int_t eNAff
number of segments selected for affine calculation
Definition: EdbPattern.h:284
Int_t eSide
emulsion side 0/1/2, defined when possible
Definition: EdbPattern.h:286
EdbID eScanID
main scanning ID for this pattern, defined when possible
Definition: EdbPattern.h:285
Float_t eStepTY
bin size for the index cell
Definition: EdbPattern.h:281
Float_t eStepTX
Definition: EdbPattern.h:281
Float_t eSigma[4]
accuracy in comparison to neibour pattern
Definition: EdbPattern.h:282
Int_t eFlag
pattern flag
Definition: EdbPattern.h:283
Float_t eStepX
Definition: EdbPattern.h:280
p
Definition: testBGReduction_AllMethods.C:8

◆ EdbPattern() [3/3]

EdbPattern::EdbPattern ( float  x0,
float  y0,
float  z0,
int  n = 0 
)
1262 : EdbSegmentsBox(x0,y0,z0,n)
1263{
1264 Set0();
1265 eCell = new TIndexCell();
1266}
brick z0
Definition: RecDispMC.C:106

◆ ~EdbPattern()

EdbPattern::~EdbPattern ( )
virtual
1285{
1286 if(eCell) { delete eCell; eCell=0; }
1287}

Member Function Documentation

◆ Cell()

TIndexCell * EdbPattern::Cell ( ) const
inline
321{return eCell;}

◆ ExtractSubPattern()

EdbPattern * EdbPattern::ExtractSubPattern ( float  min[5],
float  max[5],
int  MCevt = -1 
)

return the copy of selected segments

1471{
1475 EdbPattern *pat = new EdbPattern( X(), Y(), Z() );
1476 EdbSegP *s;
1477 int nseg = N();
1478
1479 for(int i=0; i<nseg; i++) {
1480 s = GetSegment(i);
1481
1482 if (s->MCEvt()!=MCEvt && s->MCEvt()>0 && MCEvt>=0) continue;
1483 // Do not continue not in case MCEvt was not specified at all.
1484 // This allows backward compability.
1485
1486 if(s->X() < min[0]) continue;
1487 if(s->Y() < min[1]) continue;
1488 if(s->TX() < min[2]) continue;
1489 if(s->TY() < min[3]) continue;
1490 if(s->W() < min[4]) continue;
1491
1492 if(s->X() > max[0]) continue;
1493 if(s->Y() > max[1]) continue;
1494 if(s->TX() > max[2]) continue;
1495 if(s->TY() > max[3]) continue;
1496 if(s->W() > max[4]) continue;
1497
1498 pat->AddSegment(*s);
1499 }
1500
1501 // Set ID() and PID() to have consistent values of this pattern:
1502 pat->SetID(this->ID());
1503 pat->SetPID(this->PID());
1504
1505 return pat;
1506}
float min(TClonesArray *t)
Definition: bitview.cxx:275
Definition: EdbPattern.h:273
void SetID(int id)
Definition: EdbPattern.h:309
void SetPID(int pid)
Definition: EdbPattern.h:310
int PID() const
Definition: EdbPattern.h:320
int ID() const
Definition: EdbPattern.h:319
EdbPattern()
Definition: EdbPattern.cxx:1255
Definition: EdbSegP.h:21
Float_t Z() const
Definition: EdbPattern.h:84
Int_t N() const
Definition: EdbPattern.h:86
Float_t Y() const
Definition: EdbPattern.h:83
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
Float_t X() const
Definition: EdbPattern.h:82
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72
s
Definition: check_shower.C:55
int max
Definition: check_shower.C:41

◆ FillCell()

void EdbPattern::FillCell ( float  stepx,
float  stepy,
float  steptx,
float  stepty 
)

fill cells with fixed size at z=zPat

1324{
1326
1327 TIndexCell *cell = Cell();
1328 if(cell) cell->Drop();
1329 SetStep(stepx,stepy,steptx,stepty);
1330
1331 float x,y,tx,ty,dz;
1332 Long_t val[5]; // x,y,ax,ay,i
1333 EdbSegP *p;
1334 int npat = N();
1335 for(int i=0; i<npat; i++ ) {
1336 p = GetSegment(i);
1337 dz = Z() - p->Z();
1338 tx = p->TX();
1339 ty = p->TY();
1340 x = p->X() + tx*dz;
1341 y = p->Y() + ty*dz;
1342 val[0]= (Long_t)( x / stepx );
1343 val[1]= (Long_t)( y / stepy );
1344 val[2]= (Long_t)( tx/ steptx );
1345 val[3]= (Long_t)( ty/ stepty );
1346 val[4]= (Long_t)(i);
1347 cell->Add(5,val);
1348 }
1349 cell->Sort();
1350
1351}
brick dz
Definition: RecDispMC.C:107
Int_t npat
Definition: Xi2HatStartScript.C:33
void SetStep(float stepx, float stepy, float steptx, float stepty)
Definition: EdbPattern.h:299
TIndexCell * Cell() const
Definition: EdbPattern.h:321
void Drop()
Definition: TIndexCell.cpp:237
void Sort(Int_t upto=kMaxInt)
Definition: TIndexCell.cpp:539
Int_t Add(Int_t narg, Long_t varg[])
Definition: TIndexCell.cpp:602

◆ FindCompliments()

int EdbPattern::FindCompliments ( EdbSegP s,
TObjArray &  arr,
float  nsig,
float  nsigt 
)

return the array of segments compatible with the prediction segment s with the accuracy of nsig (in sigmas)

TODO: move this cycle into iterator

1355{
1358
1359 float dz = Z()-s.Z();
1360 float x = s.X() + s.TX()*dz;
1361 float y = s.Y() + s.TY()*dz;
1362 float sx = TMath::Sqrt( s.SX() + s.STX()*dz*dz );
1363 float sy = TMath::Sqrt( s.SY() + s.STY()*dz*dz );
1364
1365 float stx = s.STX();
1366 if (stx <= 0.) stx = 1.;
1367 float sty = s.STY();
1368 if (sty <= 0.) sty = 1.;
1369
1370// printf("Step: %f %f %f %f\n",StepX(),StepY(),StepTX(),StepTY());
1371// printf("Sigm: %f %f %f %f\n",sx,sy,TMath::Sqrt(stx),TMath::Sqrt(sty));
1372
1373
1374 long vcent[4] = { (long)(x/StepX()),
1375 (long)(y/StepY()),
1376 (long)(s.TX()/StepTX()),
1377 (long)(s.TY()/StepTY()) };
1378 long vdiff[4] = { (long)(sx*nsigx/StepX()+1),
1379 (long)(sy*nsigx/StepY()+1),
1380 (long)(TMath::Sqrt(stx)*nsigt/StepTX()+1),
1381 (long)(TMath::Sqrt(sty)*nsigt/StepTY()+1) };
1382
1383 sy*=sy;
1384 sx*=sx;
1385
1386 float dx,dy,dtx,dty;
1387
1388 int nseg=0;
1389 EdbSegP *seg=0;
1390
1391
1392 //printf("find compliments in pattern %d\n",ID());
1393 //s.Print();
1394
1395 long vmin[4],vmax[4];
1396 for(int i=0; i<4; i++) {
1397 vmin[i] = vcent[i]-vdiff[i];
1398 vmax[i] = vcent[i]+vdiff[i];
1399 //printf("%f +- %f \t",(float)vcent[i],(float)vdiff[i]);
1400 }
1401 //printf("\n");
1402
1404
1405 TIndexCell *c1=0,*c2=0,*c3=0,*c4=0;
1406 for(vcent[0]=vmin[0]; vcent[0]<=vmax[0]; vcent[0]++) {
1407 c1 = eCell->Find(vcent[0]);
1408 if(!c1) continue;
1409 for(vcent[1]=vmin[1]; vcent[1]<=vmax[1]; vcent[1]++) {
1410 c2 = c1->Find(vcent[1]);
1411 if(!c2) continue;
1412 for(vcent[2]=vmin[2]; vcent[2]<=vmax[2]; vcent[2]++) {
1413 c3 = c2->Find(vcent[2]);
1414 if(!c3) continue;
1415 for(vcent[3]=vmin[3]; vcent[3]<=vmax[3]; vcent[3]++) {
1416 c4 = c3->Find(vcent[3]);
1417 if(!c4) continue;
1418
1419 for(int i=0; i<c4->N(); i++) {
1420 seg = GetSegment(c4->At(i)->Value());
1421 dtx=s.TX()-seg->TX();
1422 if( dtx*dtx > stx*nsigt*nsigt ) continue;
1423 dty=s.TY()-seg->TY();
1424 if( dty*dty > sty*nsigt*nsigt ) continue;
1425 dz = seg->Z()-s.Z();
1426 dx=s.X()+s.TX()*dz-seg->X();
1427 if( dx*dx > sx*nsigx*nsigx ) continue;
1428 dy=s.Y()+s.TY()*dz-seg->Y();
1429 if( dy*dy > sy*nsigx*nsigx ) continue;
1430 //if(!s.IsCompatible(*seg,nsigx,nsigt)) continue;
1431 arr.Add(seg);
1432 nseg++;
1433 }
1434
1435 }
1436 }
1437 }
1438 }
1439
1440
1441// TIndexCellIterV itr( eCell, 4, vcent, vdiff );
1442// const TIndexCell *c=0;
1443// while( (c=itr.Next()) )
1444// for(int i=0; i<c->N(); i++) {
1445// seg = GetSegment(c->At(i)->Value());
1446// if(!s.IsCompatible(*seg,nsigx,nsigt)) continue;
1447// arr.Add(seg);
1448// nseg++;
1449// }
1450
1451 return nseg;
1452}
float StepX() const
Definition: EdbPattern.h:315
float StepTY() const
Definition: EdbPattern.h:318
float StepTX() const
Definition: EdbPattern.h:317
float StepY() const
Definition: EdbPattern.h:316
Float_t TX() const
tangens = deltaX/deltaZ
Definition: EdbSegP.h:175
Float_t X() const
Definition: EdbSegP.h:173
Float_t Z() const
Definition: EdbSegP.h:153
Float_t Y() const
Definition: EdbSegP.h:174
Float_t TY() const
tangens = deltaY/deltaZ
Definition: EdbSegP.h:176
TIndexCell * Find(Int_t narg, Long_t varg[]) const
Definition: TIndexCell.cpp:613
TCanvas * c1
Definition: energy.C:13
TCanvas * c2
Definition: energy.C:26
MFTYPE32 long(MFTYPE MPTYPE *MOCRHOOKFCTPTR)(long HookType
Definition: milocr.h:32

◆ FindSegment()

EdbSegP * EdbPattern::FindSegment ( int  id)
1303{
1304 for(int i=0; i<N(); i++) if(GetSegment(i)->ID()==id) return GetSegment(i);
1305 return 0;
1306}
Int_t ID() const
Definition: EdbSegP.h:147
UInt_t id
Definition: tlg2couples.C:117

◆ ID()

int EdbPattern::ID ( ) const
inline
319{return eID;}

◆ NAff()

Int_t EdbPattern::NAff ( ) const
inline
314{return eNAff;}

◆ PID()

int EdbPattern::PID ( ) const
inline
320{return ePID;}

◆ Plate()

Int_t EdbPattern::Plate ( ) const
inline
327{return eScanID.ePlate;}
Int_t ePlate
Definition: EdbID.h:11

◆ Reset()

void EdbPattern::Reset ( )
1510{
1511 ((EdbSegmentsBox *)this)->Reset();
1512 Set0();
1513 if(eCell) eCell->Drop();
1514}
Definition: EdbPattern.h:27

◆ ScanID()

EdbID EdbPattern::ScanID ( ) const
inline
326{return eScanID;}

◆ Set0()

void EdbPattern::Set0 ( )
1291{
1292 eID = 0;
1293 ePID = 0;
1294 eFlag = 0;
1295 eNAff = 0;
1297 for(int i=0; i<4; i++) eSigma[i]=0;
1298 eSide=0;
1299}

◆ SetID()

void EdbPattern::SetID ( int  id)
inline
309{eID=id;}

◆ SetNAff()

void EdbPattern::SetNAff ( int  n)
inline
311{eNAff=n;}

◆ SetPID()

void EdbPattern::SetPID ( int  pid)
inline
310{ePID=pid;}
int pid[1000]
Definition: m2track.cpp:13

◆ SetScanID()

void EdbPattern::SetScanID ( EdbID  id)
inline
296{eScanID=id;}

◆ SetSegmentsPID()

void EdbPattern::SetSegmentsPID ( )
1456{
1457 int nseg = N();
1458 for(int i=0; i<nseg; i++) {
1459 GetSegment(i)->SetPID(ID()); //PID of the segment must be ID of the pattern!
1460 }
1461}
void SetPID(int pid)
Definition: EdbSegP.h:129

◆ SetSegmentsScanID()

void EdbPattern::SetSegmentsScanID ( EdbID  id)
1465{
1466 for(int i=0; i<N(); i++) GetSegment(i)->SetScanID(id);
1467}
void SetScanID(EdbID id)
Definition: EdbSegP.h:143

◆ SetSide()

void EdbPattern::SetSide ( int  side)
inline
312{eSide=side;}

◆ SetSigma()

void EdbPattern::SetSigma ( float  sx,
float  sy,
float  stx,
float  sty 
)
inline
298 {eSigma[0]=sx; eSigma[1]=sy; eSigma[2]=stx; eSigma[3]=sty;}

◆ SetStep()

void EdbPattern::SetStep ( float  stepx,
float  stepy,
float  steptx,
float  stepty 
)
inline
300 {eStepX=stepx; eStepY=stepy; eStepTX=steptx; eStepTY=stepty;}

◆ Side()

Int_t EdbPattern::Side ( ) const
inline
328{return eSide;}

◆ StepTX()

float EdbPattern::StepTX ( ) const
inline
317{return eStepTX;}

◆ StepTY()

float EdbPattern::StepTY ( ) const
inline
318{return eStepTY;}

◆ StepX()

float EdbPattern::StepX ( ) const
inline
315{return eStepX;}

◆ StepY()

float EdbPattern::StepY ( ) const
inline
316{return eStepY;}

◆ SummaryPath()

float EdbPattern::SummaryPath ( )

calculate the microscope path of the eventual prediction scan of this pattern

1310{
1312 float sum=0,dx,dy;
1313 if(N()<2) return sum;
1314 for(int i=0; i<N()-1; i++ ) {
1315 dx = GetSegment(i+1)->X()-GetSegment(i)->X();
1316 dy = GetSegment(i+1)->Y()-GetSegment(i)->Y();
1317 sum += Sqrt(dx*dx+dy*dy);
1318 }
1319 return sum;
1320}

◆ Xmean()

float EdbPattern::Xmean ( )
1518{
1519 Double_t s=0;
1520 int n=N(); if(n<1) return 0;
1521 for(int i=0; i<n; i++) s += GetSegment(i)->eX;
1522 s /= n;
1523 return (float)s;
1524}
Float_t eX
Definition: EdbSegP.h:31

◆ Ymean()

float EdbPattern::Ymean ( )
1528{
1529 Double_t s=0;
1530 int n=N(); if(n<1) return 0;
1531 for(int i=0; i<n; i++) s += GetSegment(i)->eY;
1532 s /= n;
1533 return (float)s;
1534}
Float_t eY
Definition: EdbSegP.h:31

Member Data Documentation

◆ eCell

TIndexCell* EdbPattern::eCell
private

associated with eSegments

◆ eFlag

Int_t EdbPattern::eFlag
private

pattern flag

◆ eID

Int_t EdbPattern::eID
private

pattern id in the volume

◆ eNAff

Int_t EdbPattern::eNAff
private

number of segments selected for affine calculation

◆ ePID

Int_t EdbPattern::ePID
private

correspond to the piece ID

◆ eScanID

EdbID EdbPattern::eScanID
private

main scanning ID for this pattern, defined when possible

◆ eSide

Int_t EdbPattern::eSide
private

emulsion side 0/1/2, defined when possible

◆ eSigma

Float_t EdbPattern::eSigma[4]
private

accuracy in comparison to neibour pattern

◆ eStepTX

Float_t EdbPattern::eStepTX
private

◆ eStepTY

Float_t EdbPattern::eStepTY
private

bin size for the index cell

◆ eStepX

Float_t EdbPattern::eStepX
private

◆ eStepY

Float_t EdbPattern::eStepY
private

bin size for the index cell


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