FEDRA emulsion software from the OPERA Collaboration
EdbShowAlg_RC Class Reference

#include <EdbShowAlg.h>

Inheritance diagram for EdbShowAlg_RC:
Collaboration diagram for EdbShowAlg_RC:

Public Member Functions

 ClassDef (EdbShowAlg_RC, 1)
 
 EdbShowAlg_RC ()
 
void Execute ()
 
void Finalize ()
 
Bool_t GetConeOrTubeDistanceToBTOfShowerArray (EdbSegP *sa, EdbSegP *InBT, EdbShowerP *shower, Double_t CylinderRadius, Double_t ConeAngle)
 
void Initialize ()
 
virtual ~EdbShowAlg_RC ()
 
- Public Member Functions inherited from EdbShowAlg
void AddRecoShowerArray (EdbShowerP *shower)
 
 ClassDef (EdbShowAlg, 1)
 
Double_t DeltaR_NoPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaR_WithoutPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaR_WithPropagation (EdbSegP *s, EdbSegP *stest)
 
Double_t DeltaTheta (EdbSegP *s1, EdbSegP *s2)
 
Double_t DeltaThetaComponentwise (EdbSegP *s1, EdbSegP *s2)
 
Double_t DeltaThetaSingleAngles (EdbSegP *s1, EdbSegP *s2)
 
 EdbShowAlg ()
 
 EdbShowAlg (TString AlgName, Int_t AlgValue)
 
virtual void Execute ()
 
virtual void Finalize ()
 
TString GetAlgName () const
 
Int_t GetAlgValue () const
 
Double_t GetMinimumDist (EdbSegP *seg1, EdbSegP *seg2)
 
TObjArray * GetRecoShowerArray () const
 
Int_t GetRecoShowerArrayN () const
 
EdbShowerPGetShower (Int_t i) const
 
Double_t GetSpatialDist (EdbSegP *s1, EdbSegP *s2)
 
void Help ()
 
virtual void Initialize ()
 
Bool_t IsInConeTube (EdbSegP *sTest, EdbSegP *sStart, Double_t CylinderRadius, Double_t ConeAngle)
 
void Print ()
 
void PrintAll ()
 
void PrintMore ()
 
void PrintParameters ()
 
void PrintParametersShort ()
 
void SetActualAlgParameterset (Int_t ActualAlgParametersetNr)
 
void SetEdbPVRec (EdbPVRec *Ali)
 
void SetEdbPVRecPIDNumbers (Int_t FirstPlate_eAliPID, Int_t LastPlate_eAliPID, Int_t MiddlePlate_eAliPID, Int_t NumberPlate_eAliPID)
 
void SetInBTArray (TObjArray *InBTArray)
 
void SetParameter (Int_t parNr, Float_t par)
 
void SetParameters (Float_t *par)
 
void SetRecoShowerArray (TObjArray *RecoShowerArray)
 
void SetRecoShowerArrayN (Int_t RecoShowerArrayN)
 
void SetUseAliSub (Bool_t UseAliSub)
 
void Transform_eAli (EdbSegP *InitiatorBT, Float_t ExtractSize)
 
void UpdateShowerIDs ()
 
void UpdateShowerMetaData ()
 
virtual ~EdbShowAlg ()
 

Private Member Functions

Bool_t FindPrecedingBTs (EdbSegP *s, EdbSegP *InBT, EdbPVRec *gAli, EdbShowerP *shower)
 

Additional Inherited Members

- Protected Member Functions inherited from EdbShowAlg
void Set0 ()
 
- Protected Attributes inherited from EdbShowAlg
Int_t eActualAlgParametersetNr
 
TString eAlgName
 
Int_t eAlgValue
 
EdbPVReceAli
 
EdbPVReceAli_Sub
 
Int_t eAli_SubNpat
 
Int_t eAliNpat
 
Int_t eFirstPlate_eAliPID
 
TObjArray * eInBTArray
 
Int_t eInBTArrayN
 
Int_t eLastPlate_eAliPID
 
Int_t eMiddlePlate_eAliPID
 
Int_t eNumberPlate_eAliPID
 
Int_t eParaN
 
TString eParaString [10]
 
Float_t eParaValue [10]
 
EdbShowerPeRecoShower
 
TObjArray * eRecoShowerArray
 
Int_t eRecoShowerArrayN
 
Int_t eUseAliSub
 

Constructor & Destructor Documentation

◆ EdbShowAlg_RC()

EdbShowAlg_RC::EdbShowAlg_RC ( )
1555{
1556// Default Constructor
1557 cout << "EdbShowAlg_RC::EdbShowAlg_RC() Default Constructor Recursive Cone."<<endl;
1558
1559// Reset all:
1560 Set0();
1561
1562 eAlgName="RC";
1563 eAlgValue=7; // see default.par_SHOWREC for labeling (labeling identical with ShowRec program)
1564 eParaN =4; // depends on the algorithm, how much parameters it has got.
1565
1566 eParaValue[0]=700.0;
1567 eParaString[0]="CylinderRadius"; // micrometer
1568 eParaValue[1]=0.025;
1569 eParaString[1]="ConeAngle"; // radiant
1570 eParaValue[2]=150.0;
1571 eParaString[2]="DeltaRToPreceedingBT"; // micrometer
1572 eParaValue[3]=0.13;
1573 eParaString[3]="DeltaThetaComponentwiseToPreceedingBT"; // radiant
1574}
Float_t eParaValue[10]
Definition: EdbShowAlg.h:52
TString eParaString[10]
Definition: EdbShowAlg.h:53
TString eAlgName
Definition: EdbShowAlg.h:49
void Set0()
Definition: EdbShowAlg.cxx:53
Int_t eAlgValue
Definition: EdbShowAlg.h:50
Int_t eParaN
Definition: EdbShowAlg.h:51

◆ ~EdbShowAlg_RC()

EdbShowAlg_RC::~EdbShowAlg_RC ( )
virtual
1579{
1580 cout << "EdbShowAlg_RC::~EdbShowAlg_RC()"<<endl;
1581// Default Destructor
1582}

Member Function Documentation

◆ ClassDef()

EdbShowAlg_RC::ClassDef ( EdbShowAlg_RC  ,
 
)

◆ Execute()

void EdbShowAlg_RC::Execute ( )
virtual

Reimplemented from EdbShowAlg.

1594{
1595 cout << "EdbShowAlg_RC::Execute()" << endl;
1596 cout << "EdbShowAlg_RC::Execute DOING MAIN SHOWER RECONSTRUCTION HERE" << endl;
1597
1598 EdbSegP* InBT;
1599 EdbSegP* Segment;
1600 EdbShowerP* RecoShower;
1601
1602 Bool_t StillToLoop=kTRUE;
1603 Int_t ActualPID;
1604 Int_t newActualPID;
1605 Int_t STEP=-1;
1606 Int_t NLoopedPattern=0;
1608 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_SA::Execute--- STEP for patternloop direction = " << STEP << endl;
1609
1610//--- Loop over InBTs:
1611 if (gEDBDEBUGLEVEL>2) cout << "EdbShowAlg_RC::Execute Loop over InBTs:" << endl;
1612
1613
1614
1615// Since eInBTArray is filled in ascending ordering by zpositon
1616// We use the descending loop to begin with BT with lowest z first.
1617 for (Int_t i=eInBTArrayN-1; i>=0; --i) {
1618
1619// CounterOutPut
1620 if (gEDBDEBUGLEVEL==2) if ((i%100)==0) cout << eInBTArrayN <<" InBT in total, still to do:"<<Form("%4d",i)<< "\r\r\r\r"<<flush;
1621
1622//-----------------------------------
1623// 1) Make local_gAli with cut parameters:
1624//-----------------------------------
1625
1626// Create new EdbShowerP Object for storage;
1627// See EdbShowerP why I have to call the Constructor as "unique" ideficable value
1628 RecoShower = new EdbShowerP(i,eAlgValue,eActualAlgParametersetNr);
1629
1630// Get InitiatorBT from eInBTArray
1631 InBT=(EdbSegP*)eInBTArray->At(i);
1632 if (gEDBDEBUGLEVEL>2) InBT->PrintNice();
1633 Float_t X0=InBT->X();
1634 Float_t Y0=InBT->Y();
1635
1636// Add InBT to RecoShower:
1637// This has to be done, since by definition the first BT in the RecoShower is the InBT.
1638// Otherwise, also the definition of shower axis and transversal profiles is wrong!
1639 RecoShower -> AddSegment(InBT);
1640 if (gEDBDEBUGLEVEL>4) cout << "Segment (InBT) " << Segment << " was added to RecoShower." << endl;
1641
1642// Transform (make size smaller, extract only events having same MC) the eAli object:
1643 Transform_eAli(InBT);
1644 if (gEDBDEBUGLEVEL>3) eAli_Sub->Print();
1645
1646//-----------------------------------
1647// 2) Loop over (whole) eAli, check BT for Cuts
1648// eAli_Sub
1649//-----------------------------------
1650 ActualPID= InBT->PID() ;
1651 newActualPID= InBT->PID() ;
1652
1653 while (StillToLoop) {
1654 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- --- Doing patterloop " << ActualPID << " for patterns Z position=" << eAli_Sub->GetPattern(ActualPID)->Z() << endl;
1655
1656 for (Int_t btloop_cnt=0; btloop_cnt<eAli_Sub->GetPattern(ActualPID)->GetN(); ++btloop_cnt) {
1657
1658 Segment = (EdbSegP*)eAli_Sub->GetPattern(ActualPID)->GetSegment(btloop_cnt);
1659 if (gEDBDEBUGLEVEL>4) Segment->PrintNice();
1660
1661
1662// Now apply cut conditions: RC RECURSIVE CONE Alg --------------------
1663 if (Segment->MCEvt() > 0 ) if ( Segment->MCEvt()!=InBT->MCEvt() ) continue; // MCEvtNr (>0) or BgMCNr (-999)
1664 if ( Abs(Segment->X()-X0) > 7000 ) continue;
1665 if ( Abs(Segment->Y()-Y0) > 7000 ) continue;
1666 if (!GetConeOrTubeDistanceToBTOfShowerArray(Segment, InBT, RecoShower, eParaValue[0], eParaValue[1])) continue;
1667 if (!FindPrecedingBTs(Segment, InBT, eAli_Sub, RecoShower)) continue;
1668// end of cut conditions: RC RECURSIVE CONE Alg --------------------
1669
1670
1671// If we arrive here, Basetrack Segment has passed criteria
1672// and is then added to the RecoShower:
1673// Check if its not the InBT which is already added:
1674 if (Segment->X()==InBT->X()&&Segment->Y()==InBT->Y()) {
1675 ; // is InBT, do nothing;
1676 }
1677 else {
1678 RecoShower -> AddSegment(Segment);
1679 }
1680 if (gEDBDEBUGLEVEL>4) cout << "Segment " << Segment << " was added at &Segment : " << &Segment << endl;
1681
1682 } // of btloop_cnt
1683
1684 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- --- ActualPID= " << newActualPID << " done. Reconstructed shower has up to now: " << RecoShower->N() << " Segments." << endl;
1685
1686//------------
1687 newActualPID=ActualPID+STEP;
1688 ++NLoopedPattern;
1689
1690 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- --- newActualPID= " << newActualPID << endl;
1691 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- --- NLoopedPattern= " << NLoopedPattern << endl;
1692 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- --- eNumberPlate_eAliPID= " << eNumberPlate_eAliPID << endl;
1693 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- --- StillToLoop= " << StillToLoop << endl;
1694
1695// This if holds in the case of STEP== +1
1696 if (STEP==1) {
1697 if (newActualPID>eLastPlate_eAliPID) StillToLoop=kFALSE;
1698 if (newActualPID>eLastPlate_eAliPID) cout << "EdbShowAlg_RC::Execute--- ---Stop Loop since: newActualPID>eLastPlate_eAliPID"<<endl;
1699 }
1700// This if holds in the case of STEP== -1
1701 if (STEP==-1) {
1702 if (newActualPID<eLastPlate_eAliPID) StillToLoop=kFALSE;
1703 if (newActualPID<eLastPlate_eAliPID && gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- ---Stop Loop since: newActualPID<eLastPlate_eAliPID"<<endl;
1704 }
1705// This if holds general, since eNumberPlate_eAliPID is not dependent of the structure of the gAli subject:
1706 if (NLoopedPattern>eNumberPlate_eAliPID) StillToLoop=kFALSE;
1707 if (NLoopedPattern>eNumberPlate_eAliPID && gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- ---Stop Loop since: NLoopedPattern>eNumberPlate_eAliPID"<<endl;
1708
1709 ActualPID=newActualPID;
1710 } // of // while (StillToLoop)
1711
1712// Obligatory when Shower Reconstruction is finished:
1713 RecoShower ->Update();
1714
1715 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintBasics();
1716 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintNice();
1717 if (gEDBDEBUGLEVEL>3) RecoShower ->PrintSegments();
1718
1719// if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- Before adding to array delete the histograms...by finalize() of shower."<<endl;
1720// RecoShower ->Finalize();
1721// if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::Execute--- Finalize done..."<<endl;
1722
1723
1724// Add Shower Object to Shower Reco Array.
1725// Not, if its empty or containing only one BT:
1726 if (RecoShower->N()>1) eRecoShowerArray->Add(RecoShower);
1727 if (RecoShower->N()>1) {
1728 AddRecoShowerArray(RecoShower);
1729 }
1730 else {
1731 // Delete the created shower object. Saves _lots_ of memory and speeds up reconstruction for large #InBTs!
1732 // TODO Also Include this in the other RecoAlgs ...
1733 Log(3,"EdbShowAlg_RC::Execute()","InBT # %d RecoShower->N()<=1: Dont add shower to RecoShowerArray. Delete this shower for memory safing.");
1734 delete RecoShower;
1735 RecoShower=NULL;
1736 }
1737
1738
1739// Set back loop values:
1740 StillToLoop=kTRUE;
1741 NLoopedPattern=0;
1742 } // of // for (Int_t i=eInBTArrayN-1; i>=0; --i) {
1743
1744
1745// Set new value for eRecoShowerArrayN (may now be < eInBTArrayN).
1747
1748 cout << "EdbShowAlg_RC::eRecoShowerArray() Entries: " << eRecoShowerArray->GetEntries() << endl;
1749 cout << "EdbShowAlg_RC::Execute()...done." << endl;
1750 return;
1751}
bool Log(int level, const char *location, const char *fmt,...)
Definition: EdbLog.cxx:75
brick X0
Definition: RecDispMC.C:112
void Print() const
Definition: EdbPattern.cxx:1693
EdbPattern * GetPattern(int id) const
Definition: EdbPattern.cxx:1721
Definition: EdbSegP.h:21
Float_t X() const
Definition: EdbSegP.h:173
Float_t Y() const
Definition: EdbSegP.h:174
void PrintNice() const
Definition: EdbSegP.cxx:392
Int_t PID() const
Definition: EdbSegP.h:148
Int_t MCEvt() const
Definition: EdbSegP.h:145
Float_t Z() const
Definition: EdbPattern.h:84
EdbSegP * GetSegment(int i) const
Definition: EdbPattern.h:66
Bool_t FindPrecedingBTs(EdbSegP *s, EdbSegP *InBT, EdbPVRec *gAli, EdbShowerP *shower)
Definition: EdbShowAlg.cxx:1772
Bool_t GetConeOrTubeDistanceToBTOfShowerArray(EdbSegP *sa, EdbSegP *InBT, EdbShowerP *shower, Double_t CylinderRadius, Double_t ConeAngle)
Definition: EdbShowAlg.cxx:1829
void Transform_eAli(EdbSegP *InitiatorBT, Float_t ExtractSize)
Definition: EdbShowAlg.cxx:107
void AddRecoShowerArray(EdbShowerP *shower)
Definition: EdbShowAlg.h:120
void SetRecoShowerArrayN(Int_t RecoShowerArrayN)
Definition: EdbShowAlg.h:117
TObjArray * eRecoShowerArray
Definition: EdbShowAlg.h:70
Int_t eFirstPlate_eAliPID
Definition: EdbShowAlg.h:65
Int_t eLastPlate_eAliPID
Definition: EdbShowAlg.h:66
Int_t eInBTArrayN
Definition: EdbShowAlg.h:63
TObjArray * eInBTArray
Definition: EdbShowAlg.h:62
Int_t eNumberPlate_eAliPID
Definition: EdbShowAlg.h:68
EdbPVRec * eAli_Sub
Definition: EdbShowAlg.h:74
Int_t eActualAlgParametersetNr
Definition: EdbShowAlg.h:54
Definition: EdbShowerP.h:28
Int_t N() const
Definition: EdbShowerP.h:412
void PrintNice()
Definition: EdbShowerP.cxx:2339
void PrintSegments()
Definition: EdbShowerP.cxx:2367
void PrintBasics()
Definition: EdbShowerP.cxx:2360
void Update()
Definition: EdbShowerP.cxx:975
gEDBDEBUGLEVEL
Definition: energy.C:7
#define NULL
Definition: nidaqmx.h:84

◆ Finalize()

void EdbShowAlg_RC::Finalize ( )
virtual

Reimplemented from EdbShowAlg.

1757{
1758 return;
1759}

◆ FindPrecedingBTs()

Bool_t EdbShowAlg_RC::FindPrecedingBTs ( EdbSegP s,
EdbSegP InBT,
EdbPVRec gAli,
EdbShowerP shower 
)
private
1773{
1774 if (gEDBDEBUGLEVEL>3) cout << "EdbShowAlg_RC::FindPrecedingBTs Find BTs to be connected with the Test BT:" << endl;
1775
1776 EdbSegP* s_TestBT;
1777 Int_t nentries=shower->N();
1778 Double_t dZ;
1779
1780// In case the shower is empty we do not have to search for a preceeding BT:
1781 if (nentries==0) return kTRUE;
1782
1783// We do not test BaseTracks which are before the initiator BaseTrack!
1784 if (s->Z()<InBT->Z()) {
1785 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_RC::FindPrecedingBTs(): s->Z()<InBT->Z()..: We do not test BaseTracks which are before the initiator BaseTrack! return kFALSE;" << endl;
1786 return kFALSE;
1787 }
1788
1789// For the very first Z position we do not test
1790// if testBT has Preceeders, only if it it has a BT around (case for e+e- coming from gammma):
1791// Take 50microns and 80mrad in (dR/dT) around.
1792// This does not affect the normal results, but helps for
1793// events which may have a second BT close to InBT (like in e+e-)
1794 if (s->Z()==InBT->Z()) {
1795//cout << DeltaThetaComponentwise(s, InBT) << " " << DeltaR_WithPropagation(s, InBT) << endl;
1796 if (DeltaThetaComponentwise(s, InBT) < 0.08 && DeltaR_WithoutPropagation(s, InBT) < 50.0 ) return kTRUE;
1797 }
1798
1799 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_RC::FindPrecedingBTs(): Testing " << nentries << " entries of the shower." << endl;
1800
1801 for (Int_t i=nentries-1; i>=0; --i) {
1802 s_TestBT = (EdbSegP*)( shower->GetSegment(i) );
1803
1804 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_RC::FindPrecedingBTs(): Do s_TestBT->ID() s->ID() s_TestBT->MCEvt() s_TestBT->Z() s->Z() "<< s_TestBT->ID() << " " << s->ID() << " " << s_TestBT->MCEvt() <<" " << s_TestBT->Z() << " " << s->Z() <<endl;
1805 if (gEDBDEBUGLEVEL>3) s_TestBT->PrintNice();
1806 if (gEDBDEBUGLEVEL>3) s->PrintNice();
1807
1808
1809 dZ=TMath::Abs(s_TestBT->Z()-s->Z());
1810 if (dZ<30) continue; // Exclude the case of same Zpositions...
1811 if (dZ>(3*1300.0)+30.0) continue; // Exclude the case of more than 3 plates before...
1812
1813 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_RC::FindPrecedingBTs(): Checking dT,dR and dZ for i: " << i << " " << DeltaThetaComponentwise(s, s_TestBT) << " " << DeltaR_WithPropagation(s, s_TestBT) << " "<<dZ << endl;
1814
1815 if (DeltaThetaComponentwise(s, s_TestBT) > eParaValue[3] ) continue;
1816 if (DeltaR_WithPropagation(s, s_TestBT) > eParaValue[2]) continue;
1817
1818 if (gEDBDEBUGLEVEL>3) cout << "--- --- EdbShowAlg_RC::FindPrecedingBTs(): Checking dT,dR and dZ for i: " << i << " " << DeltaThetaComponentwise(s, s_TestBT) << " " << DeltaR_WithPropagation(s, s_TestBT) << " "<<dZ << " ok!"<<endl;
1819 return kTRUE;
1820 }
1821//---------------------------------------------
1822 return kFALSE;
1823}
Int_t ID() const
Definition: EdbSegP.h:147
Float_t Z() const
Definition: EdbSegP.h:153
Double_t DeltaR_WithoutPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:404
Double_t DeltaThetaComponentwise(EdbSegP *s1, EdbSegP *s2)
Definition: EdbShowAlg.cxx:441
Double_t DeltaR_WithPropagation(EdbSegP *s, EdbSegP *stest)
Definition: EdbShowAlg.cxx:411
EdbSegP * GetSegment(int i) const
Definition: EdbShowerP.h:435
int nentries
Definition: check_shower.C:40
s
Definition: check_shower.C:55

◆ GetConeOrTubeDistanceToBTOfShowerArray()

Bool_t EdbShowAlg_RC::GetConeOrTubeDistanceToBTOfShowerArray ( EdbSegP sa,
EdbSegP InBT,
EdbShowerP shower,
Double_t  CylinderRadius,
Double_t  ConeAngle 
)
1830{
1831 Bool_t isTrueForBT=kFALSE;
1832 Int_t lastI=0;
1833 Double_t factor=1.0;
1834
1835 EdbSegP* s_TestBT;
1836 Int_t nentries=shower->GetNBT();
1837
1838// Now call IsInConeTube for every BT which was reconstructed up to now in the shower:
1839 for (Int_t i=0; i<nentries; i++) {
1840 s_TestBT = (EdbSegP*)( shower->GetSegment(i) );
1841// Dont check the BT which is before or same Z as the BaseTrack(i) position:
1842// But ---not--- for the first Basetrack. There we allow it!
1843//if (s_TestBT->Z()>=sa->Z() && i>0) continue;
1844 if (s_TestBT->Z()>sa->Z() ) continue;
1845 lastI=i;
1846 if (i==0) {
1847 factor=1.0;
1848 }
1849 else {
1850 factor=3.0;
1851 }
1852
1853// Bool_t EdbShowAlg::IsInConeTube(EdbSegP* TestingSegment, EdbSegP* StartingSegment, Double_t CylinderRadius, Double_t ConeAngle)
1854// if (GetConeOrTubeDistanceToInBT(sa, s_TestBT, factor*CylinderRadius, factor*ConeAngle)==kTRUE) {
1855 if (IsInConeTube(sa, s_TestBT, factor*CylinderRadius, factor*ConeAngle)==kTRUE) {
1856 isTrueForBT=kTRUE;
1857 break;
1858 }
1859 }
1860
1861 if (isTrueForBT) {
1862//if (lastI>0) cout <<"for i= " << lastI<< " return true with factor "<< factor << endl;
1863 return kTRUE;
1864 }
1865
1866 return kFALSE;
1867}
Bool_t IsInConeTube(EdbSegP *sTest, EdbSegP *sStart, Double_t CylinderRadius, Double_t ConeAngle)
Definition: EdbShowAlg.cxx:248
Int_t GetNBT() const
Definition: EdbShowerP.h:513

◆ Initialize()

void EdbShowAlg_RC::Initialize ( )
virtual

Reimplemented from EdbShowAlg.

1588{
1589 return;
1590}

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