CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
ElectronRegressionValueMapProducer.cc
Go to the documentation of this file.
3 
6 
8 
11 
14 
18 
19 #include "TVector2.h"
20 
21 #include <memory>
22 #include <vector>
23 #include <unordered_map>
24 
25 namespace {
26  enum reg_float_vars { k_sigmaIEtaIPhi = 0,
27  k_eMax,
28  k_e2nd,
29  k_eTop,
30  k_eBottom,
31  k_eLeft,
32  k_eRight,
33  k_clusterMaxDR,
34  k_clusterMaxDRDPhi,
35  k_clusterMaxDRDEta,
36  k_clusterMaxDRRawEnergy,
37  k_clusterRawEnergy0,
38  k_clusterRawEnergy1,
39  k_clusterRawEnergy2,
40  k_clusterDPhiToSeed0,
41  k_clusterDPhiToSeed1,
42  k_clusterDPhiToSeed2,
43  k_clusterDEtaToSeed0,
44  k_clusterDEtaToSeed1,
45  k_clusterDEtaToSeed2,
46  k_cryPhi,
47  k_cryEta,
48  k_NFloatVars };
49 
50  enum reg_int_vars { k_iPhi = 0,
51  k_iEta,
52  k_NIntVars };
53 
54  static const std::vector<std::string> float_var_names( { "sigmaIEtaIPhi",
55  "eMax",
56  "e2nd",
57  "eTop",
58  "eBottom",
59  "eLeft",
60  "eRight",
61  "clusterMaxDR",
62  "clusterMaxDRDPhi",
63  "clusterMaxDRDEta",
64  "clusterMaxDRRawEnergy",
65  "clusterRawEnergy0",
66  "clusterRawEnergy1",
67  "clusterRawEnergy2",
68  "clusterDPhiToSeed0",
69  "clusterDPhiToSeed1",
70  "clusterDPhiToSeed2",
71  "clusterDEtaToSeed0",
72  "clusterDEtaToSeed1",
73  "clusterDEtaToSeed2",
74  "cryPhi",
75  "cryEta" } );
76 
77  static const std::vector<std::string> integer_var_names( { "iPhi", "iEta" } );
78 
79  inline void set_map_val( const reg_float_vars index, const float value,
80  std::unordered_map<std::string,float>& map) {
81  map[float_var_names[index]] = value;
82  }
83  inline void set_map_val( const reg_int_vars index, const int value,
84  std::unordered_map<std::string,int>& map) {
85  map[integer_var_names[index]] = value;
86  }
87 
88  template<typename T>
89  inline void check_map(const std::unordered_map<std::string,T>& map, unsigned exp_size) {
90  if( map.size() != exp_size ) {
91  throw cms::Exception("ElectronRegressionWeirdConfig")
92  << "variable map size: " << map.size()
93  << " not equal to expected size: " << exp_size << " !"
94  << " The regression variable calculation code definitely has a bug, fix it!";
95  }
96  }
97 
98  template<typename LazyTools>
99  void calculateValues(EcalClusterLazyToolsBase* tools_tocast,
100  const edm::Ptr<reco::GsfElectron>& iEle,
101  const edm::EventSetup& iSetup,
102  std::unordered_map<std::string,float>& float_vars,
103  std::unordered_map<std::string,int>& int_vars ) {
104  LazyTools* tools = static_cast<LazyTools*>(tools_tocast);
105 
106  const auto& the_sc = iEle->superCluster();
107  const auto& theseed = the_sc->seed();
108 
109  const int numberOfClusters = the_sc->clusters().size();
110  const bool missing_clusters = !the_sc->clusters()[numberOfClusters-1].isAvailable();
111 
112  std::vector<float> vCov = tools->localCovariances( *theseed );
113 
114  const float eMax = tools->eMax( *theseed );
115  const float e2nd = tools->e2nd( *theseed );
116  const float eTop = tools->eTop( *theseed );
117  const float eLeft = tools->eLeft( *theseed );
118  const float eRight = tools->eRight( *theseed );
119  const float eBottom = tools->eBottom( *theseed );
120 
121  float dummy;
122  int iPhi;
123  int iEta;
124  float cryPhi;
125  float cryEta;
126  EcalClusterLocal _ecalLocal;
127  if (iEle->isEB())
128  _ecalLocal.localCoordsEB(*theseed, iSetup, cryEta, cryPhi, iEta, iPhi, dummy, dummy);
129  else
130  _ecalLocal.localCoordsEE(*theseed, iSetup, cryEta, cryPhi, iEta, iPhi, dummy, dummy);
131 
132  double see = (isnan(vCov[0]) ? 0. : sqrt(vCov[0]));
133  double spp = (isnan(vCov[2]) ? 0. : sqrt(vCov[2]));
134  double sep;
135  if (see*spp > 0)
136  sep = vCov[1] / (see * spp);
137  else if (vCov[1] > 0)
138  sep = 1.0;
139  else
140  sep = -1.0;
141 
142  set_map_val(k_sigmaIEtaIPhi,sep,float_vars);
143  set_map_val(k_eMax,eMax,float_vars);
144  set_map_val(k_e2nd,e2nd,float_vars);
145  set_map_val(k_eTop,eTop,float_vars);
146  set_map_val(k_eBottom,eBottom,float_vars);
147  set_map_val(k_eLeft,eLeft,float_vars);
148  set_map_val(k_eRight,eRight,float_vars);
149  set_map_val(k_cryPhi,cryPhi,float_vars);
150  set_map_val(k_cryEta,cryEta,float_vars);
151 
152  set_map_val(k_iPhi,iPhi,int_vars);
153  set_map_val(k_iEta,iEta,int_vars);
154 
155  std::vector<float> _clusterRawEnergy;
156  _clusterRawEnergy.resize(std::max(3, numberOfClusters), 0);
157  std::vector<float> _clusterDEtaToSeed;
158  _clusterDEtaToSeed.resize(std::max(3, numberOfClusters), 0);
159  std::vector<float> _clusterDPhiToSeed;
160  _clusterDPhiToSeed.resize(std::max(3, numberOfClusters), 0);
161  float _clusterMaxDR = 999.;
162  float _clusterMaxDRDPhi = 999.;
163  float _clusterMaxDRDEta = 999.;
164  float _clusterMaxDRRawEnergy = 0.;
165 
166  size_t iclus = 0;
167  float maxDR = 0;
169  if( !missing_clusters ) {
170  // loop over all clusters that aren't the seed
171  auto clusend = the_sc->clustersEnd();
172  for( auto clus = the_sc->clustersBegin(); clus != clusend; ++clus ) {
173  pclus = *clus;
174 
175  if(theseed == pclus )
176  continue;
177  _clusterRawEnergy[iclus] = pclus->energy();
178  _clusterDPhiToSeed[iclus] = reco::deltaPhi(pclus->phi(),theseed->phi());
179  _clusterDEtaToSeed[iclus] = pclus->eta() - theseed->eta();
180 
181  // find cluster with max dR
182  if(reco::deltaR(*pclus, *theseed) > maxDR) {
183  maxDR = reco::deltaR(*pclus, *theseed);
184  _clusterMaxDR = maxDR;
185  _clusterMaxDRDPhi = _clusterDPhiToSeed[iclus];
186  _clusterMaxDRDEta = _clusterDEtaToSeed[iclus];
187  _clusterMaxDRRawEnergy = _clusterRawEnergy[iclus];
188  }
189  ++iclus;
190  }
191  }
192 
193  set_map_val(k_clusterMaxDR,_clusterMaxDR,float_vars);
194  set_map_val(k_clusterMaxDRDPhi,_clusterMaxDRDPhi,float_vars);
195  set_map_val(k_clusterMaxDRDEta,_clusterMaxDRDEta,float_vars);
196  set_map_val(k_clusterMaxDRRawEnergy,_clusterMaxDRRawEnergy,float_vars);
197  set_map_val(k_clusterRawEnergy0,_clusterRawEnergy[0],float_vars);
198  set_map_val(k_clusterRawEnergy1,_clusterRawEnergy[1],float_vars);
199  set_map_val(k_clusterRawEnergy2,_clusterRawEnergy[2],float_vars);
200  set_map_val(k_clusterDPhiToSeed0,_clusterDPhiToSeed[0],float_vars);
201  set_map_val(k_clusterDPhiToSeed1,_clusterDPhiToSeed[1],float_vars);
202  set_map_val(k_clusterDPhiToSeed2,_clusterDPhiToSeed[2],float_vars);
203  set_map_val(k_clusterDEtaToSeed0,_clusterDEtaToSeed[0],float_vars);
204  set_map_val(k_clusterDEtaToSeed1,_clusterDEtaToSeed[1],float_vars);
205  set_map_val(k_clusterDEtaToSeed2,_clusterDEtaToSeed[2],float_vars);
206  }
207 }
208 
210 
211  public:
212 
215 
216  static void fillDescriptions(edm::ConfigurationDescriptions& descriptions);
217 
218  private:
219 
220  virtual void produce(edm::Event&, const edm::EventSetup&) override;
221 
222  template<typename T>
225  const std::vector<T> & values,
226  const std::string & label) const ;
227 
228  std::unique_ptr<EcalClusterLazyToolsBase> lazyTools;
229 
230  // for AOD case
235 
236  // for miniAOD case
241 
242  const bool use_full5x5_;
243 };
244 
246  use_full5x5_(iConfig.getParameter<bool>("useFull5x5")) {
247 
248  //
249  // Declare consummables, handle both AOD and miniAOD case
250  //
251  ebReducedRecHitCollection_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
252  ("ebReducedRecHitCollection"));
253  ebReducedRecHitCollectionMiniAOD_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
254  ("ebReducedRecHitCollectionMiniAOD"));
255 
256  eeReducedRecHitCollection_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
257  ("eeReducedRecHitCollection"));
258  eeReducedRecHitCollectionMiniAOD_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
259  ("eeReducedRecHitCollectionMiniAOD"));
260 
261  esReducedRecHitCollection_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
262  ("esReducedRecHitCollection"));
263  esReducedRecHitCollectionMiniAOD_ = mayConsume<EcalRecHitCollection>(iConfig.getParameter<edm::InputTag>
264  ("esReducedRecHitCollectionMiniAOD"));
265 
266  src_ = mayConsume<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>("src"));
267  srcMiniAOD_ = mayConsume<edm::View<reco::GsfElectron> >(iConfig.getParameter<edm::InputTag>("srcMiniAOD"));
268 
269  for( const std::string& name : float_var_names ) {
270  produces<edm::ValueMap<float> >(name);
271  }
272 
273  for( const std::string& name : integer_var_names ) {
274  produces<edm::ValueMap<int> >(name);
275  }
276 }
277 
279 }
280 
282 
283  using namespace edm;
284 
286 
287  // Retrieve the collection of electrons from the event.
288  // If we fail to retrieve the collection with the standard AOD
289  // name, we next look for the one with the stndard miniAOD name.
290  bool isAOD = true;
291  iEvent.getByToken(src_, src);
292 
293  if( !src.isValid() ) {
294  isAOD = false;
295  iEvent.getByToken(srcMiniAOD_,src);
296  }
297 
298  // configure lazy tools
300 
301  if( isAOD ) {
305  } else {
309  }
310 
311  if( use_full5x5_ ) {
312  lazyTools.reset( new noZS::EcalClusterLazyTools(iEvent, iSetup,
313  ebrh, eerh, esrh ) );
314  } else {
315  lazyTools.reset( new EcalClusterLazyTools(iEvent, iSetup,
316  ebrh, eerh, esrh ) );
317  }
318 
319  std::vector<std::vector<float> > float_vars(k_NFloatVars);
320  std::vector<std::vector<int> > int_vars(k_NIntVars);
321 
322  std::unordered_map<std::string,float> float_vars_map;
323  std::unordered_map<std::string,int> int_vars_map;
324 
325  // reco::GsfElectron::superCluster() is virtual so we can exploit polymorphism
326  for (size_t i = 0; i < src->size(); ++i){
327  auto iEle = src->ptrAt(i);
328 
329  if( use_full5x5_ ) {
330  calculateValues<noZS::EcalClusterLazyTools>(lazyTools.get(),
331  iEle,
332  iSetup,
333  float_vars_map,
334  int_vars_map);
335  } else {
336  calculateValues<EcalClusterLazyTools>(lazyTools.get(),
337  iEle,
338  iSetup,
339  float_vars_map,
340  int_vars_map);
341  }
342 
343  check_map(float_vars_map, k_NFloatVars);
344  check_map(int_vars_map, k_NIntVars);
345 
346  for( unsigned i = 0; i < float_vars.size(); ++i ) {
347  float_vars[i].emplace_back(float_vars_map.at(float_var_names[i]));
348  }
349 
350  for( unsigned i = 0; i < int_vars.size(); ++i ) {
351  int_vars[i].emplace_back(int_vars_map.at(integer_var_names[i]));
352  }
353  }
354 
355  for( unsigned i = 0; i < float_vars.size(); ++i ) {
356  writeValueMap(iEvent, src, float_vars[i], float_var_names[i]);
357  }
358 
359  for( unsigned i = 0; i < int_vars.size(); ++i ) {
360  writeValueMap(iEvent, src, int_vars[i], integer_var_names[i]);
361  }
362 
363  lazyTools.reset();
364 }
365 
366 template<typename T>
369  const std::vector<T> & values,
370  const std::string & label) const
371 {
372  using namespace edm;
373  using namespace std;
374  typedef ValueMap<T> TValueMap;
375 
376  auto_ptr<TValueMap> valMap(new TValueMap());
377  typename TValueMap::Filler filler(*valMap);
378  filler.insert(handle, values.begin(), values.end());
379  filler.fill();
380  iEvent.put(valMap, label);
381 }
382 
384  //The following says we do not know what parameters are allowed so do no validation
385  // Please change this to state exactly what you do use, even if it is no parameters
387  desc.setUnknown();
388  descriptions.addDefault(desc);
389 }
390 
T getParameter(std::string const &) const
edm::EDGetTokenT< EcalRecHitCollection > eeReducedRecHitCollection_
int i
Definition: DBlmapReader.cc:9
void writeValueMap(edm::Event &iEvent, const edm::Handle< edm::View< reco::GsfElectron > > &handle, const std::vector< T > &values, const std::string &label) const
EcalClusterLazyToolsT<::EcalClusterTools > EcalClusterLazyTools
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
edm::EDGetTokenT< EcalRecHitCollection > esReducedRecHitCollection_
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:449
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
bool isAvailable() const
Definition: Ptr.h:249
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
edm::EDGetTokenT< EcalRecHitCollection > esReducedRecHitCollectionMiniAOD_
int iEvent
Definition: GenABIO.cc:230
void addDefault(ParameterSetDescription const &psetDescription)
OrphanHandle< PROD > put(std::auto_ptr< PROD > product)
Put a new product.
Definition: Event.h:113
bool isnan(float x)
Definition: math.h:13
edm::EDGetTokenT< EcalRecHitCollection > ebReducedRecHitCollection_
T sqrt(T t)
Definition: SSEVec.h:48
edm::EDGetTokenT< EcalRecHitCollection > eeReducedRecHitCollectionMiniAOD_
tuple handle
Definition: patZpeak.py:22
std::unique_ptr< EcalClusterLazyToolsBase > lazyTools
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
edm::EDGetTokenT< EcalRecHitCollection > ebReducedRecHitCollectionMiniAOD_
virtual void produce(edm::Event &, const edm::EventSetup &) override
void localCoordsEB(const reco::CaloCluster &bclus, const edm::EventSetup &es, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const
void localCoordsEE(const reco::CaloCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const