CMS 3D CMS Logo

PFTrajectoryPoint.cc
Go to the documentation of this file.
2 #include <ostream>
3 #include <algorithm>
4 
5 using namespace reco;
6 
7 // To be kept in synch with the enumerator definitions in PFTrajectoryPoint.h file
8 // Don't consider "Unknown" and "NLayers"
9 std::array<std::string, PFTrajectoryPoint::NLayers> const PFTrajectoryPoint::layerTypeNames{{"ClosestApproach",
10  "BeamPipeOrEndVertex",
11  "PS1",
12  "PS2",
13  "ECALEntrance",
14  "ECALShowerMax",
15  "HCALEntrance",
16  "HCALExit",
17  "HOLayer",
18  "VFcalEntrance"}};
19 
20 PFTrajectoryPoint::PFTrajectoryPoint() : isTrackerLayer_(false), detId_(-1), layer_(-1) {}
21 
23  int layer,
24  const math::XYZPoint& posxyz,
25  const math::XYZTLorentzVector& momentum)
26  : isTrackerLayer_(false), detId_(detId), layer_(layer), posxyz_(posxyz), momentum_(momentum) {
27  if (detId)
28  isTrackerLayer_ = true;
29  posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi());
30 }
31 
33  : isTrackerLayer_(other.isTrackerLayer_),
34  detId_(other.detId_),
36  posxyz_(other.posxyz_),
37  posrep_(other.posrep_),
38  momentum_(other.momentum_) {}
39 
41 
43  auto it = std::find(layerTypeNames.begin(), layerTypeNames.end(), name);
44  if (it == layerTypeNames.end()) {
45  return Unknown; // better this or throw()?
46  }
47  int index = std::distance(layerTypeNames.begin(), it);
48  return LayerType(index);
49 }
50 
52  if (posxyz_ == other.posxyz_ && momentum_ == other.momentum_)
53  return true;
54  else
55  return false;
56 }
57 
58 std::ostream& reco::operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint) {
59  if (!out)
60  return out;
61 
62  const math::XYZPoint& posxyz = trajPoint.position();
63 
64  out << "Traj point id = " << trajPoint.detId() << ", layer = " << trajPoint.layer() << ", Eta,Phi = " << posxyz.Eta()
65  << "," << posxyz.Phi() << ", X,Y = " << posxyz.X() << "," << posxyz.Y() << ", R,Z = " << posxyz.Rho() << ","
66  << posxyz.Z() << ", E,Pt = " << trajPoint.momentum().E() << "," << trajPoint.momentum().Pt();
67 
68  return out;
69 }
bool operator==(const reco::PFTrajectoryPoint &other) const
static const char layer_[]
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
static const std::array< std::string, NLayers > layerTypeNames
int layer() const
trajectory point layer
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
constexpr std::array< uint8_t, layerIndexSize > layer
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:66
bool isTrackerLayer_
Is the measurement corresponding to a tracker layer? or was it obtained by propagating the track to a...
static LayerType layerTypeByName(const std::string &name)
const math::XYZPoint & position() const
cartesian position (x, y, z)
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
LayerType
Define the different layers where the track can be propagated.
REPPoint posrep_
position in (rho, eta, phi) base (transient)