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)
9 {
10 }
11 
13 
14 void MuonKinkFinder::init(const edm::EventSetup &iSetup) {
15  refitter_.setServices(iSetup);
16 }
17 
19  std::vector<Trajectory> traj = refitter_.transform(track);
20  if (traj.size() != 1) {
21  quality.trkKink = 999;
22  quality.tkKink_position = math::XYZPoint(0,0,0);
23  return false;
24  }
25  return fillTrkKink(quality, traj.front());
26 }
27 
29  const std::vector<TrajectoryMeasurement> &tms = trajectory.measurements();
30  quality.trkKink = -1.0;
31  quality.tkKink_position = math::XYZPoint(0,0,0);
32  bool found = false;
33  for (int itm = 3, nm = tms.size() - 3; itm < nm; ++itm) {
34  TrajectoryStateOnSurface pre = tms[itm].forwardPredictedState();
35  TrajectoryStateOnSurface post = tms[itm].backwardPredictedState();
36  if (!pre.isValid() || !post.isValid()) continue;
37  found = true;
38  double c2f = getChi2(pre,post);
39  if (c2f > quality.trkKink) {
40  quality.trkKink = c2f;
41  GlobalPoint pos = (tms[itm].updatedState().isValid() ? tms[itm].updatedState() : pre).globalPosition();
42  quality.tkKink_position = math::XYZPoint(pos.x(), pos.y(), pos.z());
43  }
44  }
45  if (!found) quality.trkKink = 999;
46  return found;
47 }
48 
50  if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
52  if (start.hasError()) cov += start.localError().matrix();
53  if (other.hasError()) cov += other.localError().matrix();
54  cropAndInvert(cov);
56  return ROOT::Math::Similarity(diff, cov);
57 }
58 
60  if (usePosition_) {
61  if (diagonalOnly_) {
62  for (size_t i = 0; i < 5; ++i) { for (size_t j = i+1; j < 5; ++j) {
63  cov(i,j) = 0;
64  } }
65  }
66  cov.Invert();
67  } else {
68  // get 3x3 covariance
69  AlgebraicSymMatrix33 momCov = cov.Sub<AlgebraicSymMatrix33>(0,0); // get 3x3 matrix
70  if (diagonalOnly_) { momCov(0,1) = 0; momCov(0,2) = 0; momCov(1,2) = 0; }
71  // invert
72  momCov.Invert();
73  // place it
74  cov.Place_at(momCov,0,0);
75  // zero the rest
76  for (size_t i = 3; i < 5; ++i) { for (size_t j = i; j < 5; ++j) {
77  cov(i,j) = 0;
78  } }
79  }
80 }
81 
Definition: start.py:1
double getChi2(const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other) const
void cropAndInvert(AlgebraicSymMatrix55 &cov) const
bool diagonalOnly_
use only on-diagonal terms of the covariance matrices
const LocalTrajectoryParameters & localParameters() const
TrackTransformer refitter_
Track Transformer.
T y() const
Definition: PV3DBase.h:63
float trkKink
value of the kink algorithm applied to the inner track stub
Definition: MuonQuality.h:11
MuonKinkFinder(const edm::ParameterSet &iConfig)
ROOT::Math::SMatrix< double, 5, 5, ROOT::Math::MatRepSym< double, 5 > > AlgebraicSymMatrix55
DataContainer const & measurements() const
Definition: Trajectory.h:196
ROOT::Math::SMatrix< double, 3, 3, ROOT::Math::MatRepSym< double, 3 > > AlgebraicSymMatrix33
void init(const edm::EventSetup &iSetup)
T z() const
Definition: PV3DBase.h:64
const AlgebraicSymMatrix55 & matrix() const
const LocalTrajectoryError & localError() const
virtual void setServices(const edm::EventSetup &)
set the services needed by the TrackTransformer
ROOT::Math::SVector< double, 5 > AlgebraicVector5
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
AlgebraicVector5 mixedFormatVector() const
math::XYZPoint tkKink_position
Kink position for the tracker stub and global track.
Definition: MuonQuality.h:32
virtual std::vector< Trajectory > transform(const reco::Track &) const
Convert a reco::Track into Trajectory.
bool usePosition_
if true, use full 5x5 track state; if false, use only the track direction
T x() const
Definition: PV3DBase.h:62
bool fillTrkKink(reco::MuonQuality &quality, const Trajectory &trajectory) const