CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_4_5_patch3/src/RecoBTag/ImpactParameter/interface/PromptTrackCountingComputer.h

Go to the documentation of this file.
00001 #ifndef ImpactParameter_PromptTrackCountingComputer_h
00002 #define ImpactParameter_PromptTrackCountingComputer_h
00003 
00004 // This returns a discriminator equal to the number of prompt tracks in the jet
00005 // It is intended for exotica physics, not b tagging.
00006 // It closely resembles the TrackCountingComputer, but with a different discrinator definition and slightly different cuts.
00007 // Author: Ian Tomalin
00008 
00009 #include "DataFormats/TrackReco/interface/Track.h"
00010 #include "DataFormats/VertexReco/interface/Vertex.h"
00011 #include "DataFormats/BTauReco/interface/TrackCountingTagInfo.h"
00012 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
00013 #include "Math/GenVector/VectorUtil.h"
00014 #include "RecoBTau/JetTagComputer/interface/JetTagComputer.h"
00015 
00016 class PromptTrackCountingComputer : public JetTagComputer
00017 {
00018  public:
00019   PromptTrackCountingComputer(const edm::ParameterSet  & parameters )
00020   {
00021      m_nthTrack         = parameters.getParameter<int>("nthTrack");
00022      m_ipType           = parameters.getParameter<int>("impactParameterType");
00023      // Maximum and minimum allowed deltaR respectively. 
00024      m_deltaR           = parameters.getParameter<double>("deltaR");
00025      m_deltaRmin      = parameters.getParameter<double>("deltaRmin");
00026      maxImpactParameter    = parameters.getParameter<double>("maxImpactParameter");
00027      maxImpactParameterSig = parameters.getParameter<double>("maxImpactParameterSig");
00028      m_cutMaxDecayLen   = parameters.getParameter<double>("maximumDecayLength"); //used
00029      m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis"); //used
00030      //
00031      // access track quality class; "any" takes everything
00032      //
00033      std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass"); //used
00034      m_trackQuality =  reco::TrackBase::qualityByName(trackQualityType);
00035      m_useAllQualities = false;
00036      if (trackQualityType == "any" || 
00037          trackQualityType == "Any" || 
00038          trackQualityType == "ANY" ) m_useAllQualities = true;
00039 
00040      uses("ipTagInfos");
00041   }
00042   
00043   float discriminator(const TagInfoHelper & ti) const 
00044    {
00045      const reco::TrackIPTagInfo & tkip = ti.get<reco::TrackIPTagInfo>();
00046      std::multiset<float> significances = orderedSignificances(tkip);
00047      std::multiset<float>::iterator sig;
00048      unsigned int nPromptTrk = 0;
00049      for(sig=significances.begin(); sig!=significances.end(); sig++) {
00050        if (fabs(*sig) < maxImpactParameterSig) nPromptTrk++;
00051        //       edm::LogDebug("") << "Track "<< nPromptTrk << " sig=" << *sig;       
00052      }
00053      return double(nPromptTrk);
00054    }
00055 
00056  protected:
00057      std::multiset<float> orderedSignificances(const reco::TrackIPTagInfo & tkip)   const  {
00058 
00059           const std::vector<reco::TrackIPTagInfo::TrackIPData> & impactParameters((tkip.impactParameterData()));
00060           const edm::RefVector<reco::TrackCollection> & tracks(tkip.selectedTracks());
00061           std::multiset<float> significances;
00062           int i=0;
00063           if(tkip.primaryVertex().isNull())  {  return std::multiset<float>();}
00064 
00065           GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z());
00066 
00067           for(std::vector<reco::TrackIPTagInfo::TrackIPData>::const_iterator it = impactParameters.begin(); it!=impactParameters.end(); ++it, i++)
00068            {
00069            if(   fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis  &&        // distance to JetAxis
00070                  (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen  &&      // max decay len
00071                  (m_useAllQualities  == true || (*tracks[i]).quality(m_trackQuality)) // use selected track qualities
00072              )
00073              {
00074                if ( ( m_deltaR    <=0  || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR ) &&
00075                     ( m_deltaRmin <=0  || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) > m_deltaRmin ) ) {
00076                  if ( fabs(((m_ipType==0)?it->ip3d:it->ip2d).value()) < maxImpactParameter ) {
00077                      significances.insert( ((m_ipType==0)?it->ip3d:it->ip2d).significance() );
00078                  }
00079                }
00080              }
00081           }
00082  
00083          return significances;    
00084    }
00085     
00086    int m_nthTrack;
00087    int m_ipType;
00088    double m_deltaR;
00089    double m_deltaRmin;
00090    double maxImpactParameter;
00091    double maxImpactParameterSig;
00092    double  m_cutMaxDecayLen;
00093    double m_cutMaxDistToAxis;
00094    reco::TrackBase::TrackQuality   m_trackQuality;
00095    bool m_useAllQualities;
00096 };
00097 
00098 #endif // ImpactParameter_PromptTrackCountingComputer_h