Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-08-21 04:46:42

0001 #ifndef RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
0002 #define RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h
0003 
0004 #include "DataFormats/BTauReco/interface/CandSoftLeptonTagInfo.h"
0005 #include "DataFormats/BTauReco/interface/TaggingVariable.h"
0006 #include "DataFormats/BTauReco/interface/TemplatedSoftLeptonTagInfo.h"
0007 #include "DataFormats/BTauReco/interface/VertexTypes.h"
0008 #include "DataFormats/TrackReco/interface/TrackFwd.h"
0009 
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 
0012 #include "RecoBTag/SecondaryVertex/interface/CombinedSVComputer.h"
0013 
0014 class CombinedSVSoftLeptonComputer : public CombinedSVComputer {
0015 public:
0016   explicit CombinedSVSoftLeptonComputer(const edm::ParameterSet &params);
0017   ~CombinedSVSoftLeptonComputer() override = default;
0018   double flipSoftLeptonValue(double value) const;
0019 
0020   template <class IPTI, class SVTI>
0021   reco::TaggingVariableList operator()(const IPTI &ipInfo,
0022                                        const SVTI &svInfo,
0023                                        const reco::CandSoftLeptonTagInfo &muonInfo,
0024                                        const reco::CandSoftLeptonTagInfo &elecInfo) const;
0025 
0026 private:
0027   bool SoftLeptonFlip;
0028   using CombinedSVComputer::operator();
0029 };
0030 
0031 inline double CombinedSVSoftLeptonComputer::flipSoftLeptonValue(double value) const {
0032   return SoftLeptonFlip ? -value : value;
0033 }
0034 
0035 template <class IPTI, class SVTI>
0036 reco::TaggingVariableList CombinedSVSoftLeptonComputer::operator()(const IPTI &ipInfo,
0037                                                                    const SVTI &svInfo,
0038                                                                    const reco::CandSoftLeptonTagInfo &muonInfo,
0039                                                                    const reco::CandSoftLeptonTagInfo &elecInfo) const {
0040   using namespace reco;
0041 
0042   // call the inherited operator()
0043   TaggingVariableList vars = CombinedSVComputer::operator()(ipInfo, svInfo);
0044 
0045   //Jets with vtxCategory 99 cause problems
0046   unsigned int vtxType =
0047       (vars.checkTag(reco::btau::vertexCategory) ? (unsigned int)(vars.get(reco::btau::vertexCategory)) : 99);
0048   if (vtxType == 99)
0049     return vars;
0050 
0051   // the following is specific to soft leptons
0052   int leptonCategory = 0;  // 0 = no lepton, 1 = muon, 2 = electron
0053 
0054   for (unsigned int i = 0; i < muonInfo.leptons();
0055        ++i)  // loop over all muons, not optimal -> find the best or use ranking from best to worst
0056   {
0057     leptonCategory = 1;  // muon category
0058     const SoftLeptonProperties &propertiesMuon = muonInfo.properties(i);
0059     vars.insert(btau::leptonPtRel, propertiesMuon.ptRel, true);
0060     vars.insert(btau::leptonSip3d, flipSoftLeptonValue(propertiesMuon.sip3d), true);
0061     vars.insert(btau::leptonDeltaR, propertiesMuon.deltaR, true);
0062     vars.insert(btau::leptonRatioRel, propertiesMuon.ratioRel, true);
0063     vars.insert(btau::leptonEtaRel, propertiesMuon.etaRel, true);
0064     vars.insert(btau::leptonRatio, propertiesMuon.ratio, true);
0065   }
0066 
0067   if (leptonCategory != 1)  // no soft muon found, try soft electron
0068   {
0069     for (unsigned int i = 0; i < elecInfo.leptons();
0070          ++i)  // loop over all electrons, not optimal -> find the best or use ranking from best to worst
0071     {
0072       leptonCategory = 2;  // electron category
0073       const SoftLeptonProperties &propertiesElec = elecInfo.properties(i);
0074       vars.insert(btau::leptonPtRel, propertiesElec.ptRel, true);
0075       vars.insert(btau::leptonSip3d, flipSoftLeptonValue(propertiesElec.sip3d), true);
0076       vars.insert(btau::leptonDeltaR, propertiesElec.deltaR, true);
0077       vars.insert(btau::leptonRatioRel, propertiesElec.ratioRel, true);
0078       vars.insert(btau::leptonEtaRel, propertiesElec.etaRel, true);
0079       vars.insert(btau::leptonRatio, propertiesElec.ratio, true);
0080     }
0081   }
0082 
0083   // set the default value for vertexLeptonCategory to 2 (= NoVertexNoSoftLepton)
0084   int vertexLepCat = 2;
0085 
0086   if (leptonCategory == 0)  // no soft lepton
0087   {
0088     if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
0089       vertexLepCat = 0;
0090     else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
0091       vertexLepCat = 1;
0092     else
0093       vertexLepCat = 2;
0094   } else if (leptonCategory == 1)  // soft muon
0095   {
0096     if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
0097       vertexLepCat = 3;
0098     else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
0099       vertexLepCat = 4;
0100     else
0101       vertexLepCat = 5;
0102   } else if (leptonCategory == 2)  // soft electron
0103   {
0104     if (vtxType == (unsigned int)(btag::Vertices::RecoVertex))
0105       vertexLepCat = 6;
0106     else if (vtxType == (unsigned int)(btag::Vertices::PseudoVertex))
0107       vertexLepCat = 7;
0108     else
0109       vertexLepCat = 8;
0110   }
0111   vars.insert(btau::vertexLeptonCategory, vertexLepCat, true);
0112 
0113   vars.finalize();
0114   return vars;
0115 }
0116 
0117 #endif  // RecoBTag_SecondaryVertex_CombinedSVSoftLeptonComputer_h