![]() |
![]() |
00001 #include <vector> 00002 #include <cstring> 00003 00004 #include "FWCore/Utilities/interface/EDMException.h" 00005 00006 #include "DataFormats/BTauReco/interface/SoftLeptonTagInfo.h" 00007 #include "DataFormats/TrackReco/interface/Track.h" 00008 #include "DataFormats/JetReco/interface/Jet.h" 00009 00010 namespace reco { 00011 00012 using namespace btau; 00013 00014 const float SoftLeptonProperties::quality::undef = -999.0; 00015 00016 unsigned int SoftLeptonProperties::quality::internalByName(const char *name) 00017 { 00018 if (std::strcmp(name, "") == 0) 00019 return 0; 00020 00021 if (std::strcmp(name, "leptonId") == 0) 00022 return leptonId; 00023 else if (std::strcmp(name, "btagLeptonCands") == 0) 00024 return btagLeptonCands; 00025 00026 if (std::strcmp(name, "pfElectronId") == 0) 00027 return pfElectronId; 00028 else if (std::strcmp(name, "btagElectronCands") == 0) 00029 return btagElectronCands; 00030 00031 if (std::strcmp(name, "muonId") == 0) 00032 return muonId; 00033 else if (std::strcmp(name, "btagMuonCands") == 0) 00034 return btagMuonCands; 00035 00036 throw edm::Exception(edm::errors::Configuration) 00037 << "Requested lepton quality \"" << name 00038 << "\" not found in SoftLeptonProperties::quality:byName" 00039 << std::endl; 00040 } 00041 00042 float SoftLeptonProperties::quality(unsigned int index, bool throwIfUndefined) const 00043 { 00044 float qual = quality::undef; 00045 if (index < qualities_.size()) 00046 qual = qualities_[index]; 00047 00048 if (qual == quality::undef && throwIfUndefined) 00049 throw edm::Exception(edm::errors::InvalidReference) 00050 << "Requested lepton quality not found in SoftLeptonProperties::quality" 00051 << std::endl; 00052 00053 return qual; 00054 } 00055 00056 void SoftLeptonProperties::setQuality(unsigned int index, float qual) 00057 { 00058 if (qualities_.size() <= index) 00059 qualities_.resize(index + 1, quality::undef); 00060 00061 qualities_[index] = qual; 00062 } 00063 00064 TaggingVariableList SoftLeptonTagInfo::taggingVariables(void) const { 00065 TaggingVariableList list; 00066 00067 const Jet & jet = *( this->jet() ); 00068 list.insert( TaggingVariable(jetEnergy, jet.energy()), true ); 00069 list.insert( TaggingVariable(jetPt, jet.et()), true ); 00070 list.insert( TaggingVariable(jetEta, jet.eta()), true ); 00071 list.insert( TaggingVariable(jetPhi, jet.phi()), true ); 00072 00073 for (unsigned int i = 0; i < m_leptons.size(); ++i) { 00074 const Track & track = *(m_leptons[i].first); 00075 list.insert( TaggingVariable(trackMomentum, track.p()), true ); 00076 list.insert( TaggingVariable(trackEta, track.eta()), true ); 00077 list.insert( TaggingVariable(trackPhi, track.phi()), true ); 00078 list.insert( TaggingVariable(trackChi2, track.normalizedChi2()), true ); 00079 const SoftLeptonProperties & data = m_leptons[i].second; 00080 list.insert( TaggingVariable(leptonQuality, data.quality(SoftLeptonProperties::quality::leptonId, false)), true ); 00081 list.insert( TaggingVariable(leptonQuality2, data.quality(SoftLeptonProperties::quality::btagLeptonCands, false)), true ); 00082 list.insert( TaggingVariable(trackSip2dSig, data.sip2d), true ); 00083 list.insert( TaggingVariable(trackSip3dSig, data.sip3d), true ); 00084 list.insert( TaggingVariable(trackPtRel, data.ptRel), true ); 00085 list.insert( TaggingVariable(trackP0Par, data.ptRel), true ); 00086 list.insert( TaggingVariable(trackEtaRel, data.etaRel), true ); 00087 list.insert( TaggingVariable(trackDeltaR, data.deltaR), true ); 00088 list.insert( TaggingVariable(trackPParRatio, data.ratioRel), true ); 00089 } 00090 00091 list.finalize(); 00092 return list; 00093 } 00094 00095 } 00096