CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateElectronExtra.h
Go to the documentation of this file.
1 #ifndef ParticleFlowCandidate_PFCandidateElectronExtra_h
2 #define ParticleFlowCandidate_PFCandidateElectronExtra_h
3 
7 
8 #include <iosfwd>
9 
10 namespace reco {
17  public:
18  enum StatusFlag {
19  X=0, // undefined
20  Selected, // selected
21  ECALDrivenPreselected, // ECAL-driven electron pre-selected
22  MVASelected, // Passed the internal particle-flow selection (mva selection)
23  Rejected // Rejected
24  };
25 
26  // if you had a variable update NMvaVariables
27  enum MvaVariable {
46  };
47 
48 
49  public:
56 
58  void setGsfTrackRef(const reco::GsfTrackRef& ref);
59 
61  void setKfTrackRef(const reco::TrackRef & ref);
62 
65 
67  reco::TrackRef kfTrackRef() const { return kfTrackRef_; }
68 
70  void setLateBrem(float val);
72  void setEarlyBrem(float val);
73 
75  void setGsfTrackPout(const math::XYZTLorentzVector& pout);
76 
78  void setClusterEnergies(const std::vector<float>& energies);
79 
81  void setSigmaEtaEta(float val);
82 
84  void setDeltaEta(float val);
85 
87  void setHadEnergy(float val);
88 
90  void setMVA(float val);
91 
93  void setStatus(StatusFlag type,bool status=true);
94 
96  bool electronStatus(StatusFlag) const ;
97 
99  int electronStatus() const {return status_;}
100 
102  bool mvaStatus(MvaVariable flag) const;
103 
105  const std::vector<float> & mvaVariables() const {return mvaVariables_;}
106 
108  float mvaVariable(MvaVariable var) const;
109 
111  float hadEnergy() const {return hadEnergy_;}
112  float sigmaEtaEta() const {return sigmaEtaEta_;}
113 
114 
115  private:
116  void setVariable(MvaVariable type,float var);
117 
118  private:
123 
125  std::vector<float> clusterEnergies_;
126 
128  std::vector<float> mvaVariables_;
129 
132 
134  int status_;
135 
138  float earlyBrem_;
139  float lateBrem_;
141  float hadEnergy_;
142  float deltaEta_;
143  };
144 
146  std::ostream& operator<<( std::ostream& out, const PFCandidateElectronExtra& c );
147 
148 }
149 #endif
type
Definition: HCALResponse.h:21
reco::GsfTrackRef gsfTrackRef() const
return a reference to the corresponding GSF track
std::vector< float > clusterEnergies_
energy of individual clusters (corrected). The first cluster is the seed
void setGsfTrackRef(const reco::GsfTrackRef &ref)
set gsftrack reference
void setMVA(float val)
set the result (mostly for debugging)
void setLateBrem(float val)
set LateBrem
void setVariable(MvaVariable type, float var)
math::XYZTLorentzVector pout_
Variables entering the MVA that should be saved.
void setSigmaEtaEta(float val)
set the sigmaetaeta
void setGsfTrackPout(const math::XYZTLorentzVector &pout)
set the pout (not trivial to get from the GSF track)
void setClusterEnergies(const std::vector< float > &energies)
set the cluster energies. the Pout should be saved first
void setEarlyBrem(float val)
set EarlyBrem
int mvaStatus_
status of mva variables
int electronStatus() const
access to the status
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
reco::TrackRef kfTrackRef_
Ref to the KF track.
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
reco::TrackRef kfTrackRef() const
return a reference to the corresponding KF track
float hadEnergy() const
access to specific variables
reco::GsfTrackRef gsfTrackRef_
Ref to the GSF track.
void setStatus(StatusFlag type, bool status=true)
set status
void setDeltaEta(float val)
set the delta eta
int status_
Status of the electron.
std::vector< float > mvaVariables_
mva variables - transient !
void setKfTrackRef(const reco::TrackRef &ref)
set kf track reference
bool mvaStatus(MvaVariable flag) const
access to mva variable status
const std::vector< float > & mvaVariables() const
access to the mva variables
float mvaVariable(MvaVariable var) const
access to any variable
tuple status
Definition: ntuplemaker.py:245
void setHadEnergy(float val)
set the had energy. The cluster energies should be entered before