CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
VertexSorting.h
Go to the documentation of this file.
1 #ifndef RecoBTag_SecondaryVertex_VertexSorting_h
2 #define RecoBTag_SecondaryVertex_VertexSorting_h
3 
4 #include <vector>
5 #include <string>
6 
8 #include <map>
14 
15 
16 
17 namespace reco {
18 
19 template <class SecondaryVertex>
21  public:
23  sortCriterium(getSortCriterium(params.getParameter<std::string>("sortCriterium"))){}
25 
26  std::vector<unsigned int>
27  operator () (const std::vector<SecondaryVertex>&svCandidates) const;
28 
29  private:
37  };
38 
39  static SortCriterium getSortCriterium(const std::string &criterium);
40 
42 };
43 
44 template <class SecondaryVertex>
47 {
48  if (criterium == "dist3dError")
49  return sortDist3dErr;
50  if (criterium == "dist3dValue")
51  return sortDist3dVal;
52  if (criterium == "dist3dSignificance")
53  return sortDist3dSig;
54  if (criterium == "dist2dError")
55  return sortDist2dErr;
56  if (criterium == "dist2dValue")
57  return sortDist2dVal;
58  if (criterium == "dist2dSignificance")
59  return sortDist2dSig;
60 
61  throw cms::Exception("InvalidArgument")
62  << "Vertex sort criterium \"" << criterium << "\" is invalid."
63  << std::endl;
64 }
65 
66 // identify most probable SV (closest to interaction point, significance-wise)
67 // FIXME: identify if this is the best strategy!
68 template <class SecondaryVertex>
70  const std::vector<SecondaryVertex> &svCandidates) const
71 {
72  Measurement1D (SecondaryVertex::*measurementFn)() const = 0;
73  switch(sortCriterium) {
74  case sortDist3dErr:
75  case sortDist3dVal:
76  case sortDist3dSig:
77  measurementFn = &SecondaryVertex::dist3d;
78  break;
79  case sortDist2dErr:
80  case sortDist2dVal:
81  case sortDist2dSig:
82  measurementFn = &SecondaryVertex::dist2d;
83  break;
84  }
85 
86  double (Measurement1D::*valueFn)() const = 0;
87  switch(sortCriterium) {
88  case sortDist3dErr:
89  case sortDist2dErr:
90  valueFn = &Measurement1D::error;
91  break;
92  case sortDist3dVal:
93  case sortDist2dVal:
94  valueFn = &Measurement1D::value;
95  break;
96  case sortDist3dSig:
97  case sortDist2dSig:
98  valueFn = &Measurement1D::significance;
99  break;
100  }
101 
102  std::multimap<double, unsigned int> sortedMap;
103  unsigned int i = 0;
104  for(typename std::vector<SecondaryVertex>::const_iterator iter =
105  svCandidates.begin(); iter != svCandidates.end(); iter++) {
106 
107  double value = std::abs((((*iter).*measurementFn)().*valueFn)());
108  sortedMap.insert(std::make_pair(value, i++));
109  }
110 
111  std::vector<unsigned int> result;
112  for(std::multimap<double, unsigned int>::const_iterator iter =
113  sortedMap.begin(); iter != sortedMap.end(); iter++)
114  result.push_back(iter->second);
115 
116  return result;
117 }
118 } // namespace reco
119 
120 #endif // RecoBTag_SecondaryVertex_VertexSorting_h
int i
Definition: DBlmapReader.cc:9
double error() const
Definition: Measurement1D.h:30
SortCriterium sortCriterium
Definition: VertexSorting.h:41
std::vector< unsigned int > operator()(const std::vector< SecondaryVertex > &svCandidates) const
Definition: VertexSorting.h:69
tuple result
Definition: query.py:137
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
VertexSorting(const edm::ParameterSet &params)
Definition: VertexSorting.h:22
double significance() const
Definition: Measurement1D.h:32
double value() const
Definition: Measurement1D.h:28
static SortCriterium getSortCriterium(const std::string &criterium)
Definition: VertexSorting.h:46