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
00027
00028 std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass");
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 useVariableJTA_ = parameters.getParameter<bool>("useVariableJTA");
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 uses("ipTagInfos");
00049 }
00050
00051 float discriminator(const TagInfoHelper & ti) const
00052 {
00053 const reco::TrackIPTagInfo & tkip = ti.get<reco::TrackIPTagInfo>();
00054
00055 const edm::RefVector<reco::TrackCollection> & tracks(tkip.selectedTracks());
00056 const std::vector<float> & allProbabilities((tkip.probabilities(m_ipType)));
00057 const std::vector<reco::TrackIPTagInfo::TrackIPData> & impactParameters((tkip.impactParameterData()));
00058
00059 if (tkip.primaryVertex().isNull()) return 0 ;
00060
00061 GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z());
00062
00063 std::vector<float> probabilities;
00064 int i=0;
00065 for(std::vector<float>::const_iterator it = allProbabilities.begin(); it!=allProbabilities.end(); ++it, i++)
00066 {
00067 if ( fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis &&
00068 (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen &&
00069 (m_useAllQualities == true || (*tracks[i]).quality(m_trackQuality))
00070 )
00071 {
00072 float p;
00073 if(m_trackSign ==0 )
00074 {
00075 if (*it >=0){p=*it/2.;}else{p=1.+*it/2.;}
00076 }
00077 else if(m_trackSign > 0)
00078 {
00079 if(*it >=0 ) p=*it; else continue;
00080 } else
00081 {
00082 if(*it <=0 ) p= -*it; else continue;
00083 }
00084 if (useVariableJTA_) {
00085 if (tkip.variableJTA( varJTApars )[i])
00086 probabilities.push_back(p);
00087 }else{
00088 if(m_deltaR <= 0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR )
00089 probabilities.push_back(p);
00090 }
00091 }
00092 }
00093 return jetProbability(probabilities);
00094 }
00095
00096 double jetProbability( const std::vector<float> & v ) const
00097 {
00098 int ngoodtracks=v.size();
00099 double SumJet=0.;
00100
00101 for(std::vector<float>::const_iterator q = v.begin(); q != v.end(); q++){
00102 SumJet+=(*q>m_minTrackProb)?log(*q):log(m_minTrackProb);
00103 }
00104
00105 double ProbJet;
00106 double Loginvlog=0;
00107
00108 if(SumJet<0.){
00109 if(ngoodtracks>=2){
00110 Loginvlog=log(-SumJet);
00111 }
00112 double Prob=1.;
00113 double lfact=1.;
00114 for(int l=1; l!=ngoodtracks; l++){
00115 lfact*=l;
00116 Prob+=exp(l*Loginvlog-log(1.*lfact));
00117 }
00118 double LogProb=log(Prob);
00119 ProbJet=
00120 std::min(exp(std::max(LogProb+SumJet,-30.)),1.);
00121 }else{
00122 ProbJet=1.;
00123 }
00124 if(ProbJet>1)
00125 std::cout << "ProbJet too high: " << ProbJet << std::endl;
00126
00127
00128
00129 return -log10(ProbJet)/4.;
00130 }
00131 private:
00132 bool useVariableJTA_;
00133 reco::TrackIPTagInfo::variableJTAParameters varJTApars;
00134 double m_minTrackProb;
00135 int m_ipType;
00136 double m_deltaR;
00137 int m_trackSign;
00138 double m_cutMaxDecayLen;
00139 double m_cutMaxDistToAxis;
00140 reco::TrackBase::TrackQuality m_trackQuality;
00141 bool m_useAllQualities;
00142
00143
00144 };
00145
00146 #endif // ImpactParameter_JetProbabilityComputer_h