File indexing completed on 2024-05-20 22:40:11
0001 #include "FWCore/Framework/interface/ConsumesCollector.h"
0002 #include "FWCore/Framework/interface/Event.h"
0003 #include "FWCore/Framework/interface/ProducesCollector.h"
0004 #include "FWCore/MessageLogger/interface/MessageLogger.h"
0005 #include "FWCore/ParameterSet/interface/ParameterSet.h"
0006 #include "SimGeneral/MixingModule/interface/PileUpEventPrincipal.h"
0007
0008 #include "DataFormats/Common/interface/Handle.h"
0009 #include "SimDataFormats/CaloAnalysis/interface/MtdSimLayerCluster.h"
0010 #include "SimDataFormats/CaloAnalysis/interface/MtdSimLayerClusterFwd.h"
0011
0012 #include "SimGeneral/PreMixingModule/interface/PreMixingWorker.h"
0013 #include "SimGeneral/PreMixingModule/interface/PreMixingWorkerFactory.h"
0014
0015 class PreMixingMtdTruthWorker : public PreMixingWorker {
0016 public:
0017 PreMixingMtdTruthWorker(const edm::ParameterSet &ps, edm::ProducesCollector, edm::ConsumesCollector &&iC);
0018 ~PreMixingMtdTruthWorker() override = default;
0019
0020 void initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
0021 void addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) override;
0022 void addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) override;
0023 void put(edm::Event &iEvent,
0024 edm::EventSetup const &iSetup,
0025 std::vector<PileupSummaryInfo> const &ps,
0026 int bunchSpacing) override;
0027
0028 private:
0029 void add(const MtdSimLayerClusterCollection &clusters);
0030
0031 edm::EDGetTokenT<MtdSimLayerClusterCollection> sigClusterToken_;
0032
0033 edm::InputTag clusterPileInputTag_;
0034 std::string mtdSimLCCollectionDM_;
0035
0036 std::unique_ptr<MtdSimLayerClusterCollection> newClusters_;
0037 };
0038
0039 PreMixingMtdTruthWorker::PreMixingMtdTruthWorker(const edm::ParameterSet &ps,
0040 edm::ProducesCollector producesCollector,
0041 edm::ConsumesCollector &&iC)
0042 : sigClusterToken_(iC.consumes<MtdSimLayerClusterCollection>(ps.getParameter<edm::InputTag>("labelSig"))),
0043 clusterPileInputTag_(ps.getParameter<edm::InputTag>("pileInputTag")),
0044 mtdSimLCCollectionDM_(ps.getParameter<std::string>("collectionDM")) {
0045 producesCollector.produces<MtdSimLayerClusterCollection>(mtdSimLCCollectionDM_);
0046 }
0047
0048 void PreMixingMtdTruthWorker::initializeEvent(edm::Event const &iEvent, edm::EventSetup const &iSetup) {
0049 newClusters_ = std::make_unique<MtdSimLayerClusterCollection>();
0050 }
0051
0052 void PreMixingMtdTruthWorker::addSignals(edm::Event const &iEvent, edm::EventSetup const &iSetup) {
0053 edm::Handle<MtdSimLayerClusterCollection> clusters;
0054 iEvent.getByToken(sigClusterToken_, clusters);
0055
0056 if (clusters.isValid()) {
0057 add(*clusters);
0058 }
0059 }
0060
0061 void PreMixingMtdTruthWorker::addPileups(PileUpEventPrincipal const &pep, edm::EventSetup const &iSetup) {
0062 edm::Handle<MtdSimLayerClusterCollection> clusters;
0063 pep.getByLabel(clusterPileInputTag_, clusters);
0064
0065 if (clusters.isValid()) {
0066 add(*clusters);
0067 }
0068 }
0069
0070 void PreMixingMtdTruthWorker::add(const MtdSimLayerClusterCollection &clusters) {
0071
0072 newClusters_->reserve(newClusters_->size() + clusters.size());
0073 std::copy(clusters.begin(), clusters.end(), std::back_inserter(*newClusters_));
0074 }
0075
0076 void PreMixingMtdTruthWorker::put(edm::Event &iEvent,
0077 edm::EventSetup const &iSetup,
0078 std::vector<PileupSummaryInfo> const &ps,
0079 int bunchSpacing) {
0080 iEvent.put(std::move(newClusters_), mtdSimLCCollectionDM_);
0081 }
0082
0083 DEFINE_PREMIXING_WORKER(PreMixingMtdTruthWorker);