File indexing completed on 2024-04-06 12:04:51
0001 #include "DataFormats/ParticleFlowCandidate/interface/PFCandidateElectronExtra.h"
0002 #include "DataFormats/GsfTrackReco/interface/GsfTrack.h"
0003 #include <ostream>
0004 #include <iomanip>
0005
0006 using namespace reco;
0007
0008 PFCandidateElectronExtra::PFCandidateElectronExtra() {
0009 status_ = 0;
0010 mvaStatus_ = 0;
0011 pout_ = math::XYZTLorentzVector(0., 0., 0., 0.);
0012 hadEnergy_ = -9999.;
0013 sigmaEtaEta_ = -9999.;
0014
0015 for (MvaVariable m = MVA_FIRST; m < MVA_LAST; m = MvaVariable(m + 1))
0016 mvaVariables_.push_back(-9999.);
0017
0018 gsfTrackRef_ = GsfTrackRef();
0019 kfTrackRef_ = TrackRef();
0020 }
0021
0022 PFCandidateElectronExtra::PFCandidateElectronExtra(const reco::GsfTrackRef& gsfTrack) {
0023 status_ = 0;
0024 mvaStatus_ = 0;
0025 pout_ = math::XYZTLorentzVector(0., 0., 0., 0.);
0026 hadEnergy_ = -9999.;
0027 sigmaEtaEta_ = -9999.;
0028
0029 for (MvaVariable m = MVA_FIRST; m < MVA_LAST; m = MvaVariable(m + 1))
0030 mvaVariables_.push_back(-9999.);
0031
0032 gsfTrackRef_ = gsfTrack;
0033 kfTrackRef_ = TrackRef();
0034
0035 setVariable(MVA_LnPtGsf, log(gsfTrackRef_->ptMode()));
0036 setVariable(MVA_EtaGsf, gsfTrackRef_->etaMode());
0037 setVariable(MVA_Chi2Gsf, gsfTrackRef_->normalizedChi2());
0038 float ptmodeerror = gsfTrackRef_->ptModeError();
0039 if (ptmodeerror > 0.)
0040 setVariable(MVA_SigmaPtOverPt, ptmodeerror / gsfTrackRef_->ptMode());
0041 else
0042 setVariable(MVA_SigmaPtOverPt, -999.);
0043
0044 setVariable(MVA_Fbrem, (gsfTrackRef_->ptMode() - pout_.pt()) / gsfTrackRef_->ptMode());
0045 }
0046
0047 void PFCandidateElectronExtra::setGsfTrackRef(const reco::GsfTrackRef& ref) { gsfTrackRef_ = ref; }
0048
0049 void PFCandidateElectronExtra::setGsfTrackPout(const math::XYZTLorentzVector& pout) { pout_ = pout; }
0050
0051 void PFCandidateElectronExtra::setKfTrackRef(const reco::TrackRef& ref) {
0052 kfTrackRef_ = ref;
0053 float nhit_kf = 0;
0054 float chi2_kf = -0.01;
0055
0056 if (kfTrackRef_.isNonnull()) {
0057 nhit_kf = (float)kfTrackRef_->hitPattern().trackerLayersWithMeasurement();
0058 chi2_kf = kfTrackRef_->normalizedChi2();
0059 }
0060 setVariable(MVA_NhitsKf, nhit_kf);
0061 setVariable(MVA_Chi2Kf, chi2_kf);
0062 }
0063
0064 void PFCandidateElectronExtra::setLateBrem(float val) {
0065 lateBrem_ = val;
0066 setVariable(MVA_LateBrem, val);
0067 }
0068
0069 void PFCandidateElectronExtra::setEarlyBrem(float val) {
0070 earlyBrem_ = val;
0071 setVariable(MVA_FirstBrem, val);
0072 }
0073
0074 void PFCandidateElectronExtra::setHadEnergy(float val) {
0075 hadEnergy_ = val;
0076 if (!clusterEnergies_.empty())
0077 setVariable(MVA_HOverHE, hadEnergy_ / (hadEnergy_ + clusterEnergies_[0]));
0078 }
0079
0080 void PFCandidateElectronExtra::setSigmaEtaEta(float val) {
0081 sigmaEtaEta_ = val;
0082 setVariable(MVA_LogSigmaEtaEta, val);
0083 }
0084
0085 void PFCandidateElectronExtra::setDeltaEta(float val) {
0086 deltaEta_ = val;
0087 setVariable(MVA_DeltaEtaTrackCluster, val);
0088 }
0089
0090 void PFCandidateElectronExtra::setClusterEnergies(const std::vector<float>& energies) {
0091 clusterEnergies_ = energies;
0092
0093 if (pout_.t() != 0.)
0094 setVariable(MVA_EseedOverPout, clusterEnergies_[0] / pout_.t());
0095
0096 const float m_el2 = 0.00051 * 0.00051;
0097 float Ein_gsf = sqrt(gsfTrackRef_->pMode() * gsfTrackRef_->pMode() + m_el2);
0098
0099 float etot = 0;
0100 unsigned size = clusterEnergies_.size();
0101
0102 float ebrem = 0.;
0103 for (unsigned ic = 0; ic < size; ++ic) {
0104 etot += clusterEnergies_[ic];
0105 if (ic > 0)
0106 ebrem += clusterEnergies_[ic];
0107 }
0108 setVariable(MVA_EtotOverPin, etot / Ein_gsf);
0109 setVariable(MVA_EbremOverDeltaP, ebrem / (Ein_gsf - pout_.t()));
0110
0111
0112 if (hadEnergy_ != -9999.)
0113 setHadEnergy(hadEnergy_);
0114 }
0115
0116 void PFCandidateElectronExtra::setMVA(float val) { setVariable(MVA_MVA, val); }
0117
0118 void PFCandidateElectronExtra::setVariable(MvaVariable type, float val) {
0119 mvaVariables_[type] = val;
0120 mvaStatus_ |= (1 << (type));
0121 }
0122
0123 void PFCandidateElectronExtra::setStatus(StatusFlag type, bool status) {
0124 if (status) {
0125 status_ |= (1 << type);
0126 } else {
0127 status_ &= ~(1 << type);
0128 }
0129 }
0130
0131 bool PFCandidateElectronExtra::electronStatus(StatusFlag flag) const { return status_ & (1 << flag); }
0132
0133 bool PFCandidateElectronExtra::mvaStatus(MvaVariable flag) const { return mvaStatus_ & (1 << (flag)); }
0134
0135 float PFCandidateElectronExtra::mvaVariable(MvaVariable var) const {
0136 return (mvaStatus(var) ? mvaVariables_[var] : -9999.);
0137 }
0138
0139 #include <string>
0140
0141 static char const* const listVar[] = {"LogPt",
0142 "Eta",
0143 "SigmaPtOverPt",
0144 "fbrem",
0145 "Chi2Gsf",
0146 "NhitsKf",
0147 "Chi2Kf",
0148 "EtotOverPin",
0149 "EseedOverPout",
0150 "EbremOverDeltaP",
0151 "DeltaEtaTrackCluster",
0152 "LogSigmaEtaEta",
0153 "H/(H+E)",
0154 "LateBrem",
0155 "FirstBrem",
0156 "MVA"};
0157
0158 std::ostream& reco::operator<<(std::ostream& out, const PFCandidateElectronExtra& extra) {
0159 if (!out)
0160 return out;
0161
0162 out << std::setiosflags(std::ios::left) << std::setw(20) << "Variable index" << std::setw(20) << "Name"
0163 << std::setw(10) << "Set(0/1)" << std::setw(8) << "value" << std::endl;
0164 for (PFCandidateElectronExtra::MvaVariable i = PFCandidateElectronExtra::MVA_FIRST;
0165 i < PFCandidateElectronExtra::MVA_LAST;
0166 i = PFCandidateElectronExtra::MvaVariable(i + 1)) {
0167 out << std::setw(20) << i << std::setw(20) << listVar[i] << std::setw(10) << extra.mvaStatus(i) << std::setw(8)
0168 << extra.mvaVariable(i) << std::endl;
0169 }
0170
0171 return out;
0172 }