CMS 3D CMS Logo

EcalRecalibRecHitProducer.cc
Go to the documentation of this file.
1 
16 
19 
21 
22 #include <iostream>
23 #include <cmath>
24 #include <vector>
25 
30 
32  : EBRecHitCollection_(ps.getParameter<edm::InputTag>("EBRecHitCollection")),
33  EERecHitCollection_(ps.getParameter<edm::InputTag>("EERecHitCollection")),
34  EBRecHitToken_((not EBRecHitCollection_.label().empty()) ? consumes<EBRecHitCollection>(EBRecHitCollection_)
35  : edm::EDGetTokenT<EBRecHitCollection>()),
36  EERecHitToken_((not EERecHitCollection_.label().empty()) ? consumes<EERecHitCollection>(EERecHitCollection_)
37  : edm::EDGetTokenT<EERecHitCollection>()),
38  EBRecalibRecHitCollection_(ps.getParameter<std::string>("EBRecalibRecHitCollection")),
39  EERecalibRecHitCollection_(ps.getParameter<std::string>("EERecalibRecHitCollection")),
40  doEnergyScale_(ps.getParameter<bool>("doEnergyScale")),
41  doIntercalib_(ps.getParameter<bool>("doIntercalib")),
42  doLaserCorrections_(ps.getParameter<bool>("doLaserCorrections")),
43 
44  doEnergyScaleInverse_(ps.getParameter<bool>("doEnergyScaleInverse")),
45  doIntercalibInverse_(ps.getParameter<bool>("doIntercalibInverse")),
46  doLaserCorrectionsInverse_(ps.getParameter<bool>("doLaserCorrectionsInverse")) {
47  produces<EBRecHitCollection>(EBRecalibRecHitCollection_);
48  produces<EERecHitCollection>(EERecalibRecHitCollection_);
49 }
50 
52  using namespace edm;
53  Handle<EBRecHitCollection> pEBRecHits;
54  Handle<EERecHitCollection> pEERecHits;
55 
56  const EBRecHitCollection* EBRecHits = nullptr;
57  const EERecHitCollection* EERecHits = nullptr;
58 
59  if (not EBRecHitCollection_.label().empty()) {
60  evt.getByToken(EBRecHitToken_, pEBRecHits);
61  EBRecHits = pEBRecHits.product(); // get a ptr to the product
62  }
63  if (not EERecHitCollection_.label().empty()) {
64  evt.getByToken(EERecHitToken_, pEERecHits);
65  EERecHits = pEERecHits.product(); // get a ptr to the product
66  }
67 
68  // collection of rechits to put in the event
69  auto EBRecalibRecHits = std::make_unique<EBRecHitCollection>();
70  auto EERecalibRecHits = std::make_unique<EERecHitCollection>();
71 
72  // now fetch all conditions we need to make rechits
73  // ADC to GeV constant
75  const EcalADCToGeVConstant* agc = nullptr;
76  float agc_eb = 1.;
77  float agc_ee = 1.;
78  if (doEnergyScale_) {
79  es.get<EcalADCToGeVConstantRcd>().get(pAgc);
80  agc = pAgc.product();
81  // use this value in the algorithm
82  agc_eb = float(agc->getEBValue());
83  agc_ee = float(agc->getEEValue());
84  }
85  // Intercalib constants
87  const EcalIntercalibConstants* ical = nullptr;
88  if (doIntercalib_) {
89  es.get<EcalIntercalibConstantsRcd>().get(pIcal);
90  ical = pIcal.product();
91  }
92  // Laser corrections
94  es.get<EcalLaserDbRecord>().get(pLaser);
95 
97  agc_eb = 1.0 / agc_eb;
98  agc_ee = 1.0 / agc_ee;
99  }
100 
101  if (EBRecHits) {
102  // loop over uncalibrated rechits to make calibrated ones
103  for (EBRecHitCollection::const_iterator it = EBRecHits->begin(); it != EBRecHits->end(); ++it) {
104  EcalIntercalibConstant icalconst = 1.;
105  if (doIntercalib_) {
106  // find intercalib constant for this xtal
107  const EcalIntercalibConstantMap& icalMap = ical->getMap();
108  EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
109  if (icalit != icalMap.end()) {
110  icalconst = (*icalit);
111  } else {
112  edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(it->id())
113  << "! something wrong with EcalIntercalibConstants in your DB? ";
114  }
115  }
116  // get laser coefficient
117  float lasercalib = 1;
118  if (doLaserCorrections_) {
119  lasercalib = pLaser->getLaserCorrection(EBDetId(it->id()), evt.time());
120  }
121 
122  // make the rechit and put in the output collection
123  // must implement op= for EcalRecHit
124 
125  if (doIntercalibInverse_) {
126  icalconst = 1.0 / icalconst;
127  }
129  lasercalib = 1.0 / lasercalib;
130  }
131 
132  EcalRecHit aHit((*it).id(), (*it).energy() * agc_eb * icalconst * lasercalib, (*it).time());
133  EBRecalibRecHits->push_back(aHit);
134  }
135  }
136 
137  if (EERecHits) {
138  // loop over uncalibrated rechits to make calibrated ones
139  for (EERecHitCollection::const_iterator it = EERecHits->begin(); it != EERecHits->end(); ++it) {
140  // find intercalib constant for this xtal
141  EcalIntercalibConstant icalconst = 1.;
142  if (doIntercalib_) {
143  const EcalIntercalibConstantMap& icalMap = ical->getMap();
144  EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
145  if (icalit != icalMap.end()) {
146  icalconst = (*icalit);
147  } else {
148  edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EEDetId(it->id())
149  << "! something wrong with EcalIntercalibConstants in your DB? ";
150  }
151  }
152  // get laser coefficient
153  float lasercalib = 1;
154  if (doLaserCorrections_) {
155  lasercalib = pLaser->getLaserCorrection(EEDetId(it->id()), evt.time());
156  }
157 
158  if (doIntercalibInverse_) {
159  icalconst = 1.0 / icalconst;
160  }
162  lasercalib = 1.0 / lasercalib;
163  }
164 
165  // make the rechit and put in the output collection
166  EcalRecHit aHit((*it).id(), (*it).energy() * agc_ee * icalconst * lasercalib, (*it).time());
167  EERecalibRecHits->push_back(aHit);
168  }
169  }
170  // put the collection of recunstructed hits in the event
171  LogInfo("EcalRecalibRecHitInfo") << "total # EB re-calibrated rechits: " << EBRecalibRecHits->size();
172  LogInfo("EcalRecalibRecHitInfo") << "total # EE re-calibrated rechits: " << EERecalibRecHits->size();
173 
174  evt.put(std::move(EBRecalibRecHits), EBRecalibRecHitCollection_);
175  evt.put(std::move(EERecalibRecHits), EERecalibRecHitCollection_);
176 }
177 
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:131
void produce(edm::StreamID sid, edm::Event &evt, const edm::EventSetup &es) const override
const edm::InputTag EERecHitCollection_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:525
const self & getMap() const
float getLaserCorrection(DetId const &xid, edm::Timestamp const &iTime) const
std::vector< EcalRecHit >::const_iterator const_iterator
char const * label
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
const_iterator end() const
EcalRecalibRecHitProducer(const edm::ParameterSet &ps)
T const * product() const
Definition: Handle.h:69
std::vector< Item >::const_iterator const_iterator
const edm::EDGetTokenT< EERecHitCollection > EERecHitToken_
const std::string EERecalibRecHitCollection_
std::string const & label() const
Definition: InputTag.h:36
const edm::EDGetTokenT< EBRecHitCollection > EBRecHitToken_
const edm::InputTag EBRecHitCollection_
HLT enums.
T get() const
Definition: EventSetup.h:73
const_iterator find(uint32_t rawId) const
const_iterator end() const
const std::string EBRecalibRecHitCollection_
edm::Timestamp time() const
Definition: EventBase.h:60
def move(src, dest)
Definition: eostools.py:511
const_iterator begin() const
float EcalIntercalibConstant