CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
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 
12  charge_(0.),
13  indexInnermost_(0),
14  indexOutermost_(0),
15  color_(1) {
16 
17  // prepare vector of trajectory points for propagated positions
19 }
20 
21 
23  charge_(charge),
24  indexInnermost_(0),
25  indexOutermost_(0),
26  color_(1) {
27 
28  // prepare vector of trajectory points for propagated positions
30 }
31 
32 
33 PFTrack::PFTrack(const PFTrack& other) :
34  charge_(other.charge_),
35  trajectoryPoints_(other.trajectoryPoints_),
36  indexInnermost_(other.indexInnermost_),
37  indexOutermost_(other.indexOutermost_),
38  color_(other.color_)
39 {}
40 
41 
42 void PFTrack::addPoint(const PFTrajectoryPoint& trajPt) {
43 
44  // cout<<"adding "<<trajPt<<endl;
45 
46  if (trajPt.isTrackerLayer()) {
47  if (!indexOutermost_) { // first time a measurement is added
49  PFTrajectoryPoint dummyPt;
50  for (unsigned iPt = trajectoryPoints_.size(); iPt < PFTrajectoryPoint::BeamPipeOrEndVertex + 1; iPt++)
51  trajectoryPoints_.push_back(dummyPt);
53  // throw an exception here
54  // edm::LogError("PFTrack")<<"trajectoryPoints_.size() is too large = "
55  // <<trajectoryPoints_.size()<<"\n";
56  }
58  } else
60  }
61  // Use push_back instead of insert in order to gain time
62  trajectoryPoints_.push_back(trajPt);
63 
64  // cout<<"adding point "<<*this<<endl;
65 }
66 
67 
69 
70  for(unsigned i=0; i<trajectoryPoints_.size(); i++) {
71  trajectoryPoints_[i].calculatePositionREP();
72  }
73 }
74 
75 
76 const reco::PFTrajectoryPoint& PFTrack::extrapolatedPoint(unsigned layerid) const {
77 
78  if( layerid >= reco::PFTrajectoryPoint::NLayers ||
79  nTrajectoryMeasurements() + layerid >= trajectoryPoints_.size() ) {
80 
81  // cout<<(*this)<<endl;
82  // cout<<"lid "<<layerid<<" "<<nTrajectoryMeasurements()<<" "<<trajectoryPoints_.size()<<endl;
83 
84  throw cms::Exception("SizeError")<<"PFRecTrack::extrapolatedPoint: cannot access "
85  <<layerid
86  <<" #traj meas = "<<nTrajectoryMeasurements()
87  <<" #traj points = "<<trajectoryPoints_.size()
88  <<endl
89  <<(*this);
90  // assert(0);
91  }
92  if (layerid < indexInnermost_)
93  return trajectoryPoints_[ layerid ];
94  else
95  return trajectoryPoints_[ nTrajectoryMeasurements() + layerid ];
96 }
97 
98 
99 ostream& reco::operator<<(ostream& out,
100  const PFTrack& track) {
101  if (!out) return out;
102 
103  const reco::PFTrajectoryPoint& closestApproach =
105 
106  out<<"Track charge = "<<track.charge()
107  <<", Pt = "<<closestApproach.momentum().Pt()
108  <<", P = "<<closestApproach.momentum().P()<<endl
109  <<"\tR0 = "<<closestApproach.position().Rho()
110  <<" Z0 = "<<closestApproach.position().Z()<<endl
111  <<"\tnumber of tracker measurements = "
112  <<track.nTrajectoryMeasurements()<<endl;
113  for(unsigned i=0; i<track.trajectoryPoints_.size(); i++)
114  out<<track.trajectoryPoints_[i]<<endl;
115 
116 
117  return out;
118 }
int i
Definition: DBlmapReader.cc:9
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
double charge(const std::vector< uint8_t > &Ampls)
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:64
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
std::vector< reco::PFTrajectoryPoint > trajectoryPoints_
vector of trajectory points
Definition: PFTrack.h:136
bool isTrackerLayer() const
is this point corresponding to an intersection with a tracker layer ?
const reco::PFTrajectoryPoint & extrapolatedPoint(unsigned layerid) const
Definition: PFTrack.cc:76
void addPoint(const reco::PFTrajectoryPoint &trajPt)
Definition: PFTrack.cc:42
tuple out
Definition: dbtoconf.py:99
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
const reco::PFTrajectoryPoint & trajectoryPoint(unsigned index) const
Definition: PFTrack.h:102
static const unsigned int nMaxTrackingLayers_
maximal number of tracking layers
Definition: PFTrack.h:130
void calculatePositionREP()
Definition: PFTrack.cc:68
unsigned int nTrajectoryMeasurements() const
Definition: PFTrack.h:94
unsigned int indexOutermost_
index outermost tracker measurement
Definition: PFTrack.h:142
unsigned int indexInnermost_
index innermost tracker measurement
Definition: PFTrack.h:139
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
double charge() const
Definition: PFTrack.h:87