CMS 3D CMS Logo

CMSSW_4_4_3_patch1/src/RecoMuon/MuonIdentification/src/MuonKinkFinder.cc

Go to the documentation of this file.
00001 #include "RecoMuon/MuonIdentification/interface/MuonKinkFinder.h"
00002 #include "FWCore/ParameterSet/interface/ParameterSet.h"
00003 #include "TrackingTools/PatternTools/interface/Trajectory.h"
00004 
00005 MuonKinkFinder::MuonKinkFinder(const edm::ParameterSet &iConfig) :
00006     diagonalOnly_(iConfig.getParameter<bool>("diagonalOnly")),
00007     usePosition_(iConfig.getParameter<bool>("usePosition")),
00008     refitter_(iConfig)
00009 { 
00010 }
00011 
00012 MuonKinkFinder::~MuonKinkFinder() { }
00013 
00014 void MuonKinkFinder::init(const edm::EventSetup &iSetup) { 
00015     refitter_.setServices(iSetup);
00016 }
00017 
00018 bool MuonKinkFinder::fillTrkKink(reco::MuonQuality & quality, const reco::Track &track) const {
00019     std::vector<Trajectory> traj  = refitter_.transform(track);
00020     if (traj.size() != 1) {
00021         quality.trkKink = -1.0;
00022         quality.tkKink_position = math::XYZPoint(0,0,0);
00023         return false;
00024     }
00025     return fillTrkKink(quality, traj.front());
00026 }
00027 
00028 bool MuonKinkFinder::fillTrkKink(reco::MuonQuality & quality, const Trajectory &trajectory) const {
00029     const std::vector<TrajectoryMeasurement> &tms = trajectory.measurements();
00030     quality.trkKink = -1.0;
00031     quality.tkKink_position = math::XYZPoint(0,0,0);
00032     bool found = false; 
00033     for (int itm = 3, nm = tms.size() - 3;  itm < nm; ++itm) {
00034         TrajectoryStateOnSurface pre  = tms[itm].forwardPredictedState();
00035         TrajectoryStateOnSurface post = tms[itm].backwardPredictedState();
00036         if (!pre.isValid() || !post.isValid()) continue;
00037         found = true;
00038         double c2f = getChi2(pre,post);
00039         if (c2f > quality.trkKink) { 
00040             quality.trkKink = c2f; 
00041             GlobalPoint pos = (tms[itm].updatedState().isValid() ? tms[itm].updatedState() : pre).globalPosition();
00042             quality.tkKink_position = math::XYZPoint(pos.x(), pos.y(), pos.z());
00043         }
00044     }
00045     return found;
00046 }
00047 
00048 double MuonKinkFinder::getChi2(const TrajectoryStateOnSurface &start, const TrajectoryStateOnSurface &other) const {
00049     if (!start.hasError() && !other.hasError()) throw cms::Exception("LogicError") << "At least one of the two states must have errors to make chi2s.\n";
00050     AlgebraicSymMatrix55 cov;
00051     if (start.hasError()) cov += start.localError().matrix();
00052     if (other.hasError()) cov += other.localError().matrix();
00053     cropAndInvert(cov);
00054     AlgebraicVector5 diff(start.localParameters().mixedFormatVector() - other.localParameters().mixedFormatVector());
00055     return ROOT::Math::Similarity(diff, cov);
00056 }
00057 
00058 void MuonKinkFinder::cropAndInvert(AlgebraicSymMatrix55 &cov) const {
00059     if (usePosition_) {
00060         if (diagonalOnly_) {
00061             for (size_t i = 0; i < 5; ++i) { for (size_t j = i+1; j < 5; ++j) {
00062                 cov(i,j) = 0;
00063             } }
00064         }
00065         cov.Invert();
00066     } else {
00067         // get 3x3 covariance
00068         AlgebraicSymMatrix33 momCov = cov.Sub<AlgebraicSymMatrix33>(0,0); // get 3x3 matrix
00069         if (diagonalOnly_) { momCov(0,1) = 0; momCov(0,2) = 0; momCov(1,2) = 0; }
00070         // invert        
00071         momCov.Invert();
00072         // place it
00073         cov.Place_at(momCov,0,0);
00074         // zero the rest
00075         for (size_t i = 3; i < 5; ++i) { for (size_t j = i; j < 5; ++j) { 
00076             cov(i,j) = 0;
00077         } }
00078     }
00079 }
00080