File indexing completed on 2024-07-03 04:18:11
0001 #include <iostream>
0002 #include <string>
0003
0004 #include "DataFormats/HcalRecHit/interface/HcalRecHitHostCollection.h"
0005 #include "DataFormats/HcalRecHit/interface/HcalRecHitCollections.h"
0006 #include "FWCore/Framework/interface/Event.h"
0007 #include "FWCore/Framework/interface/EventSetup.h"
0008 #include "FWCore/Framework/interface/MakerMacros.h"
0009 #include "FWCore/Framework/interface/stream/EDProducer.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0011 #include "FWCore/Utilities/interface/EDGetToken.h"
0012 #include "FWCore/Utilities/interface/EDPutToken.h"
0013
0014 class HcalRecHitSoAToLegacy : public edm::stream::EDProducer<> {
0015 public:
0016 explicit HcalRecHitSoAToLegacy(edm::ParameterSet const& ps);
0017 ~HcalRecHitSoAToLegacy() override = default;
0018
0019 static void fillDescriptions(edm::ConfigurationDescriptions&);
0020
0021 private:
0022 void produce(edm::Event&, edm::EventSetup const&) override;
0023
0024 private:
0025 const edm::EDGetTokenT<hcal::RecHitHostCollection> recHitsTokenIn_;
0026 const edm::EDPutTokenT<HBHERecHitCollection> recHitsLegacyTokenOut_;
0027 };
0028
0029 void HcalRecHitSoAToLegacy::fillDescriptions(edm::ConfigurationDescriptions& confDesc) {
0030 edm::ParameterSetDescription desc;
0031
0032 desc.add<edm::InputTag>("src", edm::InputTag{"hbheRecHitProducerPortable"});
0033
0034 confDesc.addWithDefaultLabel(desc);
0035 }
0036
0037 HcalRecHitSoAToLegacy::HcalRecHitSoAToLegacy(const edm::ParameterSet& ps)
0038 : recHitsTokenIn_{consumes<hcal::RecHitHostCollection>(ps.getParameter<edm::InputTag>("src"))},
0039 recHitsLegacyTokenOut_{produces<HBHERecHitCollection>()} {}
0040
0041 void HcalRecHitSoAToLegacy::produce(edm::Event& event, edm::EventSetup const& setup) {
0042
0043 auto recHitsLegacy = std::make_unique<HBHERecHitCollection>();
0044
0045
0046 auto const& hcalRechitSoAView = event.get(recHitsTokenIn_).const_view();
0047
0048 recHitsLegacy->reserve(hcalRechitSoAView.metadata().size());
0049
0050 for (auto i = 0; i < hcalRechitSoAView.metadata().size(); i++) {
0051 auto const& rechit = hcalRechitSoAView[i];
0052
0053 if (rechit.chi2() < 0)
0054 continue;
0055
0056
0057 recHitsLegacy->emplace_back(HcalDetId{rechit.detId()},
0058 rechit.energy(),
0059 0
0060 );
0061
0062 recHitsLegacy->back().setChiSquared(rechit.chi2());
0063 recHitsLegacy->back().setRawEnergy(rechit.energyM0());
0064 }
0065
0066
0067 event.put(recHitsLegacyTokenOut_, std::move(recHitsLegacy));
0068 }
0069
0070 DEFINE_FWK_MODULE(HcalRecHitSoAToLegacy);