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