CMS 3D CMS Logo

MergeClusterProducer.cc
Go to the documentation of this file.
1 // Authors: Olivie Franklova - olivie.abigail.franklova@cern.ch
2 // Date: 03/2023
3 // @file merge layer clusters which were produce by HGCalLayerClusterProducer
4 
7 
11 
14 
16 public:
24  ~MergeClusterProducer() override {}
30  static void fillDescriptions(edm::ConfigurationDescriptions &descriptions);
31 
38  void produce(edm::Event &, const edm::EventSetup &) override;
39 
40 private:
44 
49 
58  void mergeTogether(std::vector<reco::CaloCluster> &merge,
59  const std::vector<reco::CaloCluster> &EE,
60  const std::vector<reco::CaloCluster> &HSi,
61  const std::vector<reco::CaloCluster> &HSci);
62 
69  void addTo(std::vector<std::pair<float, float>> &to, const edm::ValueMap<std::pair<float, float>> &vm) {
70  size_t size = vm.size();
71  for (size_t i = 0; i < size; ++i) {
72  to.push_back(vm.get(i));
73  }
74  }
82  void mergeTime(edm::Event &evt, size_t size, std::vector<std::pair<float, float>> &times) {
84  // get values from all three part of detectors
88 
89  times.reserve(size);
90  addTo(times, *EE);
91  addTo(times, *HSi);
92  addTo(times, *HSci);
93  }
105  template <typename T>
107  const edm::EDGetTokenT<T> &EE_token,
108  const edm::EDGetTokenT<T> &HSi_token,
109  const edm::EDGetTokenT<T> &HSci_token,
110  T &merge) {
111  edm::Handle<T> EE, HSi, HSci;
112  // get values from all three part of detectors
113  evt.getByToken(EE_token, EE);
114  evt.getByToken(HSi_token, HSi);
115  evt.getByToken(HSci_token, HSci);
116  mergeTogether(merge, *EE, *HSi, *HSci);
117  }
118 };
119 
121  : timeClname_(ps.getParameter<std::string>("timeClname")),
122  clustersTimeEE_token_(
123  consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersEE"))),
124  clustersTimeHSi_token_(
125  consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersHSi"))),
126  clustersTimeHSci_token_(
127  consumes<edm::ValueMap<std::pair<float, float>>>(ps.getParameter<edm::InputTag>("time_layerclustersHSci"))) {
128  EEclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersEE"));
129  HSiclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersHSi"));
130  HSciclusters_token_ = consumes<std::vector<reco::CaloCluster>>(ps.getParameter<edm::InputTag>("layerClustersHSci"));
131 
132  produces<std::vector<float>>("InitialLayerClustersMask");
133  produces<std::vector<reco::BasicCluster>>();
134  produces<std::vector<reco::BasicCluster>>("sharing");
135  //time for layer clusters
136  produces<edm::ValueMap<std::pair<float, float>>>(timeClname_);
137 }
138 
140  // hgcalMergeLayerClusters
142  //layer clusters
143  desc.add<edm::InputTag>("layerClustersEE", edm::InputTag("hgcalLayerClustersEE"));
144  desc.add<edm::InputTag>("layerClustersHSi", edm::InputTag("hgcalLayerClustersHSi"));
145  desc.add<edm::InputTag>("layerClustersHSci", edm::InputTag("hgcalLayerClustersHSci"));
146 
147  //time
148  desc.add<edm::InputTag>("time_layerclustersEE", edm::InputTag("hgcalLayerClustersEE", "timeLayerCluster"));
149  desc.add<edm::InputTag>("time_layerclustersHSi", edm::InputTag("hgcalLayerClustersHSi", "timeLayerCluster"));
150  desc.add<edm::InputTag>("time_layerclustersHSci", edm::InputTag("hgcalLayerClustersHSci", "timeLayerCluster"));
151 
152  desc.add<std::string>("timeClname", "timeLayerCluster");
153  descriptions.add("hgcalMergeLayerClusters", desc);
154 }
155 
157  //merge clusters
158  std::unique_ptr<std::vector<reco::BasicCluster>> clusters(new std::vector<reco::BasicCluster>);
160  //put new clusters to event
161  auto clusterHandle = evt.put(std::move(clusters));
162 
163  //create layer cluster mask
164  std::unique_ptr<std::vector<float>> layerClustersMask(new std::vector<float>);
165  layerClustersMask->resize(clusterHandle->size(), 1.0);
166  //put it into event
167  evt.put(std::move(layerClustersMask), "InitialLayerClustersMask");
168 
169  //time
170  std::vector<std::pair<float, float>> times;
171  mergeTime(evt, clusterHandle->size(), times);
172 
173  auto timeCl = std::make_unique<edm::ValueMap<std::pair<float, float>>>();
175  filler.insert(clusterHandle, times.begin(), times.end());
176  filler.fill();
177  evt.put(std::move(timeCl), timeClname_);
178 }
179 
180 void MergeClusterProducer::mergeTogether(std::vector<reco::CaloCluster> &merge,
181  const std::vector<reco::CaloCluster> &EE,
182  const std::vector<reco::CaloCluster> &HSi,
183  const std::vector<reco::CaloCluster> &HSci) {
184  auto clusterSize = EE.size() + HSi.size() + HSci.size();
185  merge.reserve(clusterSize);
186 
187  merge.insert(merge.end(), EE.begin(), EE.end());
188  merge.insert(merge.end(), HSi.begin(), HSi.end());
189  merge.insert(merge.end(), HSci.begin(), HSci.end());
190 }
191 
size
Write out results.
const edm::EDGetTokenT< edm::ValueMap< std::pair< float, float > > > clustersTimeHSi_token_
Definition: merge.py:1
edm::EDGetTokenT< std::vector< reco::CaloCluster > > HSiclusters_token_
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
MergeClusterProducer(const edm::ParameterSet &)
Constructor with parameter settings - which can be changed in ...todo. Constructor will set all varia...
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Method fill description which will be used in pyhton file.
void addTo(std::vector< std::pair< float, float >> &to, const edm::ValueMap< std::pair< float, float >> &vm)
copy all values from vm to to
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:528
const edm::EDGetTokenT< edm::ValueMap< std::pair< float, float > > > clustersTimeHSci_token_
edm::EDGetTokenT< std::vector< reco::CaloCluster > > EEclusters_token_
void produce(edm::Event &, const edm::EventSetup &) override
Method will merge the producers and put them back to event.
void mergeTogether(std::vector< reco::CaloCluster > &merge, const std::vector< reco::CaloCluster > &EE, const std::vector< reco::CaloCluster > &HSi, const std::vector< reco::CaloCluster > &HSci)
method merge three vectors of reco::CaloCluster to one
void mergeTime(edm::Event &evt, size_t size, std::vector< std::pair< float, float >> &times)
Merge value map of time for all parts of detector together to vector times.
size_t size() const
Definition: Event.cc:252
edm::EDGetTokenT< std::vector< reco::CaloCluster > > HSciclusters_token_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void add(std::string const &label, ParameterSetDescription const &psetDescription)
HLT enums.
const edm::EDGetTokenT< edm::ValueMap< std::pair< float, float > > > clustersTimeEE_token_
void createMerge(edm::Event &evt, const edm::EDGetTokenT< T > &EE_token, const edm::EDGetTokenT< T > &HSi_token, const edm::EDGetTokenT< T > &HSci_token, T &merge)
get info form event and then call merge
long double T
def move(src, dest)
Definition: eostools.py:511