CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
CombinedSVSoftLeptonComputer.h
Go to the documentation of this file.
1 #ifndef RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
2 #define RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
3 
5 
7 
8 //(z)
12 
13 
15  public:
16  explicit CombinedSVSoftLeptonComputer(const edm::ParameterSet &params);
17 
18  double flipSoftLeptonValue(double value) const;
19 
20  template <class IPTI,class SVTI>
22  operator () (const IPTI &ipInfo, const SVTI &svInfo,
23  const reco::CandSoftLeptonTagInfo &muonInfo,
24  const reco::CandSoftLeptonTagInfo &elecInfo ) const;
25 
26  private:
28 };
29 
31 {
32  return SoftLeptonFlip ? -value : value;
33 }
34 
35 template <class IPTI,class SVTI>
37  const reco::CandSoftLeptonTagInfo &muonInfo,
38  const reco::CandSoftLeptonTagInfo &elecInfo) const
39 {
40  using namespace reco;
41 
42  // call the inherited operator()
44 
45  //Jets with vtxCategory 99 cause problems
46  unsigned int vtxType = ( vars.checkTag(reco::btau::vertexCategory) ? (unsigned int)(vars.get(reco::btau::vertexCategory)) : 99 );
47  if (vtxType == 99)
48  return vars;
49 
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(); ++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(); ++i) // loop over all electrons, not optimal -> find the best or use ranking from best to worst
69  {
70  leptonCategory = 2; // electron category
71  const SoftLeptonProperties & propertiesElec = elecInfo.properties(i);
72  vars.insert(btau::leptonPtRel,propertiesElec.ptRel , true);
73  vars.insert(btau::leptonSip3d,flipSoftLeptonValue(propertiesElec.sip3d) , true);
74  vars.insert(btau::leptonDeltaR,propertiesElec.deltaR , true);
75  vars.insert(btau::leptonRatioRel,propertiesElec.ratioRel , true);
76  vars.insert(btau::leptonEtaRel,propertiesElec.etaRel , true);
77  vars.insert(btau::leptonRatio,propertiesElec.ratio , true);
78  }
79  }
80 
81 
82  // set the default value for vertexLeptonCategory to 2 (= NoVertexNoSoftLepton)
83  int vertexLepCat = 2;
84 
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  }
95  else if(leptonCategory == 1) // soft muon
96  {
97  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
98  vertexLepCat = 3;
99  else if(vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
100  vertexLepCat = 4;
101  else
102  vertexLepCat = 5;
103  }
104  else if(leptonCategory == 2) // soft electron
105  {
106  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
107  vertexLepCat = 6;
108  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
109  vertexLepCat = 7;
110  else
111  vertexLepCat = 8;
112  }
113  vars.insert(btau::vertexLeptonCategory, vertexLepCat , true);
114 
115  vars.finalize();
116  return vars;
117 }
118 
119 #endif // RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
int i
Definition: DBlmapReader.cc:9
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(const edm::ParameterSet &params)
reco::TaggingVariableList operator()(const IPTI &ipInfo, const SVTI &svInfo, const reco::CandSoftLeptonTagInfo &muonInfo, const reco::CandSoftLeptonTagInfo &elecInfo) const