CMS 3D CMS Logo

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

#include <MuonAnalysis/MuonAssociators/interface/MatcherByPullsAlgorithm.cc>

Public Member Functions

void fillInvCov (const reco::Track &tk, AlgebraicSymMatrix55 &invCov) const
 Fill the inverse covariance matrix for the match(track, candidate, invCov) method. More...
 
std::pair< bool, float > match (const reco::Track &tk, const reco::Candidate &mc, const AlgebraicSymMatrix55 &invertedCovariance) const
 
std::pair< int, float > match (const reco::RecoCandidate &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good) const
 
std::pair< int, float > match (const reco::Track &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good) const
 
 MatcherByPullsAlgorithm (const edm::ParameterSet &)
 
void matchMany (const reco::RecoCandidate &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good, std::vector< std::pair< double, int > > &matchesToFill) const
 
void matchMany (const reco::Track &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good, std::vector< std::pair< double, int > > &matchesToFill) const
 
 ~MatcherByPullsAlgorithm ()
 

Private Types

enum  TrackChoice { StaTrack, TrkTrack, GlbTrack }
 Enum to define which track to use. More...
 

Private Member Functions

const reco::Tracktrack (const reco::RecoCandidate &src) const
 Get track out of Candidate, NULL if missing. More...
 

Private Attributes

double cut_
 Cut on the pull. More...
 
bool diagOnly_
 Use only the diagonal terms of the covariance matrix. More...
 
double dr2_
 DeltaR of the matching cone. More...
 
TrackChoice track_
 Track to be used in matching. More...
 
bool useVertex_
 Use also dxy / dsz in the matching. More...
 

Detailed Description

Description: Matches a RecoCandidate to a GenParticle (or any other Candidate) using the pulls of the helix parameters

Implementation: <Notes on="" implementation>="">

Definition at line 34 of file MatcherByPullsAlgorithm.h.

Member Enumeration Documentation

◆ TrackChoice

Enum to define which track to use.

Enumerator
StaTrack 
TrkTrack 
GlbTrack 

Definition at line 83 of file MatcherByPullsAlgorithm.h.

Constructor & Destructor Documentation

◆ MatcherByPullsAlgorithm()

MatcherByPullsAlgorithm::MatcherByPullsAlgorithm ( const edm::ParameterSet iConfig)
explicit

Definition at line 13 of file MatcherByPullsAlgorithm.cc.

References Exception, edm::ParameterSet::getParameter(), GlbTrack, StaTrack, AlCaHLTBitMon_QueryRunRegistry::string, track(), track_, and TrkTrack.

14  : dr2_(std::pow(iConfig.getParameter<double>("maxDeltaR"), 2)),
15  cut_(iConfig.getParameter<double>("maxPull")),
16  diagOnly_(iConfig.getParameter<bool>("diagonalElementsOnly")),
17  useVertex_(iConfig.getParameter<bool>("useVertexVariables")) {
18  std::string track = iConfig.getParameter<std::string>("track");
19  if (track == "standAloneMuon")
20  track_ = StaTrack;
21  else if (track == "combinedMuon")
22  track_ = GlbTrack;
23  else if (track == "track")
24  track_ = TrkTrack;
25  else
26  throw cms::Exception("Configuration")
27  << "MatcherByPullsAlgorithm: track '" << track << "' is not known\n"
28  << "Allowed values are: 'track', 'combinedMuon', 'standAloneMuon' (as per reco::RecoCandidate object)\n";
29 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool diagOnly_
Use only the diagonal terms of the covariance matrix.
constexpr int pow(int x)
Definition: conifer.h:24
const reco::Track * track(const reco::RecoCandidate &src) const
Get track out of Candidate, NULL if missing.
double dr2_
DeltaR of the matching cone.
bool useVertex_
Use also dxy / dsz in the matching.
double cut_
Cut on the pull.
TrackChoice track_
Track to be used in matching.

◆ ~MatcherByPullsAlgorithm()

MatcherByPullsAlgorithm::~MatcherByPullsAlgorithm ( )

Definition at line 31 of file MatcherByPullsAlgorithm.cc.

31 {}

Member Function Documentation

◆ fillInvCov()

void MatcherByPullsAlgorithm::fillInvCov ( const reco::Track tk,
AlgebraicSymMatrix55 invCov 
) const

Fill the inverse covariance matrix for the match(track, candidate, invCov) method.

Definition at line 143 of file MatcherByPullsAlgorithm.cc.

References reco::TrackBase::covariance(), diagOnly_, mps_fire::i, dqmiolumiharvest::j, and useVertex_.

Referenced by match(), and matchMany().

143  {
144  if (useVertex_) {
145  invCov = tk.covariance();
146  if (diagOnly_) {
147  for (size_t i = 0; i < 5; ++i) {
148  for (size_t j = i + 1; j < 5; ++j) {
149  invCov(i, j) = 0;
150  }
151  }
152  }
153  invCov.Invert();
154  } else {
155  AlgebraicSymMatrix33 momCov = tk.covariance().Sub<AlgebraicSymMatrix33>(0, 0); // get 3x3 matrix
156  if (diagOnly_) {
157  momCov(0, 1) = 0;
158  momCov(0, 2) = 0;
159  momCov(1, 2) = 0;
160  }
161  momCov.Invert();
162  invCov.Place_at(momCov, 0, 0);
163  }
164 }
bool diagOnly_
Use only the diagonal terms of the covariance matrix.
bool useVertex_
Use also dxy / dsz in the matching.
CovarianceMatrix covariance() const
return track covariance matrix
Definition: TrackBase.h:716
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ match() [1/3]

std::pair< bool, float > MatcherByPullsAlgorithm::match ( const reco::Track tk,
const reco::Candidate mc,
const AlgebraicSymMatrix55 invertedCovariance 
) const

Match Track to MC Candidate, using already inverted covariance matrix Return status of match and pull, or (-1,9e9)

Definition at line 40 of file MatcherByPullsAlgorithm.cc.

References HltBtagPostValidation_cff::c, reco::TrackBase::charge(), gather_cfg::cout, cut_, SiPixelRawToDigiRegional_cfi::deltaPhi, HLTMuonOfflineAnalyzer_cfi::deltaR2, change_name::diff, dr2_, reco::TrackBase::dsz(), reco::TrackBase::dxy(), reco::TrackBase::eta(), mps_fire::i, reco::TrackBase::phi(), reco::TrackBase::pt(), reco::TrackBase::qoverp(), mathSSE::sqrt(), reco::TrackBase::theta(), reco::TrackBase::vx(), reco::TrackBase::vy(), and reco::TrackBase::vz().

Referenced by match(), and matchMany().

42  {
43  if (::deltaR2(tk, c) <= dr2_) {
44  AlgebraicVector5 diff(tk.qoverp() - c.charge() / c.p(),
45  tk.theta() - c.theta(),
46  ::deltaPhi(tk.phi(), c.phi()),
47  tk.dxy(c.vertex()),
48  tk.dsz(c.vertex()));
49  double pull = ROOT::Math::Similarity(diff, invCov);
50 #if 0
51  std::cout << "Tk charge/pt/eta/phi/vx/vy/vz " << tk.charge() << "\t" << tk.pt() << "\t" << tk.eta() << "\t" << tk.phi() << "\t" << tk.vx() << "\t" << tk.vy() << "\t" << tk.vz() << std::endl;
52  std::cout << "MC charge/pt/eta/phi/vx/vy/vz " << c.charge() << "\t" << c.pt() << "\t" << c.eta() << "\t" << c.phi() << "\t" << c.vx() << "\t" << c.vy() << "\t" << c.vz() << std::endl;
53  std::cout << "Delta: " << diff << std::endl;
54  std::cout << "Sigmas: ";
55  for (size_t i = 0; i < 5; ++i) {
56  if (invCov(i,i) == 0) std::cout << "---\t";
57  else std::cout << std::sqrt(1.0/invCov(i,i)) << "\t";
58  }
59  std::cout << std::endl;
60  std::cout << "Items: ";
61  for (size_t i = 0; i < 5; ++i) {
62  if (invCov(i,i) == 0) std::cout << "---\t";
63  else std::cout << diff(i)*std::sqrt(invCov(i,i)) << "\t";
64  }
65  std::cout << std::endl;
66  std::cout << "Pull: " << pull << std::endl;
67 #endif
68  return std::pair<bool, float>(pull < cut_, pull);
69  }
70  return std::pair<bool, float>(false, 9e9);
71 }
double qoverp() const
q / p
Definition: TrackBase.h:599
double vx() const
x coordinate of the reference point on track
Definition: TrackBase.h:655
double dsz() const
dsz parameter (THIS IS NOT the SZ impact parameter to (0,0,0) if refPoint is far from (0...
Definition: TrackBase.h:614
double dr2_
DeltaR of the matching cone.
double pt() const
track transverse momentum
Definition: TrackBase.h:637
int charge() const
track electric charge
Definition: TrackBase.h:596
double vz() const
z coordinate of the reference point on track
Definition: TrackBase.h:661
T sqrt(T t)
Definition: SSEVec.h:19
ROOT::Math::SVector< double, 5 > AlgebraicVector5
double phi() const
azimuthal angle of momentum vector
Definition: TrackBase.h:649
double cut_
Cut on the pull.
double eta() const
pseudorapidity of momentum vector
Definition: TrackBase.h:652
double theta() const
polar angle
Definition: TrackBase.h:602
double vy() const
y coordinate of the reference point on track
Definition: TrackBase.h:658
double dxy() const
dxy parameter. (This is the transverse impact parameter w.r.t. to (0,0,0) ONLY if refPoint is close t...
Definition: TrackBase.h:608

◆ match() [2/3]

std::pair< int, float > MatcherByPullsAlgorithm::match ( const reco::RecoCandidate src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good 
) const

Match Reco Candidate to MC Candidates, skipping the ones which are not good Return index of matchin and pull, or (-1,9e9)

Definition at line 73 of file MatcherByPullsAlgorithm.cc.

References HLT_2023v12_cff::cands, caHitNtupletGeneratorKernels::good, match(), TrackRefitter_38T_cff::src, and track().

75  {
76  const reco::Track *tk = track(src);
77  return (tk == nullptr ? std::pair<int, float>(-1, 9e9) : match(*tk, cands, good));
78 }
const reco::Track * track(const reco::RecoCandidate &src) const
Get track out of Candidate, NULL if missing.
std::pair< bool, float > match(const reco::Track &tk, const reco::Candidate &mc, const AlgebraicSymMatrix55 &invertedCovariance) const

◆ match() [3/3]

std::pair< int, float > MatcherByPullsAlgorithm::match ( const reco::Track src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good 
) const

Match Reco Track to MC Tracks, skipping the ones which are not good Return index of matchin and pull, or (-1,9e9)

Definition at line 80 of file MatcherByPullsAlgorithm.cc.

References HLT_2023v12_cff::cands, fillInvCov(), caHitNtupletGeneratorKernels::good, mps_fire::i, visualization-live-secondInstance_cfg::m, match(), and dqmiodumpmetadata::n.

82  {
83  std::pair<int, float> best(-1, 9e9);
84 
85  AlgebraicSymMatrix55 invCov;
86  fillInvCov(tk, invCov);
87  for (int i = 0, n = cands.size(); i < n; ++i) {
88  if (!good[i])
89  continue;
90  std::pair<bool, float> m = match(tk, cands[i], invCov);
91  if (m.first && (m.second < best.second)) {
92  best.first = i;
93  best.second = m.second;
94  }
95  }
96  return best;
97 }
void fillInvCov(const reco::Track &tk, AlgebraicSymMatrix55 &invCov) const
Fill the inverse covariance matrix for the match(track, candidate, invCov) method.
std::pair< bool, float > match(const reco::Track &tk, const reco::Candidate &mc, const AlgebraicSymMatrix55 &invertedCovariance) const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55

◆ matchMany() [1/2]

void MatcherByPullsAlgorithm::matchMany ( const reco::RecoCandidate src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good,
std::vector< std::pair< double, int > > &  matchesToFill 
) const

Match Reco Candidate to MC Candidates, allowing multiple matches and skipping the ones which are not good It will fill in the vector of <double,int> with pull and index for all matching candidates, already sorted by pulls. This method assumes that matchesToFill is empty when the method is called

Definition at line 99 of file MatcherByPullsAlgorithm.cc.

References HLT_2023v12_cff::cands, caHitNtupletGeneratorKernels::good, TrackRefitter_38T_cff::src, and track().

102  {
103  const reco::Track *tk = track(src);
104  if (tk != nullptr)
105  matchMany(*tk, cands, good, matchesToFill);
106 }
const reco::Track * track(const reco::RecoCandidate &src) const
Get track out of Candidate, NULL if missing.
void matchMany(const reco::RecoCandidate &src, const std::vector< reco::GenParticle > &cands, const std::vector< uint8_t > &good, std::vector< std::pair< double, int > > &matchesToFill) const

◆ matchMany() [2/2]

void MatcherByPullsAlgorithm::matchMany ( const reco::Track src,
const std::vector< reco::GenParticle > &  cands,
const std::vector< uint8_t > &  good,
std::vector< std::pair< double, int > > &  matchesToFill 
) const

Match Reco Track to MC Tracks, allowing multiple matches and skipping the ones which are not good It will fill in the vector of <double,int> with pull and index for all matching candidates, already sorted by pulls. This method assumes that matchesToFill is empty when the method is called

Definition at line 108 of file MatcherByPullsAlgorithm.cc.

References HLT_2023v12_cff::cands, fillInvCov(), caHitNtupletGeneratorKernels::good, mps_fire::i, visualization-live-secondInstance_cfg::m, match(), dqmiodumpmetadata::n, and jetUpdater_cfi::sort.

111  {
112  AlgebraicSymMatrix55 invCov;
113  fillInvCov(tk, invCov);
114  for (int i = 0, n = cands.size(); i < n; ++i) {
115  if (!good[i])
116  continue;
117  std::pair<bool, double> m = match(tk, cands[i], invCov);
118  if (m.first)
119  matchesToFill.push_back(std::make_pair(m.second, i));
120  }
121  std::sort(matchesToFill.begin(), matchesToFill.end());
122 }
void fillInvCov(const reco::Track &tk, AlgebraicSymMatrix55 &invCov) const
Fill the inverse covariance matrix for the match(track, candidate, invCov) method.
std::pair< bool, float > match(const reco::Track &tk, const reco::Candidate &mc, const AlgebraicSymMatrix55 &invertedCovariance) const
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55

◆ track()

const reco::Track * MatcherByPullsAlgorithm::track ( const reco::RecoCandidate src) const
private

Get track out of Candidate, NULL if missing.

Definition at line 131 of file MatcherByPullsAlgorithm.cc.

References cms::cuda::assert(), GlbTrack, StaTrack, track_, and TrkTrack.

Referenced by match(), MatcherByPullsAlgorithm(), and matchMany().

131  {
132  switch (track_) {
133  case StaTrack:
134  return muon.standAloneMuon().isNonnull() ? muon.standAloneMuon().get() : nullptr;
135  case GlbTrack:
136  return muon.combinedMuon().isNonnull() ? muon.combinedMuon().get() : nullptr;
137  case TrkTrack:
138  return muon.track().isNonnull() ? muon.track().get() : nullptr;
139  }
140  assert(false);
141 }
assert(be >=bs)
TrackChoice track_
Track to be used in matching.

Member Data Documentation

◆ cut_

double MatcherByPullsAlgorithm::cut_
private

Cut on the pull.

Definition at line 92 of file MatcherByPullsAlgorithm.h.

Referenced by match().

◆ diagOnly_

bool MatcherByPullsAlgorithm::diagOnly_
private

Use only the diagonal terms of the covariance matrix.

Definition at line 95 of file MatcherByPullsAlgorithm.h.

Referenced by fillInvCov().

◆ dr2_

double MatcherByPullsAlgorithm::dr2_
private

DeltaR of the matching cone.

Definition at line 89 of file MatcherByPullsAlgorithm.h.

Referenced by match().

◆ track_

TrackChoice MatcherByPullsAlgorithm::track_
private

Track to be used in matching.

Definition at line 86 of file MatcherByPullsAlgorithm.h.

Referenced by MatcherByPullsAlgorithm(), and track().

◆ useVertex_

bool MatcherByPullsAlgorithm::useVertex_
private

Use also dxy / dsz in the matching.

Definition at line 98 of file MatcherByPullsAlgorithm.h.

Referenced by fillInvCov().