CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
PFCandidateEGammaExtra.h
Go to the documentation of this file.
1 #ifndef ParticleFlowCandidate_PFCandidateEGammaExtra_h
2 #define ParticleFlowCandidate_PFCandidateEGammaExtra_h
3 
13 
14 #include <iosfwd>
15 
16 namespace reco {
22  typedef std::pair<reco::PFBlockRef, unsigned> ElementInBlock;
23  typedef std::vector< ElementInBlock > ElementsInBlocks;
24 
26  public:
27  enum StatusFlag {
28  X=0, // undefined
29  Selected, // selected
30  ECALDrivenPreselected, // ECAL-driven electron pre-selected
31  MVASelected, // Passed the internal particle-flow selection (mva selection)
32  Rejected // Rejected
33  };
34 
35  // if you had a variable update NMvaVariables
36  enum MvaVariable {
55  };
56 
59  kKFTracksOnGSFCluster, // any number of additional tracks on GSF cluster
60  kFailsTrackAndHCALIso, // > 3 kfs on Gsf-cluster, bad h/e
61  kKillAdditionalKFs, // tracks with hcal linkbut good gsf etot/p_in
62  kItIsAPion, // bad H/P_in, H/H+E, and E_tot/P_in
63  kCrazyEoverP, // screwey track linking / weird GSFs
64  kTooLargeAngle, // angle between GSF and RSC centroid too large
66  };
67 
68  enum PhotonVetoes {
69  kFailsTrackIso, // the photon fails tracker isolation
71  };
72 
73  public:
80 
82  void setGsfTrackRef(const reco::GsfTrackRef& ref);
83 
85  void setKfTrackRef(const reco::TrackRef & ref);
86 
89  const reco::PFBlockElementCluster& ref) {
90  eleGsfCluster_ = ElementInBlock(blk,ref.index());
91  }
92 
95 
98 
101  return eleGsfCluster_;
102  }
103 
106 
109 
112 
115 
117  void addSingleLegConvTrackRef(const reco::TrackRef& trackref);
118 
120  const std::vector<reco::TrackRef>& singleLegConvTrackRef() const {return assoSingleLegRefTrack_;}
121 
123  void addSingleLegConvMva(const float& mvasingleleg);
124 
126  const std::vector<float>& singleLegConvMva() const {return assoSingleLegMva_;}
127 
129  void addConversionRef(const reco::ConversionRef& convref);
130 
133 
135  void setLateBrem(float val);
137  void setEarlyBrem(float val);
138 
140  void setGsfTrackPout(const math::XYZTLorentzVector& pout);
141 
143  void setClusterEnergies(const std::vector<float>& energies);
144 
146  void setSigmaEtaEta(float val);
147 
149  void setDeltaEta(float val);
150 
152  void setHadEnergy(float val);
153 
155  void setMVA(float val);
156 
158  void setStatus(StatusFlag type,bool status=true);
159 
161  bool electronStatus(StatusFlag) const ;
162 
164  int electronStatus() const {return status_;}
165 
167  bool mvaStatus(MvaVariable flag) const;
168 
170  const std::vector<float> & mvaVariables() const {return mvaVariables_;}
171 
173  float mvaVariable(MvaVariable var) const;
174 
176  float hadEnergy() const {return hadEnergy_;}
177  float sigmaEtaEta() const {return sigmaEtaEta_;}
178 
181  const reco::PFBlockElementTrack& tkref) {
183  assoNonConvExtraTracks_.push_back(std::make_pair(blk,tkref.index()));
184  }
185  }
188  }
189 
190  private:
191  void setVariable(MvaVariable type,float var);
192 
193  private:
200 
203 
206 
208  std::vector<reco::TrackRef> assoSingleLegRefTrack_;
209 
210  // information for track matching
212 
214  std::vector<float> assoSingleLegMva_;
215 
218 
221  std::vector<float> clusterEnergies_;
222 
224  std::vector<float> mvaVariables_;
225 
228 
230  int status_;
231 
234  float earlyBrem_;
235  float lateBrem_;
237  float hadEnergy_;
238  float deltaEta_;
239  };
240 
242  std::ostream& operator<<( std::ostream& out, const PFCandidateEGammaExtra& c );
243 
244 }
245 #endif
type
Definition: HCALResponse.h:21
void addSingleLegConvMva(const float &mvasingleleg)
add Single Leg Conversion mva
void setSuperClusterRef(reco::SuperClusterRef sc)
set reference to the corresponding supercluster
reco::TrackRef kfTrackRef_
Ref to the KF track.
void addSingleLegConvTrackRef(const reco::TrackRef &trackref)
add Single Leg Conversion TrackRef
void setStatus(StatusFlag type, bool status=true)
set status
int status_
Status of the electron.
void setGsfElectronClusterRef(const reco::PFBlockRef &blk, const reco::PFBlockElementCluster &ref)
set gsf electron cluster ref
const ElementInBlock & gsfElectronClusterRef() const
return a reference to the electron cluster ref
void setHadEnergy(float val)
set the had energy. The cluster energies should be entered before
ElementInBlock eleGsfCluster_
Ref to the electron gsf cluster;.
reco::TrackRef kfTrackRef() const
return a reference to the corresponding KF track
int mvaStatus_
status of mva variables
reco::SuperClusterRef superClusterRef() const
return a reference to the corresponding supercluster
void setGsfTrackRef(const reco::GsfTrackRef &ref)
set gsftrack reference
float hadEnergy() const
access to specific variables
reco::ConversionRefVector assoConversionsRef_
vector of ConversionRef from PF
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
void setSigmaEtaEta(float val)
set the sigmaetaeta
const std::vector< reco::TrackRef > & singleLegConvTrackRef() const
return vector of Single Leg Conversion TrackRef from
std::vector< float > mvaVariables_
mva variables - transient !
std::ostream & operator<<(std::ostream &, BeamSpot beam)
Definition: BeamSpot.cc:71
unsigned index() const
void setEarlyBrem(float val)
set EarlyBrem
void setVariable(MvaVariable type, float var)
void setGsfTrackPout(const math::XYZTLorentzVector &pout)
set the pout (not trivial to get from the GSF track)
reco::GsfTrackRef gsfTrackRef_
Ref to the GSF track.
int electronStatus() const
access to the status
const std::vector< float > & singleLegConvMva() const
return Single Leg Conversion mva
tuple out
Definition: dbtoconf.py:99
void setClusterEnergies(const std::vector< float > &energies)
set the cluster energies. the Pout should be saved first
virtual bool trackType(TrackType trType) const
std::vector< ElementInBlock > ElementsInBlocks
reco::GsfTrackRef gsfTrackRef() const
return a reference to the corresponding GSF track
reco::SuperClusterRef scRef_
Ref to (refined) supercluster.
math::XYZTLorentzVector pout_
Variables entering the MVA that should be saved.
void setDeltaEta(float val)
set the delta eta
std::pair< reco::PFBlockRef, unsigned > ElementInBlock
bool mvaStatus(MvaVariable flag) const
access to mva variable status
void setSuperClusterPFECALRef(reco::SuperClusterRef sc)
set reference to the corresponding supercluster
void addConversionRef(const reco::ConversionRef &convref)
add Conversions from PF
reco::SuperClusterRef scPFECALRef_
Ref to PF-ECAL only supercluster.
const ElementsInBlocks & extraNonConvTracks() const
void setKfTrackRef(const reco::TrackRef &ref)
set kf track reference
void setMVA(float val)
set the result (mostly for debugging)
reco::ConversionRefVector conversionRef() const
return Conversions from PF
void setLateBrem(float val)
set LateBrem
tuple status
Definition: ntuplemaker.py:245
float mvaVariable(MvaVariable var) const
access to any variable
reco::SuperClusterRef superClusterPFECALRef() const
return a reference to the corresponding box supercluster
void addExtraNonConvTrack(const reco::PFBlockRef &blk, const reco::PFBlockElementTrack &tkref)
track counting for electrons and photons
std::vector< reco::TrackRef > assoSingleLegRefTrack_
vector of TrackRef from Single Leg conversions
const std::vector< float > & mvaVariables() const
access to the mva variables
std::vector< float > assoSingleLegMva_
vector of Mvas from Single Leg conversions