CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
TemplatedSimpleSecondaryVertexComputer.h
Go to the documentation of this file.
1 #ifndef RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
2 #define RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
3 
4 #include <cmath>
5 
6 #include "Math/GenVector/VectorUtil.h"
7 
11 
13 
15 
16 template <class IPTI, class VTX>
18 public:
19  using Tokens = void;
20 
22 
24  : use2d(!parameters.getParameter<bool>("use3d")),
25  useSig(parameters.getParameter<bool>("useSignificance")),
26  unBoost(parameters.getParameter<bool>("unBoost")),
27  minTracks(parameters.getParameter<unsigned int>("minTracks")),
28  minVertices_(1) {
29  uses("svTagInfos");
30  minVertices_ =
31  parameters.existsAs<unsigned int>("minVertices") ? parameters.getParameter<unsigned int>("minVertices") : 1;
32  }
33 
34  float discriminator(const TagInfoHelper &tagInfos) const override {
35  const TagInfo &info = tagInfos.get<TagInfo>();
36  if (info.nVertices() < minVertices_)
37  return -1;
38  unsigned int idx = 0;
39  while (idx < info.nVertices()) {
40  if (info.nVertexTracks(idx) >= minTracks)
41  break;
42  idx++;
43  }
44  if (idx >= info.nVertices())
45  return -1.0;
46 
47  double gamma;
48  if (unBoost) {
49  reco::TrackKinematics kinematics(info.secondaryVertex(idx));
50  gamma = kinematics.vectorSum().Gamma();
51  } else
52  gamma = 1.0;
53 
54  double value;
55  if (useSig)
56  value = info.flightDistance(idx, use2d).significance();
57  else
58  value = info.flightDistance(idx, use2d).value();
59 
60  value /= gamma;
61 
62  if (useSig)
63  value = (value > 0) ? +std::log(1 + value) : -std::log(1 - value);
64 
65  return value;
66  }
67 
68 private:
69  bool use2d;
70  bool useSig;
71  bool unBoost;
72  unsigned int minTracks;
73  unsigned int minVertices_;
74 };
75 
76 #endif // RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
TemplatedSimpleSecondaryVertexComputer(const edm::ParameterSet &parameters)
static const TGPicture * info(bool iBackgroundIsBlack)
static std::vector< std::string > checklist log
const T & get(unsigned int index=0) const
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:171
const VTX & secondaryVertex(unsigned int index) const
reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX > TagInfo
void uses(unsigned int id, const std::string &label)
double significance() const
Definition: Measurement1D.h:29
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
double value() const
Definition: Measurement1D.h:25
float discriminator(const TagInfoHelper &tagInfos) const override
const math::XYZTLorentzVector & vectorSum() const
Measurement1D flightDistance(unsigned int index, int dim=0) const