CMS 3D CMS Logo

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:
20 
22  use2d(!parameters.getParameter<bool>("use3d")),
23  useSig(parameters.getParameter<bool>("useSignificance")),
24  unBoost(parameters.getParameter<bool>("unBoost")),
25  minTracks(parameters.getParameter<unsigned int>("minTracks")),
26  minVertices_(1)
27  {
28  uses("svTagInfos");
29  minVertices_ = parameters.existsAs<unsigned int>("minVertices") ? parameters.getParameter<unsigned int>("minVertices") : 1 ;
30  }
31 
32  float discriminator(const TagInfoHelper &tagInfos) const
33  {
34  const TagInfo &info =
35  tagInfos.get<TagInfo>();
36  if(info.nVertices() < minVertices_) return -1;
37  unsigned int idx = 0;
38  while(idx < info.nVertices()) {
39  if (info.nVertexTracks(idx) >= minTracks)
40  break;
41  idx++;
42  }
43  if (idx >= info.nVertices())
44  return -1.0;
45 
46  double gamma;
47  if (unBoost) {
48  reco::TrackKinematics kinematics(
49  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)
64  : -std::log(1 - value);
65 
66  return value;
67  }
68 
69  private:
70  bool use2d;
71  bool useSig;
72  bool unBoost;
73  unsigned int minTracks;
74  unsigned int minVertices_;
75 };
76 
77 #endif // RecoBTag_SecondaryVertex_TemplatedSimpleSecondaryVertexComputer_h
TemplatedSimpleSecondaryVertexComputer(const edm::ParameterSet &parameters)
T getParameter(std::string const &) const
static const TGPicture * info(bool iBackgroundIsBlack)
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:186
float discriminator(const TagInfoHelper &tagInfos) const
const VTX & secondaryVertex(unsigned int index) const
void uses(unsigned int id, const std::string &label)
reco::TemplatedSecondaryVertexTagInfo< IPTI, VTX > TagInfo
double significance() const
Definition: Measurement1D.h:32
double value() const
Definition: Measurement1D.h:28
const math::XYZTLorentzVector & vectorSum() const
Measurement1D flightDistance(unsigned int index, int dim=0) const