CMS 3D CMS Logo

PFTrajectoryPoint.cc
Go to the documentation of this file.
2 #include <ostream>
3 
4 using namespace reco;
5 
7  isTrackerLayer_(false),
8  detId_(-1),
9  layer_(-1) {}
10 
11 
13  int layer,
14  const math::XYZPoint& posxyz,
17  detId_(detId),
18  layer_(layer),
19  posxyz_(posxyz),
20  momentum_(momentum)
21 {
22  if (detId) isTrackerLayer_ = true;
23  posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi());
24 }
25 
26 
29  detId_(other.detId_),
30  layer_(other.layer_),
31  posxyz_(other.posxyz_),
32  posrep_(other.posrep_),
33  momentum_(other.momentum_) { }
34 
35 
37 {}
38 
39 
41  if( posxyz_ == other.posxyz_ &&
42  momentum_ == other.momentum_ ) return true;
43  else return false;
44 }
45 
46 std::ostream& reco::operator<<(std::ostream& out,
47  const reco::PFTrajectoryPoint& trajPoint) {
48  if(!out) return out;
49 
50  const math::XYZPoint& posxyz = trajPoint.position();
51 
52  out<<"Traj point id = "<<trajPoint.detId()
53  <<", layer = "<<trajPoint.layer()
54  <<", Eta,Phi = "<<posxyz.Eta()<<","<<posxyz.Phi()
55  <<", X,Y = "<<posxyz.X()<<","<<posxyz.Y()
56  <<", R,Z = "<<posxyz.Rho()<<","<<posxyz.Z()
57  <<", E,Pt = "<<trajPoint.momentum().E()<<","<<trajPoint.momentum().Pt();
58 
59  return out;
60 }
static const char layer_[]
bool operator==(const reco::PFTrajectoryPoint &other) const
int detId_
detid if measurement is corresponding to a tracker layer
const math::XYZPoint & position() const
cartesian position (x, y, z)
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
math::XYZPoint posxyz_
cartesian position (x, y, z)
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
int layer() const
trajectory point layer
bool isTrackerLayer_
Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a...
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
int layer_
propagated layer
PFTrajectoryPoint()
default constructor. Set variables at default dummy values
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
math::XYZTLorentzVector momentum_
momentum quadrivector
fixed size matrix
virtual ~PFTrajectoryPoint()
destructor
A PFTrack holds several trajectory points, which basically contain the position and momentum of a tra...
int detId() const
measurement detId
REPPoint posrep_
position in (rho, eta, phi) base (transient)