CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ReducedRecHitCollectionProducer.cc
Go to the documentation of this file.
2 
6 
7 
8 
11 
12 
13 
16 
20 
22 
23 #include <iostream>
24 
26 {
27 
28  recHitsToken_ =
29  consumes<EcalRecHitCollection>(iConfig.getParameter< edm::InputTag > ("recHitsLabel"));
30 
33  iConfig.getParameter<std::vector<edm::InputTag>>("interestingDetIdCollections"),
34  [this](edm::InputTag const & tag) {
35  return consumes<DetIdCollection>(tag);
36  }
37  );
38 
39  reducedHitsCollection_ = iConfig.getParameter<std::string>("reducedHitsCollection");
40 
41  //register your products
42  produces< EcalRecHitCollection > (reducedHitsCollection_) ;
43 
44 }
45 
46 
48 {}
49 
50 
51 // ------------ method called to produce the data ------------
52 void
54  const edm::EventSetup& iSetup)
55 {
56  using namespace edm;
57  using namespace std;
58 
59  if (interestingDetIdCollections_.size() < 1)
60  {
61  edm::LogError("ReducedRecHitCollectionProducer") << "VInputTag collections empty" ;
62  return;
63  }
64 
65 
67  iEvent.getByToken(interestingDetIdCollections_[0],detIds);
68  std::vector<DetId> xtalsToStore((*detIds).size());
69  std::copy( (*detIds).begin() , (*detIds).end() , xtalsToStore.begin() );
70 
71  //Merging DetIds from different collections
72  for( unsigned int t = 1; t < interestingDetIdCollections_.size(); ++t )
73  {
76  if(!detId.isValid())
77  {
78  Labels labels;
80  edm::LogError("MissingInput")<<"no reason to skip detid from : (" << labels.module << ", "
81  << labels.productInstance << ", "
82  << labels.process << ")" << std::endl;
83  continue;
84  }
85 
86 
87  for (unsigned int ii=0;ii<(*detId).size();ii++)
88  {
89  if (std::find(xtalsToStore.begin(),xtalsToStore.end(),(*detId)[ii]) == xtalsToStore.end())
90  xtalsToStore.push_back((*detId)[ii]);
91  }
92  }
93 
94  Handle<EcalRecHitCollection> recHitsHandle;
95  iEvent.getByToken(recHitsToken_,recHitsHandle);
96  if( !recHitsHandle.isValid() )
97  {
98  edm::LogError("ReducedRecHitCollectionProducer") << "RecHit collection not found";
99  return;
100  }
101 
102  //Create empty output collections
103  std::auto_ptr< EcalRecHitCollection > miniRecHitCollection (new EcalRecHitCollection) ;
104 
105  for (unsigned int iCry=0;iCry<xtalsToStore.size();iCry++)
106  {
107  EcalRecHitCollection::const_iterator iRecHit = recHitsHandle->find(xtalsToStore[iCry]);
108  if ( (iRecHit != recHitsHandle->end()) && (miniRecHitCollection->find(xtalsToStore[iCry]) == miniRecHitCollection->end()) )
109  miniRecHitCollection->push_back(*iRecHit);
110  }
111 
112  std::sort(xtalsToStore.begin(), xtalsToStore.end());
113  std::unique(xtalsToStore.begin(), xtalsToStore.end());
114 
115  // std::cout << "New Collection " << reducedHitsCollection_ << " size is " << miniRecHitCollection->size() << " original is " << recHitsHandle->size() << std::endl;
116  iEvent.put( miniRecHitCollection,reducedHitsCollection_ );
117 }
T getParameter(std::string const &) const
tuple t
Definition: tree.py:139
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:464
std::vector< EcalRecHit >::const_iterator const_iterator
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
int ii
Definition: cuy.py:588
void find(edm::Handle< EcalRecHitCollection > &hits, DetId thisDet, std::vector< EcalRecHitCollection::const_iterator > &hit, bool debug=false)
Definition: FindCaloHit.cc:7
char const * process
Definition: ProductLabels.h:7
int iEvent
Definition: GenABIO.cc:230
edm::EDGetTokenT< EcalRecHitCollection > recHitsToken_
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:120
bool isValid() const
Definition: HandleBase.h:75
char const * module
Definition: ProductLabels.h:5
ReducedRecHitCollectionProducer(const edm::ParameterSet &)
ctor
std::vector< edm::EDGetTokenT< DetIdCollection > > interestingDetIdCollections_
virtual void produce(edm::Event &, const edm::EventSetup &)
producer
void labelsForToken(EDGetToken iToken, Labels &oLabels) const
char const * productInstance
Definition: ProductLabels.h:6