Back to home page

Project CMSSW displayed by LXR

 
 

    


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   /**\class PFTrajectoryPoint
0017      \brief A PFTrack holds several trajectory points, which basically 
0018      contain the position and momentum of a track at a given position.
0019      
0020      \todo   detId_, layer_, isTrackerLayer_ seem to be redundant
0021      \todo   deal with origin and end vertices of PFSimParticles
0022      \todo   remove HCAL exit
0023      \author Renaud Bruneliere
0024      \date   July 2006
0025   */
0026   class PFTrajectoryPoint {
0027   public:
0028     // Next typedef uses double in ROOT 6 rather than Double32_t due to a bug in ROOT 5,
0029     // which otherwise would make ROOT5 files unreadable in ROOT6.  This does not increase
0030     // the size on disk, because due to the bug, double was actually stored on disk in ROOT 5.
0031     typedef ROOT::Math::PositionVector3D<ROOT::Math::CylindricalEta3D<double> > REPPoint;
0032 
0033     /// Define the different layers where the track can be propagated
0034     enum LayerType {
0035       /// Point of closest approach from beam axis (initial point in the case of PFSimParticle)
0036       ClosestApproach = 0,
0037       BeamPipeOrEndVertex = 1,
0038       /// Preshower layer 1
0039       PS1 = 2,
0040       /// Preshower layer 2
0041       PS2 = 3,
0042       /// ECAL front face
0043       ECALEntrance = 4,
0044       /// expected maximum of the shower in ECAL, for an e/gamma particle
0045       /// \todo: add an ECALShowerMaxHadrons
0046       ECALShowerMax = 5,
0047       /// HCAL front face
0048       HCALEntrance = 6,
0049       /// HCAL exit
0050       HCALExit = 7,
0051       /// HO layer
0052       HOLayer = 8,
0053       /// VFcal(HF) front face
0054       VFcalEntrance = 9,
0055       // Number of valid layers
0056       NLayers = 10,
0057       // Unknown
0058       Unknown = -1
0059     };
0060 
0061     static const std::array<std::string, NLayers> layerTypeNames;
0062     static LayerType layerTypeByName(const std::string& name);
0063 
0064     /// default constructor. Set variables at default dummy values
0065     PFTrajectoryPoint();
0066 
0067     /// \brief constructor from values.
0068     /// set detId to -1 if this point is not from a tracker layer
0069     PFTrajectoryPoint(int detId, int layer, const math::XYZPoint& posxyz, const math::XYZTLorentzVector& momentum);
0070 
0071     /// measurement detId
0072     int detId() const { return detId_; }
0073 
0074     /// trajectory point layer
0075     int layer() const { return layer_; }
0076 
0077     /// is this point valid ?
0078     bool isValid() const {
0079       if (layer_ == -1 && detId_ == -1)
0080         return false;
0081       else
0082         return true;
0083     }
0084 
0085     /// is this point corresponding to an intersection with a tracker layer ?
0086     bool isTrackerLayer() const {
0087       if (detId_ >= 0)
0088         return true;
0089       else
0090         return false;
0091     }
0092 
0093     /// cartesian position (x, y, z)
0094     const math::XYZPoint& position() const { return posxyz_; }
0095 
0096     /// trajectory position in (rho, eta, phi) base
0097     const REPPoint& positionREP() const { return posrep_; }
0098 
0099     /// calculate posrep_ once and for all
0100     void calculatePositionREP() { posrep_.SetCoordinates(posxyz_.Rho(), posxyz_.Eta(), posxyz_.Phi()); }
0101 
0102     /// 4-momenta quadrivector
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     /// \brief Is the measurement corresponding to a tracker layer?
0111     /// or was it obtained by propagating the track to a certain position?
0112     bool isTrackerLayer_;
0113 
0114     /// detid if measurement is corresponding to a tracker layer
0115     int detId_;
0116 
0117     /// propagated layer
0118     int layer_;
0119 
0120     /// cartesian position (x, y, z)
0121     math::XYZPoint posxyz_;
0122 
0123     /// position in (rho, eta, phi) base (transient)
0124     REPPoint posrep_;
0125 
0126     /// momentum quadrivector
0127     math::XYZTLorentzVector momentum_;
0128   };
0129 
0130   std::ostream& operator<<(std::ostream& out, const reco::PFTrajectoryPoint& trajPoint);
0131 }  // namespace reco
0132 
0133 #endif