Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-09-07 04:35:56

0001 #ifndef DataFormats_ParticleFlowReco_GsfPFRecTrack_h
0002 #define DataFormats_ParticleFlowReco_GsfPFRecTrack_h
0003 
0004 #include "DataFormats/ParticleFlowReco/interface/PFRecTrack.h"
0005 #include "DataFormats/ParticleFlowReco/interface/PFRecTrackFwd.h"
0006 /* #include "DataFormats/Common/interface/RefToBase.h" */
0007 #include "DataFormats/GsfTrackReco/interface/GsfTrackFwd.h"
0008 #include "DataFormats/ParticleFlowReco/interface/GsfPFRecTrackFwd.h"
0009 #include "DataFormats/ParticleFlowReco/interface/PFBrem.h"
0010 #include <iostream>
0011 
0012 namespace reco {
0013 
0014   /**\class PFRecTrack
0015      \brief reconstructed track used as an input to particle flow    
0016 
0017      Additional information w/r to PFTrack: 
0018      - algorithm used to reconstruct the track
0019      - track ID, soon to be replaced by a RefToBase to the corresponding Track
0020 
0021      \author Renaud Bruneliere, Michele Pioppi, Daniele Benedetti
0022      \date   July 2006
0023   */
0024   class GsfPFRecTrack : public PFRecTrack {
0025   public:
0026     GsfPFRecTrack() {}
0027     GsfPFRecTrack(double charge,
0028                   AlgoType_t algoType,
0029                   int trackId,
0030                   const reco::GsfTrackRef& gtrackref,
0031                   const edm::Ref<std::vector<PFRecTrack> >& kfpfrectrackref);
0032 
0033     /// \return reference to corresponding gsftrack
0034     const reco::GsfTrackRef& gsfTrackRef() const { return gsfTrackRef_; }
0035 
0036     /// \return reference to corresponding KF PFRecTrack  (only for GSF PFRecTrack)
0037     const edm::Ref<std::vector<PFRecTrack> >& kfPFRecTrackRef() const { return kfPFRecTrackRef_; }
0038     /// add a Bremsstrahlung photon
0039     void addBrem(const reco::PFBrem& brem);
0040 
0041     /// calculate posrep_ once and for all for each brem
0042     void calculateBremPositionREP();
0043 
0044     /// \return the vector of PFBrem
0045     const std::vector<reco::PFBrem>& PFRecBrem() const { return pfBremVec_; }
0046 
0047     /// \return id
0048     int trackId() const { return trackId_; }
0049 
0050     /// \add PFRecTrackRef from conv Brems
0051     void addConvBremPFRecTrackRef(const reco::PFRecTrackRef& pfrectracksref);
0052 
0053     /// \return vector of PFRecTrackRef from Conv Brem
0054     const std::vector<reco::PFRecTrackRef>& convBremPFRecTrackRef() const { return assoPFRecTrack_; }
0055 
0056     /// \add GsfPFRecTrackRef from duplicates
0057     void addConvBremGsfPFRecTrackRef(const reco::GsfPFRecTrackRef& gsfpfrectracksref);
0058 
0059     /// \return vector of GsfPFRecTrackRef from duplicates
0060     const std::vector<reco::GsfPFRecTrackRef>& convBremGsfPFRecTrackRef() const { return assoGsfPFRecTrack_; }
0061 
0062   private:
0063     /// reference to corresponding gsf track
0064     reco::GsfTrackRef gsfTrackRef_;
0065 
0066     ///ref to the corresponfing PfRecTrack with KF algo (only for PFRecTrack built from GSF track)
0067     reco::PFRecTrackRef kfPFRecTrackRef_;
0068 
0069     /// vector of PFBrem (empty for KF tracks)
0070     std::vector<reco::PFBrem> pfBremVec_;
0071 
0072     /// vector of PFRecTrackRef from conv Brems
0073     std::vector<reco::PFRecTrackRef> assoPFRecTrack_;
0074 
0075     /// vector of GsfPFRecTrackRef from duplicates
0076     std::vector<reco::GsfPFRecTrackRef> assoGsfPFRecTrack_;
0077 
0078     /// track id
0079     int trackId_;
0080   };
0081 
0082 }  // namespace reco
0083 
0084 #endif