CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoBTag/ImpactParameter/interface/JetBProbabilityComputer.h

Go to the documentation of this file.
00001 #ifndef ImpactParameter_JetBProbabilityComputer_h
00002 #define ImpactParameter_JetBProbabilityComputer_h
00003 
00004 #include "DataFormats/TrackReco/interface/Track.h"
00005 #include "DataFormats/BTauReco/interface/TrackProbabilityTagInfo.h"
00006 #include "DataFormats/BTauReco/interface/TrackIPTagInfo.h"
00007 #include "DataFormats/VertexReco/interface/Vertex.h"
00008 #include "Math/GenVector/VectorUtil.h"
00009 #include "RecoBTau/JetTagComputer/interface/JetTagComputer.h"
00010 #include <algorithm>
00011 #include <iostream>
00012 
00013 class JetBProbabilityComputer : public JetTagComputer
00014 {
00015  public:
00016   JetBProbabilityComputer(const edm::ParameterSet  & parameters )
00017   { 
00018      m_ipType           = parameters.getParameter<int>("impactParameterType");
00019      m_minTrackProb     = parameters.getParameter<double>("minimumProbability");
00020      m_deltaR           = parameters.getParameter<double>("deltaR");
00021      m_trackSign        = parameters.getParameter<int>("trackIpSign");
00022      m_nbTracks         = parameters.getParameter<unsigned int>("numberOfBTracks");
00023      m_cutMaxDecayLen   = parameters.getParameter<double>("maximumDecayLength");
00024      m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis");
00025  
00026      //
00027      // access track quality class; "any" takes everything
00028      //
00029      std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass"); //used
00030      m_trackQuality =  reco::TrackBase::qualityByName(trackQualityType);
00031      m_useAllQualities = false;
00032      if (trackQualityType == "any" || 
00033          trackQualityType == "Any" || 
00034          trackQualityType == "ANY" ) m_useAllQualities = true;
00035 
00036      uses("ipTagInfos");
00037   }
00038  
00039   float discriminator(const TagInfoHelper & ti) const 
00040    {
00041       const reco::TrackIPTagInfo & tkip = ti.get<reco::TrackIPTagInfo>();
00042       const edm::RefVector<reco::TrackCollection> & tracks(tkip.selectedTracks());
00043       const std::vector<float> & allProbabilities((tkip.probabilities(m_ipType)));
00044       const std::vector<reco::TrackIPTagInfo::TrackIPData> & impactParameters((tkip.impactParameterData()));
00045 
00046       if(tkip.primaryVertex().isNull()) return 0;
00047 
00048       GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z());
00049 
00050       std::vector<float> probabilities;
00051       std::vector<float> probabilitiesB;
00052       int i=0;
00053       for(std::vector<float>::const_iterator it = allProbabilities.begin(); it!=allProbabilities.end(); ++it, i++)
00054        {
00055         if (fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis  &&        // distance to JetAxis
00056             (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen  &&      // max decay len
00057             (m_useAllQualities  == true || (*tracks[i]).quality(m_trackQuality)) // use selected track qualities
00058            )
00059          {
00060     // Use only positive(or negative) tracks for B
00061            float p=fabs(*it);
00062            if( m_deltaR < 0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR)
00063              {
00064                if(m_trackSign>0 || *it >0 ) probabilities.push_back(p); //Use all tracks for positive tagger and only negative for negative tagger
00065 
00066                if(m_trackSign>0 && *it >=0){probabilitiesB.push_back(*it);} //Use only positive tracks for positive tagger
00067                if(m_trackSign<0 && *it <=0){probabilitiesB.push_back(- *it);} //Use only negative tracks for negative tagger 
00068              }
00069          }
00070        }
00071 
00072       float all = jetProbability(probabilities); 
00073       std::sort(probabilitiesB.begin(), probabilitiesB.end());
00074       if(probabilitiesB.size() > m_nbTracks )  probabilitiesB.resize(m_nbTracks);
00075       float b = jetProbability(probabilitiesB);
00076         
00077       return -log(b)/4-log(all)/4; 
00078    }
00079 
00080 double jetProbability( const std::vector<float> & v ) const
00081 {
00082    int ngoodtracks=v.size();
00083    double SumJet=0.;
00084 
00085   for(std::vector<float>::const_iterator q = v.begin(); q != v.end(); q++){
00086     SumJet+=(*q>m_minTrackProb)?log(*q):log(m_minTrackProb);
00087   }
00088 
00089   double ProbJet;
00090   double Loginvlog=0;
00091 
00092   if(SumJet<0.){
00093     if(ngoodtracks>=2){
00094       Loginvlog=log(-SumJet);
00095     }
00096     double Prob=1.;
00097     double lfact=1.;
00098     for(int l=1; l!=ngoodtracks; l++){
00099        lfact*=l;
00100       Prob+=exp(l*Loginvlog-log(1.*lfact));
00101     }
00102     double LogProb=log(Prob);
00103     ProbJet=
00104       std::min(exp(std::max(LogProb+SumJet,-30.)),1.);
00105   }else{
00106     ProbJet=1.;
00107   }
00108   if(ProbJet>1)
00109    std::cout << "ProbJet too high: "  << ProbJet << std::endl;
00110 
00111   //double LogProbJet=-log(ProbJet);
00112   //  //return 1.-ProbJet;
00113       return ProbJet;
00114   }
00115  private:
00116    double m_minTrackProb;
00117    int m_ipType;
00118    double m_deltaR;
00119    int m_trackSign;
00120    unsigned int m_nbTracks;
00121    double  m_cutMaxDecayLen;
00122    double m_cutMaxDistToAxis;
00123    reco::TrackBase::TrackQuality   m_trackQuality;
00124    bool m_useAllQualities;
00125 
00126 };
00127 
00128 #endif // ImpactParameter_JetBProbabilityComputer_h