Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:24:35

0001 #include <cmath>
0002 
0003 #include <Math/GenVector/PxPyPzM4D.h>
0004 #include "DataFormats/BTauReco/interface/ParticleMasses.h"
0005 #include "DataFormats/Math/interface/LorentzVector.h"
0006 #include "RecoBTag/SecondaryVertex/interface/V0Filter.h"
0007 
0008 using namespace reco;
0009 
0010 V0Filter::V0Filter(const edm::ParameterSet &params) : k0sMassWindow(params.getParameter<double>("k0sMassWindow")) {}
0011 
0012 bool V0Filter::operator()(const reco::Track *const *tracks, unsigned int n) const {
0013   // only check for K0s for now
0014 
0015   if (n != 2)
0016     return true;
0017 
0018   if (tracks[0]->charge() * tracks[1]->charge() > 0)
0019     return true;
0020 
0021   ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec1;
0022   ROOT::Math::LorentzVector<ROOT::Math::PxPyPzM4D<double> > vec2;
0023 
0024   vec1.SetPx(tracks[0]->px());
0025   vec1.SetPy(tracks[0]->py());
0026   vec1.SetPz(tracks[0]->pz());
0027   vec1.SetM(ParticleMasses::piPlus);
0028 
0029   vec2.SetPx(tracks[1]->px());
0030   vec2.SetPy(tracks[1]->py());
0031   vec2.SetPz(tracks[1]->pz());
0032   vec2.SetM(ParticleMasses::piPlus);
0033 
0034   double invariantMass = (vec1 + vec2).M();
0035   if (std::abs(invariantMass - ParticleMasses::k0) < k0sMassWindow)
0036     return false;
0037 
0038   return true;
0039 }
0040 
0041 bool V0Filter::operator()(const reco::TrackRef *tracks, unsigned int n) const {
0042   std::vector<const reco::Track *> trackPtrs(n);
0043   for (unsigned int i = 0; i < n; i++)
0044     trackPtrs[i] = &*tracks[i];
0045 
0046   return (*this)(&trackPtrs[0], n);
0047 }
0048 
0049 bool V0Filter::operator()(const std::vector<reco::CandidatePtr> &tracks) const {
0050   if (tracks.size() != 2)
0051     return true;
0052 
0053   if (tracks[0]->charge() * tracks[1]->charge() > 0)
0054     return true;
0055 
0056   double invariantMass = (tracks[0]->p4() + tracks[1]->p4()).M();
0057   if (std::abs(invariantMass - ParticleMasses::k0) < k0sMassWindow)
0058     return false;
0059 
0060   return true;
0061 }
0062 
0063 bool V0Filter::operator()(const std::vector<const reco::Track *> &tracks) const {
0064   return (*this)(&tracks[0], tracks.size());
0065 }
0066 
0067 bool V0Filter::operator()(const reco::Track *tracks, unsigned int n) const {
0068   std::vector<const reco::Track *> trackPtrs(n);
0069   for (unsigned int i = 0; i < n; i++)
0070     trackPtrs[i] = &tracks[i];
0071 
0072   return (*this)(&trackPtrs[0], n);
0073 }