CMS 3D CMS Logo

CalibratedPhotonProducers.cc
Go to the documentation of this file.
1 //author: Alan Smithee
2 //description:
3 // this class allows the residual scale and smearing to be applied to photon
4 // it will write out all the calibration info in the event, such as scale correction value,
5 // smearing correction value, random nr used, energy post calibration, energy pre calibration
6 // can optionally write out a new collection of photon with the energy corrected by default
7 // a port of EgammaAnalysis/ElectronTools/CalibratedPhotonProducerRun2
8 
18 
25 
26 #include "TRandom2.h"
27 
28 #include <memory>
29 
30 #include <random>
31 #include <vector>
32 
33 template <typename T>
35 public:
38  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
39  void produce(edm::Event&, const edm::EventSetup&) override;
40 
41 private:
42  void setSemiDetRandomSeed(const edm::Event& iEvent, const T& obj, size_t nrObjs, size_t objNr);
43 
46  std::unique_ptr<TRandom> semiDeterministicRng_;
50 
51  static const std::vector<int> valMapsToStore_;
52 };
53 
54 template <typename T>
55 const std::vector<int> CalibratedPhotonProducerT<T>::valMapsToStore_ = {
63 
64 namespace {
65  template <typename HandleType, typename ValType>
66  void fillAndStoreValueMap(edm::Event& iEvent,
67  HandleType objHandle,
68  const std::vector<ValType>& vals,
69  const std::string& name) {
70  auto valMap = std::make_unique<edm::ValueMap<ValType>>();
71  typename edm::ValueMap<ValType>::Filler filler(*valMap);
72  filler.insert(objHandle, vals.begin(), vals.end());
73  filler.fill();
74  iEvent.put(std::move(valMap), name);
75  }
76 } // namespace
77 
78 template <typename T>
80  : photonToken_(consumes(conf.getParameter<edm::InputTag>("src"))),
81  energyCorrector_(conf.getParameter<std::string>("correctionFile")),
82  recHitCollectionEBToken_(consumes(conf.getParameter<edm::InputTag>("recHitCollectionEB"))),
83  recHitCollectionEEToken_(consumes(conf.getParameter<edm::InputTag>("recHitCollectionEE"))),
84  produceCalibratedObjs_(conf.getParameter<bool>("produceCalibratedObjs")) {
85  energyCorrector_.setMinEt(conf.getParameter<double>("minEtToCalibrate"));
86 
87  if (conf.getParameter<bool>("semiDeterministic")) {
88  semiDeterministicRng_ = std::make_unique<TRandom2>();
90  }
91 
93  produces<std::vector<T>>();
94 
95  for (const auto& toStore : valMapsToStore_) {
96  produces<edm::ValueMap<float>>(EGEnergySysIndex::name(toStore));
97  }
98 }
99 
100 template <typename T>
103  desc.add<edm::InputTag>("src", edm::InputTag("gedGsfElectrons"));
104  desc.add<edm::InputTag>("recHitCollectionEB", edm::InputTag("reducedEcalRecHitsEB"));
105  desc.add<edm::InputTag>("recHitCollectionEE", edm::InputTag("reducedEcalRecHitsEE"));
106  desc.add<std::string>("correctionFile", std::string());
107  desc.add<double>("minEtToCalibrate", 5.0);
108  desc.add<bool>("produceCalibratedObjs", true);
109  desc.add<bool>("semiDeterministic", true);
110  std::vector<std::string> valMapsProduced;
111  valMapsProduced.reserve(valMapsToStore_.size());
112  for (auto varToStore : valMapsToStore_)
113  valMapsProduced.push_back(EGEnergySysIndex::name(varToStore));
114  desc.add<std::vector<std::string>>("valueMapsStored", valMapsProduced)
115  ->setComment(
116  "provides to python configs the list of valuemaps stored, can not be overriden in the python config");
117  descriptions.addWithDefaultLabel(desc);
118 }
119 
120 template <typename T>
122  auto inHandle = iEvent.getHandle(photonToken_);
123 
124  auto recHitCollectionEBHandle = iEvent.getHandle(recHitCollectionEBToken_);
125  auto recHitCollectionEEHandle = iEvent.getHandle(recHitCollectionEEToken_);
126 
127  std::unique_ptr<std::vector<T>> out = std::make_unique<std::vector<T>>();
128 
129  size_t nrObj = inHandle->size();
130  std::array<std::vector<float>, EGEnergySysIndex::kNrSysErrs> results;
131  for (auto& res : results)
132  res.reserve(nrObj);
133 
134  const PhotonEnergyCalibrator::EventType evtType =
136 
137  for (const auto& pho : *inHandle) {
138  out->emplace_back(pho);
139 
140  if (semiDeterministicRng_)
141  setSemiDetRandomSeed(iEvent, pho, nrObj, out->size());
142 
144  (pho.isEB()) ? recHitCollectionEBHandle.product() : recHitCollectionEEHandle.product();
145  std::array<float, EGEnergySysIndex::kNrSysErrs> uncertainties =
146  energyCorrector_.calibrate(out->back(), iEvent.id().run(), recHits, iEvent.streamID(), evtType);
147 
148  for (size_t index = 0; index < EGEnergySysIndex::kNrSysErrs; index++) {
149  results[index].push_back(uncertainties[index]);
150  }
151  }
152 
153  auto fillAndStore = [&](auto handle) {
154  for (const auto& mapToStore : valMapsToStore_) {
155  fillAndStoreValueMap(iEvent, handle, results[mapToStore], EGEnergySysIndex::name(mapToStore));
156  }
157  };
158 
159  if (produceCalibratedObjs_) {
160  fillAndStore(iEvent.put(std::move(out)));
161  } else {
162  fillAndStore(inHandle);
163  }
164 }
165 
166 //needs to be synced to CalibratedElectronProducers, want the same seed for a given SC
167 template <typename T>
169  const T& obj,
170  size_t nrObjs,
171  size_t objNr) {
172  if (obj.superCluster().isNonnull()) {
173  semiDeterministicRng_->SetSeed(egamma::getRandomSeedFromSC(iEvent, obj.superCluster()));
174  } else {
175  semiDeterministicRng_->SetSeed(egamma::getRandomSeedFromObj(iEvent, obj, nrObjs, objNr));
176  }
177 }
178 
181 
183 
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
T getParameter(std::string const &) const
Definition: ParameterSet.h:307
std::unique_ptr< TRandom > semiDeterministicRng_
void initPrivateRng(TRandom *rnd)
static const std::string & name(size_t index)
void produce(edm::Event &, const edm::EventSetup &) override
Definition: Electron.h:6
int iEvent
Definition: GenABIO.cc:224
PhotonEnergyCalibrator energyCorrector_
void setSemiDetRandomSeed(const edm::Event &iEvent, const T &obj, size_t nrObjs, size_t objNr)
CalibratedPhotonProducerT(const edm::ParameterSet &)
edm::EDGetTokenT< EcalRecHitCollection > recHitCollectionEEToken_
static constexpr size_t kNrSysErrs
uint32_t getRandomSeedFromObj(const edm::Event &iEvent, const T &obj, size_t nrObjs, size_t objNr)
edm::EDGetTokenT< EcalRecHitCollection > recHitCollectionEBToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< edm::View< T > > photonToken_
uint32_t getRandomSeedFromSC(const edm::Event &iEvent, const reco::SuperClusterRef scRef)
HLT enums.
results
Definition: mysort.py:8
static const std::vector< int > valMapsToStore_
long double T
def move(src, dest)
Definition: eostools.py:511