CMS 3D CMS Logo

CalibratedElectronProducers.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 electrons
4 // the scale and smearing is on the ecal part of the energy
5 // hence the E/p combination needs to be re-don, hence the E/p Combination Tools
6 // it re-applies the regression with the new corrected ecal energy
7 // returns a vector of calibrated energies and correction data, indexed by EGEnergySysIndex
8 // a port of EgammaAnalysis/ElectronTools/CalibratedElectronProducerRun2
9 
19 
24 
31 
32 #include "TRandom2.h"
33 
34 #include <vector>
35 
36 
37 template<typename T>
39 {
40 public:
43  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
44  void produce( edm::Event &, const edm::EventSetup & ) override ;
45 
46 private:
47  void setSemiDetRandomSeed(const edm::Event& iEvent,const T& obj,size_t nrObjs,size_t objNr);
48 
50 
53  std::unique_ptr<TRandom> semiDeterministicRng_;
57  static const std::vector<int> valMapsToStore_;
58 
59 };
60 
61 template<typename T>
62 const std::vector<int> CalibratedElectronProducerT<T>::valMapsToStore_ = {
88 };
89 
90 namespace{
91  template<typename HandleType,typename ValType>
92  void fillAndStoreValueMap(edm::Event& iEvent,HandleType objHandle,
93  const std::vector<ValType>& vals,const std::string& name)
94  {
95  auto valMap = std::make_unique<edm::ValueMap<ValType> >();
96  typename edm::ValueMap<ValType>::Filler filler(*valMap);
97  filler.insert(objHandle,vals.begin(),vals.end());
98  filler.fill();
99  iEvent.put(std::move(valMap),name);
100  }
101 }
102 
103 template<typename T>
105  electronToken_(consumes<edm::View<T>>(conf.getParameter<edm::InputTag>("src"))),
106  epCombinationTool_(conf.getParameter<edm::ParameterSet>("epCombConfig")),
107  energyCorrector_(epCombinationTool_, conf.getParameter<std::string>("correctionFile")),
108  recHitCollectionEBToken_(consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitCollectionEB"))),
109  recHitCollectionEEToken_(consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitCollectionEE"))),
110  produceCalibratedObjs_(conf.getParameter<bool>("produceCalibratedObjs"))
111 {
112  energyCorrector_.setMinEt(conf.getParameter<double>("minEtToCalibrate"));
113 
114  if (conf.getParameter<bool>("semiDeterministic")) {
115  semiDeterministicRng_.reset(new TRandom2());
117  }
118 
119  if(produceCalibratedObjs_) produces<std::vector<T>>();
120 
121  for(const auto& toStore : valMapsToStore_){
122  produces<edm::ValueMap<float>>(EGEnergySysIndex::name(toStore));
123  }
124 }
125 
126 template<typename T>
128 {
130  desc.add<edm::InputTag>("src",edm::InputTag("gedPhotons"));
132  desc.add<edm::InputTag>("recHitCollectionEB",edm::InputTag("reducedEcalRecHitsEB"));
133  desc.add<edm::InputTag>("recHitCollectionEE",edm::InputTag("reducedEcalRecHitsEE"));
134  desc.add<std::string>("correctionFile",std::string());
135  desc.add<double>("minEtToCalibrate",5.0);
136  desc.add<bool>("produceCalibratedObjs",true);
137  desc.add<bool>("semiDeterministic",true);
138  std::vector<std::string> valMapsProduced;
139  for(auto varToStore : valMapsToStore_) valMapsProduced.push_back(EGEnergySysIndex::name(varToStore));
140  desc.add<std::vector<std::string> >("valueMapsStored",valMapsProduced)->setComment("provides to python configs the list of valuemaps stored, can not be overriden in the python config");
141  descriptions.addWithDefaultLabel(desc);
142 }
143 
144 template<typename T>
145 void
147 {
148 
150 
151  edm::Handle<edm::View<T>> inHandle;
152  iEvent.getByToken(electronToken_, inHandle);
153 
154  edm::Handle<EcalRecHitCollection> recHitCollectionEBHandle;
155  edm::Handle<EcalRecHitCollection> recHitCollectionEEHandle;
156 
157  iEvent.getByToken(recHitCollectionEBToken_, recHitCollectionEBHandle);
158  iEvent.getByToken(recHitCollectionEEToken_, recHitCollectionEEHandle);
159 
160  std::unique_ptr<std::vector<T>> out = std::make_unique<std::vector<T>>();
161 
162  size_t nrObj = inHandle->size();
163  std::array<std::vector<float>,EGEnergySysIndex::kNrSysErrs> results;
164  for(auto& res : results) res.reserve(nrObj);
165 
166  const ElectronEnergyCalibrator::EventType evtType = iEvent.isRealData() ?
168 
169 
170  for (const auto& ele : *inHandle) {
171  out->push_back(ele);
172 
173  if(semiDeterministicRng_) setSemiDetRandomSeed(iEvent,ele,nrObj,out->size());
174 
175  const EcalRecHitCollection* recHits = (ele.isEB()) ? recHitCollectionEBHandle.product() : recHitCollectionEEHandle.product();
176  std::array<float,EGEnergySysIndex::kNrSysErrs> uncertainties = energyCorrector_.calibrate(out->back(), iEvent.id().run(), recHits, iEvent.streamID(), evtType);
177 
179  results[index].push_back(uncertainties[index]);
180  }
181  }
182 
183  auto fillAndStore = [&](auto handle){
184  for(const auto& mapToStore : valMapsToStore_){
185  fillAndStoreValueMap(iEvent,handle,results[mapToStore],EGEnergySysIndex::name(mapToStore));
186  }
187  };
188 
190  fillAndStore(iEvent.put(std::move(out)));
191  }else{
192  fillAndStore(inHandle);
193  }
194 
195 }
196 
197 template<typename T>
198 void CalibratedElectronProducerT<T>::setSemiDetRandomSeed(const edm::Event& iEvent,const T& obj,size_t nrObjs,size_t objNr)
199 {
200  if(obj.superCluster().isNonnull()){
201  semiDeterministicRng_->SetSeed(egamma::getRandomSeedFromSC(iEvent,obj.superCluster()));
202  }else{
203  semiDeterministicRng_->SetSeed(egamma::getRandomSeedFromObj(iEvent,obj,nrObjs,objNr));
204  }
205 }
206 
209 
211 
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
CalibratedElectronProducerT(const edm::ParameterSet &)
void addWithDefaultLabel(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:137
size_t size() const
Definition: Event.cc:278
static edm::ParameterSetDescription makePSetDescription()
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:579
static const std::string & name(size_t index)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
ElectronEnergyCalibrator energyCorrector_
void calibrate(SimpleElectron &electron, edm::StreamID const &)
bool isRealData() const
Definition: EventBase.h:64
Definition: Electron.h:6
edm::EDGetTokenT< edm::View< T > > electronToken_
edm::EDGetTokenT< EcalRecHitCollection > recHitCollectionEBToken_
int iEvent
Definition: GenABIO.cc:230
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
static constexpr size_t kNrSysErrs
uint32_t getRandomSeedFromObj(const edm::Event &iEvent, const T &obj, size_t nrObjs, size_t objNr)
void produce(edm::Event &, const edm::EventSetup &) override
ParameterDescriptionBase * add(U const &iLabel, T const &value)
void setSemiDetRandomSeed(const edm::Event &iEvent, const T &obj, size_t nrObjs, size_t objNr)
void setEventContent(const edm::EventSetup &iSetup)
T const * product() const
Definition: Handle.h:81
std::unique_ptr< TRandom > semiDeterministicRng_
edm::EDGetTokenT< EcalRecHitCollection > recHitCollectionEEToken_
edm::EventID id() const
Definition: EventBase.h:60
uint32_t getRandomSeedFromSC(const edm::Event &iEvent, const reco::SuperClusterRef scRef)
HLT enums.
static const std::vector< int > valMapsToStore_
StreamID streamID() const
Definition: Event.h:96
long double T
def move(src, dest)
Definition: eostools.py:510