CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
MatchJet Class Reference

#include <MatchJet.h>

Public Member Functions

void matchCollections (const edm::RefToBaseVector< reco::Jet > &refJets, const edm::RefToBaseVector< reco::Jet > &recJets, const reco::JetCorrector *corrector)
 match the collections More...
 
 MatchJet ()
 
 MatchJet (const edm::ParameterSet &pSet)
 
edm::RefToBase< reco::Jetoperator() (const edm::RefToBase< reco::Jet > &recJet) const
 Returns the matched "reference" jet. More...
 
void setThreshold (const double &energy)
 

Private Attributes

double maxChi2
 
CorrectJet recJetCorrector
 
edm::RefToBaseVector< reco::JetrecJets
 
std::vector< int > recToRef
 
CorrectJet refJetCorrector
 
edm::RefToBaseVector< reco::JetrefJets
 
std::vector< int > refToRec
 
double sigmaDeltaE2
 
double sigmaDeltaR2
 
double threshold
 

Detailed Description

Match jets

Definition at line 19 of file MatchJet.h.

Constructor & Destructor Documentation

MatchJet::MatchJet ( )
inline

Definition at line 22 of file MatchJet.h.

22 {}
MatchJet::MatchJet ( const edm::ParameterSet pSet)

Definition at line 20 of file MatchJet.cc.

20  :
21  maxChi2(pSet.getParameter<double>("maxChi2")),
22  sigmaDeltaR2(sqr(pSet.getParameter<double>("sigmaDeltaR"))),
23  sigmaDeltaE2(sqr(pSet.getParameter<double>("sigmaDeltaE"))),
24  threshold(1.0)
25 {
26 }
T getParameter(std::string const &) const
double sigmaDeltaE2
Definition: MatchJet.h:47
double maxChi2
Definition: MatchJet.h:45
double sigmaDeltaR2
Definition: MatchJet.h:46
Square< F >::type sqr(const F &f)
Definition: Square.h:13
double threshold
Definition: MatchJet.h:48

Member Function Documentation

void MatchJet::matchCollections ( const edm::RefToBaseVector< reco::Jet > &  refJets,
const edm::RefToBaseVector< reco::Jet > &  recJets,
const reco::JetCorrector corrector 
)

match the collections

Definition at line 28 of file MatchJet.cc.

References edm::RefToBaseVector< T >::begin(), EnergyCorrector::c, edm::RefToBaseVector< T >::clear(), edm::RefToBaseVector< T >::end(), reco::LeafCandidate::energy(), metsig::jet, MultipleCompare::Match(), patRefSel_triggerMatching_cfi::matches, PFTauMatching_cfi::matching, maxChi2, edm::RefToBaseVector< T >::push_back(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), recJetCorrector, recJets, recToRef, refJetCorrector, refJets, refToRec, CorrectJet::setCorrector(), sigmaDeltaE2, sigmaDeltaR2, edm::RefToBaseVector< T >::size(), funct::sqr(), threshold, and x.

Referenced by BTagPerformanceAnalyzerMC::getJetWithFlavour(), and setThreshold().

33 {
34  refJetCorrector.setCorrector(corrector);
35  recJetCorrector.setCorrector(corrector);
36 
37  typedef ROOT::Math::Cartesian3D<double> Vector;
38 
39  std::vector<Vector> corrRefJets;
40  refJets.clear();
42  iter != refJets_.end(); ++iter) {
43  edm::RefToBase<reco::Jet> jetRef = *iter;
44  reco::Jet jet = refJetCorrector(*jetRef);
45  if (jet.energy() < threshold)
46  continue;
47 
48  corrRefJets.push_back(Vector(jet.px(), jet.py(), jet.pz()));
49  refJets.push_back(jetRef);
50  }
51 
52  std::vector<Vector> corrRecJets;
53  recJets.clear();
55  iter != recJets_.end(); ++iter) {
56  edm::RefToBase<reco::Jet> jetRec = *iter;
57  reco::Jet jet = recJetCorrector(*jetRec);
58  if (jet.energy() < threshold)
59  continue;
60 
61  corrRecJets.push_back(Vector(jet.px(), jet.py(), jet.pz()));
62  recJets.push_back(jetRec);
63  }
64 
65  this->refJets = refJets;
66  refToRec.clear();
67  refToRec.resize(refJets.size(), -1);
68 
69  this->recJets = recJets;
70  recToRef.clear();
71  recToRef.resize(recJets.size(), -1);
72 
73  Matching<double> matching(corrRefJets, corrRecJets,
74  [this](auto& v1, auto& v2)
75  {
76  double x = ROOT::Math::VectorUtil::DeltaR2(v1, v2) / this->sigmaDeltaR2 +
77  sqr(2. * (v1.R() - v2.R()) / (v1.R() + v2.R())) / this->sigmaDeltaE2;
78  // Possible debug output
79  // std::cout << "xxx " << ROOT::Math::VectorUtil::DeltaPhi(v1, v2) << " " << (v1.Eta() - v2.Eta())
80  // << " " << (v1.R() - v2.R()) / (v1.R() + v2.R()) << " " << x << std::endl;
81  return x;
82  });
83 
85 
86  const std::vector<Match>& matches =
87  matching.match(std::less<double>(),
88  [&](auto &c){ return c < this->maxChi2;});
89  for(std::vector<Match>::const_iterator iter = matches.begin();
90  iter != matches.end(); ++iter) {
91  refToRec[iter->index1] = iter->index2;
92  recToRef[iter->index2] = iter->index1;
93  }
94 }
double sigmaDeltaE2
Definition: MatchJet.h:47
ROOT::Math::Plane3D::Vector Vector
Definition: EcalHitMaker.cc:29
Base class for all types of Jets.
Definition: Jet.h:20
double px() const final
x coordinate of momentum vector
edm::RefToBaseVector< reco::Jet > refJets
Definition: MatchJet.h:39
std::vector< int > recToRef
Definition: MatchJet.h:38
edm::RefToBaseVector< reco::Jet > recJets
Definition: MatchJet.h:40
std::vector< int > refToRec
Definition: MatchJet.h:38
double maxChi2
Definition: MatchJet.h:45
double pz() const final
z coordinate of momentum vector
double energy() const final
energy
double sigmaDeltaR2
Definition: MatchJet.h:46
size_type size() const
CorrectJet refJetCorrector
Definition: MatchJet.h:42
CorrectJet recJetCorrector
Definition: MatchJet.h:43
void setCorrector(const reco::JetCorrector *corrector)
Returns the corrected jet.
Definition: CorrectJet.h:22
double py() const final
y coordinate of momentum vector
const_iterator begin() const
def Match(required, got)
Square< F >::type sqr(const F &f)
Definition: Square.h:13
void push_back(const RefToBase< T > &)
double threshold
Definition: MatchJet.h:48
edm::RefToBase< reco::Jet > MatchJet::operator() ( const edm::RefToBase< reco::Jet > &  recJet) const

Returns the matched "reference" jet.

Definition at line 97 of file MatchJet.cc.

References mps_fire::i, edm::RefToBaseVector< T >::id(), edm::RefToBase< T >::id(), match(), recJets, recToRef, refJets, mps_fire::result, and edm::RefToBaseVector< T >::size().

Referenced by setThreshold().

98 {
100  if (recJet.id() != recJets.id())
101  return result;
102 
103  for(unsigned int i = 0; i != recJets.size(); ++i) {
104  if (recJets[i] == recJet) {
105  int match = recToRef[i];
106  if (match >= 0)
107  result = refJets[match];
108  break;
109  }
110  }
111 
112  return result;
113 }
edm::RefToBaseVector< reco::Jet > refJets
Definition: MatchJet.h:39
std::vector< int > recToRef
Definition: MatchJet.h:38
edm::RefToBaseVector< reco::Jet > recJets
Definition: MatchJet.h:40
ProductID id() const
Definition: RefToBase.h:242
size_type size() const
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10
ProductID id() const
void MatchJet::setThreshold ( const double &  energy)
inline

Member Data Documentation

double MatchJet::maxChi2
private

Definition at line 45 of file MatchJet.h.

Referenced by matchCollections().

CorrectJet MatchJet::recJetCorrector
private

Definition at line 43 of file MatchJet.h.

Referenced by matchCollections().

edm::RefToBaseVector<reco::Jet> MatchJet::recJets
private

Definition at line 40 of file MatchJet.h.

Referenced by matchCollections(), operator()(), and setThreshold().

std::vector<int> MatchJet::recToRef
private

Definition at line 38 of file MatchJet.h.

Referenced by matchCollections(), and operator()().

CorrectJet MatchJet::refJetCorrector
private

Definition at line 42 of file MatchJet.h.

Referenced by matchCollections().

edm::RefToBaseVector<reco::Jet> MatchJet::refJets
private

Definition at line 39 of file MatchJet.h.

Referenced by matchCollections(), operator()(), and setThreshold().

std::vector<int> MatchJet::refToRec
private

Definition at line 38 of file MatchJet.h.

Referenced by matchCollections().

double MatchJet::sigmaDeltaE2
private

Definition at line 47 of file MatchJet.h.

Referenced by matchCollections().

double MatchJet::sigmaDeltaR2
private

Definition at line 46 of file MatchJet.h.

Referenced by matchCollections().

double MatchJet::threshold
private

Definition at line 48 of file MatchJet.h.

Referenced by utils.StatisticalTest::get_status(), matchCollections(), and setThreshold().