Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2024-04-06 12:26:59

0001 /**  \class HLTMuonL2SelectorForL3IO
0002  * 
0003  *   L2 muon selector for L3 IO:
0004  *   finds L2 muons not previous converted into (good) L3 muons
0005  *
0006  *   \author  Benjamin Radburn-Smith, Santiago Folgueras - Purdue University
0007  */
0008 
0009 #include "RecoMuon/L3TrackFinder/interface/HLTMuonL2SelectorForL3IO.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0012 #include "DataFormats/Math/interface/deltaR.h"
0013 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
0014 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
0015 
0016 /// constructor with config
0017 HLTMuonL2SelectorForL3IO::HLTMuonL2SelectorForL3IO(const edm::ParameterSet& iConfig)
0018     : l2Src_(consumes<reco::TrackCollection>(iConfig.getParameter<edm::InputTag>("l2Src"))),
0019       l3OISrc_(consumes<reco::RecoChargedCandidateCollection>(iConfig.getParameter<edm::InputTag>("l3OISrc"))),
0020       l3linkToken_(consumes<reco::MuonTrackLinksCollection>(iConfig.getParameter<edm::InputTag>("InputLinks"))),
0021       applyL3Filters_(iConfig.getParameter<bool>("applyL3Filters")),
0022       max_NormalizedChi2_(iConfig.getParameter<double>("MaxNormalizedChi2")),
0023       max_PtDifference_(iConfig.getParameter<double>("MaxPtDifference")),
0024       min_Nhits_(iConfig.getParameter<int>("MinNhits")),
0025       min_NmuonHits_(iConfig.getParameter<int>("MinNmuonHits")) {
0026   LogTrace("Muon|RecoMuon|HLTMuonL2SelectorForL3IO") << "constructor called";
0027   produces<reco::TrackCollection>();
0028 }
0029 
0030 /// destructor
0031 HLTMuonL2SelectorForL3IO::~HLTMuonL2SelectorForL3IO() {}
0032 
0033 /// create collection of L2 muons not already reconstructed as L3 muons
0034 void HLTMuonL2SelectorForL3IO::produce(edm::Event& iEvent, const edm::EventSetup& iSetup) {
0035   const std::string metname = "Muon|RecoMuon|HLTMuonL2SelectorForL3IO";
0036 
0037   //  IN:
0038   edm::Handle<reco::TrackCollection> l2muonH;
0039   iEvent.getByToken(l2Src_, l2muonH);
0040 
0041   edm::Handle<reco::RecoChargedCandidateCollection> l3muonH;
0042   iEvent.getByToken(l3OISrc_, l3muonH);
0043 
0044   // Read Links collection:
0045   edm::Handle<reco::MuonTrackLinksCollection> links;
0046   iEvent.getByToken(l3linkToken_, links);
0047 
0048   //    OUT:
0049   std::unique_ptr<reco::TrackCollection> result(new reco::TrackCollection());
0050 
0051   for (unsigned int il2 = 0; il2 != l2muonH->size(); ++il2) {
0052     reco::TrackRef l2muRef(l2muonH, il2);
0053     bool re_do_this_L2 = true;
0054 
0055     for (unsigned int il3 = 0; il3 != l3muonH->size(); ++il3) {
0056       reco::RecoChargedCandidateRef cand(l3muonH, il3);
0057       reco::TrackRef tk = cand->track();
0058 
0059       bool useThisLink = false;
0060       for (unsigned int l(0); l < links->size() && !useThisLink; ++l) {
0061         const reco::MuonTrackLinks* link = &links->at(l);
0062 
0063         // Check if the L3 link matches the L3 candidate
0064         const reco::Track& globalTrack = *link->globalTrack();
0065         float dR2 = deltaR2(tk->eta(), tk->phi(), globalTrack.eta(), globalTrack.phi());
0066         if (dR2 < 0.02 * 0.02 and std::abs(tk->pt() - globalTrack.pt()) < 0.001 * tk->pt()) {
0067           useThisLink = true;
0068         }
0069 
0070         if (!useThisLink)
0071           continue;
0072 
0073         // Check whether the stand-alone track matches a L2, if not, we will re-use this L2
0074         const reco::TrackRef staTrack = link->standAloneTrack();
0075         if (l2muRef == staTrack)
0076           re_do_this_L2 = false;
0077 
0078         // Check the quality of the reconstructed L3, if poor quality, we will re-use this L2
0079         if (staTrack == l2muRef && applyL3Filters_) {
0080           re_do_this_L2 = true;
0081           const reco::Track& globalTrack = *link->globalTrack();
0082           if (globalTrack.numberOfValidHits() < min_Nhits_)
0083             continue;  // cut on number of hits
0084           if (globalTrack.normalizedChi2() > max_NormalizedChi2_)
0085             continue;  //normalizedChi2 cut
0086           if (globalTrack.hitPattern().numberOfValidMuonHits() < min_NmuonHits_)
0087             continue;  //min muon hits cut
0088           if (std::abs(globalTrack.pt() - l2muRef->pt()) > max_PtDifference_ * globalTrack.pt())
0089             continue;  // pt difference
0090           re_do_this_L2 = false;
0091         }
0092       }
0093     }
0094     if (re_do_this_L2)
0095       result->push_back(*l2muRef);  // used the L2 if no L3 if matched or if the matched L3 has poor quality cuts.
0096   }
0097   iEvent.put(std::move(result));
0098 }
0099 
0100 void HLTMuonL2SelectorForL3IO::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0101   edm::ParameterSetDescription desc;
0102   desc.add<edm::InputTag>("l2Src", edm::InputTag("hltL2Muons", "UpdatedAtVtx"));
0103   desc.add<edm::InputTag>("l3OISrc", edm::InputTag("hltNewOIL3MuonCandidates"));
0104   desc.add<edm::InputTag>("InputLinks", edm::InputTag("hltNewOIL3MuonsLinksCombination"));
0105   desc.add<bool>("applyL3Filters", true);
0106   desc.add<int>("MinNhits", 1);
0107   desc.add<double>("MaxNormalizedChi2", 20.0);
0108   desc.add<int>("MinNmuonHits", 0);
0109   desc.add<double>("MaxPtDifference", 999.0);  //relative difference
0110   descriptions.add("HLTMuonL2SelectorForL3IO", desc);
0111 }