CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
EcalRecalibRecHitProducer.cc
Go to the documentation of this file.
1 
19 
22 
24 
25 #include <iostream>
26 #include <cmath>
27 #include <vector>
28 
29 //#include "CondFormats/EcalObjects/interface/EcalPedestals.h"
30 //#include "CondFormats/DataRecord/interface/EcalPedestalsRcd.h"
36 
37 
38 
40 
41  EBRecHitCollection_ = ps.getParameter<edm::InputTag>("EBRecHitCollection");
42  EERecHitCollection_ = ps.getParameter<edm::InputTag>("EERecHitCollection");
43  EBRecalibRecHitCollection_ = ps.getParameter<std::string>("EBRecalibRecHitCollection");
44  EERecalibRecHitCollection_ = ps.getParameter<std::string>("EERecalibRecHitCollection");
45  doEnergyScale_ = ps.getParameter<bool>("doEnergyScale");
46  doIntercalib_ = ps.getParameter<bool>("doIntercalib");
47  doLaserCorrections_ = ps.getParameter<bool>("doLaserCorrections");
48 
49  doEnergyScaleInverse_ = ps.getParameter<bool>("doEnergyScaleInverse");
50  doIntercalibInverse_ = ps.getParameter<bool>("doIntercalibInverse");
51  doLaserCorrectionsInverse_ = ps.getParameter<bool>("doLaserCorrectionsInverse");
52 
55 
56  produces< EBRecHitCollection >(EBRecalibRecHitCollection_);
57  produces< EERecHitCollection >(EERecalibRecHitCollection_);
58 }
59 
61 
62  if (EBalgo_) delete EBalgo_;
63  if (EEalgo_) delete EEalgo_;
64 
65 }
66 
68 {
69  using namespace edm;
72 
73  const EBRecHitCollection* EBRecHits = 0;
74  const EERecHitCollection* EERecHits = 0;
75 
76  // if ( EBRecHitCollection_.label() != "" && EBRecHitCollection_.instance() != "" ) {
77  if ( EBRecHitCollection_.label() != "" ) {
78  evt.getByLabel( EBRecHitCollection_, pEBRecHits);
79  if ( pEBRecHits.isValid() ) {
80  EBRecHits = pEBRecHits.product(); // get a ptr to the product
81 #ifdef DEBUG
82  LogDebug("EcalRecHitDebug") << "total # EB rechits to be re-calibrated: " << EBRecHits->size();
83 #endif
84  } else {
85  edm::LogError("EcalRecHitError") << "Error! can't get the product " << EBRecHitCollection_.label() ;
86  }
87  }
88 
89  // if ( EERecHitCollection_.label() != "" && EERecHitCollection_.instance() != "" ) {
90  if ( EERecHitCollection_.label() != "" ) {
91  evt.getByLabel( EERecHitCollection_, pEERecHits);
92  if ( pEERecHits.isValid() ) {
93  EERecHits = pEERecHits.product(); // get a ptr to the product
94 #ifdef DEBUG
95  LogDebug("EcalRecHitDebug") << "total # EE uncalibrated rechits to be re-calibrated: " << EERecHits->size();
96 #endif
97  } else {
98  edm::LogError("EcalRecHitError") << "Error! can't get the product " << EERecHitCollection_.label() ;
99  }
100  }
101 
102  // collection of rechits to put in the event
103  std::auto_ptr< EBRecHitCollection > EBRecalibRecHits( new EBRecHitCollection );
104  std::auto_ptr< EERecHitCollection > EERecalibRecHits( new EERecHitCollection );
105 
106  // now fetch all conditions we need to make rechits
107  // ADC to GeV constant
109  const EcalADCToGeVConstant *agc = 0;
110  float agc_eb = 1.;
111  float agc_ee = 1.;
112  if (doEnergyScale_) {
113  es.get<EcalADCToGeVConstantRcd>().get(pAgc);
114  agc = pAgc.product();
115  // use this value in the algorithm
116  agc_eb = float(agc->getEBValue());
117  agc_ee = float(agc->getEEValue());
118  }
119  // Intercalib constants
121  const EcalIntercalibConstants *ical = 0;
122  if (doIntercalib_) {
123  es.get<EcalIntercalibConstantsRcd>().get(pIcal);
124  ical = pIcal.product();
125  }
126  // Laser corrections
128  es.get<EcalLaserDbRecord>().get( pLaser );
129 
130 
132  agc_eb = 1.0/agc_eb;
133  agc_ee = 1.0/agc_ee;
134  }
135 
136 
137  if (EBRecHits) {
138  // loop over uncalibrated rechits to make calibrated ones
139  for(EBRecHitCollection::const_iterator it = EBRecHits->begin(); it != EBRecHits->end(); ++it) {
140 
141  EcalIntercalibConstant icalconst = 1.;
142  if (doIntercalib_) {
143  // find intercalib constant for this xtal
144  const EcalIntercalibConstantMap &icalMap = ical->getMap();
145  EcalIntercalibConstantMap::const_iterator icalit = icalMap.find(it->id());
146  if( icalit!=icalMap.end() ){
147  icalconst = (*icalit);
148  } else {
149  edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EBDetId(it->id()) << "! something wrong with EcalIntercalibConstants in your DB? "
150  ;
151  }
152  }
153  // get laser coefficient
154  float lasercalib = 1;
155  if (doLaserCorrections_) {
156  lasercalib = pLaser->getLaserCorrection( EBDetId(it->id()), evt.time() );
157  }
158 
159  // make the rechit and put in the output collection
160  // must implement op= for EcalRecHit
161 
163  icalconst = 1.0/icalconst;
164  }
166  lasercalib = 1.0/lasercalib;
167  }
168 
169  EcalRecHit aHit( (*it).id(), (*it).energy() * agc_eb * icalconst * lasercalib, (*it).time() );
170  EBRecalibRecHits->push_back( aHit );
171  }
172  }
173 
174  if (EERecHits)
175  {
176  // loop over uncalibrated rechits to make calibrated ones
178  it != EERecHits->end(); ++it) {
179 
180  // find intercalib constant for this xtal
181  EcalIntercalibConstant icalconst = 1.;
182  if (doIntercalib_) {
183  const EcalIntercalibConstantMap &icalMap = ical->getMap();
184  EcalIntercalibConstantMap::const_iterator icalit=icalMap.find(it->id());
185  if( icalit!=icalMap.end() ) {
186  icalconst = (*icalit);
187  } else {
188  edm::LogError("EcalRecHitError") << "No intercalib const found for xtal " << EEDetId(it->id()) << "! something wrong with EcalIntercalibConstants in your DB? ";
189  }
190  }
191  // get laser coefficient
192  float lasercalib = 1;
193  if (doLaserCorrections_) {
194  lasercalib = pLaser->getLaserCorrection( EEDetId(it->id()), evt.time() );
195  }
196 
197  // make the rechit and put in the output collection
198  // must implement op= for EcalRecHit
199  //EcalRecHit aHit( EEalgo_->makeRecHit(*it, icalconst * lasercalib) );
200 
202  icalconst = 1.0/icalconst;
203  }
205  lasercalib = 1.0/lasercalib;
206  }
207 
208  EcalRecHit aHit( (*it).id(), (*it).energy() * agc_ee * icalconst * lasercalib, (*it).time() );
209  EERecalibRecHits->push_back( aHit );
210  }
211  }
212  // put the collection of recunstructed hits in the event
213  LogInfo("EcalRecalibRecHitInfo") << "total # EB re-calibrated rechits: " << EBRecalibRecHits->size();
214  LogInfo("EcalRecalibRecHitInfo") << "total # EE re-calibrated rechits: " << EERecalibRecHits->size();
215 
216  evt.put( EBRecalibRecHits, EBRecalibRecHitCollection_ );
217  evt.put( EERecalibRecHits, EERecalibRecHitCollection_ );
218 }
219 
#define LogDebug(id)
T getParameter(std::string const &) const
const self & getMap() const
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
std::vector< EcalRecHit >::const_iterator const_iterator
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:85
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
EcalRecalibRecHitProducer(const edm::ParameterSet &ps)
const T & get() const
Definition: EventSetup.h:55
std::vector< Item >::const_iterator const_iterator
virtual void produce(edm::Event &evt, const edm::EventSetup &es)
std::string const & label() const
Definition: InputTag.h:25
const_iterator find(uint32_t rawId) const
const_iterator end() const
edm::Timestamp time() const
Definition: EventBase.h:57
float EcalIntercalibConstant