File indexing completed on 2024-09-07 04:36:38
0001
0002
0003
0004
0005
0006
0007
0008
0009 #include "DataFormats/Common/interface/Handle.h"
0010
0011 #include "DataFormats/HLTReco/interface/TriggerFilterObjectWithRefs.h"
0012 #include "DataFormats/HLTReco/interface/TriggerRefsCollections.h"
0013
0014 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0015
0016 #include "DataFormats/TrackReco/interface/Track.h"
0017 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidate.h"
0018 #include "DataFormats/RecoCandidate/interface/RecoChargedCandidateFwd.h"
0019 #include "DataFormats/MuonReco/interface/MuonFwd.h"
0020 #include "HLTMuonTrimuonL3Filter.h"
0021 #include "DataFormats/BeamSpot/interface/BeamSpot.h"
0022 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeed.h"
0023 #include "DataFormats/MuonSeed/interface/L3MuonTrajectorySeedCollection.h"
0024 #include "DataFormats/TrajectorySeed/interface/TrajectorySeedCollection.h"
0025
0026 #include "TrackingTools/PatternTools/interface/ClosestApproachInRPhi.h"
0027 #include "TrackingTools/TransientTrack/interface/TransientTrack.h"
0028 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0029 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0030 #include "FWCore/Utilities/interface/InputTag.h"
0031
0032 #include "DataFormats/Math/interface/deltaR.h"
0033
0034 using namespace edm;
0035 using namespace std;
0036 using namespace reco;
0037 using namespace trigger;
0038
0039
0040
0041
0042 HLTMuonTrimuonL3Filter::HLTMuonTrimuonL3Filter(const edm::ParameterSet& iConfig)
0043 : HLTFilter(iConfig),
0044 idealMagneticFieldRecordToken_(esConsumes()),
0045 beamspotTag_(iConfig.getParameter<edm::InputTag>("BeamSpotTag")),
0046 beamspotToken_(consumes<reco::BeamSpot>(beamspotTag_)),
0047 candTag_(iConfig.getParameter<edm::InputTag>("CandTag")),
0048 candToken_(consumes<reco::RecoChargedCandidateCollection>(candTag_)),
0049 previousCandTag_(iConfig.getParameter<InputTag>("PreviousCandTag")),
0050 previousCandToken_(consumes<trigger::TriggerFilterObjectWithRefs>(previousCandTag_)),
0051 fast_Accept_(iConfig.getParameter<bool>("FastAccept")),
0052 max_Eta_(iConfig.getParameter<double>("MaxEta")),
0053 min_Nhits_(iConfig.getParameter<int>("MinNhits")),
0054 max_Dr_(iConfig.getParameter<double>("MaxDr")),
0055 max_Dz_(iConfig.getParameter<double>("MaxDz")),
0056 chargeOpt_(iConfig.getParameter<int>("ChargeOpt")),
0057 min_PtTriplet_(iConfig.getParameter<double>("MinPtTriplet")),
0058 min_PtMax_(iConfig.getParameter<double>("MinPtMax")),
0059 min_PtMin_(iConfig.getParameter<double>("MinPtMin")),
0060 min_InvMass_(iConfig.getParameter<double>("MinInvMass")),
0061 max_InvMass_(iConfig.getParameter<double>("MaxInvMass")),
0062 min_Acop_(iConfig.getParameter<double>("MinAcop")),
0063 max_Acop_(iConfig.getParameter<double>("MaxAcop")),
0064 min_PtBalance_(iConfig.getParameter<double>("MinPtBalance")),
0065 max_PtBalance_(iConfig.getParameter<double>("MaxPtBalance")),
0066 nsigma_Pt_(iConfig.getParameter<double>("NSigmaPt")),
0067 max_DCAMuMu_(iConfig.getParameter<double>("MaxDCAMuMu")),
0068 max_YTriplet_(iConfig.getParameter<double>("MaxRapidityTriplet")),
0069 theL3LinksLabel(iConfig.getParameter<InputTag>("InputLinks")),
0070 linkToken_(consumes<reco::MuonTrackLinksCollection>(theL3LinksLabel)) {
0071 LogDebug("HLTMuonTrimuonL3Filter")
0072 << " CandTag/MinN/MaxEta/MinNhits/MaxDr/MaxDz/MinPt1/MinPt2/MinInvMass/MaxInvMass/MinAcop/MaxAcop/MinPtBalance/"
0073 "MaxPtBalance/NSigmaPt/MaxDzMuMu/MaxRapidityTriplet : "
0074 << candTag_.encode() << " " << fast_Accept_ << " " << max_Eta_ << " " << min_Nhits_ << " " << max_Dr_ << " "
0075 << max_Dz_ << " " << chargeOpt_ << " " << min_PtTriplet_ << " " << min_PtMax_ << " " << min_PtMin_ << " "
0076 << min_InvMass_ << " " << max_InvMass_ << " " << min_Acop_ << " " << max_Acop_ << " " << min_PtBalance_ << " "
0077 << max_PtBalance_ << " " << nsigma_Pt_ << " " << max_DCAMuMu_ << " " << max_YTriplet_;
0078 }
0079
0080 HLTMuonTrimuonL3Filter::~HLTMuonTrimuonL3Filter() = default;
0081
0082 void HLTMuonTrimuonL3Filter::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0083 edm::ParameterSetDescription desc;
0084 makeHLTFilterDescription(desc);
0085 desc.add<edm::InputTag>("BeamSpotTag", edm::InputTag("hltOfflineBeamSpot"));
0086 desc.add<edm::InputTag>("CandTag", edm::InputTag("hltL3MuonCandidates"));
0087 desc.add<edm::InputTag>("PreviousCandTag", edm::InputTag(""));
0088 desc.add<bool>("FastAccept", false);
0089 desc.add<double>("MaxEta", 2.5);
0090 desc.add<int>("MinNhits", 0);
0091 desc.add<double>("MaxDr", 2.0);
0092 desc.add<double>("MaxDz", 9999.0);
0093 desc.add<int>("ChargeOpt", 0);
0094 desc.add<double>("MinPtTriplet", 0.0);
0095 desc.add<double>("MinPtMax", 3.0);
0096 desc.add<double>("MinPtMin", 3.0);
0097 desc.add<double>("MinInvMass", 2.8);
0098 desc.add<double>("MaxInvMass", 3.4);
0099 desc.add<double>("MinAcop", -1.0);
0100 desc.add<double>("MaxAcop", 3.15);
0101 desc.add<double>("MinPtBalance", -1.0);
0102 desc.add<double>("MaxPtBalance", 999999.0);
0103 desc.add<double>("NSigmaPt", 0.0);
0104 desc.add<double>("MaxDCAMuMu", 99999.9);
0105 desc.add<double>("MaxRapidityTriplet", 999999.0);
0106 desc.add<edm::InputTag>("InputLinks", edm::InputTag(""));
0107 descriptions.add("hltMuonTrimuonL3Filter", desc);
0108 }
0109
0110
0111
0112
0113
0114
0115 bool HLTMuonTrimuonL3Filter::hltFilter(edm::Event& iEvent,
0116 const edm::EventSetup& iSetup,
0117 trigger::TriggerFilterObjectWithRefs& filterproduct) const {
0118 double const MuMass = 0.106;
0119 double const MuMass2 = MuMass * MuMass;
0120
0121
0122
0123
0124
0125 Handle<RecoChargedCandidateCollection> mucands;
0126 if (saveTags())
0127 filterproduct.addCollectionTag(candTag_);
0128 iEvent.getByToken(candToken_, mucands);
0129
0130
0131 if (mucands->empty())
0132 return false;
0133 auto const& tk = (*mucands)[0].track();
0134 bool useL3MTS = false;
0135
0136 if (tk->seedRef().isNonnull()) {
0137 auto a = dynamic_cast<const L3MuonTrajectorySeed*>(tk->seedRef().get());
0138 useL3MTS = a != nullptr;
0139 }
0140
0141
0142 std::map<reco::TrackRef, std::vector<RecoChargedCandidateRef> > L2toL3s;
0143
0144
0145 if (useL3MTS) {
0146 unsigned int maxI = mucands->size();
0147 for (unsigned int i = 0; i != maxI; i++) {
0148 const TrackRef& tk = (*mucands)[i].track();
0149 edm::Ref<L3MuonTrajectorySeedCollection> l3seedRef =
0150 tk->seedRef().castTo<edm::Ref<L3MuonTrajectorySeedCollection> >();
0151 TrackRef staTrack = l3seedRef->l2Track();
0152 L2toL3s[staTrack].push_back(RecoChargedCandidateRef(mucands, i));
0153 }
0154 }
0155
0156 else {
0157
0158 edm::Handle<reco::MuonTrackLinksCollection> links;
0159 iEvent.getByToken(linkToken_, links);
0160
0161
0162 for (unsigned int i(0); i < mucands->size(); ++i) {
0163 RecoChargedCandidateRef cand(mucands, i);
0164 for (auto const& link : *links) {
0165 TrackRef tk = cand->track();
0166
0167
0168
0169 const reco::Track& globalTrack = *link.globalTrack();
0170 float dR2 = deltaR2(tk->eta(), tk->phi(), globalTrack.eta(), globalTrack.phi());
0171 float dPt = std::abs(tk->pt() - globalTrack.pt()) / tk->pt();
0172 const TrackRef staTrack = link.standAloneTrack();
0173 if (dR2 < 0.02 * 0.02 and dPt < 0.001) {
0174 L2toL3s[staTrack].push_back(RecoChargedCandidateRef(cand));
0175 }
0176 }
0177 }
0178 }
0179
0180 Handle<TriggerFilterObjectWithRefs> previousLevelCands;
0181 iEvent.getByToken(previousCandToken_, previousLevelCands);
0182 BeamSpot beamSpot;
0183 Handle<BeamSpot> recoBeamSpotHandle;
0184 iEvent.getByToken(beamspotToken_, recoBeamSpotHandle);
0185 beamSpot = *recoBeamSpotHandle;
0186
0187
0188 auto const& bFieldHandle = iSetup.getHandle(idealMagneticFieldRecordToken_);
0189
0190
0191 vector<RecoChargedCandidateRef> vl2cands;
0192 previousLevelCands->getObjects(TriggerMuon, vl2cands);
0193
0194
0195 int n = 0;
0196 double e1, e2, e3;
0197 Particle::LorentzVector p, p1, p2, p3;
0198
0199 auto L2toL3s_it1 = L2toL3s.begin();
0200 auto L2toL3s_end = L2toL3s.end();
0201 bool atLeastOneTriplet = false;
0202 for (; L2toL3s_it1 != L2toL3s_end; ++L2toL3s_it1) {
0203 if (!triggeredByLevel2(L2toL3s_it1->first, vl2cands))
0204 continue;
0205
0206
0207 unsigned int iTk1 = 0;
0208 unsigned int maxItk1 = L2toL3s_it1->second.size();
0209 for (; iTk1 != maxItk1; iTk1++) {
0210 bool thisL3Index1isDone = false;
0211 RecoChargedCandidateRef& cand1 = L2toL3s_it1->second[iTk1];
0212 TrackRef tk1 = cand1->get<TrackRef>();
0213
0214 LogDebug("HLTMuonTrimuonL3Filter") << " 1st muon in loop: q*pt= " << tk1->charge() * tk1->pt() << " ("
0215 << cand1->charge() * cand1->pt() << ") "
0216 << ", eta= " << tk1->eta() << " (" << cand1->eta() << ") "
0217 << ", hits= " << tk1->numberOfValidHits();
0218
0219 if (fabs(cand1->eta()) > max_Eta_)
0220 continue;
0221
0222
0223 if (tk1->numberOfValidHits() < min_Nhits_)
0224 continue;
0225
0226
0227
0228 if (fabs((-(cand1->vx() - beamSpot.x0()) * cand1->py() + (cand1->vy() - beamSpot.y0()) * cand1->px()) /
0229 cand1->pt()) > max_Dr_)
0230 continue;
0231
0232
0233 if (fabs((cand1->vz() - beamSpot.z0()) -
0234 ((cand1->vx() - beamSpot.x0()) * cand1->px() + (cand1->vy() - beamSpot.y0()) * cand1->py()) /
0235 cand1->pt() * cand1->pz() / cand1->pt()) > max_Dz_)
0236 continue;
0237
0238
0239 double pt1 = cand1->pt();
0240
0241
0242 double ptLx1 = pt1;
0243
0244 LogDebug("HLTMuonTrimuonL3Filter") << " ... 1st muon in loop, pt1= " << pt1 << ", ptLx1= " << ptLx1;
0245 auto L2toL3s_it2 = L2toL3s_it1;
0246 L2toL3s_it2++;
0247 for (; L2toL3s_it2 != L2toL3s_end; ++L2toL3s_it2) {
0248 if (!triggeredByLevel2(L2toL3s_it2->first, vl2cands))
0249 continue;
0250
0251
0252 unsigned int iTk2 = 0;
0253 unsigned int maxItk2 = L2toL3s_it2->second.size();
0254 for (; iTk2 != maxItk2; iTk2++) {
0255 RecoChargedCandidateRef& cand2 = L2toL3s_it2->second[iTk2];
0256 TrackRef tk2 = cand2->get<TrackRef>();
0257
0258
0259 LogDebug("HLTMuonTrimuonL3Filter") << " 2nd muon in loop: q*pt= " << tk2->charge() * tk2->pt() << " ("
0260 << cand2->charge() * cand2->pt() << ") "
0261 << ", eta= " << tk2->eta() << " (" << cand2->eta() << ") "
0262 << ", hits= " << tk2->numberOfValidHits() << ", d0= " << tk2->d0();
0263 if (fabs(cand2->eta()) > max_Eta_)
0264 continue;
0265
0266
0267 if (tk2->numberOfValidHits() < min_Nhits_)
0268 continue;
0269
0270
0271
0272 if (fabs((-(cand2->vx() - beamSpot.x0()) * cand2->py() + (cand2->vy() - beamSpot.y0()) * cand2->px()) /
0273 cand2->pt()) > max_Dr_)
0274 continue;
0275
0276
0277 if (fabs((cand2->vz() - beamSpot.z0()) -
0278 ((cand2->vx() - beamSpot.x0()) * cand2->px() + (cand2->vy() - beamSpot.y0()) * cand2->py()) /
0279 cand2->pt() * cand2->pz() / cand2->pt()) > max_Dz_)
0280 continue;
0281
0282
0283 double pt2 = cand2->pt();
0284
0285
0286 double ptLx2 = pt2;
0287
0288 LogDebug("HLTMuonTrimuonL3Filter") << " ... 2nd muon in loop, pt2= " << pt2 << ", ptLx2= " << ptLx2;
0289
0290 auto L2toL3s_it3 = L2toL3s_it2;
0291 L2toL3s_it3++;
0292 for (; L2toL3s_it3 != L2toL3s_end; ++L2toL3s_it3) {
0293 if (!triggeredByLevel2(L2toL3s_it3->first, vl2cands))
0294 continue;
0295
0296
0297 unsigned int iTk3 = 0;
0298 unsigned int maxItk3 = L2toL3s_it3->second.size();
0299 for (; iTk3 != maxItk3; iTk3++) {
0300 RecoChargedCandidateRef& cand3 = L2toL3s_it3->second[iTk3];
0301 TrackRef tk3 = cand3->get<TrackRef>();
0302
0303 LogDebug("HLTMuonTrimuonL3Filter") << " 3rd muon in loop: q*pt= " << tk3->charge() * tk3->pt() << " ("
0304 << cand3->charge() * cand3->pt() << ") "
0305 << ", eta= " << tk3->eta() << " (" << cand3->eta() << ") "
0306 << ", hits= " << tk3->numberOfValidHits();
0307
0308 if (fabs(cand3->eta()) > max_Eta_)
0309 continue;
0310
0311
0312 if (tk3->numberOfValidHits() < min_Nhits_)
0313 continue;
0314
0315
0316
0317 if (fabs((-(cand3->vx() - beamSpot.x0()) * cand3->py() + (cand3->vy() - beamSpot.y0()) * cand3->px()) /
0318 cand3->pt()) > max_Dr_)
0319 continue;
0320
0321
0322 if (fabs((cand3->vz() - beamSpot.z0()) -
0323 ((cand3->vx() - beamSpot.x0()) * cand3->px() + (cand3->vy() - beamSpot.y0()) * cand3->py()) /
0324 cand3->pt() * cand3->pz() / cand3->pt()) > max_Dz_)
0325 continue;
0326
0327
0328 double pt3 = cand3->pt();
0329
0330
0331 double ptLx3 = pt3;
0332
0333 LogDebug("HLTMuonTrimuonL3Filter") << " ... 3rd muon in loop, pt3= " << pt3 << ", ptLx3= " << ptLx3;
0334
0335 if (ptLx1 > ptLx2 && ptLx1 > ptLx3 && ptLx1 < min_PtMax_)
0336 continue;
0337 else if (ptLx2 > ptLx1 && ptLx2 > ptLx3 && ptLx2 < min_PtMax_)
0338 continue;
0339 else if (ptLx3 < ptLx2 && ptLx3 > ptLx1 && ptLx3 < min_PtMax_)
0340 continue;
0341
0342 if (ptLx1 < ptLx2 && ptLx1 < ptLx3 && ptLx1 < min_PtMin_)
0343 continue;
0344 else if (ptLx2 < ptLx1 && ptLx2 < ptLx3 && ptLx2 < min_PtMin_)
0345 continue;
0346 else if (ptLx3 < ptLx2 && ptLx3 < ptLx1 && ptLx3 < min_PtMin_)
0347 continue;
0348
0349 if (chargeOpt_ > 0) {
0350 if (abs(cand1->charge() + cand2->charge() + cand3->charge()) != chargeOpt_)
0351 continue;
0352 }
0353
0354
0355 double acop = fabs(cand1->phi() - cand2->phi());
0356 if (acop > M_PI)
0357 acop = 2 * M_PI - acop;
0358 acop = M_PI - acop;
0359 LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-2 acop= " << acop;
0360 if (acop < min_Acop_)
0361 continue;
0362 if (acop > max_Acop_)
0363 continue;
0364
0365 acop = fabs(cand1->phi() - cand3->phi());
0366 if (acop > M_PI)
0367 acop = 2 * M_PI - acop;
0368 acop = M_PI - acop;
0369 LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-3 acop= " << acop;
0370 if (acop < min_Acop_)
0371 continue;
0372 if (acop > max_Acop_)
0373 continue;
0374
0375 acop = fabs(cand3->phi() - cand2->phi());
0376 if (acop > M_PI)
0377 acop = 2 * M_PI - acop;
0378 acop = M_PI - acop;
0379 LogDebug("HLTMuonTrimuonL3Filter") << " ... 3-2 acop= " << acop;
0380 if (acop < min_Acop_)
0381 continue;
0382 if (acop > max_Acop_)
0383 continue;
0384
0385
0386 double ptbalance = fabs(cand1->pt() - cand2->pt());
0387 if (ptbalance < min_PtBalance_)
0388 continue;
0389 if (ptbalance > max_PtBalance_)
0390 continue;
0391 ptbalance = fabs(cand1->pt() - cand3->pt());
0392 if (ptbalance < min_PtBalance_)
0393 continue;
0394 if (ptbalance > max_PtBalance_)
0395 continue;
0396 ptbalance = fabs(cand3->pt() - cand2->pt());
0397 if (ptbalance < min_PtBalance_)
0398 continue;
0399 if (ptbalance > max_PtBalance_)
0400 continue;
0401
0402
0403 e1 = sqrt(cand1->momentum().Mag2() + MuMass2);
0404 e2 = sqrt(cand2->momentum().Mag2() + MuMass2);
0405 e3 = sqrt(cand3->momentum().Mag2() + MuMass2);
0406 p1 = Particle::LorentzVector(cand1->px(), cand1->py(), cand1->pz(), e1);
0407 p2 = Particle::LorentzVector(cand2->px(), cand2->py(), cand2->pz(), e2);
0408 p3 = Particle::LorentzVector(cand3->px(), cand3->py(), cand3->pz(), e3);
0409 p = p1 + p2 + p3;
0410
0411 double pt123 = p.pt();
0412 LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-2 pt123= " << pt123;
0413 if (pt123 < min_PtTriplet_)
0414 continue;
0415
0416 double invmass = abs(p.mass());
0417
0418 LogDebug("HLTMuonTrimuonL3Filter") << " ... 1-2 invmass= " << invmass;
0419 if (invmass < min_InvMass_)
0420 continue;
0421 if (invmass > max_InvMass_)
0422 continue;
0423
0424
0425
0426
0427
0428
0429 TransientTrack mu1TT(*tk1, &(*bFieldHandle));
0430 TransientTrack mu2TT(*tk2, &(*bFieldHandle));
0431 TransientTrack mu3TT(*tk3, &(*bFieldHandle));
0432 TrajectoryStateClosestToPoint mu1TS = mu1TT.impactPointTSCP();
0433 TrajectoryStateClosestToPoint mu2TS = mu2TT.impactPointTSCP();
0434 TrajectoryStateClosestToPoint mu3TS = mu3TT.impactPointTSCP();
0435 if (mu1TS.isValid() && mu2TS.isValid() && mu3TS.isValid()) {
0436 ClosestApproachInRPhi cApp;
0437 cApp.calculate(mu1TS.theState(), mu2TS.theState());
0438 if (!cApp.status() || cApp.distance() > max_DCAMuMu_)
0439 continue;
0440 cApp.calculate(mu1TS.theState(), mu3TS.theState());
0441 if (!cApp.status() || cApp.distance() > max_DCAMuMu_)
0442 continue;
0443 cApp.calculate(mu3TS.theState(), mu2TS.theState());
0444 if (!cApp.status() || cApp.distance() > max_DCAMuMu_)
0445 continue;
0446 }
0447
0448
0449 double rapidity = fabs(p.Rapidity());
0450 if (rapidity > max_YTriplet_)
0451 continue;
0452
0453
0454 n++;
0455 LogDebug("HLTMuonTrimuonL3Filter")
0456 << " Track1 passing filter: pt= " << cand1->pt() << ", eta: " << cand1->eta();
0457 LogDebug("HLTMuonTrimuonL3Filter")
0458 << " Track2 passing filter: pt= " << cand2->pt() << ", eta: " << cand2->eta();
0459 LogDebug("HLTMuonTrimuonL3Filter")
0460 << " Track2 passing filter: pt= " << cand3->pt() << ", eta: " << cand3->eta();
0461 LogDebug("HLTMuonTrimuonL3Filter") << " Invmass= " << invmass;
0462
0463 bool i1done = false;
0464 bool i2done = false;
0465 bool i3done = false;
0466 vector<RecoChargedCandidateRef> vref;
0467 filterproduct.getObjects(TriggerMuon, vref);
0468 for (auto& i : vref) {
0469 RecoChargedCandidateRef candref = RecoChargedCandidateRef(i);
0470 TrackRef tktmp = candref->get<TrackRef>();
0471 if (tktmp == tk1) {
0472 i1done = true;
0473 } else if (tktmp == tk2) {
0474 i2done = true;
0475 } else if (tktmp == tk3) {
0476 i3done = true;
0477 }
0478 if (i1done && i2done && i3done)
0479 break;
0480 }
0481
0482 if (!i1done) {
0483 filterproduct.addObject(TriggerMuon, cand1);
0484 }
0485 if (!i2done) {
0486 filterproduct.addObject(TriggerMuon, cand2);
0487 }
0488 if (!i3done) {
0489 filterproduct.addObject(TriggerMuon, cand3);
0490 }
0491
0492
0493 thisL3Index1isDone = true;
0494 atLeastOneTriplet = true;
0495 break;
0496 }
0497
0498 if (atLeastOneTriplet && fast_Accept_)
0499 break;
0500 }
0501 }
0502
0503 if (atLeastOneTriplet && fast_Accept_)
0504 break;
0505 }
0506
0507
0508 if (atLeastOneTriplet && fast_Accept_)
0509 break;
0510 if (thisL3Index1isDone)
0511 break;
0512 }
0513
0514 if (atLeastOneTriplet && fast_Accept_)
0515 break;
0516 }
0517
0518
0519 const bool accept(n >= 1);
0520
0521 LogDebug("HLTMuonTrimuonL3Filter") << " >>>>> Result of HLTMuonTrimuonL3Filter is " << accept
0522 << ", number of muon triplets passing thresholds= " << n;
0523
0524 return accept;
0525 }
0526
0527 bool HLTMuonTrimuonL3Filter::triggeredByLevel2(const TrackRef& staTrack, vector<RecoChargedCandidateRef>& vcands) {
0528 bool ok = false;
0529 for (auto& vcand : vcands) {
0530 if (vcand->get<TrackRef>() == staTrack) {
0531 ok = true;
0532 LogDebug("HLTMuonL3PreFilter") << "The L2 track triggered";
0533 break;
0534 }
0535 }
0536 return ok;
0537 }
0538
0539
0540 #include "FWCore/Framework/interface/MakerMacros.h"
0541 DEFINE_FWK_MODULE(HLTMuonTrimuonL3Filter);