CMS 3D CMS Logo

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"))),
64  allowMissingInputs_ = true;
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  const HBHERecHitCollection Hithbhe = *(hbhe.product());
87  for (HBHERecHitCollection::const_iterator hbheItr = Hithbhe.begin(); hbheItr != Hithbhe.end(); hbheItr++) {
88  DetId id = hbheItr->detid();
89  float recal;
90  if (jetRecalib->exists(id))
91  recal = jetRecalib->getValues(id)->getValue();
92  else
93  recal = 1.;
94  float energy = hbheItr->energy();
95  float time = hbheItr->time();
96  HBHERecHit* hbhehit = new HBHERecHit(id, recal * energy, time);
97  miniDiJetsHBHERecHitCollection->push_back(*hbhehit);
98  }
99  } catch (cms::Exception& e) { // can't find it!
100  if (!allowMissingInputs_) {
101  edm::LogError("HitCalib") << "No HBHE collection ";
102  throw e;
103  }
104  }
105 
106  try {
107  const edm::Handle<HORecHitCollection>& ho = iEvent.getHandle(tok_ho_);
108  const HORecHitCollection Hitho = *(ho.product());
109  for (HORecHitCollection::const_iterator hoItr = Hitho.begin(); hoItr != Hitho.end(); hoItr++) {
110  DetId id = hoItr->detid();
111  float recal;
112  if (jetRecalib->exists(id))
113  recal = jetRecalib->getValues(id)->getValue();
114  else
115  recal = 1.;
116  float energy = hoItr->energy();
117  float time = hoItr->time();
118  HORecHit* hohit = new HORecHit(id, recal * energy, time);
119  miniDiJetsHORecHitCollection->push_back(*hohit);
120  }
121  } catch (cms::Exception& e) { // can't find it!
122  if (!allowMissingInputs_) {
123  edm::LogError("HitCalib") << " No HO collection ";
124  throw e;
125  }
126  }
127 
128  try {
129  const edm::Handle<HFRecHitCollection>& hf = iEvent.getHandle(tok_hf_);
130  const HFRecHitCollection Hithf = *(hf.product());
131  for (HFRecHitCollection::const_iterator hfItr = Hithf.begin(); hfItr != Hithf.end(); hfItr++) {
132  DetId id = hfItr->detid();
133  float recal;
134  if (jetRecalib->exists(id))
135  recal = jetRecalib->getValues(id)->getValue();
136  else
137  recal = 1.;
138  float energy = hfItr->energy();
139  float time = hfItr->time();
140  HFRecHit* hfhit = new HFRecHit(id, recal * energy, time);
141  miniDiJetsHFRecHitCollection->push_back(*hfhit);
142  }
143  } catch (cms::Exception& e) { // can't find it!
144  if (!allowMissingInputs_)
145  throw e;
146  }
147 
148  //Put selected information in the event
149 
150  iEvent.put(std::move(miniDiJetsHBHERecHitCollection), "DiJetsHBHEReRecHitCollection");
151 
152  iEvent.put(std::move(miniDiJetsHORecHitCollection), "DiJetsHOReRecHitCollection");
153 
154  iEvent.put(std::move(miniDiJetsHFRecHitCollection), "DiJetsHFReRecHitCollection");
155  }
156 } // namespace cms
157 
159 
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
HitReCalibrator(const edm::ParameterSet &)
const edm::EDGetTokenT< HORecHitCollection > tok_ho_
std::vector< T >::const_iterator const_iterator
Log< level::Error, false > LogError
void beginJob() override
const edm::EDGetTokenT< HFRecHitCollection > tok_hf_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const edm::EDGetTokenT< HBHERecHitCollection > tok_hbhe_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const_iterator begin() const
Namespace of DDCMS conversion namespace.
const_iterator end() const
Definition: DetId.h:17
void produce(edm::Event &, const edm::EventSetup &) override
const edm::ESGetToken< HcalRespCorrs, HcalRespCorrsRcd > tok_resp_
HLT enums.
def move(src, dest)
Definition: eostools.py:511