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