|
ATLAS Offline Software
|
Go to the documentation of this file.
6 #ifndef ANALYSISUTILS_KINEMATICSELECTORCORE_H
7 #define ANALYSISUTILS_KINEMATICSELECTORCORE_H
33 static const InterfaceID IID_KinematicSelectorCore(
"KinematicSelectorCore", 1, 0);
47 static const InterfaceID&
interfaceID() {
return IID_KinematicSelectorCore; };
73 if ( (std::isinf(
pt) || std::isnan(
pt)) )
return true;
82 if ( (std::isinf(
p) || std::isnan(
p)) )
return true;
91 if ( (std::isinf(
et) || std::isnan(
et)) )
return true;
100 if ( (std::isinf(
e) || std::isnan(
e)) )
return true;
109 if ( (std::isinf(
eta) || std::isnan(
eta)) )
return true;
118 if ( (std::isinf(
eta) || std::isnan(
eta)) )
return true;
131 if ( (std::isinf(
phi) || std::isnan(
phi)) )
return true;
144 if ( (std::isinf(
mass) || std::isnan(
mass)) )
return true;
274 std::vector<double> &
maxVals );
356 #endif // ANALYSISUTILS_KINEMATICSELECTORCORE_H
bool passPt(const INavigable4Momentum *part) const
Only transverse momentum cut is applied.
std::vector< double > m_maxEtaVeto
eta veto ranges upper range boundaries
double m_maxEta
maximum eta cut value
Dual use tool (athena/ARA) for any cuts. This is the base class.
bool passE(const INavigable4Momentum *part) const
Only energy cut is applied.
Extra patterns decribing particle interation process.
double getminMass(void)
get minimum mass cut value
bool passPhi(const INavigable4Momentum *part) const
Only phi cut is applied.
double m_maxPt
maximum transverse momentum cut value
Scalar phi() const
phi method
void setmaxMass(double val)
set maximum mass cut value
std::vector< double > m_maxPhiVeto
phi veto ranges upper range boundaries
bool passPhiVetoRanges(const INavigable4Momentum *part) const
Only the phi veto ranges cut is applied.
std::vector< double > m_minEtaVeto
eta veto ranges lower range boundaries
virtual StatusCode initialize()
Gaudi Service Interface method implementations.
Scalar eta() const
pseudorapidity method
double getmaxAbsEta(void)
get maximum |eta| cut value
bool passP(const INavigable4Momentum *part) const
Only momentum cut is applied.
void setmaxE(double val)
set maximum energy cut value
void setmaxPhi(double val)
set maximum phi cut value
double getminEta(void)
get minimum eta cut value
double getminEt(void)
get minimum transverse energy cut value
void setPhiVetoRanges(std::string &range)
set phi veto ranges
void setmaxP(double val)
set maximum momentum cut value
bool passMass(const INavigable4Momentum *part) const
Only mass cut is applied.
bool passEtaVetoRanges(const INavigable4Momentum *part) const
Only the eta veto ranges cut is applied.
double getminP(void)
get minimum momentum cut value
double getmaxP(void)
get maximum momentum cut value
double getminPhi(void)
get minimum phi cut value
double m_minE
minimum energy cut value
bool accept(const IParticle *part) const
Main method for IParticle, all cuts are applied (const method)
std::string m_etaVetoRanges
eta veto ranges string
KinematicSelectorCore(PropertyMgr *pmgr=0)
Default contructor.
void convertStringRange(std::string &range, std::vector< double > &minVals, std::vector< double > &maxVals)
Helper function.
double m_maxMass
maximum mass cut value
void setminP(double val)
set minimum momentum cut value
double m_minPt
minimum transverse momentum cut value
double getmaxE(void)
get maximum energy cut value
void setminPt(double val)
set minimum transverse momentum cut value
double m_minEt
minimum transverse energy cut value
double getmaxEt(void)
get maximum transverse energy cut value
void setminEta(double val)
set minimum eta cut value
std::string getEtaVetoRanges(void)
get eta veto ranges
std::string getPhiVetoRanges(void)
get phi veto ranges
::StatusCode StatusCode
StatusCode definition for legacy code.
void setminE(double val)
set minimum energy cut value
virtual ~KinematicSelectorCore()
Default destructor.
void setmaxEt(double val)
set maximum transverse energy cut value
void setmaxPt(double val)
set maximum transverse momentum cut value
void setminMass(double val)
set minimum mass cut value
double m_maxAbsEta
maximum |eta| cut value
double getmaxMass(void)
get maximum mass cut value
void setminEt(double val)
set minimum transverse energy cut value
Dual use tool (athena/ARA) for kinematic cuts.
bool passEta(const INavigable4Momentum *part) const
Only eta cut is applied.
double m_minMass
minimum mass cut value
void setminAbsEta(double val)
set minimum |eta| cut value
double m_minAbsEta
minimum |eta| cut value
double m_minPhi
minimum phi cut value
void setminPhi(double val)
set minimum phi cut value
virtual StatusCode finalize()
Gaudi Service Interface method implementations.
void setmaxAbsEta(double val)
set maximum |eta| cut value
bool passEt(const INavigable4Momentum *part) const
Only transverse energy cut is applied.
double getminPt(void)
get minimum transverse momentum cut value
bool passAbsEta(const INavigable4Momentum *part) const
Only |eta| cut is applied.
double getminAbsEta(void)
get minimum |eta| cut value
double getminE(void)
get minimum energy cut value
double getmaxEta(void)
get maximum eta cut value
double m_minEta
minimum eta cut value
double m_maxPhi
maximum phi cut value
double m_maxE
maximum energy cut value
void setEtaVetoRanges(std::string &range)
set eta veto ranges
double m_minP
minimum momentum cut value
std::string m_phiVetoRanges
phi veto ranges string
double m_maxP
maximum momentum cut value
static const InterfaceID & interfaceID()
AlgTool interface methods.
std::vector< double > m_minPhiVeto
phi veto ranges lower range boundaries
double m_maxEt
maximum transverse energy cut value
double getmaxPt(void)
get maximum transverse momentum cut value
void setmaxEta(double val)
set maximum eta cut value
double getmaxPhi(void)
get maximum phi cut value