Back to home page

Project CMSSW displayed by LXR

 
 

    


File indexing completed on 2022-03-23 12:04:25

0001 /** \class HLTMuonL3PreFilter
0002  *
0003  * See header file for documentation
0004  *
0005  *  \author J. Alcaraz, J-R Vlimant
0006  *
0007  */
0008 
0009 #include "HLTMuonL3PreFilter.h"
0010 
0011 #include "DataFormats/Common/interface/Handle.h"
0012 #include "DataFormats/Common/interface/RefToBase.h"
0013 
0014 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0015 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
0016 
0017 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0018 #include "DataFormats/TrackReco/interface/Track.h"
0019 #include "DataFormats/MuonReco/interface/MuonTrackLinks.h"
0020 #include "DataFormats/MuonReco/interface/Muon.h"
0021 
0022 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
0023 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
0024 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0025 
0026 #include "FWCore/Utilities/interface/InputTag.h"
0027 
0028 //
0029 // constructors and destructor
0030 //
0031 using namespace std;
0032 using namespace edm;
0033 using namespace reco;
0034 using namespace trigger;
0035 
0036 HLTMuonL3PreFilter::HLTMuonL3PreFilter(const ParameterSet& iConfig)
0037     : HLTFilter(iConfig),
0038       propSetup_(iConfig, consumesCollector()),
0039       beamspotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
0040       beamspotToken_(consumes<reco::BeamSpot>(beamspotTag_)),
0041       candTag_(iConfig.getParameter<InputTag>("CandTag")),
0042       candToken_(consumes<reco::RecoChargedCandidateCollection>(candTag_)),
0043       previousCandTag_(iConfig.getParameter<InputTag>("PreviousCandTag")),
0044       previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
0045       l1CandTag_(iConfig.getParameter<InputTag>("L1CandTag")),
0046       l1CandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(l1CandTag_)),
0047       recoMuTag_(iConfig.getParameter<InputTag>("inputMuonCollection")),
0048       recoMuToken_(consumes<reco::MuonCollection>(recoMuTag_)),
0049       min_N_(iConfig.getParameter<int>("MinN")),
0050       max_Eta_(iConfig.getParameter<double>("MaxEta")),
0051       min_Nhits_(iConfig.getParameter<int>("MinNhits")),
0052       max_Dr_(iConfig.getParameter<double>("MaxDr")),
0053       min_Dr_(iConfig.getParameter<double>("MinDr")),
0054       max_Dz_(iConfig.getParameter<double>("MaxDz")),
0055       min_DxySig_(iConfig.getParameter<double>("MinDxySig")),
0056       min_Pt_(iConfig.getParameter<double>("MinPt")),
0057       nsigma_Pt_(iConfig.getParameter<double>("NSigmaPt")),
0058       max_NormalizedChi2_(iConfig.getParameter<double>("MaxNormalizedChi2")),
0059       max_DXYBeamSpot_(iConfig.getParameter<double>("MaxDXYBeamSpot")),
0060       min_DXYBeamSpot_(iConfig.getParameter<double>("MinDXYBeamSpot")),
0061       min_NmuonHits_(iConfig.getParameter<int>("MinNmuonHits")),
0062       max_PtDifference_(iConfig.getParameter<double>("MaxPtDifference")),
0063       min_TrackPt_(iConfig.getParameter<double>("MinTrackPt")),
0064       min_MuonStations_L3fromL1_(iConfig.getParameter<int>("minMuonStations")),
0065       allowedTypeMask_L3fromL1_(iConfig.getParameter<unsigned int>("allowedTypeMask")),
0066       requiredTypeMask_L3fromL1_(iConfig.getParameter<unsigned int>("requiredTypeMask")),
0067       maxNormalizedChi2_L3fromL1_(iConfig.getParameter<double>("MaxNormalizedChi2_L3FromL1")),
0068       trkMuonId_(muon::SelectionType(iConfig.getParameter<unsigned int>("trkMuonId"))),
0069       L1MatchingdR_(iConfig.getParameter<double>("L1MatchingdR")),
0070       L1MatchingdR2_(L1MatchingdR_ * L1MatchingdR_),
0071       matchPreviousCand_(iConfig.getParameter<bool>("MatchToPreviousCand")),
0072 
0073       devDebug_(false),
0074       theL3LinksLabel(iConfig.getParameter<InputTag>("InputLinks")),
0075       linkToken_(consumes<reco::MuonTrackLinksCollection>(theL3LinksLabel)) {
0076   if (L1MatchingdR_ <= 0.) {
0077     throw cms::Exception("HLTMuonL3PreFilterConfiguration")
0078         << "invalid value for parameter \"L1MatchingdR\" (must be > 0): " << L1MatchingdR_;
0079   }
0080   LogDebug("HLTMuonL3PreFilter") << " CandTag/MinN/MaxEta/MinNhits/MaxDr/MinDr/MaxDz/MinDxySig/MinPt/NSigmaPt : "
0081                                  << candTag_.encode() << " " << min_N_ << " " << max_Eta_ << " " << min_Nhits_ << " "
0082                                  << max_Dr_ << " " << min_Dr_ << " " << max_Dz_ << " " << min_DxySig_ << " " << min_Pt_
0083                                  << " " << nsigma_Pt_;
0084 }
0085 
0086 HLTMuonL3PreFilter::~HLTMuonL3PreFilter() = default;
0087 
0088 void HLTMuonL3PreFilter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0089   edm::ParameterSetDescription desc;
0090   makeHLTFilterDescription(desc);
0091   desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOfflineBeamSpot"));
0092   desc.add<edm::InputTag>("CandTag", edm::InputTag("hltL3MuonCandidates"));
0093   //  desc.add<edm::InputTag>("PreviousCandTag",edm::InputTag("hltDiMuonL2PreFiltered0"));
0094   desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
0095   desc.add<edm::InputTag>("L1CandTag", edm::InputTag(""));
0096   desc.add<edm::InputTag>("inputMuonCollection", edm::InputTag(""));
0097   desc.add<int>("MinN", 1);
0098   desc.add<double>("MaxEta", 2.5);
0099   desc.add<int>("MinNhits", 0);
0100   desc.add<double>("MaxDr", 2.0);
0101   desc.add<double>("MinDr", -1.0);
0102   desc.add<double>("MaxDz", 9999.0);
0103   desc.add<double>("MinDxySig", -1.0);
0104   desc.add<double>("MinPt", 3.0);
0105   desc.add<double>("NSigmaPt", 0.0);
0106   desc.add<double>("MaxNormalizedChi2", 9999.0);
0107   desc.add<double>("MaxDXYBeamSpot", 9999.0);
0108   desc.add<double>("MinDXYBeamSpot", -1.0);
0109   desc.add<int>("MinNmuonHits", 0);
0110   desc.add<double>("MaxPtDifference", 9999.0);
0111   desc.add<double>("MinTrackPt", 0.0);
0112   desc.add<int>("minMuonStations", -1);
0113   desc.add<int>("minTrkHits", -1);
0114   desc.add<int>("minMuonHits", -1);
0115   desc.add<unsigned int>("allowedTypeMask", 255);
0116   desc.add<unsigned int>("requiredTypeMask", 0);
0117   desc.add<double>("MaxNormalizedChi2_L3FromL1", 9999.);
0118   desc.add<unsigned int>("trkMuonId", 0);
0119   desc.add<double>("L1MatchingdR", 0.3);
0120   desc.add<bool>("MatchToPreviousCand", true);
0121   desc.add<edm::InputTag>("InputLinks", edm::InputTag(""));
0122   PropagateToMuonSetup::fillPSetDescription(desc);
0123   descriptions.add("hltMuonL3PreFilter", desc);
0124 }
0125 
0126 //
0127 // member functions
0128 //
0129 
0130 // ------------ method called to produce the data  ------------
0131 bool HLTMuonL3PreFilter::hltFilter(Event& iEvent,
0132                                    const EventSetup& iSetup,
0133                                    trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0134   // All HLT filters must create and fill an HLT filter object,
0135   // recording any reconstructed physics objects satisfying (or not)
0136   // this HLT filter, and place it in the Event.
0137 
0138   auto const prop = propSetup_.init(iSetup);
0139 
0140   if (saveTags())
0141     filterproduct.addCollectionTag(candTag_);
0142 
0143   // Read RecoChargedCandidates from L3MuonCandidateProducer:
0144   Handle<RecoChargedCandidateCollection> mucands;
0145   iEvent.getByToken(candToken_, mucands);
0146 
0147   // Read L2 triggered objects:
0148   Handle<TriggerFilterObjectWithRefs> previousLevelCands;
0149   iEvent.getByToken(previousCandToken_, previousLevelCands);
0150   vector<RecoChargedCandidateRef> vl2cands;
0151   previousLevelCands->getObjects(TriggerMuon, vl2cands);
0152 
0153   // Read BeamSpot information:
0154   Handle<BeamSpot> recoBeamSpotHandle;
0155   iEvent.getByToken(beamspotToken_, recoBeamSpotHandle);
0156   const BeamSpot& beamSpot = *recoBeamSpotHandle;
0157 
0158   // Number of objects passing the L3 Trigger:
0159   int n = 0;
0160 
0161   // sort them by L2Track
0162   std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;
0163   // map the L3 cands matched to a L1 to their position in the recoMuon collection
0164   std::map<unsigned int, RecoChargedCandidateRef> MuonToL3s;
0165 
0166   // Test to see if we can use L3MuonTrajectorySeeds:
0167   if (mucands->empty())
0168     return false;
0169   auto const& tk = (*mucands)[0].track();
0170   bool useL3MTS = false;
0171 
0172   if (tk->seedRef().isNonnull()) {
0173     auto a = dynamic_cast<const L3MuonTrajectorySeed*>(tk->seedRef().get());
0174     useL3MTS = a != nullptr;
0175   }
0176 
0177   // If we can use L3MuonTrajectory seeds run the older code:
0178   if (useL3MTS) {
0179     LogDebug("HLTMuonL3PreFilter") << "HLTMuonL3PreFilter::hltFilter is in mode: useL3MTS";
0180 
0181     unsigned int maxI = mucands->size();
0182     for (unsigned int i = 0; i != maxI; ++i) {
0183       const TrackRef& tk = (*mucands)[i].track();
0184       edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef =
0185           tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
0186       TrackRef staTrack = l3seedRef->l2Track();
0187       LogDebug("HLTMuonL3PreFilter") << "L2 from: " << iEvent.getStableProvenance(staTrack.id()).moduleLabel()
0188                                      << " index: " << staTrack.key();
0189       L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands, i));
0190     }
0191   }  //end of useL3MTS
0192 
0193   // Using normal TrajectorySeeds:
0194   else {
0195     LogDebug("HLTMuonL3PreFilter") << "HLTMuonL3PreFilter::hltFilter is in mode: not useL3MTS";
0196 
0197     // Read Links collection:
0198     edm::Handle<reco::MuonTrackLinksCollection> links;
0199     iEvent.getByToken(linkToken_, links);
0200 
0201     edm::Handle<trigger::TriggerFilterObjectWithRefs> level1Cands;
0202     std::vector<l1t::MuonRef> vl1cands;
0203 
0204     bool check_l1match = true;
0205 
0206     // Loop over RecoChargedCandidates:
0207     for (unsigned int i(0); i < mucands->size(); ++i) {
0208       RecoChargedCandidateRef cand(mucands, i);
0209       TrackRef tk = cand->track();
0210 
0211       if (!matchPreviousCand_) {
0212         MuonToL3s[i] = RecoChargedCandidateRef(cand);
0213       } else {
0214         check_l1match = true;
0215         int nlink = 0;
0216         for (auto const& link : *links) {
0217           nlink++;
0218 
0219           // Using the same method that was used to create the links between L3 and L2
0220           // ToDo: there should be a better way than dR,dPt matching
0221           const reco::Track& trackerTrack = *link.trackerTrack();
0222 
0223           float dR2 = deltaR2(tk->eta(), tk->phi(), trackerTrack.eta(), trackerTrack.phi());
0224           float dPt = std::abs(tk->pt() - trackerTrack.pt());
0225           if (tk->pt() != 0)
0226             dPt = dPt / tk->pt();
0227 
0228           if (dR2 < 0.02 * 0.02 and dPt < 0.001) {
0229             const TrackRef staTrack = link.standAloneTrack();
0230             L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
0231             check_l1match = false;
0232           }
0233         }  //MTL loop
0234 
0235         if (!l1CandTag_.label().empty() && check_l1match) {
0236           auto const propagated = prop.extrapolate(*tk);
0237           auto const etaForMatch = propagated.isValid() ? propagated.globalPosition().eta() : cand->eta();
0238           auto const phiForMatch = propagated.isValid() ? (double)propagated.globalPosition().phi() : cand->phi();
0239           iEvent.getByToken(l1CandToken_, level1Cands);
0240           level1Cands->getObjects(trigger::TriggerL1Mu, vl1cands);
0241           const unsigned int nL1Muons(vl1cands.size());
0242           for (unsigned int il1 = 0; il1 != nL1Muons; ++il1) {
0243             if (deltaR2(etaForMatch, phiForMatch, vl1cands[il1]->eta(), vl1cands[il1]->phi()) < L1MatchingdR2_) {
0244               MuonToL3s[i] = RecoChargedCandidateRef(cand);
0245             }
0246           }
0247         }
0248       }
0249     }  //RCC loop
0250   }    //end of using normal TrajectorySeeds
0251 
0252   // look at all mucands,  check cuts and add to filter object
0253   LogDebug("HLTMuonL3PreFilter") << "looking at: " << L2toL3s.size() << " L2->L3s from: " << mucands->size();
0254   for (const auto& L2toL3s_it : L2toL3s) {
0255     if (!triggeredByLevel2(L2toL3s_it.first, vl2cands))
0256       continue;
0257 
0258     //loop over the L3Tk reconstructed for this L2.
0259     unsigned int iTk = 0;
0260     unsigned int maxItk = L2toL3s_it.second.size();
0261     for (; iTk != maxItk; iTk++) {
0262       const RecoChargedCandidateRef& cand = L2toL3s_it.second[iTk];
0263       if (!applySelection(cand, beamSpot))
0264         continue;
0265 
0266       filterproduct.addObject(TriggerMuon, cand);
0267       n++;
0268       break;  // and go on with the next L2 association
0269     }
0270   }  ////loop over L2s from L3 grouping
0271 
0272   // now loop on L3 from L1
0273   edm::Handle<reco::MuonCollection> recomuons;
0274   iEvent.getByToken(recoMuToken_, recomuons);
0275 
0276   for (const auto& MuonToL3s_it : MuonToL3s) {
0277     const reco::Muon& muon(recomuons->at(MuonToL3s_it.first));
0278 
0279     // applys specific cuts for TkMu
0280     if ((muon.type() & allowedTypeMask_L3fromL1_) == 0)
0281       continue;
0282     if ((muon.type() & requiredTypeMask_L3fromL1_) != requiredTypeMask_L3fromL1_)
0283       continue;
0284     if (muon.numberOfMatchedStations() < min_MuonStations_L3fromL1_)
0285       continue;
0286     if (!muon.globalTrack().isNull()) {
0287       if (muon.globalTrack()->normalizedChi2() > maxNormalizedChi2_L3fromL1_)
0288         continue;
0289     }
0290     if (muon.isTrackerMuon() && !muon::isGoodMuon(muon, trkMuonId_))
0291       continue;
0292 
0293     const RecoChargedCandidateRef& cand = MuonToL3s_it.second;
0294     // apply common selection
0295     if (!applySelection(cand, beamSpot))
0296       continue;
0297     filterproduct.addObject(TriggerMuon, cand);
0298     n++;
0299 
0300     break;  // and go on with the next L3 from L1
0301   }
0302 
0303   vector<RecoChargedCandidateRef> vref;
0304   filterproduct.getObjects(TriggerMuon, vref);
0305   for (auto& i : vref) {
0306     RecoChargedCandidateRef candref = RecoChargedCandidateRef(i);
0307     TrackRef tk = candref->get<TrackRef>();
0308     LogDebug("HLTMuonL3PreFilter") << " Track passing filter: trackRef pt= " << tk->pt() << " (" << candref->pt()
0309                                    << ") "
0310                                    << ", eta: " << tk->eta() << " (" << candref->eta() << ") ";
0311   }
0312 
0313   // filter decision
0314   const bool accept(n >= min_N_);
0315 
0316   LogDebug("HLTMuonL3PreFilter") << " >>>>> Result of HLTMuonL3PreFilter is " << accept
0317                                  << ", number of muons passing thresholds= " << n;
0318 
0319   return accept;
0320 }
0321 
0322 bool HLTMuonL3PreFilter::triggeredByLevel2(const TrackRef& staTrack, vector<RecoChargedCandidateRef>& vcands) const {
0323   bool ok = false;
0324   for (auto& vcand : vcands) {
0325     if (vcand->get<TrackRef>() == staTrack) {
0326       ok = true;
0327       LogDebug("HLTMuonL3PreFilter") << "The L2 track triggered";
0328       break;
0329     }
0330   }
0331   return ok;
0332 }
0333 
0334 bool HLTMuonL3PreFilter::applySelection(const RecoChargedCandidateRef& cand, const BeamSpot& beamSpot) const {
0335   // eta cut
0336   if (std::abs(cand->eta()) > max_Eta_)
0337     return false;
0338 
0339   TrackRef tk = cand->track();
0340   LogDebug("HLTMuonL3PreFilter") << " Muon in loop, q*pt= " << tk->charge() * tk->pt() << " ("
0341                                  << cand->charge() * cand->pt() << ") "
0342                                  << ", eta= " << tk->eta() << " (" << cand->eta() << ") "
0343                                  << ", hits= " << tk->numberOfValidHits() << ", d0= " << tk->d0()
0344                                  << ", dz= " << tk->dz();
0345 
0346   // cut on number of hits
0347   if (tk->numberOfValidHits() < min_Nhits_)
0348     return false;
0349 
0350   //max dr cut
0351   auto dr =
0352       std::abs((-(cand->vx() - beamSpot.x0()) * cand->py() + (cand->vy() - beamSpot.y0()) * cand->px()) / cand->pt());
0353   if (dr > max_Dr_)
0354     return false;
0355 
0356   //min dr cut
0357   if (dr < min_Dr_)
0358     return false;
0359 
0360   //dz cut
0361   if (std::abs((cand->vz() - beamSpot.z0()) -
0362                ((cand->vx() - beamSpot.x0()) * cand->px() + (cand->vy() - beamSpot.y0()) * cand->py()) / cand->pt() *
0363                    cand->pz() / cand->pt()) > max_Dz_)
0364     return false;
0365 
0366   // dxy significance cut (safeguard against bizarre values)
0367   if (min_DxySig_ > 0 && (tk->dxyError() <= 0 || std::abs(tk->dxy(beamSpot.position()) / tk->dxyError()) < min_DxySig_))
0368     return false;
0369 
0370   //normalizedChi2 cut
0371   if (tk->normalizedChi2() > max_NormalizedChi2_)
0372     return false;
0373 
0374   //dxy beamspot cut
0375   float absDxy = std::abs(tk->dxy(beamSpot.position()));
0376   if (absDxy > max_DXYBeamSpot_ || absDxy < min_DXYBeamSpot_)
0377     return false;
0378 
0379   //min muon hits cut
0380   const reco::HitPattern& trackHits = tk->hitPattern();
0381   if (trackHits.numberOfValidMuonHits() < min_NmuonHits_)
0382     return false;
0383 
0384   //pt difference cut
0385   double candPt = cand->pt();
0386   double trackPt = tk->pt();
0387 
0388   if (std::abs(candPt - trackPt) > max_PtDifference_)
0389     return false;
0390 
0391   //track pt cut
0392   if (trackPt < min_TrackPt_)
0393     return false;
0394 
0395   // Pt threshold cut
0396   double pt = cand->pt();
0397   double err0 = tk->error(0);
0398   double abspar0 = std::abs(tk->parameter(0));
0399   double ptLx = pt;
0400   // convert 50% efficiency threshold to 90% efficiency threshold
0401   if (abspar0 > 0)
0402     ptLx += nsigma_Pt_ * err0 / abspar0 * pt;
0403   LogTrace("HLTMuonL3PreFilter") << " ...Muon in loop, trackkRef pt= " << tk->pt() << ", ptLx= " << ptLx << " cand pT "
0404                                  << cand->pt();
0405   if (ptLx < min_Pt_)
0406     return false;
0407 
0408   return true;
0409 }
0410 
0411 // declare this class as a framework plugin
0412 #include "FWCore/Framework/interface/MakerMacros.h"
0413 DEFINE_FWK_MODULE(HLTMuonL3PreFilter);