CMS 3D CMS Logo

CombinedSVSoftLeptonComputer.h
Go to the documentation of this file.
1 #ifndef RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
2 #define RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
3 
9 
11 
13 
15 public:
17  ~CombinedSVSoftLeptonComputer() override = default;
18  double flipSoftLeptonValue(double value) const;
19 
20  template <class IPTI, class SVTI>
21  reco::TaggingVariableList operator()(const IPTI &ipInfo,
22  const SVTI &svInfo,
23  const reco::CandSoftLeptonTagInfo &muonInfo,
24  const reco::CandSoftLeptonTagInfo &elecInfo) const;
25 
26 private:
28 };
29 
31  return SoftLeptonFlip ? -value : value;
32 }
33 
34 template <class IPTI, class SVTI>
36  const SVTI &svInfo,
37  const reco::CandSoftLeptonTagInfo &muonInfo,
38  const reco::CandSoftLeptonTagInfo &elecInfo) const {
39  using namespace reco;
40 
41  // call the inherited operator()
43 
44  //Jets with vtxCategory 99 cause problems
45  unsigned int vtxType =
46  (vars.checkTag(reco::btau::vertexCategory) ? (unsigned int)(vars.get(reco::btau::vertexCategory)) : 99);
47  if (vtxType == 99)
48  return vars;
49 
50  // the following is specific to soft leptons
51  int leptonCategory = 0; // 0 = no lepton, 1 = muon, 2 = electron
52 
53  for (unsigned int i = 0; i < muonInfo.leptons();
54  ++i) // loop over all muons, not optimal -> find the best or use ranking from best to worst
55  {
56  leptonCategory = 1; // muon category
57  const SoftLeptonProperties &propertiesMuon = muonInfo.properties(i);
58  vars.insert(btau::leptonPtRel, propertiesMuon.ptRel, true);
59  vars.insert(btau::leptonSip3d, flipSoftLeptonValue(propertiesMuon.sip3d), true);
60  vars.insert(btau::leptonDeltaR, propertiesMuon.deltaR, true);
61  vars.insert(btau::leptonRatioRel, propertiesMuon.ratioRel, true);
62  vars.insert(btau::leptonEtaRel, propertiesMuon.etaRel, true);
63  vars.insert(btau::leptonRatio, propertiesMuon.ratio, true);
64  }
65 
66  if (leptonCategory != 1) // no soft muon found, try soft electron
67  {
68  for (unsigned int i = 0; i < elecInfo.leptons();
69  ++i) // loop over all electrons, not optimal -> find the best or use ranking from best to worst
70  {
71  leptonCategory = 2; // electron category
72  const SoftLeptonProperties &propertiesElec = elecInfo.properties(i);
73  vars.insert(btau::leptonPtRel, propertiesElec.ptRel, true);
74  vars.insert(btau::leptonSip3d, flipSoftLeptonValue(propertiesElec.sip3d), true);
75  vars.insert(btau::leptonDeltaR, propertiesElec.deltaR, true);
76  vars.insert(btau::leptonRatioRel, propertiesElec.ratioRel, true);
77  vars.insert(btau::leptonEtaRel, propertiesElec.etaRel, true);
78  vars.insert(btau::leptonRatio, propertiesElec.ratio, true);
79  }
80  }
81 
82  // set the default value for vertexLeptonCategory to 2 (= NoVertexNoSoftLepton)
83  int vertexLepCat = 2;
84 
85  if (leptonCategory == 0) // no soft lepton
86  {
87  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
88  vertexLepCat = 0;
89  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
90  vertexLepCat = 1;
91  else
92  vertexLepCat = 2;
93  } else if (leptonCategory == 1) // soft muon
94  {
95  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
96  vertexLepCat = 3;
97  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
98  vertexLepCat = 4;
99  else
100  vertexLepCat = 5;
101  } else if (leptonCategory == 2) // soft electron
102  {
103  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
104  vertexLepCat = 6;
105  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
106  vertexLepCat = 7;
107  else
108  vertexLepCat = 8;
109  }
110  vars.insert(btau::vertexLeptonCategory, vertexLepCat, true);
111 
112  vars.finalize();
113  return vars;
114 }
115 
116 #endif // RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
double flipSoftLeptonValue(double value) const
const SoftLeptonProperties & properties(size_t i) const
~CombinedSVSoftLeptonComputer() override=default
CombinedSVSoftLeptonComputer(const edm::ParameterSet &params)
Definition: value.py:1
reco::TaggingVariableList operator()(const IPTI &ipInfo, const SVTI &svInfo, const reco::CandSoftLeptonTagInfo &muonInfo, const reco::CandSoftLeptonTagInfo &elecInfo) const
fixed size matrix
vars
Definition: DeepTauId.cc:166