CMS 3D CMS Logo

TemplatedJetBProbabilityComputer.h
Go to the documentation of this file.
1 #ifndef ImpactParameter_TemplatedJetBProbabilityComputer_h
2 #define ImpactParameter_TemplatedJetBProbabilityComputer_h
3 
8 #include "Math/GenVector/VectorUtil.h"
10 #include <algorithm>
11 #include <iostream>
12 
13 template <class Container, class Base>
15 {
16  public:
18 
20  {
21  m_ipType = parameters.getParameter<int>("impactParameterType");
22  m_minTrackProb = parameters.getParameter<double>("minimumProbability");
23  m_deltaR = parameters.getParameter<double>("deltaR");
24  m_trackSign = parameters.getParameter<int>("trackIpSign");
25  m_nbTracks = parameters.getParameter<unsigned int>("numberOfBTracks");
26  m_cutMaxDecayLen = parameters.getParameter<double>("maximumDecayLength");
27  m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis");
28 
29  //
30  // access track quality class; "any" takes everything
31  //
32  std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass"); //used
34  m_useAllQualities = false;
35  if (trackQualityType == "any" ||
36  trackQualityType == "Any" ||
37  trackQualityType == "ANY" ) m_useAllQualities = true;
38 
39  useVariableJTA_ = parameters.getParameter<bool>("useVariableJTA");
40  if(useVariableJTA_)
41  varJTApars = {
42  parameters.getParameter<double>("a_dR"),
43  parameters.getParameter<double>("b_dR"),
44  parameters.getParameter<double>("a_pT"),
45  parameters.getParameter<double>("b_pT"),
46  parameters.getParameter<double>("min_pT"),
47  parameters.getParameter<double>("max_pT"),
48  parameters.getParameter<double>("min_pT_dRcut"),
49  parameters.getParameter<double>("max_pT_dRcut"),
50  parameters.getParameter<double>("max_pT_trackPTcut") };
51 
52  uses("ipTagInfos");
53  }
54 
55  float discriminator(const TagInfoHelper & ti) const override
56  {
57  const TagInfo & tkip = ti.get<TagInfo>();
58  const Container & tracks(tkip.selectedTracks());
59  const std::vector<float> & allProbabilities((tkip.probabilities(m_ipType)));
60  const std::vector<reco::btag::TrackIPData> & impactParameters((tkip.impactParameterData()));
61 
62  if(tkip.primaryVertex().isNull()) return 0;
63 
64  GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z());
65 
66  std::vector<float> probabilities;
67  std::vector<float> probabilitiesB;
68  int i=0;
69  for(std::vector<float>::const_iterator it = allProbabilities.begin(); it!=allProbabilities.end(); ++it, i++)
70  {
71  if (fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis && // distance to JetAxis
72  (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen && // max decay len
73  (m_useAllQualities == true || reco::btag::toTrack(tracks[i])->quality(m_trackQuality)) // use selected track qualities
74  )
75  {
76  // Use only positive(or negative) tracks for B
77  float p=fabs(*it);
78  if (useVariableJTA_) {
79  if (tkip.variableJTA( varJTApars )[i]) {
80  if(m_trackSign>0 || *it <0 ) probabilities.push_back(p); //Use all tracks for positive tagger and only negative for negative tagger
81  if(m_trackSign>0 && *it >=0){probabilitiesB.push_back(*it);} //Use only positive tracks for positive tagger
82  if(m_trackSign<0 && *it <=0){probabilitiesB.push_back(- *it);} //Use only negative tracks for negative tagger
83  }
84  }
85  else
86  if( m_deltaR < 0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR)
87  {
88  //if(m_trackSign>0 || *it >0 ) probabilities.push_back(p); //Use all tracks for positive tagger and only negative for negative tagger
89  if(m_trackSign>0 || *it <0 ) probabilities.push_back(p); //Use all tracks for positive tagger and only negative for negative tagger
90 
91  if(m_trackSign>0 && *it >=0){probabilitiesB.push_back(*it);} //Use only positive tracks for positive tagger
92  if(m_trackSign<0 && *it <=0){probabilitiesB.push_back(- *it);} //Use only negative tracks for negative tagger
93  }
94  }
95  }
96 
97  float all = jetProbability(probabilities);
98  std::sort(probabilitiesB.begin(), probabilitiesB.end());
99  if(probabilitiesB.size() > m_nbTracks ) probabilitiesB.resize(m_nbTracks);
100  float b = jetProbability(probabilitiesB);
101 
102  return -log(b)/4-log(all)/4;
103  }
104 
105 double jetProbability( const std::vector<float> & v ) const
106 {
107  int ngoodtracks=v.size();
108  double SumJet=0.;
109 
110  for(std::vector<float>::const_iterator q = v.begin(); q != v.end(); q++){
111  SumJet+=(*q>m_minTrackProb)?log(*q):log(m_minTrackProb);
112  }
113 
114  double ProbJet;
115  double Loginvlog=0;
116 
117  if(SumJet<0.){
118  if(ngoodtracks>=2){
119  Loginvlog=log(-SumJet);
120  }
121  double Prob=1.;
122  double lfact=1.;
123  for(int l=1; l!=ngoodtracks; l++){
124  lfact*=l;
125  Prob+=exp(l*Loginvlog-log(1.*lfact));
126  }
127  double LogProb=log(Prob);
128  ProbJet=
129  std::min(exp(std::max(LogProb+SumJet,-30.)),1.);
130  }else{
131  ProbJet=1.;
132  }
133  if(ProbJet>1)
134  std::cout << "ProbJet too high: " << ProbJet << std::endl;
135 
136  //double LogProbJet=-log(ProbJet);
137  // //return 1.-ProbJet;
138  return ProbJet;
139  }
140  private:
144  int m_ipType;
145  double m_deltaR;
147  unsigned int m_nbTracks;
152 
153 };
154 
155 #endif // ImpactParameter_TemplatedJetBProbabilityComputer_h
T getParameter(std::string const &) const
const std::vector< float > & probabilities(int ip) const
Definition: IPTagInfo.h:103
double jetProbability(const std::vector< float > &v) const
const T & get(unsigned int index=0) const
std::vector< bool > variableJTA(const btag::variableJTAParameters &params) const
Definition: IPTagInfo.h:212
TrackQuality
track quality
Definition: TrackBase.h:151
float discriminator(const TagInfoHelper &ti) const override
TemplatedJetBProbabilityComputer(const edm::ParameterSet &parameters)
const Container & selectedTracks() const
Definition: IPTagInfo.h:101
const reco::Track * toTrack(const reco::TrackBaseRef &t)
Definition: IPTagInfo.h:24
const edm::Ref< VertexCollection > & primaryVertex() const
Definition: IPTagInfo.h:133
reco::btag::variableJTAParameters varJTApars
double p4[4]
Definition: TauolaWrapper.h:92
void uses(unsigned int id, const std::string &label)
def pv(vc)
Definition: MetAnalyzer.py:7
T min(T a, T b)
Definition: MathUtil.h:58
bool isNull() const
Checks for null.
Definition: Ref.h:248
const std::vector< btag::TrackIPData > & impactParameterData() const
Definition: IPTagInfo.h:91
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:134
double b
Definition: hdecay.h:120
bool quality(const TrackQuality) const
Track quality.
Definition: TrackBase.h:549
reco::IPTagInfo< Container, Base > TagInfo
edm::AssociationVector< reco::JetRefBaseProd, Values > Container