CMS 3D CMS Logo

ReducedRecHitCollectionProducer.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: ReducedRecHitCollectionProducer
4 // Class: ReducedRecHitCollectionProducer
5 //
14 #include <iostream>
15 #include <memory>
16 
35 
37 public:
42  void produce(edm::Event&, const edm::EventSetup&) override;
43 
44 private:
45  // ----------member data ---------------------------
47  std::vector<edm::EDGetTokenT<DetIdCollection>> interestingDetIdCollections_;
49 };
50 
53 
55  recHitsToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsLabel"));
56 
58  edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag>>("interestingDetIdCollections"),
59  [this](edm::InputTag const& tag) { return consumes<DetIdCollection>(tag); });
60 
61  reducedHitsCollection_ = iConfig.getParameter<std::string>("reducedHitsCollection");
62 
63  //register your products
64  produces<EcalRecHitCollection>(reducedHitsCollection_);
65 }
66 
68 
69 // ------------ method called to produce the data ------------
71  using namespace edm;
72  using namespace std;
73 
74  if (interestingDetIdCollections_.empty()) {
75  edm::LogError("ReducedRecHitCollectionProducer") << "VInputTag collections empty";
76  return;
77  }
78 
80  iEvent.getByToken(interestingDetIdCollections_[0], detIds);
81  std::vector<DetId> xtalsToStore((*detIds).size());
82  std::copy((*detIds).begin(), (*detIds).end(), xtalsToStore.begin());
83 
84  //Merging DetIds from different collections
85  for (unsigned int t = 1; t < interestingDetIdCollections_.size(); ++t) {
87  iEvent.getByToken(interestingDetIdCollections_[t], detId);
88  if (!detId.isValid()) {
89  Labels labels;
90  labelsForToken(interestingDetIdCollections_[t], labels);
91  edm::LogError("MissingInput") << "no reason to skip detid from : (" << labels.module << ", "
92  << labels.productInstance << ", " << labels.process << ")" << std::endl;
93  continue;
94  }
95 
96  for (unsigned int ii = 0; ii < (*detId).size(); ii++) {
97  if (std::find(xtalsToStore.begin(), xtalsToStore.end(), (*detId)[ii]) == xtalsToStore.end())
98  xtalsToStore.push_back((*detId)[ii]);
99  }
100  }
101 
102  Handle<EcalRecHitCollection> recHitsHandle;
103  iEvent.getByToken(recHitsToken_, recHitsHandle);
104  if (!recHitsHandle.isValid()) {
105  edm::LogError("ReducedRecHitCollectionProducer") << "RecHit collection not found";
106  return;
107  }
108 
109  //Create empty output collections
110  auto miniRecHitCollection = std::make_unique<EcalRecHitCollection>();
111 
112  for (unsigned int iCry = 0; iCry < xtalsToStore.size(); iCry++) {
113  EcalRecHitCollection::const_iterator iRecHit = recHitsHandle->find(xtalsToStore[iCry]);
114  if ((iRecHit != recHitsHandle->end()) &&
115  (miniRecHitCollection->find(xtalsToStore[iCry]) == miniRecHitCollection->end()))
116  miniRecHitCollection->push_back(*iRecHit);
117  }
118 
119  std::sort(xtalsToStore.begin(), xtalsToStore.end());
120  std::unique(xtalsToStore.begin(), xtalsToStore.end());
121 
122  // std::cout << "New Collection " << reducedHitsCollection_ << " size is " << miniRecHitCollection->size() << " original is " << recHitsHandle->size() << std::endl;
123  iEvent.put(std::move(miniRecHitCollection), reducedHitsCollection_);
124 }
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
std::vector< EcalRecHit >::const_iterator const_iterator
Log< level::Error, false > LogError
auto vector_transform(std::vector< InputType > const &input, Function predicate) -> std::vector< typename std::remove_cv< typename std::remove_reference< decltype(predicate(input.front()))>::type >::type >
Definition: transform.h:11
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:19
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< EcalRecHitCollection > recHitsToken_
def unique(seq, keepstr=True)
Definition: tier0.py:24
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
ReducedRecHitCollectionProducer(const edm::ParameterSet &)
ctor
ii
Definition: cuy.py:589
const_iterator end() const
std::vector< edm::EDGetTokenT< DetIdCollection > > interestingDetIdCollections_
bool isValid() const
Definition: HandleBase.h:70
iterator find(key_type k)
HLT enums.
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &, const edm::EventSetup &) override
producer