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 
27 
28 #include "TRandom2.h"
29 
30 #include <vector>
31 #include <random>
32 
33 template<typename T>
35 public:
36  explicit CalibratedPhotonProducerT( const edm::ParameterSet & ) ;
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 
55 template<typename T>
56 const std::vector<int> CalibratedPhotonProducerT<T>::valMapsToStore_ = {
78 };
79 
80 namespace{
81  template<typename HandleType,typename ValType>
82  void fillAndStoreValueMap(edm::Event& iEvent,HandleType objHandle,
83  const std::vector<ValType>& vals,const std::string& name)
84  {
85  auto valMap = std::make_unique<edm::ValueMap<ValType> >();
86  typename edm::ValueMap<ValType>::Filler filler(*valMap);
87  filler.insert(objHandle,vals.begin(),vals.end());
88  filler.fill();
89  iEvent.put(std::move(valMap),name);
90  }
91 }
92 
93 template<typename T>
95  photonToken_(consumes<edm::View<T> >(conf.getParameter<edm::InputTag>("src"))),
96  energyCorrector_(conf.getParameter<std::string >("correctionFile")),
97  recHitCollectionEBToken_(consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitCollectionEB"))),
98  recHitCollectionEEToken_(consumes<EcalRecHitCollection>(conf.getParameter<edm::InputTag>("recHitCollectionEE"))),
99  produceCalibratedObjs_(conf.getParameter<bool>("produceCalibratedObjs"))
100 {
101 
102  energyCorrector_.setMinEt(conf.getParameter<double>("minEtToCalibrate"));
103 
104  if (conf.getParameter<bool>("semiDeterministic")) {
105  semiDeterministicRng_.reset(new TRandom2());
107  }
108 
109  if(produceCalibratedObjs_) produces<std::vector<T>>();
110 
111  for(const auto& toStore : valMapsToStore_){
112  produces<edm::ValueMap<float>>(EGEnergySysIndex::name(toStore));
113  }
114 }
115 
116 template<typename T>
118 {
120  desc.add<edm::InputTag>("src",edm::InputTag("gedGsfElectrons"));
121  desc.add<edm::InputTag>("recHitCollectionEB",edm::InputTag("reducedEcalRecHitsEB"));
122  desc.add<edm::InputTag>("recHitCollectionEE",edm::InputTag("reducedEcalRecHitsEE"));
123  desc.add<std::string>("correctionFile",std::string());
124  desc.add<double>("minEtToCalibrate",5.0);
125  desc.add<bool>("produceCalibratedObjs",true);
126  desc.add<bool>("semiDeterministic",true);
127  std::vector<std::string> valMapsProduced;
128  for(auto varToStore : valMapsToStore_) valMapsProduced.push_back(EGEnergySysIndex::name(varToStore));
129  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");
131 }
132 
133 template<typename T>
134 void
136 
137  edm::Handle<edm::View<T>> inHandle;
138  iEvent.getByToken(photonToken_, inHandle);
139 
140  edm::Handle<EcalRecHitCollection> recHitCollectionEBHandle;
141  edm::Handle<EcalRecHitCollection> recHitCollectionEEHandle;
142 
143  iEvent.getByToken(recHitCollectionEBToken_, recHitCollectionEBHandle);
144  iEvent.getByToken(recHitCollectionEEToken_, recHitCollectionEEHandle);
145 
146  std::unique_ptr<std::vector<T>> out = std::make_unique<std::vector<T>>();
147 
148  size_t nrObj = inHandle->size();
149  std::array<std::vector<float>,EGEnergySysIndex::kNrSysErrs> results;
150  for(auto& res : results) res.reserve(nrObj);
151 
152  const PhotonEnergyCalibrator::EventType evtType = iEvent.isRealData() ?
154 
155  for (const auto& pho : *inHandle) {
156  out->emplace_back(pho);
157 
158  if(semiDeterministicRng_) setSemiDetRandomSeed(iEvent,pho,nrObj,out->size());
159 
160  const EcalRecHitCollection* recHits = (pho.isEB()) ? recHitCollectionEBHandle.product() : recHitCollectionEEHandle.product();
161  std::array<float,EGEnergySysIndex::kNrSysErrs> uncertainties = energyCorrector_.calibrate(out->back(), iEvent.id().run(), recHits, iEvent.streamID(), evtType);
162 
164  results[index].push_back(uncertainties[index]);
165  }
166  }
167 
168  auto fillAndStore = [&](auto handle){
169  for(const auto& mapToStore : valMapsToStore_){
170  fillAndStoreValueMap(iEvent,handle,results[mapToStore],EGEnergySysIndex::name(mapToStore));
171  }
172  };
173 
175  fillAndStore(iEvent.put(std::move(out)));
176  }else{
177  fillAndStore(inHandle);
178  }
179 
180 }
181 
182 //needs to be synced to CalibratedElectronProducers, want the same seed for a given SC
183 template<typename T>
184 void CalibratedPhotonProducerT<T>::setSemiDetRandomSeed(const edm::Event& iEvent,const T& obj,size_t nrObjs,size_t objNr)
185 {
186  if(obj.superCluster().isNonnull()){
187  semiDeterministicRng_->SetSeed(egamma::getRandomSeedFromSC(iEvent,obj.superCluster()));
188  }else{
189  semiDeterministicRng_->SetSeed(egamma::getRandomSeedFromObj(iEvent,obj,nrObjs,objNr));
190  }
191 }
192 
195 
197 
RunNumber_t run() const
Definition: EventID.h:39
T getParameter(std::string const &) const
std::string defaultModuleLabel()
OrphanHandle< PROD > put(std::unique_ptr< PROD > product)
Put a new product.
Definition: Event.h:127
size_t size() const
Definition: Event.cc:249
std::unique_ptr< TRandom > semiDeterministicRng_
void initPrivateRng(TRandom *rnd)
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:508
static const std::string & name(size_t index)
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void produce(edm::Event &, const edm::EventSetup &) override
bool isRealData() const
Definition: EventBase.h:64
Definition: Electron.h:4
EDGetTokenT< ProductType > consumes(edm::InputTag const &tag)
int iEvent
Definition: GenABIO.cc:230
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_
ParameterDescriptionBase * add(U const &iLabel, T const &value)
T const * product() const
Definition: Handle.h:81
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
void add(std::string const &label, ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< edm::View< T > > photonToken_
edm::EventID id() const
Definition: EventBase.h:60
uint32_t getRandomSeedFromSC(const edm::Event &iEvent, const reco::SuperClusterRef scRef)
HLT enums.
StreamID streamID() const
Definition: Event.h:86
static const std::vector< int > valMapsToStore_
long double T
def move(src, dest)
Definition: eostools.py:510
std::array< float, EGEnergySysIndex::kNrSysErrs > calibrate(reco::Photon &photon, const unsigned int runNumber, const EcalRecHitCollection *recHits, edm::StreamID const &id, const EventType eventType) const