CMS 3D CMS Logo

ReducedRecHitCollectionProducer.cc
Go to the documentation of this file.
2 
6 
9 
12 
16 
18 
19 #include <iostream>
20 
22  recHitsToken_ = consumes<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>("recHitsLabel"));
23 
25  edm::vector_transform(iConfig.getParameter<std::vector<edm::InputTag>>("interestingDetIdCollections"),
26  [this](edm::InputTag const& tag) { return consumes<DetIdCollection>(tag); });
27 
28  reducedHitsCollection_ = iConfig.getParameter<std::string>("reducedHitsCollection");
29 
30  //register your products
31  produces<EcalRecHitCollection>(reducedHitsCollection_);
32 }
33 
35 
36 // ------------ method called to produce the data ------------
38  using namespace edm;
39  using namespace std;
40 
41  if (interestingDetIdCollections_.empty()) {
42  edm::LogError("ReducedRecHitCollectionProducer") << "VInputTag collections empty";
43  return;
44  }
45 
47  iEvent.getByToken(interestingDetIdCollections_[0], detIds);
48  std::vector<DetId> xtalsToStore((*detIds).size());
49  std::copy((*detIds).begin(), (*detIds).end(), xtalsToStore.begin());
50 
51  //Merging DetIds from different collections
52  for (unsigned int t = 1; t < interestingDetIdCollections_.size(); ++t) {
55  if (!detId.isValid()) {
56  Labels labels;
57  labelsForToken(interestingDetIdCollections_[t], labels);
58  edm::LogError("MissingInput") << "no reason to skip detid from : (" << labels.module << ", "
59  << labels.productInstance << ", " << labels.process << ")" << std::endl;
60  continue;
61  }
62 
63  for (unsigned int ii = 0; ii < (*detId).size(); ii++) {
64  if (std::find(xtalsToStore.begin(), xtalsToStore.end(), (*detId)[ii]) == xtalsToStore.end())
65  xtalsToStore.push_back((*detId)[ii]);
66  }
67  }
68 
69  Handle<EcalRecHitCollection> recHitsHandle;
70  iEvent.getByToken(recHitsToken_, recHitsHandle);
71  if (!recHitsHandle.isValid()) {
72  edm::LogError("ReducedRecHitCollectionProducer") << "RecHit collection not found";
73  return;
74  }
75 
76  //Create empty output collections
77  auto miniRecHitCollection = std::make_unique<EcalRecHitCollection>();
78 
79  for (unsigned int iCry = 0; iCry < xtalsToStore.size(); iCry++) {
80  EcalRecHitCollection::const_iterator iRecHit = recHitsHandle->find(xtalsToStore[iCry]);
81  if ((iRecHit != recHitsHandle->end()) &&
82  (miniRecHitCollection->find(xtalsToStore[iCry]) == miniRecHitCollection->end()))
83  miniRecHitCollection->push_back(*iRecHit);
84  }
85 
86  std::sort(xtalsToStore.begin(), xtalsToStore.end());
87  std::unique(xtalsToStore.begin(), xtalsToStore.end());
88 
89  // std::cout << "New Collection " << reducedHitsCollection_ << " size is " << miniRecHitCollection->size() << " original is " << recHitsHandle->size() << std::endl;
90  iEvent.put(std::move(miniRecHitCollection), reducedHitsCollection_);
91 }
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
std::vector< EcalRecHit >::const_iterator const_iterator
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
bool isValid() const
Definition: HandleBase.h:70
ReducedRecHitCollectionProducer(const edm::ParameterSet &)
ctor
ii
Definition: cuy.py:590
const_iterator end() const
std::vector< edm::EDGetTokenT< DetIdCollection > > interestingDetIdCollections_
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
iterator find(key_type k)
HLT enums.
def move(src, dest)
Definition: eostools.py:511
void produce(edm::Event &, const edm::EventSetup &) override
producer