CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_2_7_hltpatch2/src/RecoBTag/ImpactParameter/interface/JetProbabilityComputer.h

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