CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
CombinedSVSoftLeptonComputer Class Reference

#include <CombinedSVSoftLeptonComputer.h>

Inheritance diagram for CombinedSVSoftLeptonComputer:
CombinedSVComputer

Public Member Functions

 CombinedSVSoftLeptonComputer (const edm::ParameterSet &params)
 
double flipSoftLeptonValue (double value) const
 
template<class IPTI , class SVTI >
reco::TaggingVariableList operator() (const IPTI &ipInfo, const SVTI &svInfo, const reco::CandSoftLeptonTagInfo &muonInfo, const reco::CandSoftLeptonTagInfo &elecInfo) const
 
 ~CombinedSVSoftLeptonComputer () override=default
 
- Public Member Functions inherited from CombinedSVComputer
 CombinedSVComputer (const edm::ParameterSet &params)
 
edm::ParameterSet dropDeltaR (const edm::ParameterSet &pset) const
 
template<class SVTI , class IPTI >
void fillCommonVariables (reco::TaggingVariableList &vars, reco::TrackKinematics &vertexKinematics, const IPTI &ipInfo, const SVTI &svInfo, double &vtx_track_ptSum, double &vtx_track_ESum) const
 
IterationRange flipIterate (int size, bool vertex) const
 
double flipValue (double value, bool vertex) const
 
virtual reco::TaggingVariableList operator() (const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
 
virtual reco::TaggingVariableList operator() (const reco::CandIPTagInfo &ipInfo, const reco::CandSecondaryVertexTagInfo &svInfo) const
 
const reco::btag::TrackIPDatathreshTrack (const reco::CandIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
 
const reco::btag::TrackIPDatathreshTrack (const reco::TrackIPTagInfo &trackIPTagInfo, const reco::btag::SortCriteria sort, const reco::Jet &jet, const GlobalPoint &pv) const
 
virtual ~CombinedSVComputer ()=default
 

Private Attributes

bool SoftLeptonFlip
 

Detailed Description

Definition at line 14 of file CombinedSVSoftLeptonComputer.h.

Constructor & Destructor Documentation

◆ CombinedSVSoftLeptonComputer()

CombinedSVSoftLeptonComputer::CombinedSVSoftLeptonComputer ( const edm::ParameterSet params)
explicit

Definition at line 6 of file CombinedSVSoftLeptonComputer.cc.

7  : CombinedSVComputer(params), SoftLeptonFlip(params.getParameter<bool>("SoftLeptonFlip")) {}
CombinedSVComputer(const edm::ParameterSet &params)

◆ ~CombinedSVSoftLeptonComputer()

CombinedSVSoftLeptonComputer::~CombinedSVSoftLeptonComputer ( )
overridedefault

Member Function Documentation

◆ flipSoftLeptonValue()

double CombinedSVSoftLeptonComputer::flipSoftLeptonValue ( double  value) const
inline

◆ operator()()

template<class IPTI , class SVTI >
reco::TaggingVariableList CombinedSVSoftLeptonComputer::operator() ( const IPTI &  ipInfo,
const SVTI &  svInfo,
const reco::CandSoftLeptonTagInfo muonInfo,
const reco::CandSoftLeptonTagInfo elecInfo 
) const

Definition at line 36 of file CombinedSVSoftLeptonComputer.h.

References reco::SoftLeptonProperties::deltaR, reco::SoftLeptonProperties::etaRel, flipSoftLeptonValue(), mps_fire::i, createfilelist::int, reco::btau::leptonDeltaR, reco::btau::leptonEtaRel, reco::btau::leptonPtRel, reco::btau::leptonRatio, reco::btau::leptonRatioRel, reco::TemplatedSoftLeptonTagInfo< REF >::leptons(), reco::btau::leptonSip3d, CombinedSVComputer::operator()(), reco::TemplatedSoftLeptonTagInfo< REF >::properties(), reco::btag::Vertices::PseudoVertex, reco::SoftLeptonProperties::ptRel, reco::SoftLeptonProperties::ratio, reco::SoftLeptonProperties::ratioRel, reco::btag::Vertices::RecoVertex, reco::SoftLeptonProperties::sip3d, reco::btau::vertexCategory, and reco::btau::vertexLeptonCategory.

39  {
40  using namespace reco;
41 
42  // call the inherited operator()
44 
45  //Jets with vtxCategory 99 cause problems
46  unsigned int vtxType =
47  (vars.checkTag(reco::btau::vertexCategory) ? (unsigned int)(vars.get(reco::btau::vertexCategory)) : 99);
48  if (vtxType == 99)
49  return vars;
50 
51  // the following is specific to soft leptons
52  int leptonCategory = 0; // 0 = no lepton, 1 = muon, 2 = electron
53 
54  for (unsigned int i = 0; i < muonInfo.leptons();
55  ++i) // loop over all muons, not optimal -> find the best or use ranking from best to worst
56  {
57  leptonCategory = 1; // muon category
58  const SoftLeptonProperties &propertiesMuon = muonInfo.properties(i);
59  vars.insert(btau::leptonPtRel, propertiesMuon.ptRel, true);
60  vars.insert(btau::leptonSip3d, flipSoftLeptonValue(propertiesMuon.sip3d), true);
61  vars.insert(btau::leptonDeltaR, propertiesMuon.deltaR, true);
62  vars.insert(btau::leptonRatioRel, propertiesMuon.ratioRel, true);
63  vars.insert(btau::leptonEtaRel, propertiesMuon.etaRel, true);
64  vars.insert(btau::leptonRatio, propertiesMuon.ratio, true);
65  }
66 
67  if (leptonCategory != 1) // no soft muon found, try soft electron
68  {
69  for (unsigned int i = 0; i < elecInfo.leptons();
70  ++i) // loop over all electrons, not optimal -> find the best or use ranking from best to worst
71  {
72  leptonCategory = 2; // electron category
73  const SoftLeptonProperties &propertiesElec = elecInfo.properties(i);
74  vars.insert(btau::leptonPtRel, propertiesElec.ptRel, true);
75  vars.insert(btau::leptonSip3d, flipSoftLeptonValue(propertiesElec.sip3d), true);
76  vars.insert(btau::leptonDeltaR, propertiesElec.deltaR, true);
77  vars.insert(btau::leptonRatioRel, propertiesElec.ratioRel, true);
78  vars.insert(btau::leptonEtaRel, propertiesElec.etaRel, true);
79  vars.insert(btau::leptonRatio, propertiesElec.ratio, true);
80  }
81  }
82 
83  // set the default value for vertexLeptonCategory to 2 (= NoVertexNoSoftLepton)
84  int vertexLepCat = 2;
85 
86  if (leptonCategory == 0) // no soft lepton
87  {
88  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
89  vertexLepCat = 0;
90  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
91  vertexLepCat = 1;
92  else
93  vertexLepCat = 2;
94  } else if (leptonCategory == 1) // soft muon
95  {
96  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
97  vertexLepCat = 3;
98  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
99  vertexLepCat = 4;
100  else
101  vertexLepCat = 5;
102  } else if (leptonCategory == 2) // soft electron
103  {
104  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
105  vertexLepCat = 6;
106  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
107  vertexLepCat = 7;
108  else
109  vertexLepCat = 8;
110  }
111  vars.insert(btau::vertexLeptonCategory, vertexLepCat, true);
112 
113  vars.finalize();
114  return vars;
115 }
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
vars
Definition: DeepTauIdBase.h:60
double flipSoftLeptonValue(double value) const
const SoftLeptonProperties & properties(size_t i) const
fixed size matrix

Member Data Documentation

◆ SoftLeptonFlip

bool CombinedSVSoftLeptonComputer::SoftLeptonFlip
private

Definition at line 27 of file CombinedSVSoftLeptonComputer.h.

Referenced by flipSoftLeptonValue().