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:
29 };
30 
32  return SoftLeptonFlip ? -value : value;
33 }
34 
35 template <class IPTI, class SVTI>
37  const SVTI &svInfo,
38  const reco::CandSoftLeptonTagInfo &muonInfo,
39  const reco::CandSoftLeptonTagInfo &elecInfo) const {
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 }
116 
117 #endif // RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
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
~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