CMS 3D CMS Logo

EcalRecHitRecalib.cc
Go to the documentation of this file.
3 
8 
10  : ecalHitsProducer_(iConfig.getParameter<std::string>("ecalRecHitsProducer")),
11  barrelHits_(iConfig.getParameter<std::string>("barrelHitCollection")),
12  endcapHits_(iConfig.getParameter<std::string>("endcapHitCollection")),
13  recalibBarrelHits_(iConfig.getParameter<std::string>("RecalibBarrelHitCollection")),
14  recalibEndcapHits_(iConfig.getParameter<std::string>("RecalibEndcapHitCollection")),
15  refactor_(iConfig.getUntrackedParameter<double>("Refactor", (double)1)),
16  refactor_mean_(iConfig.getUntrackedParameter<double>("Refactor_mean", (double)1)),
17  ebRecHitToken_(consumes<EBRecHitCollection>(edm::InputTag(ecalHitsProducer_, barrelHits_))),
18  eeRecHitToken_(consumes<EERecHitCollection>(edm::InputTag(ecalHitsProducer_, endcapHits_))),
19  intercalibConstsToken_(esConsumes()),
20  barrelHitsToken_(produces<EBRecHitCollection>(recalibBarrelHits_)),
21  endcapHitsToken_(produces<EERecHitCollection>(recalibEndcapHits_)) {}
22 
24 
25 // ------------ method called to produce the data ------------
27  using namespace edm;
28  using namespace std;
29 
30  Handle<EBRecHitCollection> barrelRecHitsHandle;
31  Handle<EERecHitCollection> endcapRecHitsHandle;
32 
33  const EBRecHitCollection* EBRecHits = nullptr;
34  const EERecHitCollection* EERecHits = nullptr;
35 
36  iEvent.getByToken(ebRecHitToken_, barrelRecHitsHandle);
37  if (!barrelRecHitsHandle.isValid()) {
38  LogDebug("") << "EcalREcHitMiscalib: Error! can't get product!" << std::endl;
39  } else {
40  EBRecHits = barrelRecHitsHandle.product(); // get a ptr to the product
41  }
42 
43  iEvent.getByToken(eeRecHitToken_, endcapRecHitsHandle);
44  if (!endcapRecHitsHandle.isValid()) {
45  LogDebug("") << "EcalREcHitMiscalib: Error! can't get product!" << std::endl;
46  } else {
47  EERecHits = endcapRecHitsHandle.product(); // get a ptr to the product
48  }
49 
50  //Create empty output collections
51  auto RecalibEBRecHitCollection = std::make_unique<EBRecHitCollection>();
52  auto RecalibEERecHitCollection = std::make_unique<EERecHitCollection>();
53 
54  // Intercalib constants
56 
57  if (EBRecHits) {
58  //loop on all EcalRecHits (barrel)
60  for (itb = EBRecHits->begin(); itb != EBRecHits->end(); ++itb) {
61  // find intercalib constant for this xtal
62  EcalIntercalibConstantMap::const_iterator icalit = ical.getMap().find(itb->id().rawId());
63  EcalIntercalibConstant icalconst = -1;
64 
65  if (icalit != ical.getMap().end()) {
66  icalconst = (*icalit);
67 
68  } else {
69  edm::LogError("EcalRecHitRecalib") << "No intercalib const found for xtal " << EBDetId(itb->id())
70  << "! something wrong with EcalIntercalibConstants in your DB? ";
71  }
72 
73  // make the rechit with rescaled energy and put in the output collection
74  icalconst = refactor_mean_ +
75  (icalconst - refactor_mean_) * refactor_; //apply additional scaling factor (works if gaussian)
76  EcalRecHit aHit(itb->id(), itb->energy() * icalconst, itb->time());
77 
78  RecalibEBRecHitCollection->push_back(aHit);
79  }
80  }
81 
82  if (EERecHits) {
83  //loop on all EcalRecHits (barrel)
85  for (ite = EERecHits->begin(); ite != EERecHits->end(); ++ite) {
86  // find intercalib constant for this xtal
87  EcalIntercalibConstantMap::const_iterator icalit = ical.getMap().find(ite->id().rawId());
88  EcalIntercalibConstant icalconst = -1;
89 
90  if (icalit != ical.getMap().end()) {
91  icalconst = (*icalit);
92  } else {
93  edm::LogError("EcalRecHitRecalib") << "No intercalib const found for xtal " << EEDetId(ite->id())
94  << "! something wrong with EcalIntercalibConstants in your DB? ";
95  }
96 
97  // make the rechit with rescaled energy and put in the output collection
98 
99  icalconst = refactor_mean_ +
100  (icalconst - refactor_mean_) * refactor_; //apply additional scaling factor (works if gaussian)
101  EcalRecHit aHit(ite->id(), ite->energy() * icalconst, ite->time());
102 
103  RecalibEERecHitCollection->push_back(aHit);
104  }
105  }
106 
107  //Put Recalibrated rechit in the event
108  iEvent.put(barrelHitsToken_, std::move(RecalibEBRecHitCollection));
109  iEvent.put(endcapHitsToken_, std::move(RecalibEERecHitCollection));
110 }
ESGetTokenH3DDVariant esConsumes(std::string const &Record, edm::ConsumesCollector &)
Definition: DeDxTools.cc:283
const edm::EDPutTokenT< EBRecHitCollection > barrelHitsToken_
~EcalRecHitRecalib() override
T const * product() const
Definition: Handle.h:70
std::vector< EcalRecHit >::const_iterator const_iterator
EcalRecHitRecalib(const edm::ParameterSet &)
Log< level::Error, false > LogError
const edm::ESGetToken< EcalIntercalibConstants, EcalIntercalibConstantsRcd > intercalibConstsToken_
int iEvent
Definition: GenABIO.cc:224
const edm::EDPutTokenT< EERecHitCollection > endcapHitsToken_
bool getData(T &iHolder) const
Definition: EventSetup.h:122
const double refactor_mean_
std::vector< Item >::const_iterator const_iterator
const edm::EDGetTokenT< EERecHitCollection > eeRecHitToken_
const double refactor_
bool isValid() const
Definition: HandleBase.h:70
const edm::EDGetTokenT< EBRecHitCollection > ebRecHitToken_
HLT enums.
void produce(edm::Event &, const edm::EventSetup &) override
def move(src, dest)
Definition: eostools.py:511
float EcalIntercalibConstant
#define LogDebug(id)