CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros 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  {
26 
27  public:
29  GsfPFRecTrack(double charge,
31  int trackId,
32  const reco::GsfTrackRef& gtrackref,
33  const edm::Ref<std::vector<PFRecTrack> >& kfpfrectrackref);
34 
35 
36 
38  const reco::GsfTrackRef&
39  gsfTrackRef() const {return gsfTrackRef_;}
40 
43  kfPFRecTrackRef() const {return kfPFRecTrackRef_;}
45  void addBrem( const reco::PFBrem& brem);
46 
49 
51  const std::vector<reco::PFBrem>& PFRecBrem()const {return pfBremVec_;}
52 
54  int trackId() const {return trackId_;}
55 
56 
58  void addConvBremPFRecTrackRef(const reco::PFRecTrackRef& pfrectracksref);
59 
61  const std::vector<reco::PFRecTrackRef>& convBremPFRecTrackRef() const {return assoPFRecTrack_;}
62 
64  void addConvBremGsfPFRecTrackRef(const reco::GsfPFRecTrackRef& gsfpfrectracksref);
65 
67  const std::vector<reco::GsfPFRecTrackRef>& convBremGsfPFRecTrackRef() const {return assoGsfPFRecTrack_;}
68 
69  private:
72 
75 
77  std::vector<reco::PFBrem> pfBremVec_;
78 
80  std::vector<reco::PFRecTrackRef> assoPFRecTrack_;
81 
83  std::vector<reco::GsfPFRecTrackRef> assoGsfPFRecTrack_;
84 
86  int trackId_;
87  };
88 
89 
90 }
91 
92 #endif
const std::vector< reco::PFRecTrackRef > & convBremPFRecTrackRef() const
Definition: GsfPFRecTrack.h:61
reconstructed track used as an input to particle flow
Definition: PFRecTrack.h:22
const std::vector< reco::PFBrem > & PFRecBrem() const
Definition: GsfPFRecTrack.h:51
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:71
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:80
const reco::GsfTrackRef & gsfTrackRef() const
Definition: GsfPFRecTrack.h:39
int trackId_
track id
Definition: GsfPFRecTrack.h:86
int trackId() const
Definition: GsfPFRecTrack.h:54
std::vector< reco::GsfPFRecTrackRef > assoGsfPFRecTrack_
vector of GsfPFRecTrackRef from duplicates
Definition: GsfPFRecTrack.h:83
const edm::Ref< std::vector< PFRecTrack > > & kfPFRecTrackRef() const
Definition: GsfPFRecTrack.h:43
void addConvBremGsfPFRecTrackRef(const reco::GsfPFRecTrackRef &gsfpfrectracksref)
GsfPFRecTrackRef from duplicates
const std::vector< reco::GsfPFRecTrackRef > & convBremGsfPFRecTrackRef() const
Definition: GsfPFRecTrack.h:67
AlgoType_t
different types of fitting algorithms
Definition: PFRecTrack.h:27
unsigned int algoType() const
Definition: PFRecTrack.h:47
std::vector< reco::PFBrem > pfBremVec_
vector of PFBrem (empty for KF tracks)
Definition: GsfPFRecTrack.h:77
reco::PFRecTrackRef kfPFRecTrackRef_
ref to the corresponfing PfRecTrack with KF algo (only for PFRecTrack built from GSF track) ...
Definition: GsfPFRecTrack.h:74
double charge() const
Definition: PFTrack.h:87