CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC2/src/RecoBTag/ImpactParameter/interface/TrackCountingComputer.h

Go to the documentation of this file.
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 //#include "RecoBTag/BTagTools/interface/VariableJTA.h"
00011 
00012 class TrackCountingComputer : public JetTagComputer
00013 {
00014  
00015  public:
00016   TrackCountingComputer(const edm::ParameterSet  & parameters )
00017   {
00018     m_nthTrack         = parameters.getParameter<int>("nthTrack");
00019     m_ipType           = parameters.getParameter<int>("impactParameterType");
00020     m_deltaR           = parameters.getParameter<double>("deltaR");
00021     m_cutMaxDecayLen   = parameters.getParameter<double>("maximumDecayLength"); //used
00022     m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis"); //used
00023     //
00024     // access track quality class; "any" takes everything
00025     //
00026     std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass"); //used
00027     m_trackQuality =  reco::TrackBase::qualityByName(trackQualityType);
00028     m_useAllQualities = false;
00029     if (trackQualityType == "any" || 
00030         trackQualityType == "Any" || 
00031         trackQualityType == "ANY" ) m_useAllQualities = true;
00032     
00033     uses("ipTagInfos");
00034   
00035     useVariableJTA_    = parameters.existsAs<bool>("useVariableJTA") ?  parameters.getParameter<bool>("useVariableJTA") : false ;
00036     if (useVariableJTA_){ 
00037       varJTApars = {
00038         parameters.getParameter<double>("a_dR"),
00039         parameters.getParameter<double>("b_dR"),
00040         parameters.getParameter<double>("a_pT"),
00041         parameters.getParameter<double>("b_pT"),
00042         parameters.getParameter<double>("min_pT"),  
00043         parameters.getParameter<double>("max_pT"),
00044         parameters.getParameter<double>("min_pT_dRcut"),  
00045         parameters.getParameter<double>("max_pT_dRcut"),
00046         parameters.getParameter<double>("max_pT_trackPTcut") };
00047     }
00048   }
00049   
00050   
00051   float discriminator(const TagInfoHelper & ti) const 
00052   {           
00053     const reco::TrackIPTagInfo & tkip = ti.get<reco::TrackIPTagInfo>();
00054     std::multiset<float> significances = orderedSignificances(tkip);
00055     std::multiset<float>::reverse_iterator nth=significances.rbegin();
00056     for(int i=0;i<m_nthTrack-1 && nth!=significances.rend();i++) nth++;  
00057     if(nth!=significances.rend()) return *nth; else return -100.;
00058   }
00059   
00060  protected:
00061 
00062   std::multiset<float> orderedSignificances(const reco::TrackIPTagInfo & tkip)   const  {
00063     
00064     const std::vector<reco::TrackIPTagInfo::TrackIPData> & impactParameters((tkip.impactParameterData()));
00065     const edm::RefVector<reco::TrackCollection> & tracks(tkip.selectedTracks());
00066     std::multiset<float> significances;
00067     int i=0;
00068     if(tkip.primaryVertex().isNull())  {  return std::multiset<float>();}
00069     
00070     GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z());
00071     
00072     for(std::vector<reco::TrackIPTagInfo::TrackIPData>::const_iterator it = impactParameters.begin(); it!=impactParameters.end(); ++it, i++)  {
00073       if( fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis  &&        // distance to JetAxis
00074           (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen  &&      // max decay len
00075           (m_useAllQualities  == true || (*tracks[i]).quality(m_trackQuality)) // use selected track qualities
00076           ) {    
00077         if (useVariableJTA_) {
00078           if (tkip.variableJTA( varJTApars )[i])  significances.insert( ((m_ipType==0)?it->ip3d:it->ip2d).significance() );
00079         }
00080         else // no using variable JTA, use the default method 
00081           if(m_deltaR <=0  || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR)
00082             significances.insert( ((m_ipType==0)?it->ip3d:it->ip2d).significance() );
00083       }
00084     }     
00085     
00086     return significances;    
00087   }
00088   
00089   
00090   bool useVariableJTA_;
00091   reco::TrackIPTagInfo::variableJTAParameters varJTApars;
00092   
00093   int m_nthTrack;
00094   int m_ipType;
00095   double m_deltaR;
00096   double  m_cutMaxDecayLen;
00097   double m_cutMaxDistToAxis;
00098   reco::TrackBase::TrackQuality   m_trackQuality;
00099   bool m_useAllQualities;
00100 };
00101 
00102 #endif // ImpactParameter_TrackCountingComputer_h