CMS 3D CMS Logo

PFTrajectoryPoint.cc
Go to the documentation of this file.
2 #include <array>
3 #include <ostream>
4 #include <algorithm>
5 
6 using namespace reco;
7 
8 // To be kept in synch with the enumerator definitions in PFTrajectoryPoint.h file
9 // Don't consider "Unknown" and "NLayers"
10 std::array<std::string, PFTrajectoryPoint::NLayers> const PFTrajectoryPoint::layerTypeNames{{"ClosestApproach",
11  "BeamPipeOrEndVertex",
12  "PS1",
13  "PS2",
14  "ECALEntrance",
15  "ECALShowerMax",
16  "HCALEntrance",
17  "HCALExit",
18  "HOLayer",
19  "VFcalEntrance"}};
20 
21 PFTrajectoryPoint::PFTrajectoryPoint() : isTrackerLayer_(false), detId_(-1), layer_(-1) {}
22 
24  int layer,
25  const math::XYZPoint& posxyz,
26  const math::XYZTLorentzVector& momentum)
27  : isTrackerLayer_(false), detId_(detId), layer_(layer), posxyz_(posxyz), momentum_(momentum) {
28  if (detId)
29  isTrackerLayer_ = true;
30  posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi());
31 }
32 
34  auto it = std::find(layerTypeNames.begin(), layerTypeNames.end(), name);
35  if (it == layerTypeNames.end()) {
36  return Unknown; // better this or throw()?
37  }
38  int index = std::distance(layerTypeNames.begin(), it);
39  return LayerType(index);
40 }
41 
42 std::ostream& reco::operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint) {
43  if (!out)
44  return out;
45 
46  const math::XYZPoint& posxyz = trajPoint.position();
47 
48  out << "Traj point id = " << trajPoint.detId() << ", layer = " << trajPoint.layer() << ", Eta,Phi = " << posxyz.Eta()
49  << "," << posxyz.Phi() << ", X,Y = " << posxyz.X() << "," << posxyz.Y() << ", R,Z = " << posxyz.Rho() << ","
50  << posxyz.Z() << ", E,Pt = " << trajPoint.momentum().E() << "," << trajPoint.momentum().Pt();
51 
52  return out;
53 }
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
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
fixed size matrix
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)