CMS 3D CMS Logo

PFTrajectoryPoint.h
Go to the documentation of this file.
1 #ifndef DataFormats_ParticleFlowReco_PFTrajectoryPoint_h
2 #define DataFormats_ParticleFlowReco_PFTrajectoryPoint_h
3 
5 #include <vector>
6 #include <map>
7 #include <iosfwd>
8 
10 #include "Rtypes.h"
12 #include "Math/GenVector/PositionVector3D.h"
13 
14 namespace reco {
15 
27  public:
28  // Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
29  // which otherwise would make ROOT5 files unreadable in ROOT6. This does not increase
30  // the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
31  typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > REPPoint;
32 
34  enum LayerType {
39  PS1 = 2,
41  PS2 = 3,
50  HCALExit = 7,
52  HOLayer = 8,
55  // Number of valid layers
56  NLayers = 10,
57  // Unknown
58  Unknown = -1
59  };
60 
61  static const std::array<std::string, NLayers> layerTypeNames;
63 
66 
70 
73 
75  virtual ~PFTrajectoryPoint();
76 
78  int detId() const { return detId_; }
79 
81  int layer() const { return layer_; }
82 
84  bool isValid() const {
85  if (layer_ == -1 && detId_ == -1)
86  return false;
87  else
88  return true;
89  }
90 
92  bool isTrackerLayer() const {
93  if (detId_ >= 0)
94  return true;
95  else
96  return false;
97  }
98 
100  const math::XYZPoint& position() const { return posxyz_; }
101 
103  const REPPoint& positionREP() const { return posrep_; }
104 
106  void calculatePositionREP() { posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi()); }
107 
109  const math::XYZTLorentzVector& momentum() const { return momentum_; }
110 
111  bool operator==(const reco::PFTrajectoryPoint& other) const;
112 
113  friend std::ostream& operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint);
114 
115  private:
119 
121  int detId_;
122 
124  int layer_;
125 
128 
131 
134  };
135 
136  std::ostream& operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint);
137 } // namespace reco
138 
139 #endif
bool operator==(const reco::PFTrajectoryPoint &other) const
const REPPoint & positionREP() const
trajectory position in (rho, eta, phi) base
int detId_
detid if measurement is corresponding to a tracker layer
const math::XYZTLorentzVector & momentum() const
4-momenta quadrivector
static const std::array< std::string, NLayers > layerTypeNames
friend std::ostream & operator<<(std::ostream &out, const reco::PFTrajectoryPoint &trajPoint)
bool isTrackerLayer() const
is this point corresponding to an intersection with a tracker layer ?
int layer() const
trajectory point layer
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
bool isValid() const
is this point valid ?
math::XYZPoint posxyz_
cartesian position (x, y, z)
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:66
Point of closest approach from beam axis (initial point in the case of PFSimParticle) ...
ROOT::Math::PositionVector3D< ROOT::Math::CylindricalEta3D< double > > REPPoint
void calculatePositionREP()
calculate posrep_ once and for all
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)
int layer_
propagated layer
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)