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
00028
00029 std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass");
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 useVariableJTA_ = parameters.getParameter<bool>("useVariableJTA");
00037 if(useVariableJTA_)
00038 varJTApars = {
00039 parameters.getParameter<double>("a_dR"),
00040 parameters.getParameter<double>("b_dR"),
00041 parameters.getParameter<double>("a_pT"),
00042 parameters.getParameter<double>("b_pT"),
00043 parameters.getParameter<double>("min_pT"),
00044 parameters.getParameter<double>("max_pT"),
00045 parameters.getParameter<double>("min_pT_dRcut"),
00046 parameters.getParameter<double>("max_pT_dRcut"),
00047 parameters.getParameter<double>("max_pT_trackPTcut") };
00048
00049 uses("ipTagInfos");
00050 }
00051
00052 float discriminator(const TagInfoHelper & ti) const
00053 {
00054 const reco::TrackIPTagInfo & tkip = ti.get<reco::TrackIPTagInfo>();
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 std::vector<float> probabilitiesB;
00065 int i=0;
00066 for(std::vector<float>::const_iterator it = allProbabilities.begin(); it!=allProbabilities.end(); ++it, i++)
00067 {
00068 if (fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis &&
00069 (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen &&
00070 (m_useAllQualities == true || (*tracks[i]).quality(m_trackQuality))
00071 )
00072 {
00073
00074 float p=fabs(*it);
00075 if (useVariableJTA_) {
00076 if (tkip.variableJTA( varJTApars )[i]) {
00077 if(m_trackSign>0 || *it <0 ) probabilities.push_back(p);
00078 if(m_trackSign>0 && *it >=0){probabilitiesB.push_back(*it);}
00079 if(m_trackSign<0 && *it <=0){probabilitiesB.push_back(- *it);}
00080 }
00081 }
00082 else
00083 if( m_deltaR < 0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR)
00084 {
00085
00086 if(m_trackSign>0 || *it <0 ) probabilities.push_back(p);
00087
00088 if(m_trackSign>0 && *it >=0){probabilitiesB.push_back(*it);}
00089 if(m_trackSign<0 && *it <=0){probabilitiesB.push_back(- *it);}
00090 }
00091 }
00092 }
00093
00094 float all = jetProbability(probabilities);
00095 std::sort(probabilitiesB.begin(), probabilitiesB.end());
00096 if(probabilitiesB.size() > m_nbTracks ) probabilitiesB.resize(m_nbTracks);
00097 float b = jetProbability(probabilitiesB);
00098
00099 return -log(b)/4-log(all)/4;
00100 }
00101
00102 double jetProbability( const std::vector<float> & v ) const
00103 {
00104 int ngoodtracks=v.size();
00105 double SumJet=0.;
00106
00107 for(std::vector<float>::const_iterator q = v.begin(); q != v.end(); q++){
00108 SumJet+=(*q>m_minTrackProb)?log(*q):log(m_minTrackProb);
00109 }
00110
00111 double ProbJet;
00112 double Loginvlog=0;
00113
00114 if(SumJet<0.){
00115 if(ngoodtracks>=2){
00116 Loginvlog=log(-SumJet);
00117 }
00118 double Prob=1.;
00119 double lfact=1.;
00120 for(int l=1; l!=ngoodtracks; l++){
00121 lfact*=l;
00122 Prob+=exp(l*Loginvlog-log(1.*lfact));
00123 }
00124 double LogProb=log(Prob);
00125 ProbJet=
00126 std::min(exp(std::max(LogProb+SumJet,-30.)),1.);
00127 }else{
00128 ProbJet=1.;
00129 }
00130 if(ProbJet>1)
00131 std::cout << "ProbJet too high: " << ProbJet << std::endl;
00132
00133
00134
00135 return ProbJet;
00136 }
00137 private:
00138 bool useVariableJTA_;
00139 reco::TrackIPTagInfo::variableJTAParameters varJTApars;
00140 double m_minTrackProb;
00141 int m_ipType;
00142 double m_deltaR;
00143 int m_trackSign;
00144 unsigned int m_nbTracks;
00145 double m_cutMaxDecayLen;
00146 double m_cutMaxDistToAxis;
00147 reco::TrackBase::TrackQuality m_trackQuality;
00148 bool m_useAllQualities;
00149
00150 };
00151
00152 #endif // ImpactParameter_JetBProbabilityComputer_h