Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2025-01-31 02:19:52

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