Base class for particle flow input reconstructed tracks and simulated particles. More...
#include <PFTrack.h>
Public Member Functions | |
void | addPoint (const reco::PFTrajectoryPoint &trajPt) |
void | calculatePositionREP () |
double | charge () const |
int | color () const |
const reco::PFTrajectoryPoint & | extrapolatedPoint (unsigned layerid) const |
std::vector < reco::PFTrajectoryPoint > ::const_iterator | innermostMeasurement () const |
iterator on innermost tracker measurement | |
unsigned int | nTrajectoryMeasurements () const |
unsigned int | nTrajectoryPoints () const |
std::vector < reco::PFTrajectoryPoint > ::const_iterator | outermostMeasurement () const |
iterator on outermost tracker measurement | |
PFTrack (double charge) | |
PFTrack () | |
PFTrack (const PFTrack &other) | |
void | setColor (int color) |
void | setPoint (unsigned int index, const reco::PFTrajectoryPoint &measurement) |
set a trajectory point | |
const reco::PFTrajectoryPoint & | trajectoryPoint (unsigned index) const |
const std::vector < reco::PFTrajectoryPoint > & | trajectoryPoints () const |
Protected Attributes | |
double | charge_ |
charge | |
int | color_ |
color (transient) | |
unsigned int | indexInnermost_ |
index innermost tracker measurement | |
unsigned int | indexOutermost_ |
index outermost tracker measurement | |
std::vector < reco::PFTrajectoryPoint > | trajectoryPoints_ |
vector of trajectory points | |
Static Protected Attributes | |
static const unsigned int | nMaxTrackingLayers_ = 17 |
maximal number of tracking layers | |
Friends | |
std::ostream & | operator<< (std::ostream &out, const PFTrack &track) |
Base class for particle flow input reconstructed tracks and simulated particles.
A PFTrack contains a vector of PFTrajectoryPoint objects. These points are stored in a vector to benefit from the random access. One must take care of storing the points in the right order, and it might even be necessary to insert dummy points.
For a PFRecTrack, the ordering of the points is the following:
For a PFSimParticle, the ordering of the points is the following.
PFRecTracks and PFSimParticles are created in the PFTrackProducer module.
PFTrack::PFTrack | ( | ) |
Definition at line 11 of file PFTrack.cc.
References reco::PFTrajectoryPoint::NLayers, nMaxTrackingLayers_, and trajectoryPoints_.
: charge_(0.), indexInnermost_(0), indexOutermost_(0), color_(1) { // prepare vector of trajectory points for propagated positions trajectoryPoints_.reserve(PFTrajectoryPoint::NLayers + nMaxTrackingLayers_); }
PFTrack::PFTrack | ( | double | charge | ) |
Definition at line 22 of file PFTrack.cc.
References reco::PFTrajectoryPoint::NLayers, nMaxTrackingLayers_, and trajectoryPoints_.
: charge_(charge), indexInnermost_(0), indexOutermost_(0), color_(1) { // prepare vector of trajectory points for propagated positions trajectoryPoints_.reserve(PFTrajectoryPoint::NLayers + nMaxTrackingLayers_); }
PFTrack::PFTrack | ( | const PFTrack & | other | ) |
Definition at line 33 of file PFTrack.cc.
: charge_(other.charge_), trajectoryPoints_(other.trajectoryPoints_), indexInnermost_(other.indexInnermost_), indexOutermost_(other.indexOutermost_), color_(other.color_) {}
void PFTrack::addPoint | ( | const reco::PFTrajectoryPoint & | trajPt | ) |
add a trajectory measurement
Definition at line 42 of file PFTrack.cc.
References reco::PFTrajectoryPoint::BeamPipeOrEndVertex, indexInnermost_, indexOutermost_, reco::PFTrajectoryPoint::isTrackerLayer(), and trajectoryPoints_.
Referenced by PFTrackTransformer::addPoints(), and PFTrackTransformer::addPointsAndBrems().
{ // cout<<"adding "<<trajPt<<endl; if (trajPt.isTrackerLayer()) { if (!indexOutermost_) { // first time a measurement is added if (trajectoryPoints_.size() < PFTrajectoryPoint::BeamPipeOrEndVertex + 1) { PFTrajectoryPoint dummyPt; for (unsigned iPt = trajectoryPoints_.size(); iPt < PFTrajectoryPoint::BeamPipeOrEndVertex + 1; iPt++) trajectoryPoints_.push_back(dummyPt); } else if (trajectoryPoints_.size() > PFTrajectoryPoint::BeamPipeOrEndVertex + 1) { // throw an exception here // edm::LogError("PFTrack")<<"trajectoryPoints_.size() is too large = " // <<trajectoryPoints_.size()<<"\n"; } indexOutermost_ = indexInnermost_ = PFTrajectoryPoint::BeamPipeOrEndVertex + 1; } else indexOutermost_++; } // Use push_back instead of insert in order to gain time trajectoryPoints_.push_back(trajPt); // cout<<"adding point "<<*this<<endl; }
void PFTrack::calculatePositionREP | ( | ) |
calculate posrep_ once and for all for each point
Definition at line 68 of file PFTrack.cc.
References i, and trajectoryPoints_.
Referenced by PFTrackTransformer::addPointsAndBrems(), reco::GsfPFRecTrack::calculateBremPositionREP(), and ConvBremPFTrackFinder::runConvBremFinder().
{ for(unsigned i=0; i<trajectoryPoints_.size(); i++) { trajectoryPoints_[i].calculatePositionREP(); } }
double reco::PFTrack::charge | ( | ) | const [inline] |
Definition at line 87 of file PFTrack.h.
References charge_.
Referenced by PFRootEventManager::countChargedAndPhotons(), and PFRootEventManagerColin::processHIGH_E_TAUS().
{ return charge_; }
int reco::PFTrack::color | ( | ) | const [inline] |
Definition at line 122 of file PFTrack.h.
References color_.
Referenced by setColor().
{ return color_; }
const reco::PFTrajectoryPoint & PFTrack::extrapolatedPoint | ( | unsigned | layerid | ) | const |
Definition at line 76 of file PFTrack.cc.
References Exception, indexInnermost_, reco::PFTrajectoryPoint::NLayers, nTrajectoryMeasurements(), and trajectoryPoints_.
Referenced by PFRootEventManager::closestParticle(), PFRootEventManager::fillOutEventWithSimParticles(), PFAlgo::isSatelliteCluster(), DisplayManager::loadGSimParticles(), LinkByRecHit::testTrackAndClusterByRecHit(), and PFBlockAlgo::testTrackAndPS().
{ if( layerid >= reco::PFTrajectoryPoint::NLayers || nTrajectoryMeasurements() + layerid >= trajectoryPoints_.size() ) { // cout<<(*this)<<endl; // cout<<"lid "<<layerid<<" "<<nTrajectoryMeasurements()<<" "<<trajectoryPoints_.size()<<endl; throw cms::Exception("SizeError")<<"PFRecTrack::extrapolatedPoint: cannot access " <<layerid <<" #traj meas = "<<nTrajectoryMeasurements() <<" #traj points = "<<trajectoryPoints_.size() <<endl <<(*this); // assert(0); } if (layerid < indexInnermost_) return trajectoryPoints_[ layerid ]; else return trajectoryPoints_[ nTrajectoryMeasurements() + layerid ]; }
std::vector< reco::PFTrajectoryPoint >::const_iterator reco::PFTrack::innermostMeasurement | ( | ) | const [inline] |
iterator on innermost tracker measurement
Definition at line 111 of file PFTrack.h.
References indexInnermost_, and trajectoryPoints_.
{ return trajectoryPoints_.begin() + indexInnermost_; }
unsigned int reco::PFTrack::nTrajectoryMeasurements | ( | ) | const [inline] |
Definition at line 94 of file PFTrack.h.
References indexInnermost_, and indexOutermost_.
Referenced by PFRootEventManager::closestParticle(), and extrapolatedPoint().
{ return (indexOutermost_ ? indexOutermost_ - indexInnermost_ + 1 : 0); }
unsigned int reco::PFTrack::nTrajectoryPoints | ( | ) | const [inline] |
Definition at line 90 of file PFTrack.h.
References trajectoryPoints_.
Referenced by CalibratableTest::analyze(), PFRootEventManager::closestParticle(), and PFRootEventManager::fillOutEventWithSimParticles().
{ return trajectoryPoints_.size(); }
std::vector< reco::PFTrajectoryPoint >::const_iterator reco::PFTrack::outermostMeasurement | ( | ) | const [inline] |
iterator on outermost tracker measurement
Definition at line 116 of file PFTrack.h.
References indexOutermost_, and trajectoryPoints_.
{ return trajectoryPoints_.begin() + indexOutermost_; }
void reco::PFTrack::setColor | ( | int | color | ) | [inline] |
void reco::PFTrack::setPoint | ( | unsigned int | index, |
const reco::PFTrajectoryPoint & | measurement | ||
) | [inline] |
set a trajectory point
Definition at line 78 of file PFTrack.h.
References getHLTprescales::index, and trajectoryPoints_.
{ trajectoryPoints_[index] = measurement; }
const reco::PFTrajectoryPoint& reco::PFTrack::trajectoryPoint | ( | unsigned | index | ) | const [inline] |
Definition at line 102 of file PFTrack.h.
References getHLTprescales::index, and trajectoryPoints_.
Referenced by CalibratableTest::analyze(), CalibratableTest::findCandidatesInDeltaR(), and DisplayManager::loadGSimParticles().
{ return trajectoryPoints_[index]; }
const std::vector< reco::PFTrajectoryPoint >& reco::PFTrack::trajectoryPoints | ( | ) | const [inline] |
Definition at line 98 of file PFTrack.h.
References trajectoryPoints_.
Referenced by DisplayManager::loadGSimParticles(), and PFRootEventManager::trackInsideGCut().
{ return trajectoryPoints_; }
std::ostream& operator<< | ( | std::ostream & | out, |
const PFTrack & | track | ||
) | [friend] |
double reco::PFTrack::charge_ [protected] |
int reco::PFTrack::color_ [protected] |
unsigned int reco::PFTrack::indexInnermost_ [protected] |
index innermost tracker measurement
Definition at line 139 of file PFTrack.h.
Referenced by addPoint(), extrapolatedPoint(), innermostMeasurement(), and nTrajectoryMeasurements().
unsigned int reco::PFTrack::indexOutermost_ [protected] |
index outermost tracker measurement
Definition at line 142 of file PFTrack.h.
Referenced by addPoint(), nTrajectoryMeasurements(), and outermostMeasurement().
const unsigned PFTrack::nMaxTrackingLayers_ = 17 [static, protected] |
std::vector< reco::PFTrajectoryPoint > reco::PFTrack::trajectoryPoints_ [protected] |
vector of trajectory points
Definition at line 136 of file PFTrack.h.
Referenced by addPoint(), calculatePositionREP(), extrapolatedPoint(), innermostMeasurement(), nTrajectoryPoints(), outermostMeasurement(), PFTrack(), setPoint(), trajectoryPoint(), and trajectoryPoints().