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 
10  public:
11  explicit CombinedSVSoftLeptonComputer(const edm::ParameterSet &params);
12 
13  template <class IPTI,class SVTI>
15  operator () (const IPTI &ipInfo, const SVTI &svInfo,
16  const reco::SoftLeptonTagInfo &muonInfo,
17  const reco::SoftLeptonTagInfo &elecInfo ) const;
18 };
19 
20 template <class IPTI,class SVTI>
22  const reco::SoftLeptonTagInfo &muonInfo,
23  const reco::SoftLeptonTagInfo &elecInfo) const
24 {
25  using namespace reco;
26 
27  // call the inherited operator()
29 
30  // the following is specific to soft leptons
31  int leptonCategory = 0; // 0 = no lepton, 1 = muon, 2 = electron
32 
33  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
34  {
35  leptonCategory = 1; // muon category
36  const SoftLeptonProperties & propertiesMuon = muonInfo.properties(i);
37  vars.insert(btau::leptonPtRel,propertiesMuon.ptRel , true);
38  vars.insert(btau::leptonSip3d,propertiesMuon.sip3d , true);
39  vars.insert(btau::leptonDeltaR,propertiesMuon.deltaR , true);
40  vars.insert(btau::leptonRatioRel,propertiesMuon.ratioRel , true);
41  vars.insert(btau::leptonEtaRel,propertiesMuon.etaRel , true);
42  vars.insert(btau::leptonRatio,propertiesMuon.ratio , true);
43  }
44 
45  if(leptonCategory != 1) // no soft muon found, try soft electron
46  {
47  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
48  {
49  leptonCategory = 2; // electron category
50  const SoftLeptonProperties & propertiesElec = elecInfo.properties(i);
51  vars.insert(btau::leptonPtRel,propertiesElec.ptRel , true);
52  vars.insert(btau::leptonSip3d,propertiesElec.sip3d , true);
53  vars.insert(btau::leptonDeltaR,propertiesElec.deltaR , true);
54  vars.insert(btau::leptonRatioRel,propertiesElec.ratioRel , true);
55  vars.insert(btau::leptonEtaRel,propertiesElec.etaRel , true);
56  vars.insert(btau::leptonRatio,propertiesElec.ratio , true);
57  }
58  }
59 
60 
61  // set the default value for vertexLeptonCategory to 2 (= NoVertexNoSoftLepton)
62  int vertexLepCat = 2;
63 
64  unsigned int vtxType = ( vars.checkTag(reco::btau::vertexCategory) ? (unsigned int)(vars.get(reco::btau::vertexCategory)) : 99 );
65 
66  if(leptonCategory == 0) // no soft lepton
67  {
68  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
69  vertexLepCat = 0;
70  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
71  vertexLepCat = 1;
72  else
73  vertexLepCat = 2;
74  }
75  else if(leptonCategory == 1) // soft muon
76  {
77  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
78  vertexLepCat = 3;
79  else if(vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
80  vertexLepCat = 4;
81  else
82  vertexLepCat = 5;
83  }
84  else if(leptonCategory == 2) // soft electron
85  {
86  if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
87  vertexLepCat = 6;
88  else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
89  vertexLepCat = 7;
90  else
91  vertexLepCat = 8;
92  }
93  vars.insert(btau::vertexLeptonCategory, vertexLepCat , true);
94 
95  vars.finalize();
96  return vars;
97 }
98 
99 #endif // RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
int i
Definition: DBlmapReader.cc:9
virtual reco::TaggingVariableList operator()(const reco::TrackIPTagInfo &ipInfo, const reco::SecondaryVertexTagInfo &svInfo) const
const SoftLeptonProperties & properties(size_t i) const
CombinedSVSoftLeptonComputer(const edm::ParameterSet &params)
reco::TaggingVariableList operator()(const IPTI &ipInfo, const SVTI &svInfo, const reco::SoftLeptonTagInfo &muonInfo, const reco::SoftLeptonTagInfo &elecInfo) const