File indexing completed on 2024-04-06 12:18:35
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017 #include "HLTL1TMuonSelector.h"
0018
0019
0020 #include "FWCore/Framework/interface/EventSetup.h"
0021 #include "FWCore/Framework/interface/Event.h"
0022 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0023 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0024 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0025
0026 using namespace std;
0027 using namespace edm;
0028 using namespace l1t;
0029
0030
0031 HLTL1TMuonSelector::HLTL1TMuonSelector(const edm::ParameterSet& iConfig)
0032 : theSource(iConfig.getParameter<InputTag>("InputObjects")),
0033 theL1MinPt(iConfig.getParameter<double>("L1MinPt")),
0034 theL1MaxEta(iConfig.getParameter<double>("L1MaxEta")),
0035 theL1MinQuality(iConfig.getParameter<unsigned int>("L1MinQuality")),
0036 centralBxOnly_(iConfig.getParameter<bool>("CentralBxOnly")) {
0037 muCollToken_ = consumes<MuonBxCollection>(theSource);
0038
0039 produces<MuonBxCollection>();
0040 }
0041
0042
0043 HLTL1TMuonSelector::~HLTL1TMuonSelector() = default;
0044
0045 void HLTL1TMuonSelector::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0046 edm::ParameterSetDescription desc;
0047 desc.add<edm::InputTag>("InputObjects", edm::InputTag("hltGmtStage2Digis"));
0048 desc.add<double>("L1MinPt", -1.);
0049 desc.add<double>("L1MaxEta", 5.0);
0050 desc.add<unsigned int>("L1MinQuality", 0);
0051 desc.add<bool>("CentralBxOnly", true);
0052 descriptions.add("hltL1TMuonSelector", desc);
0053 }
0054
0055 void HLTL1TMuonSelector::produce(edm::StreamID, edm::Event& iEvent, const edm::EventSetup& iSetup) const {
0056 const std::string metname = "Muon|RecoMuon|HLTL1TMuonSelector";
0057
0058 unique_ptr<MuonBxCollection> output(new MuonBxCollection());
0059
0060
0061 edm::Handle<MuonBxCollection> muColl;
0062 iEvent.getByToken(muCollToken_, muColl);
0063 LogTrace(metname) << "Number of muons " << muColl->size() << endl;
0064
0065 for (int ibx = muColl->getFirstBX(); ibx <= muColl->getLastBX(); ++ibx) {
0066 if (centralBxOnly_ && (ibx != 0))
0067 continue;
0068 for (auto it = muColl->begin(ibx); it != muColl->end(ibx); it++) {
0069 unsigned int quality = it->hwQual();
0070 int valid_charge = it->hwChargeValid();
0071
0072 float pt = it->pt();
0073 float eta = it->eta();
0074 float theta = 2 * atan(exp(-eta));
0075 float phi = it->phi();
0076 int charge = it->charge();
0077
0078 if (!valid_charge)
0079 charge = 0;
0080
0081 if (pt < theL1MinPt || fabs(eta) > theL1MaxEta)
0082 continue;
0083
0084 LogTrace(metname) << "L1 Muon Found";
0085 LogTrace(metname) << "Pt = " << pt << " GeV/c";
0086 LogTrace(metname) << "eta = " << eta;
0087 LogTrace(metname) << "theta = " << theta << " rad";
0088 LogTrace(metname) << "phi = " << phi << " rad";
0089 LogTrace(metname) << "charge = " << charge;
0090
0091 if (quality <= theL1MinQuality)
0092 continue;
0093 LogTrace(metname) << "quality = " << quality;
0094
0095 output->push_back(ibx, *it);
0096 }
0097 }
0098
0099 iEvent.put(std::move(output));
0100 }
0101
0102
0103 #include "FWCore/Framework/interface/MakerMacros.h"
0104 DEFINE_FWK_MODULE(HLTL1TMuonSelector);