00001 #ifndef ImpactParameter_TrackCountingComputer_h 00002 #define ImpactParameter_TrackCountingComputer_h 00003 00004 #include "DataFormats/TrackReco/interface/Track.h" 00005 #include "DataFormats/VertexReco/interface/Vertex.h" 00006 #include "DataFormats/BTauReco/interface/TrackCountingTagInfo.h" 00007 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h" 00008 #include "Math/GenVector/VectorUtil.h" 00009 #include "RecoBTau/JetTagComputer/interface/JetTagComputer.h" 00010 00011 class TrackCountingComputer : public JetTagComputer 00012 { 00013 public: 00014 TrackCountingComputer(const edm::ParameterSet & parameters ) 00015 { 00016 m_nthTrack = parameters.getParameter<int>("nthTrack"); 00017 m_ipType = parameters.getParameter<int>("impactParameterType"); 00018 m_deltaR = parameters.getParameter<double>("deltaR"); 00019 m_cutMaxDecayLen = parameters.getParameter<double>("maximumDecayLength"); //used 00020 m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis"); //used 00021 // 00022 // access track quality class; "any" takes everything 00023 // 00024 std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass"); //used 00025 m_trackQuality = reco::TrackBase::qualityByName(trackQualityType); 00026 m_useAllQualities = false; 00027 if (trackQualityType == "any" || 00028 trackQualityType == "Any" || 00029 trackQualityType == "ANY" ) m_useAllQualities = true; 00030 00031 uses("ipTagInfos"); 00032 } 00033 00034 00035 float discriminator(const TagInfoHelper & ti) const 00036 { 00037 const reco::TrackIPTagInfo & tkip = ti.get<reco::TrackIPTagInfo>(); 00038 std::multiset<float> significances = orderedSignificances(tkip); 00039 std::multiset<float>::reverse_iterator nth=significances.rbegin(); 00040 for(int i=0;i<m_nthTrack-1 && nth!=significances.rend();i++) nth++; 00041 if(nth!=significances.rend()) return *nth; else return -100.; 00042 } 00043 00044 protected: 00045 std::multiset<float> orderedSignificances(const reco::TrackIPTagInfo & tkip) const { 00046 00047 const std::vector<reco::TrackIPTagInfo::TrackIPData> & impactParameters((tkip.impactParameterData())); 00048 const edm::RefVector<reco::TrackCollection> & tracks(tkip.selectedTracks()); 00049 std::multiset<float> significances; 00050 int i=0; 00051 if(tkip.primaryVertex().isNull()) { return std::multiset<float>();} 00052 00053 GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z()); 00054 00055 for(std::vector<reco::TrackIPTagInfo::TrackIPData>::const_iterator it = impactParameters.begin(); it!=impactParameters.end(); ++it, i++) 00056 { 00057 if( fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis && // distance to JetAxis 00058 (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen && // max decay len 00059 (m_useAllQualities == true || (*tracks[i]).quality(m_trackQuality)) // use selected track qualities 00060 ) 00061 { 00062 if(m_deltaR <=0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR) 00063 significances.insert( ((m_ipType==0)?it->ip3d:it->ip2d).significance() ); 00064 } 00065 } 00066 00067 return significances; 00068 } 00069 00070 00071 00072 int m_nthTrack; 00073 int m_ipType; 00074 double m_deltaR; 00075 double m_cutMaxDecayLen; 00076 double m_cutMaxDistToAxis; 00077 reco::TrackBase::TrackQuality m_trackQuality; 00078 bool m_useAllQualities; 00079 }; 00080 00081 #endif // ImpactParameter_TrackCountingComputer_h