File indexing completed on 2025-07-24 02:01:44
0001 #include <memory>
0002
0003 #include "FWCore/Framework/interface/Frameworkfwd.h"
0004 #include "FWCore/Framework/interface/global/EDProducer.h"
0005 #include "FWCore/Framework/interface/Event.h"
0006 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0007 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0008 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0009
0010 #include "DataFormats/Math/interface/libminifloat.h"
0011 #include "DataFormats/ParticleFlowReco/interface/PFLayer.h"
0012 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
0013 #include "DataFormats/ParticleFlowReco/interface/PFRecHitFwd.h"
0014 #include "DataFormats/Scouting/interface/Run3ScoutingEBRecHit.h"
0015 #include "DataFormats/Scouting/interface/Run3ScoutingEERecHit.h"
0016 #include "DataFormats/Scouting/interface/Run3ScoutingHBHERecHit.h"
0017
0018 class HLTScoutingRecHitProducer : public edm::global::EDProducer<> {
0019 public:
0020 explicit HLTScoutingRecHitProducer(const edm::ParameterSet&);
0021 ~HLTScoutingRecHitProducer() override = default;
0022
0023 static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
0024
0025 private:
0026 void produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const final;
0027 static void produceEcal(edm::Event& iEvent,
0028 const reco::PFRecHitCollection& inputRecHits,
0029 double minEnergyEB,
0030 double minEnergyEE,
0031 int mantissaPrecision,
0032 const std::string& tag = "");
0033 static void produceHcal(edm::Event& iEvent,
0034 const reco::PFRecHitCollection& inputRecHits,
0035 double minEnergyHBHE,
0036 int mantissaPrecision,
0037 const std::string& tag = "");
0038
0039 template <typename T>
0040 void setToken(edm::EDGetTokenT<T>& token, const edm::ParameterSet& iConfig, std::string name) {
0041 const auto inputTag = iConfig.getParameter<edm::InputTag>(name);
0042 if (!inputTag.encode().empty()) {
0043 token = consumes(inputTag);
0044 }
0045 }
0046
0047 edm::EDGetTokenT<reco::PFRecHitCollection> recoPFRecHitsTokenECAL_;
0048 edm::EDGetTokenT<reco::PFRecHitCollection> recoPFRecHitsTokenECALCleaned_;
0049 edm::EDGetTokenT<reco::PFRecHitCollection> recoPFRecHitsTokenHBHE_;
0050 const double minEnergyEB_;
0051 const double minEnergyEE_;
0052 const double minEnergyCleanedEB_;
0053 const double minEnergyCleanedEE_;
0054 const double minEnergyHBHE_;
0055 const int mantissaPrecision_;
0056 };
0057
0058 HLTScoutingRecHitProducer::HLTScoutingRecHitProducer(const edm::ParameterSet& iConfig)
0059 : minEnergyEB_(iConfig.getParameter<double>("minEnergyEB")),
0060 minEnergyEE_(iConfig.getParameter<double>("minEnergyEE")),
0061 minEnergyCleanedEB_(iConfig.getParameter<double>("minEnergyCleanedEB")),
0062 minEnergyCleanedEE_(iConfig.getParameter<double>("minEnergyCleanedEE")),
0063 minEnergyHBHE_(iConfig.getParameter<double>("minEnergyHBHE")),
0064 mantissaPrecision_(iConfig.getParameter<int>("mantissaPrecision")) {
0065
0066
0067 setToken(recoPFRecHitsTokenECAL_, iConfig, "pfRecHitsECAL");
0068 setToken(recoPFRecHitsTokenECALCleaned_, iConfig, "pfRecHitsECALCleaned");
0069 setToken(recoPFRecHitsTokenHBHE_, iConfig, "pfRecHitsHBHE");
0070 produces<Run3ScoutingEBRecHitCollection>("EB");
0071 produces<Run3ScoutingEERecHitCollection>("EE");
0072 produces<Run3ScoutingEBRecHitCollection>("EBCleaned");
0073 produces<Run3ScoutingEERecHitCollection>("EECleaned");
0074 produces<Run3ScoutingHBHERecHitCollection>("HBHE");
0075 }
0076
0077 void HLTScoutingRecHitProducer::produceEcal(edm::Event& iEvent,
0078 const reco::PFRecHitCollection& inputRecHits,
0079 double minEnergyEB,
0080 double minEnergyEE,
0081 int mantissaPrecision,
0082 const std::string& tag) {
0083 auto run3ScoutEBRecHits = std::make_unique<Run3ScoutingEBRecHitCollection>();
0084 run3ScoutEBRecHits->reserve(inputRecHits.size());
0085
0086 auto run3ScoutEERecHits = std::make_unique<Run3ScoutingEERecHitCollection>();
0087 run3ScoutEERecHits->reserve(inputRecHits.size());
0088
0089 for (auto const& rh : inputRecHits) {
0090 if (rh.layer() == PFLayer::ECAL_BARREL) {
0091 if (rh.energy() < minEnergyEB) {
0092 continue;
0093 }
0094
0095 run3ScoutEBRecHits->emplace_back(
0096 MiniFloatConverter::reduceMantissaToNbitsRounding(rh.energy(), mantissaPrecision),
0097 MiniFloatConverter::reduceMantissaToNbitsRounding(rh.time(), mantissaPrecision),
0098 rh.detId(),
0099 rh.flags());
0100 } else if (rh.layer() == PFLayer::ECAL_ENDCAP) {
0101 if (rh.energy() < minEnergyEE) {
0102 continue;
0103 }
0104
0105 run3ScoutEERecHits->emplace_back(
0106 MiniFloatConverter::reduceMantissaToNbitsRounding(rh.energy(), mantissaPrecision),
0107 MiniFloatConverter::reduceMantissaToNbitsRounding(rh.time(), mantissaPrecision),
0108 rh.detId());
0109 } else {
0110 edm::LogWarning("HLTScoutingRecHitProducer")
0111 << "Skipping PFRecHit because of unexpected PFLayer value (" << rh.layer() << ").";
0112 }
0113 }
0114 iEvent.put(std::move(run3ScoutEBRecHits), "EB" + tag);
0115 iEvent.put(std::move(run3ScoutEERecHits), "EE" + tag);
0116 }
0117
0118 void HLTScoutingRecHitProducer::produceHcal(edm::Event& iEvent,
0119 const reco::PFRecHitCollection& inputRecHits,
0120 double minEnergyHBHE,
0121 int mantissaPrecision,
0122 const std::string& tag) {
0123 auto run3ScoutHBHERecHits = std::make_unique<Run3ScoutingHBHERecHitCollection>();
0124 run3ScoutHBHERecHits->reserve(inputRecHits.size());
0125
0126 for (auto const& rh : inputRecHits) {
0127 if (rh.energy() < minEnergyHBHE) {
0128 continue;
0129 }
0130
0131 run3ScoutHBHERecHits->emplace_back(
0132 MiniFloatConverter::reduceMantissaToNbitsRounding(rh.energy(), mantissaPrecision), rh.detId());
0133 }
0134
0135 iEvent.put(std::move(run3ScoutHBHERecHits), "HBHE" + tag);
0136 }
0137
0138 void HLTScoutingRecHitProducer::produce(edm::StreamID, edm::Event& iEvent, edm::EventSetup const&) const {
0139
0140 if (!recoPFRecHitsTokenECAL_.isUninitialized()) {
0141 auto const& recoPFRecHitsECAL = iEvent.get(recoPFRecHitsTokenECAL_);
0142 produceEcal(iEvent, recoPFRecHitsECAL, minEnergyEB_, minEnergyEE_, mantissaPrecision_);
0143 }
0144
0145 if (!recoPFRecHitsTokenECALCleaned_.isUninitialized()) {
0146 auto const& recoPFRecHitsECALCleaned = iEvent.get(recoPFRecHitsTokenECALCleaned_);
0147 produceEcal(
0148 iEvent, recoPFRecHitsECALCleaned, minEnergyCleanedEB_, minEnergyCleanedEE_, mantissaPrecision_, "Cleaned");
0149 }
0150
0151 if (!recoPFRecHitsTokenHBHE_.isUninitialized()) {
0152 auto const& recoPFRecHitsHBHE = iEvent.get(recoPFRecHitsTokenHBHE_);
0153 produceHcal(iEvent, recoPFRecHitsHBHE, minEnergyHBHE_, mantissaPrecision_);
0154 }
0155 }
0156
0157 void HLTScoutingRecHitProducer::fillDescriptions(edm::ConfigurationDescriptions& descriptions) {
0158 edm::ParameterSetDescription desc;
0159 desc.add<edm::InputTag>("pfRecHitsECAL", edm::InputTag("hltParticleFlowRecHitECALUnseeded"));
0160 desc.add<edm::InputTag>("pfRecHitsECALCleaned", edm::InputTag("hltParticleFlowRecHitECALUnseeded", "Cleaned"));
0161 desc.add<edm::InputTag>("pfRecHitsHBHE", edm::InputTag("hltParticleFlowRecHitHBHE"));
0162 desc.add<double>("minEnergyEB", -1)->setComment("Minimum energy of the EcalBarrel PFRecHit in GeV");
0163 desc.add<double>("minEnergyEE", -1)->setComment("Minimum energy of the EcalEndcap PFRecHit in GeV");
0164 desc.add<double>("minEnergyCleanedEB", -1)->setComment("Minimum energy of the cleaned EcalBarrel PFRecHit in GeV");
0165 desc.add<double>("minEnergyCleanedEE", -1)->setComment("Minimum energy of the cleaned EcalEndcap PFRecHit in GeV");
0166 desc.add<double>("minEnergyHBHE", -1)->setComment("Minimum energy of the HBHE PFRecHit in GeV");
0167 desc.add<int>("mantissaPrecision", 10)->setComment("default of 10 corresponds to float16, change to 23 for float32");
0168 descriptions.add("hltScoutingRecHitProducer", desc);
0169 }
0170
0171 #include "FWCore/Framework/interface/MakerMacros.h"
0172 DEFINE_FWK_MODULE(HLTScoutingRecHitProducer);