CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
MatchJet.cc
Go to the documentation of this file.
1 #include <functional>
2 #include <algorithm>
3 #include <vector>
4 #include <memory>
5 #include <string>
6 
7 #include <Math/GenVector/Cartesian3D.h>
8 #include <Math/GenVector/VectorUtil.h>
9 
11 
12 #include "Matching.h"
13 
14 using namespace btag;
15 
16 namespace {
17  template<typename T>
18  static inline T sqr(const T & val) { return val * val; }
19 
20  template<typename T1, typename T2, typename R>
21  struct JetDistance : public std::binary_function<T1, T2, R> {
22  JetDistance(double sigmaDeltaR, double sigmaDeltaE) :
23  sigmaDeltaR2(sqr(sigmaDeltaR)),
24  sigmaDeltaE2(sqr(sigmaDeltaE)) {}
25 
26  double operator () (const T1 &v1, const T2 &v2) const
27  {
28  using namespace ROOT::Math;
29 // return VectorUtil::DeltaR2(v1, v2) / sigmaDeltaR2 +
30 // sqr(2. * (v1.R() - v2.R()) /
31 // (v1.R() + v2.R())) / sigmaDeltaE2;
32  double x = VectorUtil::DeltaR2(v1, v2) / sigmaDeltaR2 +
33  sqr(2. * (v1.R() - v2.R()) /
34  (v1.R() + v2.R())) / sigmaDeltaE2;
35 // std::cout << "xxx " << VectorUtil::DeltaPhi(v1, v2) << " " << (v1.Eta() - v2.Eta()) << " " << (v1.R() - v2.R()) / (v1.R() + v2.R()) << " " << x << std::endl;
36  return x;
37  }
38 
39  double sigmaDeltaR2, sigmaDeltaE2;
40  };
41 }
42 
44  refJetCorrector(pSet.getParameter<std::string>("refJetCorrection")),
45  recJetCorrector(pSet.getParameter<std::string>("recJetCorrection")),
46  maxChi2(pSet.getParameter<double>("maxChi2")),
47  sigmaDeltaR(pSet.getParameter<double>("sigmaDeltaR")),
48  sigmaDeltaE(pSet.getParameter<double>("sigmaDeltaE")),
49  threshold(1.0)
50 {
51 }
52 
54  const edm::RefToBaseVector<reco::Jet> & refJets_,
55  const edm::RefToBaseVector<reco::Jet> & recJets_,
56  const edm::EventSetup & es)
57 {
60 
61  typedef ROOT::Math::Cartesian3D<double> Vector;
62 
63  std::vector<Vector> corrRefJets;
64  refJets.clear();
66  iter != refJets_.end(); ++iter) {
68  reco::Jet jet = refJetCorrector(*jetRef);
69  if (jet.energy() < threshold)
70  continue;
71 
72  corrRefJets.push_back(Vector(jet.px(), jet.py(), jet.pz()));
73  refJets.push_back(jetRef);
74  }
75 
76  std::vector<Vector> corrRecJets;
77  recJets.clear();
79  iter != recJets_.end(); ++iter) {
81  reco::Jet jet = recJetCorrector(*jetRec);
82  if (jet.energy() < threshold)
83  continue;
84 
85  corrRecJets.push_back(Vector(jet.px(), jet.py(), jet.pz()));
86  recJets.push_back(jetRec);
87  }
88 
89  this->refJets = refJets;
90  refToRec.clear();
91  refToRec.resize(refJets.size(), -1);
92 
93  this->recJets = recJets;
94  recToRef.clear();
95  recToRef.resize(recJets.size(), -1);
96 
97  Matching<double> matching(corrRefJets, corrRecJets,
101 
102  const std::vector<Match>& matches =
103  matching.match(std::less<double>(),
104  std::bind2nd(std::less<double>(), maxChi2));
105  for(std::vector<Match>::const_iterator iter = matches.begin();
106  iter != matches.end(); ++iter) {
107  refToRec[iter->index1] = iter->index2;
108  recToRef[iter->index2] = iter->index1;
109  }
110 }
111 
114 {
116  if (recJet.id() != recJets.id())
117  return result;
118 
119  for(unsigned int i = 0; i != recJets.size(); ++i) {
120  if (recJets[i] == recJet) {
121  int match = recToRef[i];
122  if (match >= 0)
123  result = refJets[match];
124  break;
125  }
126  }
127 
128  return result;
129 }
virtual double energy() const GCC11_FINAL
energy
int i
Definition: DBlmapReader.cc:9
void matchCollections(const edm::RefToBaseVector< reco::Jet > &refJets, const edm::RefToBaseVector< reco::Jet > &recJets, const edm::EventSetup &es)
match the collections
Definition: MatchJet.cc:53
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:29
Base class for all types of Jets.
Definition: Jet.h:20
edm::RefToBaseVector< reco::Jet > refJets
Definition: MatchJet.h:38
const_iterator end() const
std::vector< int > recToRef
Definition: MatchJet.h:37
void setEventSetup(const edm::EventSetup &es)
Definition: CorrectJet.h:23
edm::RefToBaseVector< reco::Jet > recJets
Definition: MatchJet.h:39
virtual double pz() const GCC11_FINAL
z coordinate of momentum vector
virtual double py() const GCC11_FINAL
y coordinate of momentum vector
ProductID id() const
Definition: RefToBase.h:221
double sigmaDeltaE
Definition: MatchJet.h:46
double sigmaDeltaR
Definition: MatchJet.h:45
std::vector< int > refToRec
Definition: MatchJet.h:37
double maxChi2
Definition: MatchJet.h:44
tuple result
Definition: query.py:137
virtual double px() const GCC11_FINAL
x coordinate of momentum vector
size_type size() const
edm::RefToBase< reco::Jet > operator()(const edm::RefToBase< reco::Jet > &recJet) const
Returns the matched &quot;reference&quot; jet.
Definition: MatchJet.cc:113
CorrectJet refJetCorrector
Definition: MatchJet.h:41
CorrectJet recJetCorrector
Definition: MatchJet.h:42
const_iterator begin() const
std::vector< Match > match(SortComparator sortComparator=SortComparator(), CutCriterion cutCriterion=CutCriterion())
Definition: Matching.h:72
Square< F >::type sqr(const F &f)
Definition: Square.h:13
void push_back(const RefToBase< T > &)
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:6
double threshold
Definition: MatchJet.h:47
MatchJet()
Definition: MatchJet.h:22
Definition: DDAxes.h:10
long double T
ProductID id() const