Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:28:39

0001 #include "RecoTracker/PixelTrackFitting/interface/PixelTrackFilterByKinematics.h"
0002 #include "DataFormats/TrackReco/interface/Track.h"
0003 #include "DataFormats/TrackReco/interface/TrackBase.h"
0004 
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 
0007 template <class T>
0008 T sqr(T t) {
0009   return t * t;
0010 }
0011 
0012 #include <iostream>
0013 
0014 PixelTrackFilterByKinematics::PixelTrackFilterByKinematics(double ptmin, double tipmax, double chi2max)
0015     : theoPtMin(1 / ptmin),
0016       theNSigmaInvPtTolerance(0.),
0017       theTIPMax(tipmax),
0018       theNSigmaTipMaxTolerance(0.),
0019       theChi2Max(chi2max) {}
0020 
0021 PixelTrackFilterByKinematics::PixelTrackFilterByKinematics(
0022     float optmin, float invPtTolerance, float tipmax, float tipmaxTolerance, float chi2max)
0023     : theoPtMin(optmin),
0024       theNSigmaInvPtTolerance(invPtTolerance),
0025       theTIPMax(tipmax),
0026       theNSigmaTipMaxTolerance(tipmaxTolerance),
0027       theChi2Max(chi2max) {}
0028 
0029 PixelTrackFilterByKinematics::~PixelTrackFilterByKinematics() {}
0030 
0031 bool PixelTrackFilterByKinematics::operator()(const reco::Track* ptrack, const PixelTrackFilterBase::Hits& hits) const {
0032   if (!ptrack)
0033     return false;
0034   auto const& track = *ptrack;
0035   if (track.chi2() > theChi2Max)
0036     return false;
0037   if ((std::abs(track.d0()) - theTIPMax) > theNSigmaTipMaxTolerance * track.d0Error())
0038     return false;
0039 
0040   float pt_v = float(track.pt());
0041   float opt_v = 1.f / pt_v;
0042   float pz_v = track.pz();
0043   float p_v = float(track.p());
0044   float op_v = 1.f / p_v;
0045   float cosTheta = pz_v * op_v;
0046   float osinTheta = p_v * opt_v;
0047   float errLambda2 = track.covariance(reco::TrackBase::i_lambda, reco::TrackBase::i_lambda);
0048   float errInvP2 = track.covariance(reco::TrackBase::i_qoverp, reco::TrackBase::i_qoverp);
0049   float covIPtTheta = track.covariance(reco::TrackBase::i_qoverp, reco::TrackBase::i_lambda);
0050   float errInvPt2 =
0051       (errInvP2 + sqr(cosTheta * opt_v) * errLambda2 + 2.f * (cosTheta * opt_v) * covIPtTheta) * sqr(osinTheta);
0052 
0053   return (opt_v - theoPtMin) < theNSigmaInvPtTolerance * std::sqrt(errInvPt2);
0054 }