Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

L3 Muon Filter Fix (Master Branch) #18602

Merged
merged 4 commits into from May 8, 2017
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Jump to
Jump to file
Failed to load files.
Diff view
Diff view
3 changes: 3 additions & 0 deletions HLTrigger/Muon/interface/HLTMuonDimuonL3Filter.h
Expand Up @@ -14,6 +14,7 @@
#include "HLTrigger/HLTcore/interface/HLTFilter.h"
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "DataFormats/MuonReco/interface/MuonTrackLinks.h"

namespace edm {
class ConfigurationDescriptions;
Expand Down Expand Up @@ -59,6 +60,8 @@ class HLTMuonDimuonL3Filter : public HLTFilter {
double max_DCAMuMu_; // DCA between the two muons
double max_YPair_; // |rapidity| of pair
bool cutCowboys_; ///< if true, reject muon-track pairs that bend towards each other
const edm::InputTag theL3LinksLabel; //Needed to iterL3
const edm::EDGetTokenT<reco::MuonTrackLinksCollection> linkToken_; //Needed to iterL3
};

#endif //HLTMuonDimuonFilter_h
3 changes: 3 additions & 0 deletions HLTrigger/Muon/interface/HLTMuonTrimuonL3Filter.h
Expand Up @@ -15,6 +15,7 @@
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonReco/interface/MuonTrackLinks.h"

namespace edm {
class ConfigurationDescriptions;
Expand Down Expand Up @@ -56,6 +57,8 @@ class HLTMuonTrimuonL3Filter : public HLTFilter {
double nsigma_Pt_; // pt uncertainty margin (in number of sigmas)
double max_DCAMuMu_; // DCA between the three muons
double max_YTriplet_; // |rapidity| of triplet
const edm::InputTag theL3LinksLabel; //Needed to iterL3
const edm::EDGetTokenT<reco::MuonTrackLinksCollection> linkToken_; //Needed to iterL3
};

#endif //HLTMuonDimuonFilter_h
67 changes: 57 additions & 10 deletions HLTrigger/Muon/src/HLTMuonDimuonL3Filter.cc
Expand Up @@ -22,6 +22,7 @@
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"

#include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
Expand Down Expand Up @@ -70,7 +71,9 @@ HLTMuonDimuonL3Filter::HLTMuonDimuonL3Filter(const edm::ParameterSet& iConfig) :
nsigma_Pt_ (iConfig.getParameter<double> ("NSigmaPt")),
max_DCAMuMu_ (iConfig.getParameter<double>("MaxDCAMuMu")),
max_YPair_ (iConfig.getParameter<double>("MaxRapidityPair")),
cutCowboys_(iConfig.getParameter<bool>("CutCowboys"))
cutCowboys_(iConfig.getParameter<bool>("CutCowboys")),
theL3LinksLabel (iConfig.getParameter<InputTag>("InputLinks")),
linkToken_ (consumes<reco::MuonTrackLinksCollection>(theL3LinksLabel))
{

LogDebug("HLTMuonDimuonL3Filter")
Expand Down Expand Up @@ -132,6 +135,7 @@ HLTMuonDimuonL3Filter::fillDescriptions(edm::ConfigurationDescriptions& descript
desc.add<double>("MaxDCAMuMu",99999.9);
desc.add<double>("MaxRapidityPair",999999.0);
desc.add<bool>("CutCowboys",false);
desc.add<edm::InputTag>("InputLinks",edm::InputTag(""));
descriptions.add("hltMuonDimuonL3Filter",desc);
}

Expand Down Expand Up @@ -161,19 +165,62 @@ HLTMuonDimuonL3Filter::hltFilter(edm::Event& iEvent, const edm::EventSetup& iSet
Handle<RecoChargedCandidateCollection> mucands;
if (saveTags()) filterproduct.addCollectionTag(candTag_);
iEvent.getByToken(candToken_,mucands);

// Test to see if we can use L3MuonTrajectorySeeds:
if (mucands->empty()) return false;
auto const &tk = (*mucands)[0].track();
bool useL3MTS=false;

if (tk->seedRef().isNonnull()){
auto a = dynamic_cast<const L3MuonTrajectorySeed*>(tk->seedRef().get());
useL3MTS = a != nullptr;
}

// sort them by L2Track
std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;
unsigned int maxI = mucands->size();
for (unsigned int i=0;i!=maxI;i++){
TrackRef tk = (*mucands)[i].track();
if (previousCandIsL2_) {
edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
TrackRef staTrack = l3seedRef->l2Track();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands,i));
} else {
L2toL3s[tk].push_back(RecoChargedCandidateRef(mucands,i));

// If we can use L3MuonTrajectory seeds run the older code:
if (useL3MTS){
unsigned int maxI = mucands->size();
for (unsigned int i=0;i!=maxI;i++){
const TrackRef &tk = (*mucands)[i].track();
if (previousCandIsL2_) {
edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
TrackRef staTrack = l3seedRef->l2Track();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands,i));
} else {
L2toL3s[tk].push_back(RecoChargedCandidateRef(mucands,i));
}
}
}
// Using normal TrajectorySeeds:
else{
// Read Links collection:
edm::Handle<reco::MuonTrackLinksCollection> links;
iEvent.getByToken(linkToken_, links);

// Loop over RecoChargedCandidates:
for(unsigned int i(0); i < mucands->size(); ++i){
RecoChargedCandidateRef cand(mucands,i);
for(auto const & link : *links){
TrackRef tk = cand->track();

// Using the same method that was used to create the links between L3 and L2
// ToDo: there should be a better way than dR,dPt matching
const reco::Track& globalTrack = *link.globalTrack();
float dR2 = deltaR2(tk->eta(),tk->phi(),globalTrack.eta(),globalTrack.phi());
float dPt = std::abs(tk->pt() - globalTrack.pt())/tk->pt();
const TrackRef staTrack = link.standAloneTrack();
if (dR2 < 0.02*0.02 and dPt < 0.001 and previousCandIsL2_) {
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
}
else if (not previousCandIsL2_){
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
}
} //MTL loop
} //RCC loop
} //end of using normal TrajectorySeeds


Handle<TriggerFilterObjectWithRefs> previousLevelCands;
iEvent.getByToken(previousCandToken_,previousLevelCands);
Expand Down
129 changes: 31 additions & 98 deletions HLTrigger/Muon/src/HLTMuonL3PreFilter.cc
Expand Up @@ -138,8 +138,8 @@ bool HLTMuonL3PreFilter::hltFilter(Event& iEvent, const EventSetup& iSetup, trig
std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;

// Test to see if we can use L3MuonTrajectorySeeds:
if (mucands->size()<1) return false;
auto tk = (*mucands)[0].track();
if (mucands->empty()) return false;
auto const &tk = (*mucands)[0].track();
bool useL3MTS=false;

if (tk->seedRef().isNonnull()){
Expand All @@ -153,15 +153,42 @@ bool HLTMuonL3PreFilter::hltFilter(Event& iEvent, const EventSetup& iSetup, trig

unsigned int maxI = mucands->size();
for (unsigned int i=0;i!=maxI;++i){
TrackRef tk = (*mucands)[i].track();
const TrackRef &tk = (*mucands)[i].track();
edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
TrackRef staTrack = l3seedRef->l2Track();
LogDebug("HLTMuonL3PreFilter") <<"L2 from: "<<iEvent.getProvenance(staTrack.id()).moduleLabel() <<" index: "<<staTrack.key();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands,i));
}
} //end of useL3MTS

// Using normal TrajectorySeeds:
else{
LogDebug("HLTMuonL3PreFilter") << "HLTMuonL3PreFilter::hltFilter is in mode: not useL3MTS";

// Read Links collection:
edm::Handle<reco::MuonTrackLinksCollection> links;
iEvent.getByToken(linkToken_, links);

// Loop over RecoChargedCandidates:
for(unsigned int i(0); i < mucands->size(); ++i){
RecoChargedCandidateRef cand(mucands,i);
for(auto const & link : *links){
TrackRef tk = cand->track();

// Using the same method that was used to create the links between L3 and L2
// ToDo: there should be a better way than dR,dPt matching
const reco::Track& globalTrack = *link.globalTrack();
float dR2 = deltaR2(tk->eta(),tk->phi(),globalTrack.eta(),globalTrack.phi());
float dPt = std::abs(tk->pt() - globalTrack.pt())/tk->pt();
if (dR2 < 0.02*0.02 and dPt < 0.001) {
const TrackRef staTrack = link.standAloneTrack();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
}
} //MTL loop
} //RCC loop
} //end of using normal TrajectorySeeds

// look at all mucands, check cuts and add to filter object
int n = 0;
auto L2toL3s_it = L2toL3s.begin();
auto L2toL3s_end = L2toL3s.end();
LogDebug("HLTMuonL3PreFilter")<<"looking at: "<<L2toL3s.size()<<" L2->L3s from: "<<mucands->size();
Expand Down Expand Up @@ -251,101 +278,7 @@ bool HLTMuonL3PreFilter::hltFilter(Event& iEvent, const EventSetup& iSetup, trig
LogDebug("HLTMuonL3PreFilter") << " >>>>> Result of HLTMuonL3PreFilter is " << accept << ", number of muons passing thresholds= " << n;

return accept;
} //end of useL3MTS

// Using normal TrajectorySeeds:
else{
LogDebug("HLTMuonL3PreFilter") << "HLTMuonL3PreFilter::hltFilter is in mode: not useL3MTS";

// Read Links collection:
edm::Handle<reco::MuonTrackLinksCollection> links;
iEvent.getByToken(linkToken_, links);

// Loop over RecoChargedCandidates:
for(unsigned int i(0); i < mucands->size(); ++i){
RecoChargedCandidateRef cand(mucands,i);
for(auto const & l : *links){
const reco::MuonTrackLinks* link = &l;
bool useThisLink=false;
TrackRef tk = cand->track();
reco::TrackRef trkTrack = link->trackerTrack();

// Using the same method that was used to create the links
// ToDo: there should be a better way than dR,dPt matching
const reco::Track& globalTrack = *link->globalTrack();
float dR2 = deltaR2(tk->eta(),tk->phi(),globalTrack.eta(),globalTrack.phi());
float dPt = std::abs(tk->pt() - globalTrack.pt())/tk->pt();
if (dR2 < 0.02*0.02 and dPt < 0.001) {
useThisLink=true;
}

if (useThisLink){
const TrackRef staTrack = link->standAloneTrack();
if (!triggeredByLevel2(staTrack,vl2cands)) continue;

// eta cut
if (std::abs(cand->eta())>max_Eta_) continue;

// cut on number of hits
if (tk->numberOfValidHits()<min_Nhits_) continue;

//max dr cut
//if (std::abs(tk->d0())>max_Dr_) continue;
auto dr = std::abs( (- (cand->vx()-beamSpot.x0()) * cand->py() + (cand->vy()-beamSpot.y0()) * cand->px() ) / cand->pt() );
if (dr >max_Dr_) continue;

//min dr cut
if (dr <min_Dr_) continue;

//dz cut
if (std::abs((cand->vz()-beamSpot.z0()) - ((cand->vx()-beamSpot.x0())*cand->px()+(cand->vy()-beamSpot.y0())*cand->py())/cand->pt() * cand->pz()/cand->pt())>max_Dz_) continue;

// dxy significance cut (safeguard against bizarre values)
if (min_DxySig_ > 0 && (tk->dxyError() <= 0 || std::abs(tk->dxy(beamSpot.position())/tk->dxyError()) < min_DxySig_)) continue;

//normalizedChi2 cut
if (tk->normalizedChi2() > max_NormalizedChi2_ ) continue;

//dxy beamspot cut
float absDxy = std::abs(tk->dxy(beamSpot.position()));
if (absDxy > max_DXYBeamSpot_ || absDxy < min_DXYBeamSpot_ ) continue;

//min muon hits cut
const reco::HitPattern& trackHits = tk->hitPattern();
if (trackHits.numberOfValidMuonHits() < min_NmuonHits_ ) continue;

//pt difference cut
double candPt = cand->pt();
double trackPt = tk->pt();

if (std::abs(candPt - trackPt) > max_PtDifference_ ) continue;

//track pt cut
if (trackPt < min_TrackPt_ ) continue;

// Pt threshold cut
double pt = cand->pt();
double err0 = tk->error(0);
double abspar0 = std::abs(tk->parameter(0));
double ptLx = pt;
// convert 50% efficiency threshold to 90% efficiency threshold
if (abspar0>0) ptLx += nsigma_Pt_*err0/abspar0*pt;
LogTrace("HLTMuonL3PreFilter") << " ...Muon in loop, trackkRef pt= "
<< tk->pt() << ", ptLx= " << ptLx
<< " cand pT " << cand->pt();
if (ptLx<min_Pt_) continue;

filterproduct.addObject(TriggerMuon,cand);
n++;
break; // and go on with the next L2 association
} //end of useThisLink
} //end of muons in links collection
} //end of RecoCand collection

// filter decision:
const bool accept (n >= min_N_);
return accept;
} //not useL3MTS
}

bool
Expand Down
57 changes: 50 additions & 7 deletions HLTrigger/Muon/src/HLTMuonTrimuonL3Filter.cc
Expand Up @@ -22,6 +22,7 @@
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"

#include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
#include "TrackingTools/TransientTrack/interface/TransientTrack.h"
Expand Down Expand Up @@ -65,7 +66,9 @@ HLTMuonTrimuonL3Filter::HLTMuonTrimuonL3Filter(const edm::ParameterSet& iConfig)
max_PtBalance_ (iConfig.getParameter<double> ("MaxPtBalance")),
nsigma_Pt_ (iConfig.getParameter<double> ("NSigmaPt")),
max_DCAMuMu_ (iConfig.getParameter<double>("MaxDCAMuMu")),
max_YTriplet_ (iConfig.getParameter<double>("MaxRapidityTriplet"))
max_YTriplet_ (iConfig.getParameter<double>("MaxRapidityTriplet")),
theL3LinksLabel (iConfig.getParameter<InputTag>("InputLinks")),
linkToken_ (consumes<reco::MuonTrackLinksCollection>(theL3LinksLabel))
{

LogDebug("HLTMuonTrimuonL3Filter")
Expand Down Expand Up @@ -113,6 +116,7 @@ HLTMuonTrimuonL3Filter::fillDescriptions(edm::ConfigurationDescriptions& descrip
desc.add<double>("NSigmaPt",0.0);
desc.add<double>("MaxDCAMuMu",99999.9);
desc.add<double>("MaxRapidityTriplet",999999.0);
desc.add<edm::InputTag>("InputLinks",edm::InputTag(""));
descriptions.add("hltMuonTrimuonL3Filter",desc);
}

Expand All @@ -135,15 +139,54 @@ HLTMuonTrimuonL3Filter::hltFilter(edm::Event& iEvent, const edm::EventSetup& iSe
Handle<RecoChargedCandidateCollection> mucands;
if (saveTags()) filterproduct.addCollectionTag(candTag_);
iEvent.getByToken(candToken_,mucands);

// Test to see if we can use L3MuonTrajectorySeeds:
if (mucands->empty()) return false;
auto const &tk = (*mucands)[0].track();
bool useL3MTS=false;

if (tk->seedRef().isNonnull()){
auto a = dynamic_cast<const L3MuonTrajectorySeed*>(tk->seedRef().get());
useL3MTS = a != nullptr;
}

// sort them by L2Track
std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;
unsigned int maxI = mucands->size();
for (unsigned int i=0;i!=maxI;i++){
TrackRef tk = (*mucands)[i].track();
edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
TrackRef staTrack = l3seedRef->l2Track();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands,i));

// If we can use L3MuonTrajectory seeds run the older code:
if (useL3MTS){
unsigned int maxI = mucands->size();
for (unsigned int i=0;i!=maxI;i++){
const TrackRef &tk = (*mucands)[i].track();
edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef = tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
TrackRef staTrack = l3seedRef->l2Track();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands,i));
}
}
// Using normal TrajectorySeeds:
else{
// Read Links collection:
edm::Handle<reco::MuonTrackLinksCollection> links;
iEvent.getByToken(linkToken_, links);

// Loop over RecoChargedCandidates:
for(unsigned int i(0); i < mucands->size(); ++i){
RecoChargedCandidateRef cand(mucands,i);
for(auto const & link : *links){
TrackRef tk = cand->track();

// Using the same method that was used to create the links between L3 and L2
// ToDo: there should be a better way than dR,dPt matching
const reco::Track& globalTrack = *link.globalTrack();
float dR2 = deltaR2(tk->eta(),tk->phi(),globalTrack.eta(),globalTrack.phi());
float dPt = std::abs(tk->pt() - globalTrack.pt())/tk->pt();
const TrackRef staTrack = link.standAloneTrack();
if (dR2 < 0.02*0.02 and dPt < 0.001) {
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
}
} //MTL loop
} //RCC loop
} //end of using normal TrajectorySeeds

Handle<TriggerFilterObjectWithRefs> previousLevelCands;
iEvent.getByToken(previousCandToken_,previousLevelCands);
Expand Down