CMS 3D CMS Logo

PFTrack.cc
Go to the documentation of this file.
2 #include "Math/GenVector/PositionVector3D.h"
5 
6 using namespace reco;
7 using namespace std;
8 
9 const unsigned PFTrack::nMaxTrackingLayers_ = 17;
10 
11 PFTrack::PFTrack() : charge_(0.), indexInnermost_(0), indexOutermost_(0), color_(1) {
12  // prepare vector of trajectory points for propagated positions
14 }
15 
17  // prepare vector of trajectory points for propagated positions
19 }
20 
22  : charge_(other.charge_),
26  color_(other.color_) {}
27 
28 void PFTrack::addPoint(const PFTrajectoryPoint& trajPt) {
29  // cout<<"adding "<<trajPt<<endl;
30 
31  if (trajPt.isTrackerLayer()) {
32  if (!indexOutermost_) { // first time a measurement is added
34  PFTrajectoryPoint dummyPt;
35  for (unsigned iPt = trajectoryPoints_.size(); iPt < PFTrajectoryPoint::BeamPipeOrEndVertex + 1; iPt++)
36  trajectoryPoints_.push_back(dummyPt);
38  // throw an exception here
39  // edm::LogError("PFTrack")<<"trajectoryPoints_.size() is too large = "
40  // <<trajectoryPoints_.size()<<"\n";
41  }
43  } else
45  }
46  // Use push_back instead of insert in order to gain time
47  trajectoryPoints_.push_back(trajPt);
48 
49  // cout<<"adding point "<<*this<<endl;
50 }
51 
53  //for(unsigned i=0; i<trajectoryPoints_.size(); i++) {
54  // trajectoryPoints_[i].calculatePositionREP();
55  //}
56 }
57 
58 const reco::PFTrajectoryPoint& PFTrack::extrapolatedPoint(unsigned layerid) const {
59  const unsigned offset_layerid = nTrajectoryMeasurements() + layerid;
60  if (layerid >= reco::PFTrajectoryPoint::NLayers || offset_layerid >= trajectoryPoints_.size()) {
61  throw cms::Exception("SizeError") << "PFRecTrack::extrapolatedPoint: cannot access " << layerid
62  << " #traj meas = " << nTrajectoryMeasurements()
63  << " #traj points = " << trajectoryPoints_.size() << endl
64  << (*this);
65  // assert(0);
66  }
67  if (layerid < indexInnermost_)
68  return trajectoryPoints_[layerid];
69  else
70  return trajectoryPoints_[offset_layerid];
71 }
72 
73 ostream& reco::operator<<(ostream& out, const PFTrack& track) {
74  if (!out)
75  return out;
76 
78 
79  out << "Track charge = " << track.charge() << ", Pt = " << closestApproach.momentum().Pt()
80  << ", P = " << closestApproach.momentum().P() << endl
81  << "\tR0 = " << closestApproach.position().Rho() << " Z0 = " << closestApproach.position().Z() << endl
82  << "\tnumber of tracker measurements = " << track.nTrajectoryMeasurements() << endl;
83  for (unsigned i = 0; i < track.trajectoryPoints().size(); i++)
84  out << track.trajectoryPoints()[i] << endl;
85 
86  return out;
87 }
int color_
color (transient)
Definition: PFTrack.h:136
const math::XYZPoint & position() const
cartesian position (x, y, z)
Base class for particle flow input reconstructed tracks and simulated particles.
Definition: PFTrack.h:63
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
double charge_
charge
Definition: PFTrack.h:124
bool isTrackerLayer() const
is this point corresponding to an intersection with a tracker layer ?
const std::vector< reco::PFTrajectoryPoint > & trajectoryPoints() const
Definition: PFTrack.h:96
const reco::PFTrajectoryPoint & extrapolatedPoint(unsigned layerid) const
Definition: PFTrack.cc:58
void addPoint(const reco::PFTrajectoryPoint &trajPt)
Definition: PFTrack.cc:28
std::vector< reco::PFTrajectoryPoint > trajectoryPoints_
vector of trajectory points
Definition: PFTrack.h:127
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
const reco::PFTrajectoryPoint & trajectoryPoint(unsigned index) const
Definition: PFTrack.h:99
static const unsigned int nMaxTrackingLayers_
maximal number of tracking layers
Definition: PFTrack.h:121
void calculatePositionREP()
Definition: PFTrack.cc:52
unsigned int nTrajectoryMeasurements() const
Definition: PFTrack.h:91
unsigned int indexOutermost_
index outermost tracker measurement
Definition: PFTrack.h:133
fixed size matrix
unsigned int indexInnermost_
index innermost tracker measurement
Definition: PFTrack.h:130
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
double charge() const
Definition: PFTrack.h:85