CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
GsfPFRecTrack.h
Go to the documentation of this file.
1 #ifndef DataFormats_ParticleFlowReco_GsfPFRecTrack_h
2 #define DataFormats_ParticleFlowReco_GsfPFRecTrack_h
3 
6 /* #include "DataFormats/Common/interface/RefToBase.h" */
10 #include <iostream>
11 
12 namespace reco {
13 
24  class GsfPFRecTrack : public PFRecTrack {
25  public:
27  GsfPFRecTrack(double charge,
29  int trackId,
30  const reco::GsfTrackRef& gtrackref,
31  const edm::Ref<std::vector<PFRecTrack> >& kfpfrectrackref);
32 
34  const reco::GsfTrackRef& gsfTrackRef() const { return gsfTrackRef_; }
35 
39  void addBrem(const reco::PFBrem& brem);
40 
43 
45  const std::vector<reco::PFBrem>& PFRecBrem() const { return pfBremVec_; }
46 
48  int trackId() const { return trackId_; }
49 
51  void addConvBremPFRecTrackRef(const reco::PFRecTrackRef& pfrectracksref);
52 
54  const std::vector<reco::PFRecTrackRef>& convBremPFRecTrackRef() const { return assoPFRecTrack_; }
55 
57  void addConvBremGsfPFRecTrackRef(const reco::GsfPFRecTrackRef& gsfpfrectracksref);
58 
60  const std::vector<reco::GsfPFRecTrackRef>& convBremGsfPFRecTrackRef() const { return assoGsfPFRecTrack_; }
61 
62  private:
65 
68 
70  std::vector<reco::PFBrem> pfBremVec_;
71 
73  std::vector<reco::PFRecTrackRef> assoPFRecTrack_;
74 
76  std::vector<reco::GsfPFRecTrackRef> assoGsfPFRecTrack_;
77 
79  int trackId_;
80  };
81 
82 } // namespace reco
83 
84 #endif
const std::vector< reco::PFRecTrackRef > & convBremPFRecTrackRef() const
Definition: GsfPFRecTrack.h:54
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:20
const std::vector< reco::PFBrem > & PFRecBrem() const
Definition: GsfPFRecTrack.h:45
void addConvBremPFRecTrackRef(const reco::PFRecTrackRef &pfrectracksref)
PFRecTrackRef from conv Brems
void addBrem(const reco::PFBrem &brem)
add a Bremsstrahlung photon
reco::GsfTrackRef gsfTrackRef_
reference to corresponding gsf track
Definition: GsfPFRecTrack.h:64
void calculateBremPositionREP()
calculate posrep_ once and for all for each brem
std::vector< reco::PFRecTrackRef > assoPFRecTrack_
vector of PFRecTrackRef from conv Brems
Definition: GsfPFRecTrack.h:73
const reco::GsfTrackRef & gsfTrackRef() const
Definition: GsfPFRecTrack.h:34
int trackId_
track id
Definition: GsfPFRecTrack.h:79
int trackId() const
Definition: GsfPFRecTrack.h:48
std::vector< reco::GsfPFRecTrackRef > assoGsfPFRecTrack_
vector of GsfPFRecTrackRef from duplicates
Definition: GsfPFRecTrack.h:76
const edm::Ref< std::vector< PFRecTrack > > & kfPFRecTrackRef() const
Definition: GsfPFRecTrack.h:37
void addConvBremGsfPFRecTrackRef(const reco::GsfPFRecTrackRef &gsfpfrectracksref)
GsfPFRecTrackRef from duplicates
const std::vector< reco::GsfPFRecTrackRef > & convBremGsfPFRecTrackRef() const
Definition: GsfPFRecTrack.h:60
AlgoType_t
different types of fitting algorithms
Definition: PFRecTrack.h:23
unsigned int algoType() const
Definition: PFRecTrack.h:39
std::vector< reco::PFBrem > pfBremVec_
vector of PFBrem (empty for KF tracks)
Definition: GsfPFRecTrack.h:70
reco::PFRecTrackRef kfPFRecTrackRef_
ref to the corresponfing PfRecTrack with KF algo (only for PFRecTrack built from GSF track) ...
Definition: GsfPFRecTrack.h:67
double charge() const
Definition: PFTrack.h:81