![Logo](../../ATLAS-Logo-Square-Blue-RGB.png) |
ATLAS Offline Software
|
Go to the documentation of this file.
60 "Eta ranges in which the object should NOT be. Example: [[-1.52, -1.37], [1.37, 1.52]]");
66 "Phi ranges in which the object should NOT be. Example: [[-1.4, -1.2], [1.1, 1.5]]");
141 std::cerr <<
"ERROR! Couldn't cast to INavigable4Momentum!"
142 <<
" But this should always work from an IParticle!"
178 bool isInVetoRange =
false;
181 if ( (std::isinf(
eta) || std::isnan(
eta)) )
return true;
193 isInVetoRange =
true;
198 return !isInVetoRange;
207 bool isInVetoRange =
false;
210 if ( (std::isinf(
phi) || std::isnan(
phi)) )
return true;
222 isInVetoRange =
true;
227 return !isInVetoRange;
238 std::string rangeCopy(
range);
241 for (
unsigned int i=0;
i<rangeCopy.length();
i++ )
243 if ( rangeCopy[
i]==
' ' )
245 rangeCopy.erase(
i,1);
251 std::vector<size_t> openBracketPositions;
252 std::vector<size_t> closedBracketPositions;
253 for (
size_t i=0;
i < rangeCopy.size(); ++
i )
255 if ( rangeCopy.compare(
i, 1,
"[") == 0 )
257 openBracketPositions.push_back(
i );
259 if ( rangeCopy.compare(
i, 1,
"]") == 0 )
261 closedBracketPositions.push_back(
i );
266 if ( openBracketPositions.size() != closedBracketPositions.size() )
272 <<
"Inconsistent number of open and closed brackets! "
273 <<
" numberOpen=" << openBracketPositions.size()
274 <<
" numberClosed=" << closedBracketPositions.size()
275 <<
" and cut string=" << rangeCopy
282 for (
unsigned int i=0;
i<openBracketPositions.size(); ++
i )
284 std::string newRange;
285 newRange = rangeCopy.substr( openBracketPositions[
i]+1,
286 closedBracketPositions[
i] );
287 size_t comma = newRange.find(
",");
288 minVals.push_back( (
double)strtod( (newRange.substr( 0, comma )).c_str(),
290 maxVals.push_back( (
double)strtod( (newRange.substr( comma+1 )).c_str(),
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.
bool passPhi(const INavigable4Momentum *part) const
Only phi cut is applied.
double m_maxPt
maximum transverse momentum cut value
Scalar phi() const
phi method
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
bool passP(const INavigable4Momentum *part) const
Only momentum cut is applied.
void setPhiVetoRanges(std::string &range)
set phi veto ranges
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 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
double m_minPt
minimum transverse momentum cut value
double m_minEt
minimum transverse energy cut value
::StatusCode StatusCode
StatusCode definition for legacy code.
double m_maxAbsEta
maximum |eta| cut value
bool passEta(const INavigable4Momentum *part) const
Only eta cut is applied.
double m_minMass
minimum mass cut value
double m_minAbsEta
minimum |eta| cut value
double m_minPhi
minimum phi cut value
virtual StatusCode finalize()
Gaudi Service Interface method implementations.
bool passEt(const INavigable4Momentum *part) const
Only transverse energy cut is applied.
bool passAbsEta(const INavigable4Momentum *part) const
Only |eta| cut is applied.
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
std::vector< double > m_minPhiVeto
phi veto ranges lower range boundaries
double m_maxEt
maximum transverse energy cut value