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

update L3 filter, for 2017 online muon reconstruction #19465

Merged
merged 4 commits into from Jul 5, 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
19 changes: 19 additions & 0 deletions HLTrigger/Muon/interface/HLTMuonL3PreFilter.h
Expand Up @@ -17,8 +17,11 @@
#include "DataFormats/TrackReco/interface/TrackFwd.h"
#include "DataFormats/MuonReco/interface/MuonFwd.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonReco/interface/MuonSelectors.h"
#include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
#include "DataFormats/Math/interface/deltaR.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
#include "FWCore/ParameterSet/interface/ParameterSetDescription.h"

Expand All @@ -32,12 +35,18 @@ class HLTMuonL3PreFilter : public HLTFilter {

private:
bool triggeredByLevel2(const reco::TrackRef& track,std::vector<reco::RecoChargedCandidateRef>& vcands) const;
bool applySelection(const reco::RecoChargedCandidateRef&, const reco::BeamSpot&) const;

const edm::InputTag beamspotTag_ ;
const edm::EDGetTokenT<reco::BeamSpot> beamspotToken_ ;
const edm::InputTag candTag_; // input tag identifying product contains muons
const edm::EDGetTokenT<reco::RecoChargedCandidateCollection> candToken_; // token identifying product contains muons
const edm::InputTag previousCandTag_; // input tag identifying product contains muons passing the previous level
const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> previousCandToken_; // token identifying product contains muons passing the previous level
const edm::InputTag l1CandTag_; // input tag identifying product contains muons passing the L1 level
const edm::EDGetTokenT<trigger::TriggerFilterObjectWithRefs> l1CandToken_; // token identifying product contains muons passing the L1 level
const edm::InputTag recoMuTag_; // input tag identifying reco muons
const edm::EDGetTokenT<reco::MuonCollection> recoMuToken_; // token identifying product contains reco muons
const int min_N_; // minimum number of muons to fire the trigger
const double max_Eta_; // Eta cut
const int min_Nhits_; // threshold on number of hits on muon
Expand All @@ -53,6 +62,16 @@ class HLTMuonL3PreFilter : public HLTFilter {
const int min_NmuonHits_; // cutoff in minumum number of chi2 hits
const double max_PtDifference_; // cutoff in maximum different between global track and tracker track
const double min_TrackPt_; // cutoff in tracker track pt

// cuts for L3FromL1
const int min_MuonStations_L3fromL1_; // cut on min number of stations - for L3FromL1
const unsigned int allowedTypeMask_L3fromL1_;
const unsigned int requiredTypeMask_L3fromL1_;
double maxNormalizedChi2_L3fromL1_;
muon::SelectionType trkMuonId_;
const double L1MatchingdR_;
const bool matchPreviousCand_;

const bool devDebug_;
const edm::InputTag theL3LinksLabel;
const edm::EDGetTokenT<reco::MuonTrackLinksCollection> linkToken_;
Expand Down
254 changes: 168 additions & 86 deletions HLTrigger/Muon/src/HLTMuonL3PreFilter.cc
Expand Up @@ -16,10 +16,9 @@

#include "FWCore/MessageLogger/interface/MessageLogger.h"
#include "DataFormats/TrackReco/interface/Track.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
#include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
#include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
#include "DataFormats/BeamSpot/interface/BeamSpot.h"
#include "DataFormats/MuonReco/interface/Muon.h"

#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
#include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
#include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
Expand All @@ -41,6 +40,10 @@ HLTMuonL3PreFilter::HLTMuonL3PreFilter(const ParameterSet& iConfig) : HLTFilter(
candToken_ (consumes<reco::RecoChargedCandidateCollection>(candTag_)),
previousCandTag_ (iConfig.getParameter<InputTag > ("PreviousCandTag")),
previousCandToken_ (consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
l1CandTag_ (iConfig.getParameter<InputTag > ("L1CandTag")),
l1CandToken_ (consumes<trigger::TriggerFilterObjectWithRefs>(l1CandTag_)),
recoMuTag_ (iConfig.getParameter<InputTag > ("inputMuonCollection")),
recoMuToken_ (consumes<reco::MuonCollection>(recoMuTag_)),
min_N_ (iConfig.getParameter<int> ("MinN")),
max_Eta_ (iConfig.getParameter<double> ("MaxEta")),
min_Nhits_ (iConfig.getParameter<int> ("MinNhits")),
Expand All @@ -56,6 +59,14 @@ HLTMuonL3PreFilter::HLTMuonL3PreFilter(const ParameterSet& iConfig) : HLTFilter(
min_NmuonHits_ (iConfig.getParameter<int> ("MinNmuonHits")),
max_PtDifference_ (iConfig.getParameter<double> ("MaxPtDifference")),
min_TrackPt_ (iConfig.getParameter<double> ("MinTrackPt")),
min_MuonStations_L3fromL1_ (iConfig.getParameter<int>("minMuonStations")),
allowedTypeMask_L3fromL1_ (iConfig.getParameter<unsigned int>("allowedTypeMask")),
requiredTypeMask_L3fromL1_ (iConfig.getParameter<unsigned int>("requiredTypeMask")),
maxNormalizedChi2_L3fromL1_(iConfig.getParameter<double>("MaxNormalizedChi2_L3FromL1")),
trkMuonId_ ( muon::SelectionType(iConfig.getParameter<unsigned int>("trkMuonId"))),
L1MatchingdR_ (iConfig.getParameter<double> ("L1MatchingdR")),
matchPreviousCand_( iConfig.getParameter<bool>("MatchToPreviousCand") ),

devDebug_ (false),
theL3LinksLabel (iConfig.getParameter<InputTag>("InputLinks")),
linkToken_ (consumes<reco::MuonTrackLinksCollection>(theL3LinksLabel))
Expand Down Expand Up @@ -84,6 +95,8 @@ HLTMuonL3PreFilter::fillDescriptions(edm::ConfigurationDescriptions& description
desc.add<edm::InputTag>("CandTag",edm::InputTag("hltL3MuonCandidates"));
// desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag("hltDiMuonL2PreFiltered0"));
desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag(""));
desc.add<edm::InputTag>("L1CandTag",edm::InputTag(""));
desc.add<edm::InputTag>("inputMuonCollection",edm::InputTag(""));
desc.add<int>("MinN",1);
desc.add<double>("MaxEta",2.5);
desc.add<int>("MinNhits",0);
Expand All @@ -99,6 +112,15 @@ HLTMuonL3PreFilter::fillDescriptions(edm::ConfigurationDescriptions& description
desc.add<int>("MinNmuonHits",0);
desc.add<double>("MaxPtDifference",9999.0);
desc.add<double>("MinTrackPt",0.0);
desc.add<int>("minMuonStations",-1);
desc.add<int>("minTrkHits",-1);
desc.add<int>("minMuonHits",-1);
desc.add<unsigned int>("allowedTypeMask",255);
desc.add<unsigned int>("requiredTypeMask",0);
desc.add<double>("MaxNormalizedChi2_L3FromL1",0.0);
desc.add<unsigned int>("trkMuonId",0);
desc.add<double>("L1MatchingdR",0.3);
desc.add<bool>("MatchToPreviousCand", true);
desc.add<edm::InputTag>("InputLinks",edm::InputTag(""));
descriptions.add("hltMuonL3PreFilter",desc);
}
Expand Down Expand Up @@ -134,9 +156,11 @@ bool HLTMuonL3PreFilter::hltFilter(Event& iEvent, const EventSetup& iSetup, trig
// Number of objects passing the L3 Trigger:
int n = 0;

// sort them by L2Track
// sort them by L2Track
std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;

// map the L3 cands matched to a L1 to their position in the recoMuon collection
std::map<unsigned int, RecoChargedCandidateRef > MuonToL3s;

// Test to see if we can use L3MuonTrajectorySeeds:
if (mucands->empty()) return false;
auto const &tk = (*mucands)[0].track();
Expand Down Expand Up @@ -168,101 +192,99 @@ bool HLTMuonL3PreFilter::hltFilter(Event& iEvent, const EventSetup& iSetup, trig
// Read Links collection:
edm::Handle<reco::MuonTrackLinksCollection> links;
iEvent.getByToken(linkToken_, links);


edm::Handle<trigger::TriggerFilterObjectWithRefs> level1Cands;
std::vector<l1t::MuonRef> vl1cands;
std::vector<l1t::MuonRef>::iterator vl1cands_begin;
std::vector<l1t::MuonRef>::iterator vl1cands_end;

bool check_l1match = true;

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

// 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));
}
const reco::Track& trackerTrack = *link.trackerTrack();

float dR2 = deltaR2(tk->eta(),tk->phi(),trackerTrack.eta(),trackerTrack.phi());
float dPt = std::abs(tk->pt() - trackerTrack.pt());
if (tk->pt() != 0) dPt = dPt/tk->pt();

if (matchPreviousCand_ && dR2 < 0.02*0.02 and dPt < 0.001) {
const TrackRef staTrack = link.standAloneTrack();
L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
check_l1match = false;
}
} //MTL loop


if ( !l1CandTag_.label().empty() && check_l1match && matchPreviousCand_){
iEvent.getByToken(l1CandToken_,level1Cands);
level1Cands->getObjects(trigger::TriggerL1Mu,vl1cands);
const unsigned int nL1Muons(vl1cands.size());
for (unsigned int il1=0; il1!=nL1Muons; ++il1) {
if (deltaR(cand->eta(), cand->phi(), vl1cands[il1]->eta(), vl1cands[il1]->phi()) < L1MatchingdR_) {
MuonToL3s[i] = RecoChargedCandidateRef(cand);
}
}
}
} //RCC loop
} //end of using normal TrajectorySeeds

// look at all mucands, check cuts and add to filter object
auto L2toL3s_it = L2toL3s.begin();
auto L2toL3s_end = L2toL3s.end();
LogDebug("HLTMuonL3PreFilter")<<"looking at: "<<L2toL3s.size()<<" L2->L3s from: "<<mucands->size();
for (; L2toL3s_it!=L2toL3s_end; ++L2toL3s_it){

if (!triggeredByLevel2(L2toL3s_it->first,vl2cands)) continue;

//loop over the L3Tk reconstructed for this L2.
unsigned int iTk=0;
unsigned int maxItk=L2toL3s_it->second.size();
for (; iTk!=maxItk; iTk++){

RecoChargedCandidateRef & cand=L2toL3s_it->second[iTk];
TrackRef tk = cand->track();

LogDebug("HLTMuonL3PreFilter") << " Muon in loop, q*pt= " << tk->charge()*tk->pt() <<" (" << cand->charge()*cand->pt() << ") " << ", eta= " << tk->eta() << " (" << cand->eta() << ") " << ", hits= " << tk->numberOfValidHits() << ", d0= " << tk->d0() << ", dz= " << tk->dz();

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

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

//max dr cut
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;
// look at all mucands, check cuts and add to filter object
LogDebug("HLTMuonL3PreFilter")<<"looking at: "<<L2toL3s.size()<<" L2->L3s from: "<<mucands->size();
for ( const auto& L2toL3s_it : L2toL3s ) {

// 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;
if (matchPreviousCand_ && !( triggeredByLevel2(L2toL3s_it.first,vl2cands))) continue;

filterproduct.addObject(TriggerMuon,cand);
n++;
break; // and go on with the next L2 association
}
}////loop over L2s from L3 grouping
//loop over the L3Tk reconstructed for this L2.
unsigned int iTk=0;
unsigned int maxItk=L2toL3s_it.second.size();
for (; iTk!=maxItk; iTk++){

const RecoChargedCandidateRef & cand=L2toL3s_it.second[iTk];
if (! applySelection(cand, beamSpot)) continue;

filterproduct.addObject(TriggerMuon,cand);
n++;
break; // and go on with the next L2 association
}
}////loop over L2s from L3 grouping

// now loop on L3 from L1
edm::Handle<reco::MuonCollection> recomuons;
iEvent.getByToken(recoMuToken_,recomuons);

for ( const auto& MuonToL3s_it : MuonToL3s ) {

const reco::Muon& muon(recomuons->at(MuonToL3s_it.first));

// applys specific cuts for TkMu
if ( (muon.type() & allowedTypeMask_L3fromL1_ ) == 0 ) continue;
if ( (muon.type() & requiredTypeMask_L3fromL1_) != requiredTypeMask_L3fromL1_ ) continue;
if ( muon.numberOfMatchedStations() < min_MuonStations_L3fromL1_ ) continue;
if ( !muon.globalTrack().isNull() ){
if (muon.globalTrack()->normalizedChi2() > maxNormalizedChi2_L3fromL1_) continue;
}
if ( muon.isTrackerMuon() && !muon::isGoodMuon(muon,trkMuonId_) ) continue;

const RecoChargedCandidateRef & cand=MuonToL3s_it.second;
// apply common selection
if (! applySelection(cand, beamSpot)) continue;
filterproduct.addObject(TriggerMuon,cand);
n++;

break; // and go on with the next L3 from L1

}

vector<RecoChargedCandidateRef> vref;
filterproduct.getObjects(TriggerMuon,vref);
for (auto & i : vref) {
Expand Down Expand Up @@ -295,3 +317,63 @@ HLTMuonL3PreFilter::triggeredByLevel2(const TrackRef& staTrack,vector<RecoCharge
return ok;
}

bool
HLTMuonL3PreFilter::applySelection(const RecoChargedCandidateRef& cand, const BeamSpot& beamSpot) const
{
// eta cut
if (std::abs(cand->eta())>max_Eta_) return false;

TrackRef tk = cand->track();
LogDebug("HLTMuonL3PreFilter") << " Muon in loop, q*pt= " << tk->charge()*tk->pt() <<" (" << cand->charge()*cand->pt() << ") " << ", eta= " << tk->eta() << " (" << cand->eta() << ") " << ", hits= " << tk->numberOfValidHits() << ", d0= " << tk->d0() << ", dz= " << tk->dz();

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

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

//min dr cut
if (dr <min_Dr_) return false;

//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_) return false;

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

//normalizedChi2 cut
if (tk->normalizedChi2() > max_NormalizedChi2_ ) return false;

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

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

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

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

//track pt cut
if (trackPt < min_TrackPt_ ) return false;

// 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_) return false;

return true;

}