CMS 3D CMS Logo

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

#include <MuonKinkFinder.h>

Public Member Functions

bool fillTrkKink (reco::MuonQuality &quality, const Trajectory &trajectory) const
 
bool fillTrkKink (reco::MuonQuality &quality, const reco::Track &track) const
 
void init (const edm::EventSetup &iSetup)
 
 MuonKinkFinder (const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
 
 ~MuonKinkFinder ()
 

Private Member Functions

void cropAndInvert (AlgebraicSymMatrix55 &cov) const
 
double getChi2 (const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other) const
 

Private Attributes

bool diagonalOnly_
 use only on-diagonal terms of the covariance matrices More...
 
TrackTransformer refitter_
 Track Transformer. More...
 
bool usePosition_
 if true, use full 5x5 track state; if false, use only the track direction More...
 

Detailed Description

Definition at line 8 of file MuonKinkFinder.h.

Constructor & Destructor Documentation

◆ MuonKinkFinder()

MuonKinkFinder::MuonKinkFinder ( const edm::ParameterSet iConfig,
edm::ConsumesCollector iC 
)

Definition at line 5 of file MuonKinkFinder.cc.

6  : diagonalOnly_(iConfig.getParameter<bool>("diagonalOnly")),
7  usePosition_(iConfig.getParameter<bool>("usePosition")),
8  refitter_(iConfig, iC) {}
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
bool diagonalOnly_
use only on-diagonal terms of the covariance matrices
TrackTransformer refitter_
Track Transformer.
bool usePosition_
if true, use full 5x5 track state; if false, use only the track direction

◆ ~MuonKinkFinder()

MuonKinkFinder::~MuonKinkFinder ( )

Definition at line 10 of file MuonKinkFinder.cc.

10 {}

Member Function Documentation

◆ cropAndInvert()

void MuonKinkFinder::cropAndInvert ( AlgebraicSymMatrix55 cov) const
private

Definition at line 60 of file MuonKinkFinder.cc.

References diagonalOnly_, mps_fire::i, dqmiolumiharvest::j, and usePosition_.

Referenced by getChi2().

60  {
61  if (usePosition_) {
62  if (diagonalOnly_) {
63  for (size_t i = 0; i < 5; ++i) {
64  for (size_t j = i + 1; j < 5; ++j) {
65  cov(i, j) = 0;
66  }
67  }
68  }
69  cov.Invert();
70  } else {
71  // get 3x3 covariance
72  AlgebraicSymMatrix33 momCov = cov.Sub<AlgebraicSymMatrix33>(0, 0); // get 3x3 matrix
73  if (diagonalOnly_) {
74  momCov(0, 1) = 0;
75  momCov(0, 2) = 0;
76  momCov(1, 2) = 0;
77  }
78  // invert
79  momCov.Invert();
80  // place it
81  cov.Place_at(momCov, 0, 0);
82  // zero the rest
83  for (size_t i = 3; i < 5; ++i) {
84  for (size_t j = i; j < 5; ++j) {
85  cov(i, j) = 0;
86  }
87  }
88  }
89 }
bool diagonalOnly_
use only on-diagonal terms of the covariance matrices
bool usePosition_
if true, use full 5x5 track state; if false, use only the track direction
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33

◆ fillTrkKink() [1/2]

bool MuonKinkFinder::fillTrkKink ( reco::MuonQuality quality,
const Trajectory trajectory 
) const

Definition at line 24 of file MuonKinkFinder.cc.

References newFWLiteAna::found, getChi2(), Trajectory::measurements(), findAndChange::post, findAndChange::pre, and quality.

Referenced by fillTrkKink().

24  {
25  const std::vector<TrajectoryMeasurement> &tms = trajectory.measurements();
26  quality.trkKink = -1.0;
27  quality.tkKink_position = math::XYZPoint(0, 0, 0);
28  bool found = false;
29  for (int itm = 3, nm = tms.size() - 3; itm < nm; ++itm) {
30  TrajectoryStateOnSurface pre = tms[itm].forwardPredictedState();
31  TrajectoryStateOnSurface post = tms[itm].backwardPredictedState();
32  if (!pre.isValid() || !post.isValid())
33  continue;
34  found = true;
35  double c2f = getChi2(pre, post);
36  if (c2f > quality.trkKink) {
37  quality.trkKink = c2f;
38  GlobalPoint pos = (tms[itm].updatedState().isValid() ? tms[itm].updatedState() : pre).globalPosition();
39  quality.tkKink_position = math::XYZPoint(pos.x(), pos.y(), pos.z());
40  }
41  }
42  if (!found)
43  quality.trkKink = 999;
44  return found;
45 }
double getChi2(const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other) const
DataContainer const & measurements() const
Definition: Trajectory.h:178
string quality
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12

◆ fillTrkKink() [2/2]

bool MuonKinkFinder::fillTrkKink ( reco::MuonQuality quality,
const reco::Track track 
) const

Definition at line 14 of file MuonKinkFinder.cc.

References fillTrkKink(), quality, refitter_, HLT_2024v13_cff::track, and TrackTransformer::transform().

14  {
15  std::vector<Trajectory> traj = refitter_.transform(track);
16  if (traj.size() != 1) {
17  quality.trkKink = 999;
18  quality.tkKink_position = math::XYZPoint(0, 0, 0);
19  return false;
20  }
21  return fillTrkKink(quality, traj.front());
22 }
TrackTransformer refitter_
Track Transformer.
string quality
bool fillTrkKink(reco::MuonQuality &quality, const Trajectory &trajectory) const
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
std::vector< Trajectory > transform(const reco::Track &) const override
Convert a reco::Track into Trajectory.

◆ getChi2()

double MuonKinkFinder::getChi2 ( const TrajectoryStateOnSurface start,
const TrajectoryStateOnSurface other 
) const
private

Definition at line 47 of file MuonKinkFinder.cc.

References cropAndInvert(), change_name::diff, and trackingPlots::other.

Referenced by fillTrkKink().

47  {
48  if (!start.hasError() && !other.hasError())
49  throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
51  if (start.hasError())
52  cov += start.localError().matrix();
53  if (other.hasError())
54  cov += other.localError().matrix();
55  cropAndInvert(cov);
56  AlgebraicVector5 diff(start.localParameters().mixedFormatVector() - other.localParameters().mixedFormatVector());
57  return ROOT::Math::Similarity(diff, cov);
58 }
Definition: start.py:1
ROOT::Math::SVector< double, 5 > AlgebraicVector5
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
void cropAndInvert(AlgebraicSymMatrix55 &cov) const

◆ init()

void MuonKinkFinder::init ( const edm::EventSetup iSetup)

Definition at line 12 of file MuonKinkFinder.cc.

References refitter_, and TrackTransformer::setServices().

12 { refitter_.setServices(iSetup); }
void setServices(const edm::EventSetup &) override
set the services needed by the TrackTransformer
TrackTransformer refitter_
Track Transformer.

Member Data Documentation

◆ diagonalOnly_

bool MuonKinkFinder::diagonalOnly_
private

use only on-diagonal terms of the covariance matrices

Definition at line 24 of file MuonKinkFinder.h.

Referenced by cropAndInvert().

◆ refitter_

TrackTransformer MuonKinkFinder::refitter_
private

Track Transformer.

Definition at line 29 of file MuonKinkFinder.h.

Referenced by fillTrkKink(), and init().

◆ usePosition_

bool MuonKinkFinder::usePosition_
private

if true, use full 5x5 track state; if false, use only the track direction

Definition at line 26 of file MuonKinkFinder.h.

Referenced by cropAndInvert().