Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:04:53

0001 #include "DataFormats/ParticleFlowReco/interface/PFTrajectoryPoint.h"
0002 #include <array>
0003 #include <ostream>
0004 #include <algorithm>
0005 
0006 using namespace reco;
0007 
0008 // To be kept in synch with the enumerator definitions in PFTrajectoryPoint.h file
0009 // Don't consider "Unknown" and "NLayers"
0010 std::array<std::string, PFTrajectoryPoint::NLayers> const PFTrajectoryPoint::layerTypeNames{{"ClosestApproach",
0011                                                                                              "BeamPipeOrEndVertex",
0012                                                                                              "PS1",
0013                                                                                              "PS2",
0014                                                                                              "ECALEntrance",
0015                                                                                              "ECALShowerMax",
0016                                                                                              "HCALEntrance",
0017                                                                                              "HCALExit",
0018                                                                                              "HOLayer",
0019                                                                                              "VFcalEntrance"}};
0020 
0021 PFTrajectoryPoint::PFTrajectoryPoint() : isTrackerLayer_(false), detId_(-1), layer_(-1) {}
0022 
0023 PFTrajectoryPoint::PFTrajectoryPoint(int detId,
0024                                      int layer,
0025                                      const math::XYZPoint& posxyz,
0026                                      const math::XYZTLorentzVector& momentum)
0027     : isTrackerLayer_(false), detId_(detId), layer_(layer), posxyz_(posxyz), momentum_(momentum) {
0028   if (detId)
0029     isTrackerLayer_ = true;
0030   posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi());
0031 }
0032 
0033 PFTrajectoryPoint::LayerType PFTrajectoryPoint::layerTypeByName(const std::string& name) {
0034   auto it = std::find(layerTypeNames.begin(), layerTypeNames.end(), name);
0035   if (it == layerTypeNames.end()) {
0036     return Unknown;  // better this or throw()?
0037   }
0038   int index = std::distance(layerTypeNames.begin(), it);
0039   return LayerType(index);
0040 }
0041 
0042 std::ostream& reco::operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint) {
0043   if (!out)
0044     return out;
0045 
0046   const math::XYZPoint& posxyz = trajPoint.position();
0047 
0048   out << "Traj point id = " << trajPoint.detId() << ", layer = " << trajPoint.layer() << ", Eta,Phi = " << posxyz.Eta()
0049       << "," << posxyz.Phi() << ", X,Y = " << posxyz.X() << "," << posxyz.Y() << ", R,Z = " << posxyz.Rho() << ","
0050       << posxyz.Z() << ", E,Pt = " << trajPoint.momentum().E() << "," << trajPoint.momentum().Pt();
0051 
0052   return out;
0053 }