FEDRA emulsion software from the OPERA Collaboration
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
EdbPatternsVolume Class Reference

#include <EdbPattern.h>

Inheritance diagram for EdbPatternsVolume:
Collaboration diagram for EdbPatternsVolume:

Public Member Functions

void AddPattern (EdbPattern *pat)
 
void AddPatternAt (EdbPattern *pat, int id)
 
void Centralize ()
 
void Centralize (float xc, float yc)
 
void DropCell ()
 
int DropCouples ()
 
 EdbPatternsVolume ()
 
 EdbPatternsVolume (EdbPatternsVolume &pvol)
 
int FindComplimentsVol (EdbSegP &s, TObjArray &arr, float nsig, float nsigt, int dpat)
 
EdbPatternGetPattern (int id) const
 
EdbPatternGetPatternByPID (int pid) const
 
EdbPatternGetPatternByPlate (int plate, int side)
 
EdbPatternGetPatternByZ (float z, float tolerance=5) const
 EdbPattern* GetPatternByZ(float z) const;. More...
 
EdbPatternGetPatternNext (float z, int dir) const
 
EdbPatternGetPatternPreceding (EdbPattern *pat) const
 
EdbPatternGetPatternSucceding (EdbPattern *pat) const
 
EdbPatternGetPatternZLowestHighest (Bool_t lowestZ=kTRUE) const
 
EdbSegPGetSegment (Long_t vid) const
 
EdbPatternInsertPattern (EdbPattern *pat, Bool_t descendingZ=0)
 
EdbPatternNextPattern (float z, int dir) const
 
Int_t Npatterns () const
 
void PassProperties (EdbPatternsVolume &pvol)
 
Int_t Pid (Long_t vid) const
 
void Print () const
 
void PrintAff () const
 
void PrintStat (EdbPattern &pat) const
 
void PrintStat (Option_t *opt="") const
 
void Set0 ()
 
void SetPatternsID ()
 
void SetXYZ (float x, float y, float z)
 
void Shift (float x, float y)
 
Int_t Sid (Long_t vid) const
 
void SortPatternsByZ (Bool_t descendingZ=0)
 
void Transform (const EdbAffine2D *aff)
 
Long_t Vid (int pid, int sid) const
 
Float_t X () const
 
Float_t Xmean ()
 
Float_t Y () const
 
Float_t Ymean ()
 
Float_t Z () const
 
virtual ~EdbPatternsVolume ()
 

Public Attributes

Bool_t eDescendingZ
 if =0 - z increase in the pattrens array; if =1 - decrease More...
 
TObjArray * ePatterns
 collection of patterns More...
 
TIndexCellePatternsCell
 "pid:id1:chi2:id2" - all found couples More...
 
TIndexCelleTracksCell
 "vidt:vids" - connected segments cell More...
 

Private Attributes

Float_t eX
 
Float_t eY
 
Float_t eZ
 central point More...
 

Constructor & Destructor Documentation

◆ EdbPatternsVolume() [1/2]

EdbPatternsVolume::EdbPatternsVolume ( )
1546{
1547 ePatterns = new TObjArray();
1548 eTracksCell = 0;
1549 ePatternsCell = 0;
1550 Set0();
1551}
TIndexCell * eTracksCell
"vidt:vids" - connected segments cell
Definition: EdbPattern.h:344
TIndexCell * ePatternsCell
"pid:id1:chi2:id2" - all found couples
Definition: EdbPattern.h:345
void Set0()
Definition: EdbPattern.cxx:1585
TObjArray * ePatterns
collection of patterns
Definition: EdbPattern.h:342

◆ EdbPatternsVolume() [2/2]

EdbPatternsVolume::EdbPatternsVolume ( EdbPatternsVolume pvol)
1555{
1556 ePatterns = new TObjArray();
1557 eTracksCell = 0;
1558 ePatternsCell = 0;
1559 Set0();
1560
1561 pvol.PassProperties(*this);
1562 int npat,nseg;
1563 npat = Npatterns();
1564 for(int j=0; j<npat; j++) {
1565 nseg = GetPattern(j)->N();
1566 for(int i=0; i<nseg; i++ ) {
1567 pvol.GetPattern(j)->AddSegment( *(GetPattern(j)->GetSegment(i)) );
1568 }
1569 }
1570}
Int_t npat
Definition: Xi2HatStartScript.C:33
void PassProperties(EdbPatternsVolume &pvol)
Definition: EdbPattern.cxx:1729
EdbSegP * GetSegment(Long_t vid) const
Definition: EdbPattern.h:397
Int_t Npatterns() const
Definition: EdbPattern.h:366
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Int_t N() const
Definition: EdbPattern.h:86
EdbSegP * AddSegment(int i, EdbSegP &s)
Definition: EdbPattern.cxx:72

◆ ~EdbPatternsVolume()

EdbPatternsVolume::~EdbPatternsVolume ( )
virtual
1574{
1575 if(ePatterns) {
1576 ePatterns->Delete();
1577 delete ePatterns;
1578 ePatterns=0;
1579 }
1580 if(eTracksCell) { delete eTracksCell; eTracksCell=0; }
1581 if(ePatternsCell) { delete ePatternsCell; ePatternsCell=0; }
1582}

Member Function Documentation

◆ AddPattern()

void EdbPatternsVolume::AddPattern ( EdbPattern pat)
1708{
1709 ePatterns->Add(pat);
1710}

◆ AddPatternAt()

void EdbPatternsVolume::AddPatternAt ( EdbPattern pat,
int  id 
)
1714{
1715 if(ePatterns->GetSize()<id+1) ePatterns->Expand(id+1);
1716 pat->SetID(id);
1717 ePatterns->AddAt(pat,id);
1718}
void SetID(int id)
Definition: EdbPattern.h:309

◆ Centralize() [1/2]

void EdbPatternsVolume::Centralize ( )

find geometrical center (XY) of all patterns and set it as the center of coordinates to simplify transformations To be used before any operations on patterns

1622{
1626
1627 float xc=0;
1628 float yc=0;
1629 int npat = Npatterns();
1630 for(int i=0; i<npat; i++ ) {
1631 xc += GetPattern(i)->Xmax() + GetPattern(i)->Xmin();
1632 yc += GetPattern(i)->Ymax() + GetPattern(i)->Ymin();
1633 }
1634 xc = xc/Npatterns()/2;
1635 yc = yc/Npatterns()/2;
1636
1637 Centralize(xc,yc);
1638}
void Centralize()
Definition: EdbPattern.cxx:1621
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

◆ Centralize() [2/2]

void EdbPatternsVolume::Centralize ( float  xc,
float  yc 
)
1642{
1643 eX = xc; eY=yc;
1644 Shift(-xc,-yc);
1645 float npat = Npatterns();
1646 for(int i=0; i<npat; i++ )
1647 GetPattern(i)->SetKeep(1,0,0,1,0,0);
1648}
Float_t eX
Definition: EdbPattern.h:338
void Shift(float x, float y)
Definition: EdbPattern.cxx:1748
Float_t eY
Definition: EdbPattern.h:338
virtual void SetKeep(float a11, float a12, float a21, float a22, float b1, float b2)
Definition: EdbVirtual.cxx:143

◆ DropCell()

void EdbPatternsVolume::DropCell ( )
1761{
1762 int npat = Npatterns();
1763 for(int i=0; i<npat; i++ ) GetPattern(i)->Cell()->Drop();
1764}
TIndexCell * Cell() const
Definition: EdbPattern.h:321
void Drop()
Definition: TIndexCell.cpp:237

◆ DropCouples()

int EdbPatternsVolume::DropCouples ( )
1593{
1594 int count=0;
1595 int npat=Npatterns();
1596 for(int i=0; i<npat; i++ )
1597 count += GetPattern(i)->Cell()->DropCouples(4);
1598 if(count) printf("%d couples are dropped in volume cells\n",count);
1599 return count;
1600}
int DropCouples(int level)
Definition: TIndexCell.cpp:206

◆ FindComplimentsVol()

int EdbPatternsVolume::FindComplimentsVol ( EdbSegP s,
TObjArray &  arr,
float  nsig,
float  nsigt,
int  dpat 
)
1939{
1940 EdbPattern *pat = 0;
1941 int npat = Npatterns();
1942 bool p_inverse_z = (GetPattern(1)->Z() - GetPattern(0)->Z()) > 0. ? false : true ;
1943 int ii = 0;
1944 int dpat = 1;
1945 int pid = 0;
1946 if (p_inverse_z) dpat = -1;
1947 for (int i = 0; i<npat; i++)
1948 {
1949 ii = i;
1950 if (p_inverse_z) ii = npat-1-i;
1951 pat = GetPattern(ii);
1952 if (ss.Z() < pat->Z())
1953 {
1954 pid = ii - dpat;
1955 if ( pid >= npat ) ss.SetPID(npat-1);
1956 if ( pid < 0 ) ss.SetPID(0);
1957 else ss.SetPID(pid);
1958 break;
1959 }
1960 else if (ss.Z() == pat->Z())
1961 {
1962 ss.SetPID(ii);
1963 break;
1964 }
1965 }
1966
1967 int p0 = ss.PID();
1968 int pend = p0 + Dpat;
1969 if (pend >= npat) pend = npat - 1;
1970 int pstart = p0 - Dpat;
1971 if (pstart < 0) pstart = 0;
1972
1973 //arr.Clear();
1974 int nseg = 0;
1975 int n = 0;
1976 for(int i=pstart; i <=pend; i++ ) {
1977 pat = GetPattern(i);
1978 if(!pat) continue;
1979 n = pat->FindCompliments(ss,arr,nsig,nsigt);
1980 nseg += n;
1981 }
1982
1983 return nseg;
1984}
Definition: EdbPattern.h:273
int FindCompliments(EdbSegP &s, TObjArray &arr, float nsig, float nsigt)
Definition: EdbPattern.cxx:1354
Float_t Z() const
Definition: EdbPattern.h:84
int pid[1000]
Definition: m2track.cpp:13
ss
Definition: energy.C:62

◆ GetPattern()

EdbPattern * EdbPatternsVolume::GetPattern ( int  id) const
1722{
1723 if(Npatterns()>id) return (EdbPattern*)ePatterns->At(id);
1724 else return 0;
1725}

◆ GetPatternByPID()

EdbPattern * EdbPatternsVolume::GetPatternByPID ( int  pid) const

Return pattern having PID() == pid. If there is no such pattern, return a NULL pattern.

This is important because some codings (libShower for example). rely on getting Patterns by the PIDs.

1880{
1886
1887 EdbPattern *pat=0;
1888 for(int i=0; i<Npatterns(); i++) {
1889 pat = GetPattern(i);
1890 if (GetPattern(i)->PID() == pid ) return pat;
1891 }
1892
1893 if (gEDBDEBUGLEVEL>2)
1894 if (pat == 0) cout << "EdbPattern* EdbPatternsVolume::GetPatternNext() WARNING: pattern is NULL ! Check your input!" << endl;
1895 return NULL;
1896}
gEDBDEBUGLEVEL
Definition: energy.C:7
#define NULL
Definition: nidaqmx.h:84
Definition: Flexmotn.h:102

◆ GetPatternByPlate()

EdbPattern * EdbPatternsVolume::GetPatternByPlate ( int  plate,
int  side 
)
1768{
1769 EdbPattern *p=0;
1770 for(int i=0; i<Npatterns(); i++) {
1771 p = GetPattern(i);
1772 if(p) if( p->Plate()==plate && p->Side()==side ) return p;
1773 }
1774 return 0;
1775}
Int_t plate
Definition: merge_Energy_SytematicSources_Electron.C:1
p
Definition: testBGReduction_AllMethods.C:8

◆ GetPatternByZ()

EdbPattern * EdbPatternsVolume::GetPatternByZ ( float  z,
float  tolerance = 5 
) const

EdbPattern* GetPatternByZ(float z) const;.

Return pattern having Z() == z +- tolerance(microns) due to roundings. If there is no such pattern, return a NULL pattern.

This is important because some codings (libShower for example). rely on getting Patterns by the Z. Warning: If tolerance is to big, there is not a unique solution anymore

1900{
1907
1908 EdbPattern *pat=0;
1909 for(int i=0; i<Npatterns(); i++) {
1910 pat = GetPattern(i);
1911 if ( TMath::Abs(GetPattern(i)->Z() - z)<tolerance ) return pat;
1912 }
1913
1914 if (gEDBDEBUGLEVEL>2)
1915 if (pat == 0) cout << "EdbPattern* EdbPatternsVolume::GetPatternNext() WARNING: pattern is NULL ! Check your input!" << endl;
1916 return NULL;
1917}
Float_t Z() const
Definition: EdbPattern.h:365

◆ GetPatternNext()

EdbPattern * EdbPatternsVolume::GetPatternNext ( float  z,
int  dir 
) const

Return next pattern in either downstream (dir=+1) (i.e. Z gets bigger) or in either upstream (dir=-1) (i.e. Z gets smaller) direction. If the last/first pattern is already given and one ask therefor for a pattern that does not exists, then return a NULL pattern.

1858{
1864
1865 EdbPattern *pat=0;
1866 float dz=dir*99999999.;
1867 float zpat;
1868 for(int i=0; i<Npatterns(); i++) {
1869 zpat=GetPattern(i)->Z();
1870 if(dir>0) if(zpat>z) if(zpat-z<dz) {dz=zpat-z; pat=GetPattern(i);}
1871 if(dir<0) if(zpat<z) if(zpat-z>dz) {dz=zpat-z; pat=GetPattern(i);}
1872 }
1873 if (gEDBDEBUGLEVEL>2)
1874 if (pat == 0) cout << "EdbPattern* EdbPatternsVolume::GetPatternNext() WARNING: pattern is NULL ! Check your input!" << endl;
1875 return pat;
1876}
brick dz
Definition: RecDispMC.C:107

◆ GetPatternPreceding()

EdbPattern * EdbPatternsVolume::GetPatternPreceding ( EdbPattern pat) const
1852{
1853 return GetPatternNext(pat->Z(),-1);
1854}
EdbPattern * GetPatternNext(float z, int dir) const
Definition: EdbPattern.cxx:1857

◆ GetPatternSucceding()

EdbPattern * EdbPatternsVolume::GetPatternSucceding ( EdbPattern pat) const
1846{
1847 return GetPatternNext(pat->Z(),1);
1848}

◆ GetPatternZLowestHighest()

EdbPattern * EdbPatternsVolume::GetPatternZLowestHighest ( Bool_t  lowestZ = kTRUE) const

By default it returns the pattern with the lowest Z position. This is NOT necessarily the first pattern of the patterns volume.

1817{
1820
1821 if (Npatterns()==0) {
1822 Log(1,"EdbPatternsVolume::GetPatternZLowestHighest","ERROR: Attempt to get Pattern, but NO patterns are in EdbPatternsVolume!");
1823 return 0;
1824 }
1825 if (gEDBDEBUGLEVEL>3) cout << "EdbPatternsVolume::GetPatternZLowestHighest lowestZ = " << lowestZ << endl;
1826
1827 EdbPattern *pat=0;
1828 pat=GetPattern(0);
1829 if (Npatterns()==1) return pat;
1830 float zpat0,zpat1;
1831 zpat0=GetPattern(0)->Z();
1832 zpat1=GetPattern(Npatterns()-1)->Z();
1833 if (gEDBDEBUGLEVEL>2) {
1834 cout << "EdbPatternsVolume Pattern At Position 0 has Z = " << zpat0 << endl;
1835 cout << "EdbPatternsVolume Pattern At Position last has Z = " << zpat1 << endl;
1836 }
1837 if (zpat0>zpat1) pat=GetPattern(Npatterns()-1);
1838 if (kTRUE==lowestZ) return pat;
1839 if (kFALSE==lowestZ && zpat0<zpat1) return GetPattern(Npatterns()-1);
1840 return GetPattern(0);
1841}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75

◆ GetSegment()

EdbSegP * EdbPatternsVolume::GetSegment ( Long_t  vid) const
inline
398 {return GetPattern(Pid(vid))->GetSegment( Sid(vid) );}
Int_t Sid(Long_t vid) const
Definition: EdbPattern.h:395
Int_t Pid(Long_t vid) const
Definition: EdbPattern.h:394
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66

◆ InsertPattern()

EdbPattern * EdbPatternsVolume::InsertPattern ( EdbPattern pat,
Bool_t  descendingZ = 0 
)

insert new pattern using it's Z as a main criteria for the positioning

1779{
1781 EdbPattern *p = GetPatternByPlate( pat->Plate(), pat->Side() );
1782 if(p) {
1783 Log(1,"EdbPatternsVolume::InsertPattern",
1784 "ERROR: Attempt to insert new pattern for already existing plate/side z = %f znew = %f",
1785 pat->Plate(), pat->Side(), p->Z(), pat->Z());
1786 return p;
1787 }
1788 AddPattern(pat);
1789 SortPatternsByZ(descendingZ);
1790 return pat;
1791}
Int_t Side() const
Definition: EdbPattern.h:328
Int_t Plate() const
Definition: EdbPattern.h:327
void AddPattern(EdbPattern *pat)
Definition: EdbPattern.cxx:1707
void SortPatternsByZ(Bool_t descendingZ=0)
Definition: EdbPattern.cxx:1794
EdbPattern * GetPatternByPlate(int plate, int side)
Definition: EdbPattern.cxx:1767

◆ NextPattern()

EdbPattern * EdbPatternsVolume::NextPattern ( float  z,
int  dir 
) const

Return next pattern in either downstream (dir=+1) or in either upstream (dir=-1) direction.

1921{
1925
1926 EdbPattern *pat=0;
1927 float dz=dir*99999999.;
1928 float zpat;
1929 for(int i=0; i<Npatterns(); i++) {
1930 zpat=GetPattern(i)->Z();
1931 if(dir>0) if(zpat>z) if(zpat-z<dz) {dz=zpat-z; pat=GetPattern(i);}
1932 if(dir<0) if(zpat<z) if(zpat-z>dz) {dz=zpat-z; pat=GetPattern(i);}
1933 }
1934 return pat;
1935}

◆ Npatterns()

Int_t EdbPatternsVolume::Npatterns ( ) const
inline
366 { if(ePatterns)
367 return ePatterns->GetEntriesFast();
368 else return 0; }

◆ PassProperties()

void EdbPatternsVolume::PassProperties ( EdbPatternsVolume pvol)
1730{
1731 pvol.SetXYZ(eX,eY,eZ);
1732 EdbAffine2D a;
1733 EdbPattern *p=0;
1734
1735 int npat = Npatterns();
1736 for(int i=0; i<npat; i++ ) {
1737 p = GetPattern(i);
1738 p->GetKeep(a);
1739 EdbPattern *psel = new EdbPattern( p->X(),p->Y(),p->Z() );
1740 psel->SetKeep( a.A11(), a.A12(), a.A21(), a.A22(), a.B1(),a.B2() );
1741 pvol.AddPattern(psel);
1742 }
1743 pvol.SetPatternsID();
1744}
void a()
Definition: check_aligned.C:59
Definition: EdbAffine.h:17
void SetXYZ(float x, float y, float z)
Definition: EdbPattern.h:361
void SetPatternsID()
Definition: EdbPattern.cxx:1603
Float_t eZ
central point
Definition: EdbPattern.h:338

◆ Pid()

Int_t EdbPatternsVolume::Pid ( Long_t  vid) const
inline
394{ return vid/1000000; }

◆ Print()

void EdbPatternsVolume::Print ( ) const
1694{
1695 int npat = Npatterns();
1696 printf("\nEdbPatternsVolume with %d patterns\n",npat);
1697 EdbPattern *pat=0;
1698 for(int i=0; i<npat; i++ ) {
1699 pat = GetPattern(i);
1700 printf(" id=%3d pid=%3d x:y:z = %12.3f %12.3f %12.3f n= %8d ScanID = %s \tside=%d\n",
1701 pat->ID(), pat->PID(), pat->X(),pat->Y(),pat->Z(),pat->N(), pat->ScanID().AsString(), pat->Side() );
1702 }
1703 printf("\n");
1704}
char * AsString() const
Definition: EdbID.cxx:26
int PID() const
Definition: EdbPattern.h:320
EdbID ScanID() const
Definition: EdbPattern.h:326
int ID() const
Definition: EdbPattern.h:319
Float_t Y() const
Definition: EdbPattern.h:83
Float_t X() const
Definition: EdbPattern.h:82

◆ PrintAff()

void EdbPatternsVolume::PrintAff ( ) const
1652{
1653 EdbAffine2D a;
1654 int npat = Npatterns();
1655 for(int i=0; i<npat; i++ ) {
1656 GetPattern(i)->GetKeep(a);
1657 printf(" %3d (%5d) Z=%13.3f ",i,GetPattern(i)->NAff(), GetPattern(i)->Z() ); a.Print();
1658 }
1659}
virtual const EdbAffine2D * GetKeep() const
Definition: EdbVirtual.h:176

◆ PrintStat() [1/2]

void EdbPatternsVolume::PrintStat ( EdbPattern pat) const
1689{
1690}

◆ PrintStat() [2/2]

void EdbPatternsVolume::PrintStat ( Option_t *  opt = "") const
1663{
1664 int npat = Npatterns();
1665 printf("\nVolume statistics for %d patterns\n",npat);
1666
1667 float dx,dy;
1668 EdbPattern *pat=0;
1669 printf("pat# \t segments \t dX \t\tdY \t meanDist \n");
1670 for(int i=0; i<npat; i++ ) {
1671 pat = GetPattern(i);
1672 if(pat) {
1673 dx = pat->Xmax() - pat->Xmin();
1674 dy = pat->Ymax()- pat->Ymin();
1675 printf(" %d\t %d\t %10.2f \t %10.2f \t %10.4f \n",
1676 i, pat->GetN(),dx,dy, TMath::Sqrt(dx*dy/pat->GetN()) );
1677 }
1678 }
1679
1680 npat=Npatterns();
1681 for(int i=0; i<npat; i++ ) {
1682 pat = GetPattern(i);
1683 if(pat) if(pat->Cell()) pat->Cell()->PrintStat();
1684 }
1685}
Int_t GetN() const
Definition: EdbPattern.h:65
void PrintStat() const
Definition: TIndexCell.cpp:264

◆ Set0()

void EdbPatternsVolume::Set0 ( )
1586{
1587 eX=0; eY=0; eZ=0;
1588 eDescendingZ=0;
1589}
Bool_t eDescendingZ
if =0 - z increase in the pattrens array; if =1 - decrease
Definition: EdbPattern.h:347

◆ SetPatternsID()

void EdbPatternsVolume::SetPatternsID ( )
1604{
1605 int npat = Npatterns();
1606 for(int i=0; i<npat; i++ )
1607 GetPattern(i)->SetID(i);
1608}

◆ SetXYZ()

void EdbPatternsVolume::SetXYZ ( float  x,
float  y,
float  z 
)
inline
361{ eX=x; eY=y; eZ=z; }

◆ Shift()

void EdbPatternsVolume::Shift ( float  x,
float  y 
)
1749{
1750 EdbAffine2D aff;
1751 aff.ShiftX(x);
1752 aff.ShiftY(y);
1753
1754 int npat = Npatterns();
1755 for(int i=0; i<npat; i++)
1756 GetPattern(i)->Transform(&aff);
1757}
void ShiftX(float d)
Definition: EdbAffine.h:64
void ShiftY(float d)
Definition: EdbAffine.h:65
virtual void Transform(const EdbAffine2D *a)
Definition: EdbVirtual.cxx:155

◆ Sid()

Int_t EdbPatternsVolume::Sid ( Long_t  vid) const
inline
395{ return vid%1000000; }

◆ SortPatternsByZ()

void EdbPatternsVolume::SortPatternsByZ ( Bool_t  descendingZ = 0)
1795{
1796 int n = Npatterns();
1797 int *ind = new int[n];
1798 float *z = new float[n];
1799 for(int i=0; i<Npatterns(); i++) {
1800 ind[i] = i;
1801 z[i] = GetPattern(i)->Z();
1802 }
1803 TMath::Sort(n,z,ind,descendingZ);
1804 TObjArray *pnew = new TObjArray(ePatterns->GetSize());
1805 for(int i=0; i<n; i++) pnew->Add(GetPattern(ind[i]));
1806 ePatterns->Clear();
1807 delete ePatterns;
1808 ePatterns = pnew;
1809 delete ind;
1810 delete z;
1811}

◆ Transform()

void EdbPatternsVolume::Transform ( const EdbAffine2D aff)
1612{
1613 int npat = Npatterns();
1614 for(int i=0; i<npat; i++ ) {
1615 GetPattern(i)->Transform(aff);
1616 GetPattern(i)->TransformARot(aff);
1617 }
1618}
void TransformARot(const EdbAffine2D *affA)
Definition: EdbPattern.cxx:358

◆ Vid()

Long_t EdbPatternsVolume::Vid ( int  pid,
int  sid 
) const
inline
393{ return pid*1000000+sid; }

◆ X()

Float_t EdbPatternsVolume::X ( ) const
inline
363{return eX;}

◆ Xmean()

Float_t EdbPatternsVolume::Xmean ( )
1988{
1989 if (Npatterns()==0) return 0;
1990 EdbPattern* pat = GetPattern(0);
1991 Double_t s=0;
1992 int n=pat->N(); if(n<1) return 0;
1993 for(int i=0; i<n; i++) s += pat->GetSegment(i)->eX;
1994 s /= n;
1995 return (float)s;
1996}
Float_t eX
Definition: EdbSegP.h:31
s
Definition: check_shower.C:55

◆ Y()

Float_t EdbPatternsVolume::Y ( ) const
inline
364{return eY;}

◆ Ymean()

Float_t EdbPatternsVolume::Ymean ( )
2000{
2001 if (Npatterns()==0) return 0;
2002 EdbPattern* pat = GetPattern(0);
2003 Double_t s=0;
2004 int n=pat->N(); if(n<1) return 0;
2005 for(int i=0; i<n; i++) s += pat->GetSegment(i)->eY;
2006 s /= n;
2007 return (float)s;
2008}
Float_t eY
Definition: EdbSegP.h:31

◆ Z()

Float_t EdbPatternsVolume::Z ( ) const
inline
365{return eZ;}

Member Data Documentation

◆ eDescendingZ

Bool_t EdbPatternsVolume::eDescendingZ

if =0 - z increase in the pattrens array; if =1 - decrease

◆ ePatterns

TObjArray* EdbPatternsVolume::ePatterns

collection of patterns

◆ ePatternsCell

TIndexCell* EdbPatternsVolume::ePatternsCell

"pid:id1:chi2:id2" - all found couples

◆ eTracksCell

TIndexCell* EdbPatternsVolume::eTracksCell

"vidt:vids" - connected segments cell

◆ eX

Float_t EdbPatternsVolume::eX
private

◆ eY

Float_t EdbPatternsVolume::eY
private

◆ eZ

Float_t EdbPatternsVolume::eZ
private

central point


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