Skip to content

Commit

Permalink
Merge pull request #21333 from cms-bph-trigger/robustVertexFilters-fr…
Browse files Browse the repository at this point in the history
…om-CMSSW_9_2_X_2017-11-15-2300

Make the HLT vertex filters robust to L3Muon duplicates - 92X backport
  • Loading branch information
cmsbuild committed Nov 21, 2017
2 parents 94488b5 + 480e414 commit afe6014
Show file tree
Hide file tree
Showing 2 changed files with 30 additions and 29 deletions.
47 changes: 24 additions & 23 deletions HLTrigger/btau/src/HLTmumutkFilter.cc
Expand Up @@ -131,43 +131,44 @@ bool HLTmumutkFilter::hltFilter(edm::Event& iEvent, const edm::EventSetup& iSetu
reco::RecoChargedCandidateCollection::const_iterator tkcand ;

int iFoundRefs = 0;
bool threeMuons = false;
bool track1Matched = false;
bool track2Matched = false;
bool track3Matched = false;
for (auto cand=mucands->begin(); cand!=mucands->end(); cand++) {
reco::TrackRef tkRef = cand->get<reco::TrackRef>();
if (tkRef == vertextkRef1 && iFoundRefs==0) {mucand1 = cand; iFoundRefs++;}
else if(tkRef == vertextkRef1 && iFoundRefs==1) {mucand2 = cand; iFoundRefs++;}
else if(tkRef == vertextkRef1 && iFoundRefs==2) {threeMuons = true;}
if (tkRef == vertextkRef2 && iFoundRefs==0) {mucand1 = cand; iFoundRefs++;}
else if(tkRef == vertextkRef2 && iFoundRefs==1) {mucand2 = cand; iFoundRefs++;}
else if(tkRef == vertextkRef2 && iFoundRefs==2) {threeMuons = true;}
if (tkRef == vertextkRef3 && iFoundRefs==0) {mucand1 = cand; iFoundRefs++;}
else if(tkRef == vertextkRef3 && iFoundRefs==1) {mucand2 = cand; iFoundRefs++;}
else if(tkRef == vertextkRef3 && iFoundRefs==2) {threeMuons = true;}
if (!track1Matched) {
if (tkRef == vertextkRef1 && iFoundRefs==0) {mucand1 = cand; iFoundRefs++; track1Matched = true;}
else if(tkRef == vertextkRef1 && iFoundRefs==1) {mucand2 = cand; iFoundRefs++; track1Matched = true;}
}
if (!track2Matched) {
if (tkRef == vertextkRef2 && iFoundRefs==0) {mucand1 = cand; iFoundRefs++; track2Matched = true;}
else if(tkRef == vertextkRef2 && iFoundRefs==1) {mucand2 = cand; iFoundRefs++; track2Matched = true;}
}
if (!track3Matched) {
if (tkRef == vertextkRef3 && iFoundRefs==0) {mucand1 = cand; iFoundRefs++; track3Matched = true;}
else if(tkRef == vertextkRef3 && iFoundRefs==1) {mucand2 = cand; iFoundRefs++; track3Matched = true;}
}
}
if(threeMuons) throw cms::Exception("BadLogic") << "HLTmumutkFilterr: ERROR: the vertex must have "
<< " exactly two muons by definition." << std::endl;
if(iFoundRefs < 2) throw cms::Exception("BadLogic") << "HLTmumutkFilterr: ERROR: the vertex must have "
<< " at least two muons by definition." << std::endl;

bool twoTrks = false;
int iTrkFoundRefs = 0;
for (auto cand=trkcands->begin(); cand!=trkcands->end(); cand++) {
reco::TrackRef tkRef = cand->get<reco::TrackRef>();
if (tkRef == vertextkRef1 && iTrkFoundRefs==0) {tkcand = cand; iTrkFoundRefs++;}
else if(tkRef == vertextkRef1 && iTrkFoundRefs==1) {twoTrks = true;}
if (tkRef == vertextkRef2 && iTrkFoundRefs==0) {tkcand = cand; iTrkFoundRefs++;}
else if(tkRef == vertextkRef2 && iTrkFoundRefs==1) {twoTrks = true;}
if (tkRef == vertextkRef3 && iTrkFoundRefs==0) {tkcand = cand; iTrkFoundRefs++;}
else if(tkRef == vertextkRef3 && iTrkFoundRefs==1) {twoTrks = true;}
if (tkRef == vertextkRef1 && iTrkFoundRefs==0 && !track1Matched) {tkcand = cand; iTrkFoundRefs++; track1Matched = true; break;}
if (tkRef == vertextkRef2 && iTrkFoundRefs==0 && !track2Matched) {tkcand = cand; iTrkFoundRefs++; track2Matched = true; break;}
if (tkRef == vertextkRef3 && iTrkFoundRefs==0 && !track3Matched) {tkcand = cand; iTrkFoundRefs++; track3Matched = true; break;}
}
if(twoTrks) throw cms::Exception("BadLogic") << "HLTmumutkFilterr: ERROR: the vertex must have "
<< " exactly one track by definition." << std::endl;
if(iTrkFoundRefs == 0) throw cms::Exception("BadLogic") << "HLTmumutkFilterr: ERROR: the vertex must have "
<< " at least one track by definition." << std::endl;

// calculate three-track transverse momentum
math::XYZVector pperp(mucand1->px() + mucand2->px() + tkcand->px(),
mucand1->py() + mucand2->py() + tkcand->py(),
0.);

// get vertex position and error to calculate the decay length significance
reco::Vertex::Point vpoint=displacedVertex.position();
const reco::Vertex::Point& vpoint=displacedVertex.position();
reco::Vertex::Error verr = displacedVertex.error();
GlobalPoint secondaryVertex (vpoint.x(), vpoint.y(), vpoint.z());
GlobalError err(verr.At(0,0), verr.At(1,0), verr.At(1,1), verr.At(2,0), verr.At(2,1), verr.At(2,2) );
Expand Down Expand Up @@ -209,4 +210,4 @@ bool HLTmumutkFilter::triggerdByPreviousLevel(const reco::RecoChargedCandidateRe
if (candref == vcands[i]) return true;
}
return false;
}
}
12 changes: 6 additions & 6 deletions HLTrigger/btau/src/HLTmumutktkFilter.cc
Expand Up @@ -134,25 +134,25 @@ bool HLTmumutktkFilter::hltFilter(edm::Event& iEvent, const edm::EventSetup& iSe
if (tkRef == iVec) {mucandVec.push_back(cand); break;}
}
}
if(mucandVec.size()!= 2) throw cms::Exception("BadLogic") << "HLTmumutktkFilterr: ERROR: the vertex must have "
<< " exactly two muons by definition." << std::endl;
if(mucandVec.size() < 2) throw cms::Exception("BadLogic") << "HLTmumutktkFilterr: ERROR: the vertex must have "
<< " at least two muons by definition." << std::endl;

for (auto cand=trkcands->begin(); cand!=trkcands->end(); cand++) {
reco::TrackRef tkRef = cand->get<reco::TrackRef>();
for (auto & iVec : vertexTrksRefVec){
if (tkRef == iVec) {trkcandVec.push_back(cand); break;}
}
}
if(trkcandVec.size()!= 2 ) throw cms::Exception("BadLogic") << "HLTmumutktkFilterr: ERROR: the vertex must have "
<< " exactly two tracks by definition." << std::endl;
if(trkcandVec.size() < 2 ) throw cms::Exception("BadLogic") << "HLTmumutktkFilterr: ERROR: the vertex must have "
<< " at least two tracks by definition." << std::endl;

// calculate four-track transverse momentum
math::XYZVector pperp(mucandVec.at(0)->px() + mucandVec.at(1)->px() + trkcandVec.at(0)->px() + trkcandVec.at(1)->px(),
mucandVec.at(0)->py() + mucandVec.at(1)->py() + trkcandVec.at(0)->px() + trkcandVec.at(1)->py(),
0.);

// get vertex position and error to calculate the decay length significance
reco::Vertex::Point vpoint=displacedVertex.position();
const reco::Vertex::Point& vpoint=displacedVertex.position();
reco::Vertex::Error verr = displacedVertex.error();
GlobalPoint secondaryVertex (vpoint.x(), vpoint.y(), vpoint.z());
GlobalError err(verr.At(0,0), verr.At(1,0), verr.At(1,1), verr.At(2,0), verr.At(2,1), verr.At(2,2) );
Expand Down Expand Up @@ -198,4 +198,4 @@ bool HLTmumutktkFilter::triggerdByPreviousLevel(const reco::RecoChargedCandidate
if (candref == vcands[i]) return true;
}
return false;
}
}

0 comments on commit afe6014

Please sign in to comment.