Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:27

0001 #ifndef ImpactParameter_PromptTrackCountingComputer_h
0002 #define ImpactParameter_PromptTrackCountingComputer_h
0003 
0004 // This returns a discriminator equal to the number of prompt tracks in the jet
0005 // It is intended for exotica physics, not b tagging.
0006 // It closely resembles the TrackCountingComputer, but with a different discrinator definition and slightly different cuts.
0007 // Author: Ian Tomalin
0008 
0009 #include "DataFormats/TrackReco/interface/Track.h"
0010 #include "DataFormats/VertexReco/interface/Vertex.h"
0011 #include "DataFormats/BTauReco/interface/TrackCountingTagInfo.h"
0012 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
0013 #include "Math/GenVector/VectorUtil.h"
0014 #include "RecoBTau/JetTagComputer/interface/JetTagComputer.h"
0015 
0016 class PromptTrackCountingComputer : public JetTagComputer {
0017 public:
0018   using Tokens = void;
0019 
0020   PromptTrackCountingComputer(const edm::ParameterSet& parameters) {
0021     m_nthTrack = parameters.getParameter<int>("nthTrack");
0022     m_ipType = parameters.getParameter<int>("impactParameterType");
0023     // Maximum and minimum allowed deltaR respectively.
0024     m_deltaR = parameters.getParameter<double>("deltaR");
0025     m_deltaRmin = parameters.getParameter<double>("deltaRmin");
0026     maxImpactParameter = parameters.getParameter<double>("maxImpactParameter");
0027     maxImpactParameterSig = parameters.getParameter<double>("maxImpactParameterSig");
0028     m_cutMaxDecayLen = parameters.getParameter<double>("maximumDecayLength");          //used
0029     m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis");  //used
0030     //
0031     // access track quality class; "any" takes everything
0032     //
0033     std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass");  //used
0034     m_trackQuality = reco::TrackBase::qualityByName(trackQualityType);
0035     m_useAllQualities = false;
0036     if (trackQualityType == "any" || trackQualityType == "Any" || trackQualityType == "ANY")
0037       m_useAllQualities = true;
0038 
0039     uses("ipTagInfos");
0040   }
0041 
0042   float discriminator(const TagInfoHelper& ti) const override {
0043     const reco::TrackIPTagInfo& tkip = ti.get<reco::TrackIPTagInfo>();
0044     std::multiset<float> significances = orderedSignificances(tkip);
0045     std::multiset<float>::iterator sig;
0046     unsigned int nPromptTrk = 0;
0047     for (sig = significances.begin(); sig != significances.end(); sig++) {
0048       if (fabs(*sig) < maxImpactParameterSig)
0049         nPromptTrk++;
0050       //       edm::LogDebug("") << "Track "<< nPromptTrk << " sig=" << *sig;
0051     }
0052     return double(nPromptTrk);
0053   }
0054 
0055 protected:
0056   std::multiset<float> orderedSignificances(const reco::TrackIPTagInfo& tkip) const {
0057     const std::vector<reco::btag::TrackIPData>& impactParameters((tkip.impactParameterData()));
0058     const edm::RefVector<reco::TrackCollection>& tracks(tkip.selectedTracks());
0059     std::multiset<float> significances;
0060     int i = 0;
0061     if (tkip.primaryVertex().isNull()) {
0062       return std::multiset<float>();
0063     }
0064 
0065     GlobalPoint pv(tkip.primaryVertex()->position().x(),
0066                    tkip.primaryVertex()->position().y(),
0067                    tkip.primaryVertex()->position().z());
0068 
0069     for (std::vector<reco::btag::TrackIPData>::const_iterator it = impactParameters.begin();
0070          it != impactParameters.end();
0071          ++it, i++) {
0072       if (fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis &&  // distance to JetAxis
0073           (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen &&      // max decay len
0074           (m_useAllQualities == true || (*tracks[i]).quality(m_trackQuality))          // use selected track qualities
0075       ) {
0076         if ((m_deltaR <= 0 ||
0077              ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR) &&
0078             (m_deltaRmin <= 0 ||
0079              ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) > m_deltaRmin)) {
0080           if (fabs(((m_ipType == 0) ? it->ip3d : it->ip2d).value()) < maxImpactParameter) {
0081             significances.insert(((m_ipType == 0) ? it->ip3d : it->ip2d).significance());
0082           }
0083         }
0084       }
0085     }
0086 
0087     return significances;
0088   }
0089 
0090   int m_nthTrack;
0091   int m_ipType;
0092   double m_deltaR;
0093   double m_deltaRmin;
0094   double maxImpactParameter;
0095   double maxImpactParameterSig;
0096   double m_cutMaxDecayLen;
0097   double m_cutMaxDistToAxis;
0098   reco::TrackBase::TrackQuality m_trackQuality;
0099   bool m_useAllQualities;
0100 };
0101 
0102 #endif  // ImpactParameter_PromptTrackCountingComputer_h