Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:32

0001 #ifndef RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
0002 #define RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
0003 
0004 #include <cmath>
0005 
0006 #include "Math/GenVector/VectorUtil.h"
0007 
0008 #include "DataFormats/TrackReco/interface/Track.h"
0009 #include "DataFormats/VertexReco/interface/Vertex.h"
0010 #include "DataFormats/BTauReco/interface/TemplatedSecondaryVertexTagInfo.h"
0011 
0012 #include "RecoBTau/JetTagComputer/interface/JetTagComputer.h"
0013 
0014 #include "RecoBTag/SecondaryVertex/interface/TrackKinematics.h"
0015 
0016 template <class IPTI, class VTX>
0017 class TemplatedSimpleSecondaryVertexComputer : public JetTagComputer {
0018 public:
0019   using Tokens = void;
0020 
0021   typedef reco::TemplatedSecondaryVertexTagInfo<IPTI, VTX> TagInfo;
0022 
0023   TemplatedSimpleSecondaryVertexComputer(const edm::ParameterSet &parameters)
0024       : use2d(!parameters.getParameter<bool>("use3d")),
0025         useSig(parameters.getParameter<bool>("useSignificance")),
0026         unBoost(parameters.getParameter<bool>("unBoost")),
0027         minTracks(parameters.getParameter<unsigned int>("minTracks")),
0028         minVertices_(1) {
0029     uses("svTagInfos");
0030     minVertices_ =
0031         parameters.existsAs<unsigned int>("minVertices") ? parameters.getParameter<unsigned int>("minVertices") : 1;
0032   }
0033 
0034   float discriminator(const TagInfoHelper &tagInfos) const override {
0035     const TagInfo &info = tagInfos.get<TagInfo>();
0036     if (info.nVertices() < minVertices_)
0037       return -1;
0038     unsigned int idx = 0;
0039     while (idx < info.nVertices()) {
0040       if (info.nVertexTracks(idx) >= minTracks)
0041         break;
0042       idx++;
0043     }
0044     if (idx >= info.nVertices())
0045       return -1.0;
0046 
0047     double gamma;
0048     if (unBoost) {
0049       reco::TrackKinematics kinematics(info.secondaryVertex(idx));
0050       gamma = kinematics.vectorSum().Gamma();
0051     } else
0052       gamma = 1.0;
0053 
0054     double value;
0055     if (useSig)
0056       value = info.flightDistance(idx, use2d).significance();
0057     else
0058       value = info.flightDistance(idx, use2d).value();
0059 
0060     value /= gamma;
0061 
0062     if (useSig)
0063       value = (value > 0) ? +std::log(1 + value) : -std::log(1 - value);
0064 
0065     return value;
0066   }
0067 
0068 private:
0069   bool use2d;
0070   bool useSig;
0071   bool unBoost;
0072   unsigned int minTracks;
0073   unsigned int minVertices_;
0074 };
0075 
0076 #endif  // RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h