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() [1/2]

MatchJet::MatchJet ( )
inline

Definition at line 21 of file MatchJet.h.

21 {}

◆ MatchJet() [2/2]

MatchJet::MatchJet ( const edm::ParameterSet pSet)

Definition at line 20 of file MatchJet.cc.

21  : maxChi2(pSet.getParameter<double>("maxChi2")),
22  sigmaDeltaR2(sqr(pSet.getParameter<double>("sigmaDeltaR"))),
23  sigmaDeltaE2(sqr(pSet.getParameter<double>("sigmaDeltaE"))),
24  threshold(1.0) {}
double sigmaDeltaE2
Definition: MatchJet.h:44
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
double maxChi2
Definition: MatchJet.h:42
double sigmaDeltaR2
Definition: MatchJet.h:43
Square< F >::type sqr(const F &f)
Definition: Square.h:14
double threshold
Definition: MatchJet.h:45

Member Function Documentation

◆ matchCollections()

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 26 of file MatchJet.cc.

References edm::RefToBaseVector< T >::begin(), HltBtagPostValidation_cff::c, edm::RefToBaseVector< T >::clear(), pfClustersFromHGC3DClusters_cfi::corrector, edm::RefToBaseVector< T >::end(), metsig::jet, MultipleCompare::Match(), oniaPATMuonsWithTrigger_cff::matches, maxChi2, edm::RefToBaseVector< T >::push_back(), recJetCorrector, recJets, recToRef, refJetCorrector, refJets, refToRec, CorrectJet::setCorrector(), sigmaDeltaE2, sigmaDeltaR2, edm::RefToBaseVector< T >::size(), funct::sqr(), threshold, and x.

Referenced by BTagPerformanceAnalyzerMC::getJetWithFlavour().

28  {
31 
32  typedef ROOT::Math::Cartesian3D<double> Vector;
33 
34  std::vector<Vector> corrRefJets;
35  refJets.clear();
36  for (edm::RefToBaseVector<reco::Jet>::const_iterator iter = refJets.begin(); iter != refJets_.end(); ++iter) {
37  edm::RefToBase<reco::Jet> jetRef = *iter;
38  reco::Jet jet = refJetCorrector(*jetRef);
39  if (jet.energy() < threshold)
40  continue;
41 
42  corrRefJets.push_back(Vector(jet.px(), jet.py(), jet.pz()));
43  refJets.push_back(jetRef);
44  }
45 
46  std::vector<Vector> corrRecJets;
47  recJets.clear();
48  for (edm::RefToBaseVector<reco::Jet>::const_iterator iter = recJets.begin(); iter != recJets_.end(); ++iter) {
49  edm::RefToBase<reco::Jet> jetRec = *iter;
50  reco::Jet jet = recJetCorrector(*jetRec);
51  if (jet.energy() < threshold)
52  continue;
53 
54  corrRecJets.push_back(Vector(jet.px(), jet.py(), jet.pz()));
55  recJets.push_back(jetRec);
56  }
57 
58  refToRec.clear();
59  refToRec.resize(refJets.size(), -1);
60 
61  recToRef.clear();
62  recToRef.resize(recJets.size(), -1);
63 
64  Matching<double> matching(corrRefJets, corrRecJets, [this](auto& v1, auto& v2) {
65  double x = ROOT::Math::VectorUtil::DeltaR2(v1, v2) / this->sigmaDeltaR2 +
66  sqr(2. * (v1.R() - v2.R()) / (v1.R() + v2.R())) / this->sigmaDeltaE2;
67  // Possible debug output
68  // std::cout << "xxx " << ROOT::Math::VectorUtil::DeltaPhi(v1, v2) << " " << (v1.Eta() - v2.Eta())
69  // << " " << (v1.R() - v2.R()) / (v1.R() + v2.R()) << " " << x << std::endl;
70  return x;
71  });
72 
74 
75  const std::vector<Match>& matches = matching.match(std::less<double>(), [&](auto& c) { return c < this->maxChi2; });
76  for (std::vector<Match>::const_iterator iter = matches.begin(); iter != matches.end(); ++iter) {
77  refToRec[iter->index1] = iter->index2;
78  recToRef[iter->index2] = iter->index1;
79  }
80 }
double sigmaDeltaE2
Definition: MatchJet.h:44
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:36
std::vector< int > recToRef
Definition: MatchJet.h:35
edm::RefToBaseVector< reco::Jet > recJets
Definition: MatchJet.h:37
size_type size() const
const_iterator begin() const
std::vector< int > refToRec
Definition: MatchJet.h:35
double maxChi2
Definition: MatchJet.h:42
double sigmaDeltaR2
Definition: MatchJet.h:43
CorrectJet refJetCorrector
Definition: MatchJet.h:39
CorrectJet recJetCorrector
Definition: MatchJet.h:40
void setCorrector(const reco::JetCorrector *corrector)
Returns the corrected jet.
Definition: CorrectJet.h:21
def Match(required, got)
Square< F >::type sqr(const F &f)
Definition: Square.h:14
void push_back(const RefToBase< T > &)
double threshold
Definition: MatchJet.h:45

◆ operator()()

edm::RefToBase< reco::Jet > MatchJet::operator() ( const edm::RefToBase< reco::Jet > &  recJet) const

Returns the matched "reference" jet.

Definition at line 82 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().

82  {
84  if (recJet.id() != recJets.id())
85  return result;
86 
87  for (unsigned int i = 0; i != recJets.size(); ++i) {
88  if (recJets[i] == recJet) {
89  int match = recToRef[i];
90  if (match >= 0)
91  result = refJets[match];
92  break;
93  }
94  }
95 
96  return result;
97 }
edm::RefToBaseVector< reco::Jet > refJets
Definition: MatchJet.h:36
std::vector< int > recToRef
Definition: MatchJet.h:35
edm::RefToBaseVector< reco::Jet > recJets
Definition: MatchJet.h:37
size_type size() const
ProductID id() const
Definition: RefToBase.h:216
ProductID id() const
std::pair< typename Association::data_type::first_type, double > match(Reference key, Association association, bool bestMatchByMaxValue)
Generic matching function.
Definition: Utils.h:10

◆ setThreshold()

void MatchJet::setThreshold ( const double &  energy)
inline

Member Data Documentation

◆ maxChi2

double MatchJet::maxChi2
private

Definition at line 42 of file MatchJet.h.

Referenced by matchCollections().

◆ recJetCorrector

CorrectJet MatchJet::recJetCorrector
private

Definition at line 40 of file MatchJet.h.

Referenced by matchCollections().

◆ recJets

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

Definition at line 37 of file MatchJet.h.

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

◆ recToRef

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

Definition at line 35 of file MatchJet.h.

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

◆ refJetCorrector

CorrectJet MatchJet::refJetCorrector
private

Definition at line 39 of file MatchJet.h.

Referenced by matchCollections().

◆ refJets

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

Definition at line 36 of file MatchJet.h.

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

◆ refToRec

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

Definition at line 35 of file MatchJet.h.

Referenced by matchCollections().

◆ sigmaDeltaE2

double MatchJet::sigmaDeltaE2
private

Definition at line 44 of file MatchJet.h.

Referenced by matchCollections().

◆ sigmaDeltaR2

double MatchJet::sigmaDeltaR2
private

Definition at line 43 of file MatchJet.h.

Referenced by matchCollections().

◆ threshold

double MatchJet::threshold
private

Definition at line 45 of file MatchJet.h.

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