CMS 3D CMS Logo

MuonKinkFinder.cc
Go to the documentation of this file.
4 
6  : diagonalOnly_(iConfig.getParameter<bool>("diagonalOnly")),
7  usePosition_(iConfig.getParameter<bool>("usePosition")),
8  refitter_(iConfig, iC) {}
9 
11 
13 
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 }
23 
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 }
46 
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 }
59 
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 }
change_name.diff
diff
Definition: change_name.py:13
MuonKinkFinder::diagonalOnly_
bool diagonalOnly_
use only on-diagonal terms of the covariance matrices
Definition: MuonKinkFinder.h:24
AlgebraicSymMatrix33
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
Definition: AlgebraicROOTObjects.h:21
electrons_cff.bool
bool
Definition: electrons_cff.py:393
mps_fire.i
i
Definition: mps_fire.py:428
MuonKinkFinder::cropAndInvert
void cropAndInvert(AlgebraicSymMatrix55 &cov) const
Definition: MuonKinkFinder.cc:60
reco::MuonQuality
Definition: MuonQuality.h:6
start
Definition: start.py:1
HLT_FULL_cff.track
track
Definition: HLT_FULL_cff.py:11779
TrackTransformer::setServices
void setServices(const edm::EventSetup &) override
set the services needed by the TrackTransformer
Definition: TrackTransformer.cc:79
pos
Definition: PixelAliasList.h:18
TrackTransformer::transform
std::vector< Trajectory > transform(const reco::Track &) const override
Convert a reco::Track into Trajectory.
Definition: TrackTransformer.cc:182
MuonKinkFinder::init
void init(const edm::EventSetup &iSetup)
Definition: MuonKinkFinder.cc:12
newFWLiteAna.found
found
Definition: newFWLiteAna.py:118
TrajectoryStateOnSurface
Definition: TrajectoryStateOnSurface.h:16
MuonKinkFinder::MuonKinkFinder
MuonKinkFinder(const edm::ParameterSet &iConfig, edm::ConsumesCollector &iC)
Definition: MuonKinkFinder.cc:5
reco::Track
Definition: Track.h:27
MuonKinkFinder::refitter_
TrackTransformer refitter_
Track Transformer.
Definition: MuonKinkFinder.h:29
trackingPlots.other
other
Definition: trackingPlots.py:1467
Point3DBase< float, GlobalTag >
MuonKinkFinder::fillTrkKink
bool fillTrkKink(reco::MuonQuality &quality, const Trajectory &trajectory) const
Definition: MuonKinkFinder.cc:24
edm::ParameterSet
Definition: ParameterSet.h:47
math::XYZPoint
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
edm::EventSetup
Definition: EventSetup.h:57
AlgebraicVector5
ROOT::Math::SVector< double, 5 > AlgebraicVector5
Definition: AlgebraicROOTObjects.h:14
Trajectory::measurements
DataContainer const & measurements() const
Definition: Trajectory.h:178
Trajectory.h
qcdUeDQM_cfi.quality
quality
Definition: qcdUeDQM_cfi.py:31
Trajectory
Definition: Trajectory.h:38
MuonKinkFinder::getChi2
double getChi2(const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other) const
Definition: MuonKinkFinder.cc:47
MuonKinkFinder.h
MuonKinkFinder::~MuonKinkFinder
~MuonKinkFinder()
Definition: MuonKinkFinder.cc:10
cms::Exception
Definition: Exception.h:70
ParameterSet.h
dqmiolumiharvest.j
j
Definition: dqmiolumiharvest.py:66
MuonKinkFinder::usePosition_
bool usePosition_
if true, use full 5x5 track state; if false, use only the track direction
Definition: MuonKinkFinder.h:26
AlgebraicSymMatrix55
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
Definition: AlgebraicROOTObjects.h:23
edm::ConsumesCollector
Definition: ConsumesCollector.h:45
TrajectoryStateOnSurface::isValid
bool isValid() const
Definition: TrajectoryStateOnSurface.h:54