CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
HitReCalibrator.cc
Go to the documentation of this file.
1 // system include files
2 #include <memory>
3 #include <string>
4 
5 // user include files
13 
16 
28 
29 //
30 // class declaration
31 //
32 
33 namespace cms {
34 
36  public:
37  explicit HitReCalibrator(const edm::ParameterSet&);
38  ~HitReCalibrator() override;
39 
40  void beginJob() override;
41 
42  void produce(edm::Event&, const edm::EventSetup&) override;
43 
44  private:
45  // ----------member data ---------------------------
46 
48 
52 
54  };
55 } // end namespace cms
56 
57 namespace cms {
58 
60  tok_hbhe_ = consumes<HBHERecHitCollection>(iConfig.getParameter<edm::InputTag>("hbheInput"));
61  tok_ho_ = consumes<HORecHitCollection>(iConfig.getParameter<edm::InputTag>("hoInput"));
62  tok_hf_ = consumes<HFRecHitCollection>(iConfig.getParameter<edm::InputTag>("hfInput"));
63  allowMissingInputs_ = true;
64  tok_resp_ = esConsumes<HcalRespCorrs, HcalRespCorrsRcd>();
65 
66  //register your products
67 
68  produces<HBHERecHitCollection>("DiJetsHBHEReRecHitCollection");
69  produces<HORecHitCollection>("DiJetsHOReRecHitCollection");
70  produces<HFRecHitCollection>("DiJetsHFReRecHitCollection");
71  }
73 
75 
76  // ------------ method called to produce the data ------------
78  auto miniDiJetsHBHERecHitCollection = std::make_unique<HBHERecHitCollection>();
79  auto miniDiJetsHORecHitCollection = std::make_unique<HORecHitCollection>();
80  auto miniDiJetsHFRecHitCollection = std::make_unique<HFRecHitCollection>();
81 
82  const HcalRespCorrs* jetRecalib = &iSetup.getData(tok_resp_);
83 
84  try {
86  iEvent.getByToken(tok_hbhe_, hbhe);
87  const HBHERecHitCollection Hithbhe = *(hbhe.product());
88  for (HBHERecHitCollection::const_iterator hbheItr = Hithbhe.begin(); hbheItr != Hithbhe.end(); hbheItr++) {
89  DetId id = hbheItr->detid();
90  float recal;
91  if (jetRecalib->exists(id))
92  recal = jetRecalib->getValues(id)->getValue();
93  else
94  recal = 1.;
95  float energy = hbheItr->energy();
96  float time = hbheItr->time();
97  HBHERecHit* hbhehit = new HBHERecHit(id, recal * energy, time);
98  miniDiJetsHBHERecHitCollection->push_back(*hbhehit);
99  }
100  } catch (cms::Exception& e) { // can't find it!
101  if (!allowMissingInputs_) {
102  edm::LogError("HitCalib") << "No HBHE collection ";
103  throw e;
104  }
105  }
106 
107  try {
109  iEvent.getByToken(tok_ho_, ho);
110  const HORecHitCollection Hitho = *(ho.product());
111  for (HORecHitCollection::const_iterator hoItr = Hitho.begin(); hoItr != Hitho.end(); hoItr++) {
112  DetId id = hoItr->detid();
113  float recal;
114  if (jetRecalib->exists(id))
115  recal = jetRecalib->getValues(id)->getValue();
116  else
117  recal = 1.;
118  float energy = hoItr->energy();
119  float time = hoItr->time();
120  HORecHit* hohit = new HORecHit(id, recal * energy, time);
121  miniDiJetsHORecHitCollection->push_back(*hohit);
122  }
123  } catch (cms::Exception& e) { // can't find it!
124  if (!allowMissingInputs_) {
125  edm::LogError("HitCalib") << " No HO collection ";
126  throw e;
127  }
128  }
129 
130  try {
132  iEvent.getByToken(tok_hf_, hf);
133  const HFRecHitCollection Hithf = *(hf.product());
134  for (HFRecHitCollection::const_iterator hfItr = Hithf.begin(); hfItr != Hithf.end(); hfItr++) {
135  DetId id = hfItr->detid();
136  float recal;
137  if (jetRecalib->exists(id))
138  recal = jetRecalib->getValues(id)->getValue();
139  else
140  recal = 1.;
141  float energy = hfItr->energy();
142  float time = hfItr->time();
143  HFRecHit* hfhit = new HFRecHit(id, recal * energy, time);
144  miniDiJetsHFRecHitCollection->push_back(*hfhit);
145  }
146  } catch (cms::Exception& e) { // can't find it!
147  if (!allowMissingInputs_)
148  throw e;
149  }
150 
151  //Put selected information in the event
152 
153  iEvent.put(std::move(miniDiJetsHBHERecHitCollection), "DiJetsHBHEReRecHitCollection");
154 
155  iEvent.put(std::move(miniDiJetsHORecHitCollection), "DiJetsHOReRecHitCollection");
156 
157  iEvent.put(std::move(miniDiJetsHFRecHitCollection), "DiJetsHFReRecHitCollection");
158  }
159 } // namespace cms
160 
162 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:133
HitReCalibrator(const edm::ParameterSet &)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
std::vector< T >::const_iterator const_iterator
Log< level::Error, false > LogError
void beginJob() override
bool getData(T &iHolder) const
Definition: EventSetup.h:128
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< HORecHitCollection > tok_ho_
def move
Definition: eostools.py:511
const_iterator end() const
edm::EDGetTokenT< HFRecHitCollection > tok_hf_
Definition: DetId.h:17
void produce(edm::Event &, const edm::EventSetup &) override
T const * product() const
Definition: Handle.h:70
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
edm::ESGetToken< HcalRespCorrs, HcalRespCorrsRcd > tok_resp_
const_iterator begin() const
edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_