File indexing completed on 2024-07-03 04:18:13
0001
0002
0003
0004
0005 #include "FWCore/Framework/interface/Frameworkfwd.h"
0006 #include "FWCore/Framework/interface/stream/EDProducer.h"
0007
0008 #include "FWCore/Framework/interface/Event.h"
0009 #include "FWCore/ParameterSet/interface/ConfigurationDescriptions.h"
0010 #include "FWCore/ParameterSet/interface/ParameterSetDescription.h"
0011
0012 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
0013 #include "DataFormats/Common/interface/ValueMap.h"
0014
0015 class MergeClusterProducer : public edm::stream::EDProducer<> {
0016 public:
0017
0018
0019
0020
0021
0022
0023 MergeClusterProducer(const edm::ParameterSet &);
0024 ~MergeClusterProducer() override {}
0025
0026
0027
0028
0029
0030 static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
0031
0032
0033
0034
0035
0036
0037
0038 void produce(edm::Event &, const edm::EventSetup &) override;
0039
0040 private:
0041 edm::EDGetTokenT<std::vector<reco::CaloCluster>> EEclusters_token_;
0042 edm::EDGetTokenT<std::vector<reco::CaloCluster>> HSiclusters_token_;
0043 edm::EDGetTokenT<std::vector<reco::CaloCluster>> HSciclusters_token_;
0044
0045 std::string timeClname_;
0046 const edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> clustersTimeEE_token_;
0047 const edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> clustersTimeHSi_token_;
0048 const edm::EDGetTokenT<edm::ValueMap<std::pair<float, float>>> clustersTimeHSci_token_;
0049
0050
0051
0052
0053
0054
0055
0056
0057
0058 void mergeTogether(std::vector<reco::CaloCluster> &merge,
0059 const std::vector<reco::CaloCluster> &EE,
0060 const std::vector<reco::CaloCluster> &HSi,
0061 const std::vector<reco::CaloCluster> &HSci);
0062
0063
0064
0065
0066
0067
0068
0069 void addTo(std::vector<std::pair<float, float>> &to, const edm::ValueMap<std::pair<float, float>> &vm) {
0070 size_t size = vm.size();
0071 for (size_t i = 0; i < size; ++i) {
0072 to.push_back(vm.get(i));
0073 }
0074 }
0075
0076
0077
0078
0079
0080
0081
0082 void mergeTime(edm::Event &evt, size_t size, std::vector<std::pair<float, float>> ×) {
0083 edm::Handle<edm::ValueMap<std::pair<float, float>>> EE, HSi, HSci;
0084
0085 evt.getByToken(clustersTimeEE_token_, EE);
0086 evt.getByToken(clustersTimeHSi_token_, HSi);
0087 evt.getByToken(clustersTimeHSci_token_, HSci);
0088
0089 times.reserve(size);
0090 addTo(times, *EE);
0091 addTo(times, *HSi);
0092 addTo(times, *HSci);
0093 }
0094
0095
0096
0097
0098
0099
0100
0101
0102
0103
0104
0105 template <typename T>
0106 void createMerge(edm::Event &evt,
0107 const edm::EDGetTokenT<T> &EE_token,
0108 const edm::EDGetTokenT<T> &HSi_token,
0109 const edm::EDGetTokenT<T> &HSci_token,
0110 T &merge) {
0111 edm::Handle<T> EE, HSi, HSci;
0112
0113 evt.getByToken(EE_token, EE);
0114 evt.getByToken(HSi_token, HSi);
0115 evt.getByToken(HSci_token, HSci);
0116 mergeTogether(merge, *EE, *HSi, *HSci);
0117 }
0118 };
0119
0120 MergeClusterProducer::MergeClusterProducer(const edm::ParameterSet &ps)
0121 : timeClname_(ps.getParameter<std::string>("timeClname")),
0122 clustersTimeEE_token_(
0123 consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersEE"))),
0124 clustersTimeHSi_token_(
0125 consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersHSi"))),
0126 clustersTimeHSci_token_(
0127 consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersHSci"))) {
0128 EEclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersEE"));
0129 HSiclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersHSi"));
0130 HSciclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersHSci"));
0131
0132 produces<std::vector<float>>("InitialLayerClustersMask");
0133 produces<std::vector<reco::BasicCluster>>();
0134 produces<std::vector<reco::BasicCluster>>("sharing");
0135
0136 produces<edm::ValueMap<std::pair<float, float>>>(timeClname_);
0137 }
0138
0139 void MergeClusterProducer::fillDescriptions(edm::ConfigurationDescriptions &descriptions) {
0140
0141 edm::ParameterSetDescription desc;
0142
0143 desc.add<edm::InputTag>("layerClustersEE", edm::InputTag("hgcalLayerClustersEE"));
0144 desc.add<edm::InputTag>("layerClustersHSi", edm::InputTag("hgcalLayerClustersHSi"));
0145 desc.add<edm::InputTag>("layerClustersHSci", edm::InputTag("hgcalLayerClustersHSci"));
0146
0147
0148 desc.add<edm::InputTag>("time_layerclustersEE", edm::InputTag("hgcalLayerClustersEE", "timeLayerCluster"));
0149 desc.add<edm::InputTag>("time_layerclustersHSi", edm::InputTag("hgcalLayerClustersHSi", "timeLayerCluster"));
0150 desc.add<edm::InputTag>("time_layerclustersHSci", edm::InputTag("hgcalLayerClustersHSci", "timeLayerCluster"));
0151
0152 desc.add<std::string>("timeClname", "timeLayerCluster");
0153 descriptions.add("hgcalMergeLayerClusters", desc);
0154 }
0155
0156 void MergeClusterProducer::produce(edm::Event &evt, const edm::EventSetup &es) {
0157
0158 std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>);
0159 createMerge(evt, EEclusters_token_, HSiclusters_token_, HSciclusters_token_, *clusters);
0160
0161 auto clusterHandle = evt.put(std::move(clusters));
0162
0163
0164 std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
0165 layerClustersMask->resize(clusterHandle->size(), 1.0);
0166
0167 evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
0168
0169
0170 std::vector<std::pair<float, float>> times;
0171 mergeTime(evt, clusterHandle->size(), times);
0172
0173 auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
0174 edm::ValueMap<std::pair<float, float>>::Filler filler(*timeCl);
0175 filler.insert(clusterHandle, times.begin(), times.end());
0176 filler.fill();
0177 evt.put(std::move(timeCl), timeClname_);
0178 }
0179
0180 void MergeClusterProducer::mergeTogether(std::vector<reco::CaloCluster> &merge,
0181 const std::vector<reco::CaloCluster> &EE,
0182 const std::vector<reco::CaloCluster> &HSi,
0183 const std::vector<reco::CaloCluster> &HSci) {
0184 auto clusterSize = EE.size() + HSi.size() + HSci.size();
0185 merge.reserve(clusterSize);
0186
0187 merge.insert(merge.end(), EE.begin(), EE.end());
0188 merge.insert(merge.end(), HSi.begin(), HSi.end());
0189 merge.insert(merge.end(), HSci.begin(), HSci.end());
0190 }
0191
0192 #include "FWCore/Framework/interface/MakerMacros.h"
0193 DEFINE_FWK_MODULE(MergeClusterProducer);