File indexing completed on 2024-04-06 12:04:52
0001 #ifndef DataFormats_ParticleFlowReco_PFTrajectoryPoint_h
0002 #define DataFormats_ParticleFlowReco_PFTrajectoryPoint_h
0003
0004
0005 #include <vector>
0006 #include <map>
0007 #include <iosfwd>
0008
0009 #include "DataFormats/Math/interface/Point3D.h"
0010 #include "Rtypes.h"
0011 #include "DataFormats/Math/interface/LorentzVector.h"
0012 #include "Math/GenVector/PositionVector3D.h"
0013
0014 namespace reco {
0015
0016
0017
0018
0019
0020
0021
0022
0023
0024
0025
0026 class PFTrajectoryPoint {
0027 public:
0028
0029
0030
0031 typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > REPPoint;
0032
0033
0034 enum LayerType {
0035
0036 ClosestApproach = 0,
0037 BeamPipeOrEndVertex = 1,
0038
0039 PS1 = 2,
0040
0041 PS2 = 3,
0042
0043 ECALEntrance = 4,
0044
0045
0046 ECALShowerMax = 5,
0047
0048 HCALEntrance = 6,
0049
0050 HCALExit = 7,
0051
0052 HOLayer = 8,
0053
0054 VFcalEntrance = 9,
0055
0056 NLayers = 10,
0057
0058 Unknown = -1
0059 };
0060
0061 static const std::array<std::string, NLayers> layerTypeNames;
0062 static LayerType layerTypeByName(const std::string& name);
0063
0064
0065 PFTrajectoryPoint();
0066
0067
0068
0069 PFTrajectoryPoint(int detId, int layer, const math::XYZPoint& posxyz, const math::XYZTLorentzVector& momentum);
0070
0071
0072 int detId() const { return detId_; }
0073
0074
0075 int layer() const { return layer_; }
0076
0077
0078 bool isValid() const {
0079 if (layer_ == -1 && detId_ == -1)
0080 return false;
0081 else
0082 return true;
0083 }
0084
0085
0086 bool isTrackerLayer() const {
0087 if (detId_ >= 0)
0088 return true;
0089 else
0090 return false;
0091 }
0092
0093
0094 const math::XYZPoint& position() const { return posxyz_; }
0095
0096
0097 const REPPoint& positionREP() const { return posrep_; }
0098
0099
0100 void calculatePositionREP() { posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi()); }
0101
0102
0103 const math::XYZTLorentzVector& momentum() const { return momentum_; }
0104
0105 bool operator==(const reco::PFTrajectoryPoint& other) const;
0106
0107 friend std::ostream& operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint);
0108
0109 private:
0110
0111
0112 bool isTrackerLayer_;
0113
0114
0115 int detId_;
0116
0117
0118 int layer_;
0119
0120
0121 math::XYZPoint posxyz_;
0122
0123
0124 REPPoint posrep_;
0125
0126
0127 math::XYZTLorentzVector momentum_;
0128 };
0129
0130 std::ostream& operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint);
0131 }
0132
0133 #endif