A PFTrack holds several trajectory points, which basically contain the position and momentum of a track at a given position. More...
#include <PFTrajectoryPoint.h>
Public Types | |
enum | LayerType { ClosestApproach = 0, BeamPipeOrEndVertex = 1, PS1 = 2, PS2 = 3, ECALEntrance = 4, ECALShowerMax = 5, HCALEntrance = 6, HCALExit = 7, HOLayer = 8, NLayers = 9 } |
Define the different layers where the track can be propagated. More... | |
typedef ROOT::Math::PositionVector3D < ROOT::Math::CylindricalEta3D < Double32_t > > | REPPoint |
Public Member Functions | |
void | calculatePositionREP () |
calculate posrep_ once and for all | |
int | detId () const |
measurement detId | |
bool | isTrackerLayer () const |
is this point corresponding to an intersection with a tracker layer ? | |
bool | isValid () const |
is this point valid ? | |
int | layer () const |
trajectory point layer | |
const math::XYZTLorentzVector & | momentum () const |
4-momenta quadrivector | |
bool | operator== (const reco::PFTrajectoryPoint &other) const |
PFTrajectoryPoint () | |
default constructor. Set variables at default dummy values | |
PFTrajectoryPoint (const PFTrajectoryPoint &other) | |
copy | |
PFTrajectoryPoint (int detId, int layer, const math::XYZPoint &posxyz, const math::XYZTLorentzVector &momentum) | |
constructor from values. set detId to -1 if this point is not from a tracker layer | |
const math::XYZPoint & | position () const |
cartesian position (x, y, z) | |
const REPPoint & | positionREP () const |
trajectory position in (rho, eta, phi) base | |
virtual | ~PFTrajectoryPoint () |
destructor | |
Private Attributes | |
int | detId_ |
detid if measurement is corresponding to a tracker layer | |
bool | isTrackerLayer_ |
Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a certain position? | |
int | layer_ |
propagated layer | |
math::XYZTLorentzVector | momentum_ |
momentum quadrivector | |
REPPoint | posrep_ |
position in (rho, eta, phi) base (transient) | |
math::XYZPoint | posxyz_ |
cartesian position (x, y, z) | |
Friends | |
std::ostream & | operator<< (std::ostream &out, const reco::PFTrajectoryPoint &trajPoint) |
A PFTrack holds several trajectory points, which basically contain the position and momentum of a track at a given position.
Definition at line 26 of file PFTrajectoryPoint.h.
typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<Double32_t> > reco::PFTrajectoryPoint::REPPoint |
Definition at line 29 of file PFTrajectoryPoint.h.
Define the different layers where the track can be propagated.
ClosestApproach |
Point of closest approach from beam axis (initial point in the case of PFSimParticle) |
BeamPipeOrEndVertex | |
PS1 |
Preshower layer 1. |
PS2 |
Preshower layer 2. |
ECALEntrance |
ECAL front face. |
ECALShowerMax |
expected maximum of the shower in ECAL, for an e/gamma particle |
HCALEntrance |
HCAL front face. |
HCALExit |
HCAL exit. |
HOLayer |
HO layer. |
NLayers |
Definition at line 32 of file PFTrajectoryPoint.h.
{ ClosestApproach = 0, BeamPipeOrEndVertex = 1, PS1 = 2, PS2 = 3, ECALEntrance = 4, ECALShowerMax = 5, HCALEntrance = 6, HCALExit = 7, HOLayer = 8, NLayers = 9 };
PFTrajectoryPoint::PFTrajectoryPoint | ( | ) |
default constructor. Set variables at default dummy values
Definition at line 6 of file PFTrajectoryPoint.cc.
: isTrackerLayer_(false), detId_(-1), layer_(-1) {}
PFTrajectoryPoint::PFTrajectoryPoint | ( | int | detId, |
int | layer, | ||
const math::XYZPoint & | posxyz, | ||
const math::XYZTLorentzVector & | momentum | ||
) |
constructor from values. set detId to -1 if this point is not from a tracker layer
Definition at line 12 of file PFTrajectoryPoint.cc.
References isTrackerLayer_, posrep_, and posxyz_.
PFTrajectoryPoint::PFTrajectoryPoint | ( | const PFTrajectoryPoint & | other | ) |
copy
Definition at line 27 of file PFTrajectoryPoint.cc.
PFTrajectoryPoint::~PFTrajectoryPoint | ( | ) | [virtual] |
void reco::PFTrajectoryPoint::calculatePositionREP | ( | ) | [inline] |
int reco::PFTrajectoryPoint::detId | ( | ) | const [inline] |
measurement detId
Definition at line 73 of file PFTrajectoryPoint.h.
References detId_.
Referenced by reco::operator<<().
{ return detId_; }
bool reco::PFTrajectoryPoint::isTrackerLayer | ( | ) | const [inline] |
is this point corresponding to an intersection with a tracker layer ?
Definition at line 85 of file PFTrajectoryPoint.h.
References detId_.
Referenced by reco::PFTrack::addPoint().
{ if(detId_ >= 0 ) return true; else return false; }
bool reco::PFTrajectoryPoint::isValid | ( | void | ) | const [inline] |
is this point valid ?
Definition at line 79 of file PFTrajectoryPoint.h.
References detId_, and layer_.
Referenced by PFElecTkProducer::isSharingEcalEnergyWithEgSC(), PFBlockAlgo::link(), PFElecTkProducer::minTangDist(), reco::PFBlockElementBrem::PFBlockElementBrem(), reco::PFBlockElementGsfTrack::PFBlockElementGsfTrack(), reco::PFBlockElementTrack::PFBlockElementTrack(), PFConversionProducer::produce(), PFTrackProducer::produce(), ConvBremPFTrackFinder::runConvBremFinder(), KDTreeLinkerTrackHcal::searchLinks(), KDTreeLinkerTrackEcal::searchLinks(), LinkByRecHit::testTrackAndClusterByRecHit(), and PFBlockAlgo::testTrackAndPS().
int reco::PFTrajectoryPoint::layer | ( | ) | const [inline] |
trajectory point layer
Definition at line 76 of file PFTrajectoryPoint.h.
References layer_.
Referenced by reco::operator<<().
{ return layer_; }
const math::XYZTLorentzVector& reco::PFTrajectoryPoint::momentum | ( | ) | const [inline] |
4-momenta quadrivector
Definition at line 102 of file PFTrajectoryPoint.h.
References momentum_.
Referenced by CalibratableTest::analyze(), PFRootEventManager::closestParticle(), PFRootEventManager::fillOutEventWithSimParticles(), EFilter::filter(), DisplayManager::loadGGsfRecTracks(), DisplayManager::loadGRecTracks(), DisplayManager::loadGSimParticles(), reco::operator<<(), ConvBremPFTrackFinder::runConvBremFinder(), KDTreeLinkerTrackEcal::searchLinks(), PFRootEventManager::tauBenchmark(), and LinkByRecHit::testTrackAndClusterByRecHit().
{ return momentum_; }
bool PFTrajectoryPoint::operator== | ( | const reco::PFTrajectoryPoint & | other | ) | const |
const math::XYZPoint& reco::PFTrajectoryPoint::position | ( | ) | const [inline] |
cartesian position (x, y, z)
Definition at line 91 of file PFTrajectoryPoint.h.
References posxyz_.
Referenced by PFRootEventManager::closestParticle(), PFRootEventManager::fillOutEventWithSimParticles(), PFBlockAlgo::link(), DisplayManager::loadGSimParticles(), reco::operator<<(), reco::PFBlockElementBrem::PFBlockElementBrem(), reco::PFBlockElementGsfTrack::PFBlockElementGsfTrack(), reco::PFBlockElementTrack::PFBlockElementTrack(), KDTreeLinkerTrackEcal::searchLinks(), LinkByRecHit::testTrackAndClusterByRecHit(), and PFBlockAlgo::testTrackAndPS().
{ return posxyz_; }
const REPPoint& reco::PFTrajectoryPoint::positionREP | ( | ) | const [inline] |
trajectory position in (rho, eta, phi) base
Definition at line 94 of file PFTrajectoryPoint.h.
References posrep_.
Referenced by CalibratableTest::analyze(), CalibratableTest::findCandidatesInDeltaR(), PFBlockAlgo::link(), PFElecTkProducer::minTangDist(), KDTreeLinkerTrackHcal::searchLinks(), KDTreeLinkerTrackEcal::searchLinks(), LinkByRecHit::testTrackAndClusterByRecHit(), and KDTreeLinkerTrackHcal::updatePFBlockEltWithLinks().
{ return posrep_; }
std::ostream& operator<< | ( | std::ostream & | out, |
const reco::PFTrajectoryPoint & | trajPoint | ||
) | [friend] |
int reco::PFTrajectoryPoint::detId_ [private] |
detid if measurement is corresponding to a tracker layer
Definition at line 115 of file PFTrajectoryPoint.h.
Referenced by detId(), isTrackerLayer(), and isValid().
bool reco::PFTrajectoryPoint::isTrackerLayer_ [private] |
Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a certain position?
Definition at line 112 of file PFTrajectoryPoint.h.
Referenced by PFTrajectoryPoint().
int reco::PFTrajectoryPoint::layer_ [private] |
propagated layer
Definition at line 118 of file PFTrajectoryPoint.h.
momentum quadrivector
Definition at line 127 of file PFTrajectoryPoint.h.
Referenced by momentum(), and operator==().
REPPoint reco::PFTrajectoryPoint::posrep_ [private] |
position in (rho, eta, phi) base (transient)
Definition at line 124 of file PFTrajectoryPoint.h.
Referenced by calculatePositionREP(), PFTrajectoryPoint(), and positionREP().
cartesian position (x, y, z)
Definition at line 121 of file PFTrajectoryPoint.h.
Referenced by calculatePositionREP(), operator==(), PFTrajectoryPoint(), and position().