CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
GsfElectronGSCrysFixer.cc
Go to the documentation of this file.
1 #ifndef RecoEgamma_EgammaElectronProducers_GsfElectronGSCrysFixer_h
2 #define RecoEgamma_EgammaElectronProducers_GsfElectronGSCrysFixer_h
3 
7 
11 
12 
17 
27 
32 
34 
36 
38 
39 #include <iostream>
40 #include <string>
41 
43 public:
44  explicit GsfElectronGSCrysFixer(const edm::ParameterSet& );
46 
47  void produce(edm::Event&, const edm::EventSetup& ) override;
49  edm::EventSetup const&) override;
50 
51 
52  template<typename T>
54  auto tag(pset.getParameter<edm::InputTag>(label));
55  if (!instance.empty())
56  tag = edm::InputTag(tag.label(), instance, tag.process());
57 
58  token = consumes<T>(tag);
59  }
60 private:
63 
68 
71  std::unique_ptr<ModifyObjectValueBase> gedRegression_;
72 };
73 
74 
75 namespace {
76  template<typename T> edm::Handle<T> getHandle(const edm::Event& iEvent,const edm::EDGetTokenT<T>& token){
78  iEvent.getByToken(token,handle);
79  return handle;
80  }
81 }
82 
83 
84 
86 
88 {
89  getToken(oldGsfElesToken_,pset,"oldEles");
90  getToken(newCoresToken_,pset,"newCores");
91  getToken(newCoresToOldSCMapToken_, pset, "newCores");
92  getToken(ebRecHitsToken_,pset,"ebRecHits");
93 
94  //ripped wholesale from GEDGsfElectronFinalizer
95  if( pset.existsAs<edm::ParameterSet>("regressionConfig") ) {
96  const edm::ParameterSet& iconf = pset.getParameterSet("regressionConfig");
97  const std::string& mname = iconf.getParameter<std::string>("modifierName");
99  ModifyObjectValueFactory::get()->create(mname,iconf);
100  gedRegression_.reset(plugin);
102  gedRegression_->setConsumes(sumes);
103  } else {
104  gedRegression_.reset(nullptr);
105  }
106 
107  produces<reco::GsfElectronCollection >();
108  // new -> old
109  produces<ElectronRefMap>();
110 }
111 
112 namespace {
113 
114  reco::GsfElectronCoreRef getNewCore(reco::GsfElectron const& oldEle,
116  edm::ValueMap<reco::SuperClusterRef> const& newToOldCoresMap)
117  {
118  for(size_t coreNr=0;coreNr<newCores->size();coreNr++){
119  reco::GsfElectronCoreRef coreRef(newCores,coreNr);
120  auto oldRef = newToOldCoresMap[coreRef];
121  if( oldRef.isNonnull() && oldRef==oldEle.superCluster()){
122  return coreRef;
123  }
124  }
125  return reco::GsfElectronCoreRef(newCores.id());
126  }
127 }
128 
130 {
131  auto outEles = std::make_unique<reco::GsfElectronCollection>();
132 
133  if( gedRegression_ ) {
134  gedRegression_->setEvent(iEvent);
135  gedRegression_->setEventContent(iSetup);
136  }
137 
138  auto elesHandle = getHandle(iEvent,oldGsfElesToken_);
139  auto& ebRecHits = *getHandle(iEvent,ebRecHitsToken_);
140  auto newCoresHandle = getHandle(iEvent,newCoresToken_);
141  auto newCoresToOldSCMap = *getHandle(iEvent, newCoresToOldSCMapToken_);
142 
143  std::vector<reco::GsfElectronRef> oldElectrons;
144 
145  unsigned iE(0);
146  for (auto& oldEle : *elesHandle) {
147  auto&& refToBase(elesHandle->refAt(iE++));
148  oldElectrons.emplace_back(refToBase.id(), &oldEle, refToBase.key());
149 
150  reco::GsfElectronCoreRef newCoreRef = getNewCore(oldEle,newCoresHandle,newCoresToOldSCMap);
151 
152  outEles->emplace_back(oldEle,newCoreRef);
153  auto& newEle(outEles->back());
154 
155  if (GainSwitchTools::hasEBGainSwitchIn5x5(*oldEle.superCluster(), &ebRecHits, topology_)) {
156 
157  reco::GsfElectron::ShowerShape full5x5ShowerShape = GainSwitchTools::redoEcalShowerShape<true>(newEle.full5x5_showerShape(),newEle.superCluster(),&ebRecHits,topology_,geometry_);
158  reco::GsfElectron::ShowerShape showerShape = GainSwitchTools::redoEcalShowerShape<false>(newEle.showerShape(),newEle.superCluster(),&ebRecHits,topology_,geometry_);
159 
160  float eNewSCOverEOldSC = newEle.superCluster()->energy()/oldEle.superCluster()->energy();
163 
164  newEle.full5x5_setShowerShape(full5x5ShowerShape);
165  newEle.setShowerShape(showerShape);
166 
167 
168  if( gedRegression_ )
169  gedRegression_->modifyObject(newEle);
170  }
171  }
172 
173  auto&& newElectronsHandle(iEvent.put(std::move(outEles)));
174 
175  std::auto_ptr<ElectronRefMap> pRefMap(new ElectronRefMap);
176  ElectronRefMap::Filler refMapFiller(*pRefMap);
177  refMapFiller.insert(newElectronsHandle, oldElectrons.begin(), oldElectrons.end());
178  refMapFiller.fill();
179  iEvent.put(pRefMap);
180 }
181 
183  edm::EventSetup const& es) {
184  edm::ESHandle<CaloGeometry> caloGeom ;
185  edm::ESHandle<CaloTopology> caloTopo ;
186  es.get<CaloGeometryRecord>().get(caloGeom);
187  es.get<CaloTopologyRecord>().get(caloTopo);
188  geometry_ = caloGeom.product();
189  topology_ = caloTopo.product();
190 }
191 
192 
193 
195 #endif
196 
T getParameter(std::string const &) const
edm::Ref< GsfElectronCoreCollection > GsfElectronCoreRef
bool existsAs(std::string const &parameterName, bool trackiness=true) const
checks if a parameter exists as a given type
Definition: ParameterSet.h:186
auto_ptr< JetDefinition::Plugin > plugin
ProductID id() const
Definition: HandleBase.cc:15
static PFTauRenderPlugin instance
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
edm::EDGetTokenT< SCRefMap > newCoresToOldSCMapToken_
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
void getToken(edm::EDGetTokenT< T > &token, const edm::ParameterSet &pset, const std::string &label, const std::string &instance="")
void insert(const H &h, I begin, I end)
Definition: ValueMap.h:52
edm::View< reco::GsfElectron > GsfElectronView
const CaloTopology * topology_
std::unique_ptr< ModifyObjectValueBase > gedRegression_
const CaloGeometry * geometry_
int iEvent
Definition: GenABIO.cc:230
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:121
ConsumesCollector consumesCollector()
Use a ConsumesCollector to gather consumes information from helper functions.
virtual SuperClusterRef superCluster() const
reference to a SuperCluster
Definition: GsfElectron.h:182
def move
Definition: eostools.py:510
tuple handle
Definition: patZpeak.py:22
edm::EDGetTokenT< GsfElectronView > oldGsfElesToken_
GsfElectronGSCrysFixer(const edm::ParameterSet &)
static void correctHadem(reco::GsfElectron::ShowerShape &showerShape, float eNewOverEOld, const GainSwitchTools::ShowerShapeType ssType)
void beginLuminosityBlock(edm::LuminosityBlock const &, edm::EventSetup const &) override
edm::ValueMap< reco::SuperClusterRef > SCRefMap
edm::ValueMap< reco::GsfElectronRef > ElectronRefMap
edm::EDGetTokenT< EcalRecHitCollection > ebRecHitsToken_
ParameterSet const & getParameterSet(std::string const &) const
void produce(edm::Event &, const edm::EventSetup &) override
const T & get() const
Definition: EventSetup.h:56
T const * product() const
Definition: ESHandle.h:86
edm::EDGetTokenT< reco::GsfElectronCoreCollection > newCoresToken_
static bool hasEBGainSwitchIn5x5(const reco::SuperCluster &superClus, const EcalRecHitCollection *recHits, const CaloTopology *topology)
T get(const Candidate &c)
Definition: component.h:55