CMS 3D CMS Logo

List of all members | Public Types | Public Member Functions | Protected Member Functions | Protected Attributes
PromptTrackCountingComputer Class Reference

#include <PromptTrackCountingComputer.h>

Inheritance diagram for PromptTrackCountingComputer:
JetTagComputer

Public Types

using Tokens = void
 

Public Member Functions

float discriminator (const TagInfoHelper &ti) const override
 
 PromptTrackCountingComputer (const edm::ParameterSet &parameters)
 
- Public Member Functions inherited from JetTagComputer
const std::vector< std::string > & getInputLabels () const
 
virtual void initialize (const JetTagComputerRecord &)
 
 JetTagComputer ()
 
 JetTagComputer (const edm::ParameterSet &configuration)
 
float operator() (const reco::BaseTagInfo &info) const
 
float operator() (const TagInfoHelper &helper) const
 
void setupDone ()
 
virtual ~JetTagComputer ()
 

Protected Member Functions

std::multiset< float > orderedSignificances (const reco::TrackIPTagInfo &tkip) const
 
- Protected Member Functions inherited from JetTagComputer
virtual float discriminator (const reco::BaseTagInfo &) const
 
void uses (unsigned int id, const std::string &label)
 
void uses (const std::string &label)
 

Protected Attributes

double m_cutMaxDecayLen
 
double m_cutMaxDistToAxis
 
double m_deltaR
 
double m_deltaRmin
 
int m_ipType
 
int m_nthTrack
 
reco::TrackBase::TrackQuality m_trackQuality
 
bool m_useAllQualities
 
double maxImpactParameter
 
double maxImpactParameterSig
 

Detailed Description

Definition at line 16 of file PromptTrackCountingComputer.h.

Member Typedef Documentation

◆ Tokens

Definition at line 18 of file PromptTrackCountingComputer.h.

Constructor & Destructor Documentation

◆ PromptTrackCountingComputer()

PromptTrackCountingComputer::PromptTrackCountingComputer ( const edm::ParameterSet parameters)
inline

Definition at line 20 of file PromptTrackCountingComputer.h.

References m_cutMaxDecayLen, m_cutMaxDistToAxis, m_deltaR, m_deltaRmin, m_ipType, m_nthTrack, m_trackQuality, m_useAllQualities, maxImpactParameter, maxImpactParameterSig, reco::TrackBase::qualityByName(), AlCaHLTBitMon_QueryRunRegistry::string, and JetTagComputer::uses().

20  {
21  m_nthTrack = parameters.getParameter<int>("nthTrack");
22  m_ipType = parameters.getParameter<int>("impactParameterType");
23  // Maximum and minimum allowed deltaR respectively.
24  m_deltaR = parameters.getParameter<double>("deltaR");
25  m_deltaRmin = parameters.getParameter<double>("deltaRmin");
26  maxImpactParameter = parameters.getParameter<double>("maxImpactParameter");
27  maxImpactParameterSig = parameters.getParameter<double>("maxImpactParameterSig");
28  m_cutMaxDecayLen = parameters.getParameter<double>("maximumDecayLength"); //used
29  m_cutMaxDistToAxis = parameters.getParameter<double>("maximumDistanceToJetAxis"); //used
30  //
31  // access track quality class; "any" takes everything
32  //
33  std::string trackQualityType = parameters.getParameter<std::string>("trackQualityClass"); //used
35  m_useAllQualities = false;
36  if (trackQualityType == "any" || trackQualityType == "Any" || trackQualityType == "ANY")
37  m_useAllQualities = true;
38 
39  uses("ipTagInfos");
40  }
reco::TrackBase::TrackQuality m_trackQuality
void uses(unsigned int id, const std::string &label)
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:126

Member Function Documentation

◆ discriminator()

float PromptTrackCountingComputer::discriminator ( const TagInfoHelper ti) const
inlineoverridevirtual

Reimplemented from JetTagComputer.

Definition at line 42 of file PromptTrackCountingComputer.h.

References JetTagComputer::TagInfoHelper::get(), maxImpactParameterSig, and orderedSignificances().

42  {
43  const reco::TrackIPTagInfo& tkip = ti.get<reco::TrackIPTagInfo>();
44  std::multiset<float> significances = orderedSignificances(tkip);
45  std::multiset<float>::iterator sig;
46  unsigned int nPromptTrk = 0;
47  for (sig = significances.begin(); sig != significances.end(); sig++) {
48  if (fabs(*sig) < maxImpactParameterSig)
49  nPromptTrk++;
50  // edm::LogDebug("") << "Track "<< nPromptTrk << " sig=" << *sig;
51  }
52  return double(nPromptTrk);
53  }
std::multiset< float > orderedSignificances(const reco::TrackIPTagInfo &tkip) const

◆ orderedSignificances()

std::multiset<float> PromptTrackCountingComputer::orderedSignificances ( const reco::TrackIPTagInfo tkip) const
inlineprotected

Definition at line 56 of file PromptTrackCountingComputer.h.

References electronAnalyzer_cfi::DeltaR, mps_fire::i, reco::IPTagInfo< Container, Base >::impactParameterData(), edm::Ref< C, T, F >::isNull(), m_cutMaxDecayLen, m_cutMaxDistToAxis, m_deltaR, m_deltaRmin, m_ipType, m_trackQuality, m_useAllQualities, maxImpactParameter, reco::IPTagInfo< Container, Base >::primaryVertex(), MetAnalyzer::pv(), quality, reco::IPTagInfo< Container, Base >::selectedTracks(), tracks, and relativeConstraints::value.

Referenced by discriminator().

56  {
57  const std::vector<reco::btag::TrackIPData>& impactParameters((tkip.impactParameterData()));
59  std::multiset<float> significances;
60  int i = 0;
61  if (tkip.primaryVertex().isNull()) {
62  return std::multiset<float>();
63  }
64 
65  GlobalPoint pv(tkip.primaryVertex()->position().x(),
66  tkip.primaryVertex()->position().y(),
67  tkip.primaryVertex()->position().z());
68 
69  for (std::vector<reco::btag::TrackIPData>::const_iterator it = impactParameters.begin();
70  it != impactParameters.end();
71  ++it, i++) {
72  if (fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis && // distance to JetAxis
73  (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen && // max decay len
74  (m_useAllQualities == true || (*tracks[i]).quality(m_trackQuality)) // use selected track qualities
75  ) {
76  if ((m_deltaR <= 0 ||
77  ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR) &&
78  (m_deltaRmin <= 0 ||
79  ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) > m_deltaRmin)) {
80  if (fabs(((m_ipType == 0) ? it->ip3d : it->ip2d).value()) < maxImpactParameter) {
81  significances.insert(((m_ipType == 0) ? it->ip3d : it->ip2d).significance());
82  }
83  }
84  }
85  }
86 
87  return significances;
88  }
reco::TrackBase::TrackQuality m_trackQuality
const Container & selectedTracks() const
Definition: IPTagInfo.h:99
def pv(vc)
Definition: MetAnalyzer.py:7
const std::vector< btag::TrackIPData > & impactParameterData() const
Definition: IPTagInfo.h:90
bool isNull() const
Checks for null.
Definition: Ref.h:235
auto const & tracks
cannot be loose
const edm::Ref< VertexCollection > & primaryVertex() const
Definition: IPTagInfo.h:133
string quality

Member Data Documentation

◆ m_cutMaxDecayLen

double PromptTrackCountingComputer::m_cutMaxDecayLen
protected

◆ m_cutMaxDistToAxis

double PromptTrackCountingComputer::m_cutMaxDistToAxis
protected

◆ m_deltaR

double PromptTrackCountingComputer::m_deltaR
protected

◆ m_deltaRmin

double PromptTrackCountingComputer::m_deltaRmin
protected

◆ m_ipType

int PromptTrackCountingComputer::m_ipType
protected

◆ m_nthTrack

int PromptTrackCountingComputer::m_nthTrack
protected

Definition at line 90 of file PromptTrackCountingComputer.h.

Referenced by PromptTrackCountingComputer().

◆ m_trackQuality

reco::TrackBase::TrackQuality PromptTrackCountingComputer::m_trackQuality
protected

◆ m_useAllQualities

bool PromptTrackCountingComputer::m_useAllQualities
protected

◆ maxImpactParameter

double PromptTrackCountingComputer::maxImpactParameter
protected

◆ maxImpactParameterSig

double PromptTrackCountingComputer::maxImpactParameterSig
protected

Definition at line 95 of file PromptTrackCountingComputer.h.

Referenced by discriminator(), and PromptTrackCountingComputer().