Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-06-06 07:38:23

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::PFTrajectoryPoint(const PFTrajectoryPoint& other)
0034     : isTrackerLayer_(other.isTrackerLayer_),
0035       detId_(other.detId_),
0036       layer_(other.layer_),
0037       posxyz_(other.posxyz_),
0038       posrep_(other.posrep_),
0039       momentum_(other.momentum_) {}
0040 
0041 PFTrajectoryPoint::~PFTrajectoryPoint() {}
0042 
0043 PFTrajectoryPoint::LayerType PFTrajectoryPoint::layerTypeByName(const std::string& name) {
0044   auto it = std::find(layerTypeNames.begin(), layerTypeNames.end(), name);
0045   if (it == layerTypeNames.end()) {
0046     return Unknown;  // better this or throw()?
0047   }
0048   int index = std::distance(layerTypeNames.begin(), it);
0049   return LayerType(index);
0050 }
0051 
0052 bool PFTrajectoryPoint::operator==(const reco::PFTrajectoryPoint& other) const {
0053   if (posxyz_ == other.posxyz_ && momentum_ == other.momentum_)
0054     return true;
0055   else
0056     return false;
0057 }
0058 
0059 std::ostream& reco::operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint) {
0060   if (!out)
0061     return out;
0062 
0063   const math::XYZPoint& posxyz = trajPoint.position();
0064 
0065   out << "Traj point id = " << trajPoint.detId() << ", layer = " << trajPoint.layer() << ", Eta,Phi = " << posxyz.Eta()
0066       << "," << posxyz.Phi() << ", X,Y = " << posxyz.X() << "," << posxyz.Y() << ", R,Z = " << posxyz.Rho() << ","
0067       << posxyz.Z() << ", E,Pt = " << trajPoint.momentum().E() << "," << trajPoint.momentum().Pt();
0068 
0069   return out;
0070 }