CMS 3D CMS Logo

CaloRecHitsBeamHaloCleaned.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: RecoMET/METProducers
4 // Class: CaloRecHitsBeamHaloCleaned
5 //
13 //
14 // Original Author: Thomas Laurent
15 // Created: Tue, 09 Feb 2016 13:09:37 GMT
16 //
17 //
18 
19 
20 // system include files
21 #include <memory>
22 
23 // user include files
26 
29 
32 
35 #include <vector>
36 #include <iostream>
46 
47 using namespace reco;
48 
50 public:
52  ~CaloRecHitsBeamHaloCleaned() override;
53 
54  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
55 
56 private:
57  void produce(edm::Event&, const edm::EventSetup&) override;
58 
59 
64 
65  //Input tags
70 
71 
72  bool ishlt;
73 };
74 
75 //
76 // constructors and destructor
77 //
79 {
80 
81  ishlt = iConfig.getParameter< bool> ("IsHLT");
82 
83 
84 
85  produces<EcalRecHitCollection>("EcalRecHitsEB");
86  produces<EcalRecHitCollection>("EcalRecHitsEE");
87  produces<HBHERecHitCollection>();
88 
89 
90  it_EBRecHits = iConfig.getParameter<edm::InputTag>("EBRecHitsLabel");
91  it_EERecHits = iConfig.getParameter<edm::InputTag>("EERecHitsLabel");
92  it_HBHERecHits = iConfig.getParameter<edm::InputTag>("HBHERecHitsLabel");
93  it_GlobalHaloData = iConfig.getParameter<edm::InputTag>("GlobalHaloDataLabel");
94 
95 
96  ecalebhits_token= consumes<EcalRecHitCollection>(it_EBRecHits);
97  ecaleehits_token= consumes<EcalRecHitCollection>(it_EERecHits);
98  hbhehits_token= consumes<HBHERecHitCollection>(it_HBHERecHits);
99 
100  globalhalo_token=consumes<GlobalHaloData>(it_GlobalHaloData);
101 }
102 
103 
105 {
106 
107  // do anything here that needs to be done at destruction time
108  // (e.g. close files, deallocate resources etc.)
109 
110 }
111 
112 
113 //
114 // member functions
115 //
116 
117 // ------------ method called to produce the data ------------
118 void
120 {
121  using namespace edm;
122  using namespace reco;
123  using namespace std;
124 
125 
126  Handle<EcalRecHitCollection> ebrhitsuncleaned;
127  iEvent.getByToken(ecalebhits_token, ebrhitsuncleaned );
128 
129  Handle<EcalRecHitCollection> eerhitsuncleaned;
130  iEvent.getByToken(ecaleehits_token, eerhitsuncleaned );
131 
132  Handle<HBHERecHitCollection> hbherhitsuncleaned;
133  iEvent.getByToken(hbhehits_token, hbherhitsuncleaned );
134 
135  // Get GlobalHaloData
136  edm::Handle<reco::GlobalHaloData> TheGlobalHaloData;
137  iEvent.getByToken(globalhalo_token, TheGlobalHaloData);
138 
139 
140 
141  const GlobalHaloData TheSummaryHalo = (*TheGlobalHaloData );
142 
143 
144  //Cleaning of the various rechits collections:
145 
146  // EcalRecHit EB
147  auto ebrhitscleaned = std::make_unique<EcalRecHitCollection>();
148  for(unsigned int i = 0; i < ebrhitsuncleaned->size(); i++){
149  const EcalRecHit & rhit = (*ebrhitsuncleaned)[i];
150  bool isclean(true);
151  const edm::RefVector<EcalRecHitCollection>& refbeamhalorechits = TheSummaryHalo.GetEBRechits();
152  for(unsigned int j = 0; j <refbeamhalorechits.size() ; j++){
153  const EcalRecHit &rhitbeamhalo = *(refbeamhalorechits)[j];
154  if( rhit.detid() == rhitbeamhalo.detid() ) {
155  isclean = false;
156  break;
157  }
158  }
159  if(isclean) ebrhitscleaned->push_back(rhit);
160  }
161 
162  // EcalRecHit EE
163  auto eerhitscleaned = std::make_unique<EcalRecHitCollection>();
164  for(unsigned int i = 0; i < eerhitsuncleaned->size(); i++){
165  const EcalRecHit & rhit = (*eerhitsuncleaned)[i];
166  bool isclean(true);
167  const edm::RefVector<EcalRecHitCollection>& refbeamhalorechits = TheSummaryHalo.GetEERechits();
168  for(unsigned int j = 0; j <refbeamhalorechits.size() ; j++){
169  const EcalRecHit &rhitbeamhalo = *(refbeamhalorechits)[j];
170  if( rhit.detid() == rhitbeamhalo.detid() ) {
171  isclean = false;
172  break;
173  }
174  }
175  if(isclean) eerhitscleaned->push_back(rhit);
176  }
177 
178  // HBHERecHit
179  auto hbherhitscleaned = std::make_unique<HBHERecHitCollection>();
180  for(unsigned int i = 0; i < hbherhitsuncleaned->size(); i++){
181  const HBHERecHit & rhit = (*hbherhitsuncleaned)[i];
182  bool isclean(true);
183  const edm::RefVector<HBHERecHitCollection>& refbeamhalorechits = TheSummaryHalo.GetHBHERechits();
184  for(unsigned int j = 0; j <refbeamhalorechits.size() ; j++){
185  const HBHERecHit &rhitbeamhalo = *(refbeamhalorechits)[j];
186  if( rhit.detid() == rhitbeamhalo.detid() ) {
187  isclean = false;
188  break;
189  }
190  }
191  if(isclean) hbherhitscleaned->push_back(rhit);
192  }
193 
194 
195 
196  iEvent.put(std::move(ebrhitscleaned),"EcalRecHitsEB");
197  iEvent.put(std::move(eerhitscleaned),"EcalRecHitsEE");
198  iEvent.put(std::move(hbherhitscleaned));
199 
200 }
201 
202 
203 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
204 void
207  desc.add<edm::InputTag>("EBRecHitsLabel", edm::InputTag("EcalRecHit","EcalRecHitsEB"));
208  desc.add<edm::InputTag>("EERecHitsLabel", edm::InputTag("EcalRecHit","EcalRecHitsEE"));
209  desc.add<edm::InputTag>("HBHERecHitsLabel", edm::InputTag("hbhereco"));
210  desc.add<edm::InputTag>("GlobalHaloDataLabel", edm::InputTag("GlobalHaloData"));
211  desc.add<bool>("IsHLT",false);
212  descriptions.add("caloRecHitsBeamHaloCleaned",desc);
213 
214 }
215 
216 //define this as a plug-in
T getParameter(std::string const &) const
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
const DetId & detid() const
Definition: CaloRecHit.h:21
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
CaloRecHitsBeamHaloCleaned(const edm::ParameterSet &)
const DetId & detid() const
Definition: EcalRecHit.h:72
int iEvent
Definition: GenABIO.cc:230
edm::RefVector< HBHERecHitCollection > & GetHBHERechits()
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
ParameterDescriptionBase * add(U const &iLabel, T const &value)
edm::EDGetTokenT< GlobalHaloData > globalhalo_token
edm::EDGetTokenT< EcalRecHitCollection > ecalebhits_token
void produce(edm::Event &, const edm::EventSetup &) override
void add(std::string const &label, ParameterSetDescription const &psetDescription)
fixed size matrix
HLT enums.
size_type size() const
size_type size() const
Size of the RefVector.
Definition: RefVector.h:107
edm::RefVector< EcalRecHitCollection > & GetEBRechits()
edm::EDGetTokenT< HBHERecHitCollection > hbhehits_token
edm::RefVector< EcalRecHitCollection > & GetEERechits()
def move(src, dest)
Definition: eostools.py:510
edm::EDGetTokenT< EcalRecHitCollection > ecaleehits_token