CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Protected Member Functions | Protected Attributes
PromptTrackCountingComputer Class Reference

#include <PromptTrackCountingComputer.h>

Inheritance diagram for PromptTrackCountingComputer:
JetTagComputer

Public Member Functions

float discriminator (const TagInfoHelper &ti) const
 
 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.

Constructor & Destructor Documentation

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

Definition at line 19 of file PromptTrackCountingComputer.h.

References edm::ParameterSet::getParameter(), 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" ||
37  trackQualityType == "Any" ||
38  trackQualityType == "ANY" ) m_useAllQualities = true;
39 
40  uses("ipTagInfos");
41  }
T getParameter(std::string const &) const
reco::TrackBase::TrackQuality m_trackQuality
void uses(unsigned int id, const std::string &label)
static TrackQuality qualityByName(const std::string &name)
Definition: TrackBase.cc:46

Member Function Documentation

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

Reimplemented from JetTagComputer.

Definition at line 43 of file PromptTrackCountingComputer.h.

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

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

Definition at line 57 of file PromptTrackCountingComputer.h.

References i, reco::TrackIPTagInfo::impactParameterData(), PyquenDefaultSettings_cff::impactParameters, edm::Ref< C, T, F >::isNull(), reco::JTATagInfo::jet(), m_cutMaxDecayLen, m_cutMaxDistToAxis, m_deltaR, m_deltaRmin, m_ipType, m_trackQuality, m_useAllQualities, maxImpactParameter, p4, reco::TrackIPTagInfo::primaryVertex(), reco::TrackIPTagInfo::selectedTracks(), testEve_cfg::tracks, and relativeConstraints::value.

Referenced by discriminator().

57  {
58 
59  const std::vector<reco::TrackIPTagInfo::TrackIPData> & impactParameters((tkip.impactParameterData()));
61  std::multiset<float> significances;
62  int i=0;
63  if(tkip.primaryVertex().isNull()) { return std::multiset<float>();}
64 
65  GlobalPoint pv(tkip.primaryVertex()->position().x(),tkip.primaryVertex()->position().y(),tkip.primaryVertex()->position().z());
66 
67  for(std::vector<reco::TrackIPTagInfo::TrackIPData>::const_iterator it = impactParameters.begin(); it!=impactParameters.end(); ++it, i++)
68  {
69  if( fabs(impactParameters[i].distanceToJetAxis.value()) < m_cutMaxDistToAxis && // distance to JetAxis
70  (impactParameters[i].closestToJetAxis - pv).mag() < m_cutMaxDecayLen && // max decay len
71  (m_useAllQualities == true || (*tracks[i]).quality(m_trackQuality)) // use selected track qualities
72  )
73  {
74  if ( ( m_deltaR <=0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) < m_deltaR ) &&
75  ( m_deltaRmin <=0 || ROOT::Math::VectorUtil::DeltaR((*tkip.jet()).p4().Vect(), (*tracks[i]).momentum()) > m_deltaRmin ) ) {
76  if ( fabs(((m_ipType==0)?it->ip3d:it->ip2d).value()) < maxImpactParameter ) {
77  significances.insert( ((m_ipType==0)?it->ip3d:it->ip2d).significance() );
78  }
79  }
80  }
81  }
82 
83  return significances;
84  }
int i
Definition: DBlmapReader.cc:9
reco::TrackBase::TrackQuality m_trackQuality
const edm::Ref< VertexCollection > & primaryVertex() const
bool isNull() const
Checks for null.
Definition: Ref.h:247
double p4[4]
Definition: TauolaWrapper.h:92
virtual edm::RefToBase< Jet > jet(void) const
returns a polymorphic reference to the tagged jet
Definition: JTATagInfo.h:20
tuple tracks
Definition: testEve_cfg.py:39
const std::vector< TrackIPData > & impactParameterData() const
const edm::RefVector< TrackCollection > & selectedTracks() const

Member Data Documentation

double PromptTrackCountingComputer::m_cutMaxDecayLen
protected
double PromptTrackCountingComputer::m_cutMaxDistToAxis
protected
double PromptTrackCountingComputer::m_deltaR
protected
double PromptTrackCountingComputer::m_deltaRmin
protected
int PromptTrackCountingComputer::m_ipType
protected
int PromptTrackCountingComputer::m_nthTrack
protected

Definition at line 86 of file PromptTrackCountingComputer.h.

Referenced by PromptTrackCountingComputer().

reco::TrackBase::TrackQuality PromptTrackCountingComputer::m_trackQuality
protected
bool PromptTrackCountingComputer::m_useAllQualities
protected
double PromptTrackCountingComputer::maxImpactParameter
protected
double PromptTrackCountingComputer::maxImpactParameterSig
protected

Definition at line 91 of file PromptTrackCountingComputer.h.

Referenced by discriminator(), and PromptTrackCountingComputer().