CMS 3D CMS Logo

L1TEGMultiMerger.cc
Go to the documentation of this file.
6 
12 #include <vector>
13 #include <iostream>
14 
16 public:
17  explicit L1TEGMultiMerger(const edm::ParameterSet&);
18  ~L1TEGMultiMerger() override;
19 
20 private:
21  void produce(edm::StreamID, edm::Event&, const edm::EventSetup&) const override;
22 
23  struct RefRemapper {
25  std::map<edm::Ref<BXVector<l1t::EGamma>>, edm::Ref<BXVector<l1t::EGamma>>> old2newRefMap;
26  };
27 
28  template <class T>
30  public:
32  : instanceLabel_(conf.getParameter<std::string>("instance")) {
33  for (const auto& producer_tag : conf.getParameter<std::vector<edm::InputTag>>("pfProducers")) {
34  tokens_.push_back(
35  prod->consumes<T>(edm::InputTag(producer_tag.label(), producer_tag.instance(), producer_tag.process())));
36  }
37  // FIXME: move this outside
38  prod->produces<T>(instanceLabel_);
39  }
40 
41  void produce(edm::Event& iEvent, RefRemapper& refRemapper) const {
43  auto out = std::make_unique<T>();
44  for (const auto& token : tokens_) {
45  iEvent.getByToken(token, handle);
46  populate(out, handle, refRemapper);
47  }
48  remapRefs(iEvent, out, refRemapper);
50  }
51 
52  private:
53  template <class TT>
54  void remapRefs(edm::Event& iEvent, std::unique_ptr<TT>& out, RefRemapper& refRemapper) const {
55  for (auto& egobj : *out) {
56  auto newref = refRemapper.old2newRefMap.find(egobj.EGRef());
57  if (newref != refRemapper.old2newRefMap.end()) {
58  egobj.setEGRef(newref->second);
59  }
60  }
61  }
62 
65  edm::Ref<BXVector<l1t::EGamma>>::key_type idx = 0;
66  for (std::size_t ix = 0; ix < out->size(); ix++) {
67  refRemapper.old2newRefMap[refRemapper.oldRefs[ix]] = edm::Ref<BXVector<l1t::EGamma>>(ref_egs, idx++);
68  }
69  }
70 
71  template <class TT>
72  void populate(std::unique_ptr<TT>& out, const edm::Handle<TT>& in, RefRemapper& refRemapper) const {
73  out->insert(out->end(), in->begin(), in->end());
74  }
75 
78  RefRemapper& refRemapper) const {
79  edm::Ref<BXVector<l1t::EGamma>>::key_type idx = 0;
80  for (int bx = in->getFirstBX(); bx <= in->getLastBX(); bx++) {
81  for (auto egee_itr = in->begin(bx); egee_itr != in->end(bx); egee_itr++) {
82  out->push_back(bx, *egee_itr);
83  // this to ensure that the old ref and the new object have the same index in the BXVector collection so that we can still match them
84  // no matter which BX we will insert next
86  }
87  }
88  }
89 
90  std::vector<edm::EDGetTokenT<T>> tokens_;
92  };
93 
94  std::vector<InstanceMerger<l1t::TkElectronCollection>> tkEleMerger;
95  std::vector<InstanceMerger<l1t::TkEmCollection>> tkEmMerger;
96  std::vector<InstanceMerger<BXVector<l1t::EGamma>>> tkEGMerger;
97 };
98 
100  for (const auto& config : conf.getParameter<std::vector<edm::ParameterSet>>("tkEgs")) {
102  }
103  for (const auto& config : conf.getParameter<std::vector<edm::ParameterSet>>("tkElectrons")) {
105  }
106  for (const auto& config : conf.getParameter<std::vector<edm::ParameterSet>>("tkEms")) {
108  }
109 }
110 
112 
114  RefRemapper refmapper;
115  for (const auto& egMerger : tkEGMerger)
116  egMerger.produce(iEvent, refmapper);
117  for (const auto& eleMerger : tkEleMerger)
118  eleMerger.produce(iEvent, refmapper);
119  for (const auto& emMerger : tkEmMerger)
120  emMerger.produce(iEvent, refmapper);
121 }
122 
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< InstanceMerger< BXVector< l1t::EGamma > > > tkEGMerger
void produce(edm::Event &iEvent, RefRemapper &refRemapper) const
InstanceMerger(L1TEGMultiMerger *prod, const edm::ParameterSet &conf)
std::vector< InstanceMerger< l1t::TkEmCollection > > tkEmMerger
Definition: config.py:1
std::map< edm::Ref< BXVector< l1t::EGamma > >, edm::Ref< BXVector< l1t::EGamma > > > old2newRefMap
std::unique_ptr< T, impl::DeviceDeleter > unique_ptr
void populate(std::unique_ptr< TT > &out, const edm::Handle< TT > &in, RefRemapper &refRemapper) const
int iEvent
Definition: GenABIO.cc:224
std::vector< InstanceMerger< l1t::TkElectronCollection > > tkEleMerger
L1TEGMultiMerger(const edm::ParameterSet &)
void produce(edm::StreamID, edm::Event &, const edm::EventSetup &) const override
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< edm::EDGetTokenT< T > > tokens_
void remapRefs(edm::Event &iEvent, std::unique_ptr< TT > &out, RefRemapper &refRemapper) const
void populate(std::unique_ptr< BXVector< l1t::EGamma >> &out, const edm::Handle< BXVector< l1t::EGamma >> &in, RefRemapper &refRemapper) const
long double T
BXVector< edm::Ref< BXVector< l1t::EGamma > > > oldRefs
~L1TEGMultiMerger() override
def move(src, dest)
Definition: eostools.py:511
void push_back(int bx, T object)
void remapRefs(edm::Event &iEvent, std::unique_ptr< BXVector< l1t::EGamma >> &out, RefRemapper &refRemapper) const