CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
Classes | Functions
ElectronRegressionValueMapProducer.cc File Reference
#include "FWCore/Framework/interface/Frameworkfwd.h"
#include "FWCore/Framework/interface/stream/EDProducer.h"
#include "FWCore/Framework/interface/Event.h"
#include "FWCore/Framework/interface/MakerMacros.h"
#include "FWCore/ParameterSet/interface/ParameterSet.h"
#include "DataFormats/Common/interface/ValueMap.h"
#include "DataFormats/Common/interface/View.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectron.h"
#include "DataFormats/EgammaCandidates/interface/GsfElectronFwd.h"
#include "DataFormats/PatCandidates/interface/Electron.h"
#include "RecoEcal/EgammaCoreTools/interface/EcalClusterLazyTools.h"
#include "RecoEgamma/EgammaTools/interface/EcalClusterLocal.h"
#include "TVector2.h"
#include <memory>
#include <vector>

Go to the source code of this file.

Classes

class  ElectronRegressionValueMapProducer
 

Functions

template<typename LazyTools >
void calculateValues (EcalClusterLazyToolsBase *tools_tocast, const edm::Ptr< reco::GsfElectron > &iEle, const edm::EventSetup &iSetup, std::vector< float > &vsigmaIEtaIPhi, std::vector< float > &veMax, std::vector< float > &ve2nd, std::vector< float > &veTop, std::vector< float > &veBottom, std::vector< float > &veLeft, std::vector< float > &veRight, std::vector< float > &vclusterMaxDR, std::vector< float > &vclusterMaxDRDPhi, std::vector< float > &vclusterMaxDRDEta, std::vector< float > &vclusterMaxDRRawEnergy, std::vector< float > &vclusterRawEnergy0, std::vector< float > &vclusterRawEnergy1, std::vector< float > &vclusterRawEnergy2, std::vector< float > &vclusterDPhiToSeed0, std::vector< float > &vclusterDPhiToSeed1, std::vector< float > &vclusterDPhiToSeed2, std::vector< float > &vclusterDEtaToSeed0, std::vector< float > &vclusterDEtaToSeed1, std::vector< float > &vclusterDEtaToSeed2, std::vector< int > &veleIPhi, std::vector< int > &veleIEta, std::vector< float > &veleCryPhi, std::vector< float > &veleCryEta)
 
static const
edm::ParameterSetDescriptionFillerPluginFactory::PMaker
< edm::ParameterSetDescriptionFiller
< ElectronRegressionValueMapProducer > > 
s_filler__LINE__ ("ElectronRegressionValueMapProducer")
 
static const
edm::MakerPluginFactory::PMaker
< edm::WorkerMaker
< ElectronRegressionValueMapProducer > > 
s_maker__LINE__ ("ElectronRegressionValueMapProducer")
 

Function Documentation

template<typename LazyTools >
void calculateValues ( EcalClusterLazyToolsBase tools_tocast,
const edm::Ptr< reco::GsfElectron > &  iEle,
const edm::EventSetup iSetup,
std::vector< float > &  vsigmaIEtaIPhi,
std::vector< float > &  veMax,
std::vector< float > &  ve2nd,
std::vector< float > &  veTop,
std::vector< float > &  veBottom,
std::vector< float > &  veLeft,
std::vector< float > &  veRight,
std::vector< float > &  vclusterMaxDR,
std::vector< float > &  vclusterMaxDRDPhi,
std::vector< float > &  vclusterMaxDRDEta,
std::vector< float > &  vclusterMaxDRRawEnergy,
std::vector< float > &  vclusterRawEnergy0,
std::vector< float > &  vclusterRawEnergy1,
std::vector< float > &  vclusterRawEnergy2,
std::vector< float > &  vclusterDPhiToSeed0,
std::vector< float > &  vclusterDPhiToSeed1,
std::vector< float > &  vclusterDPhiToSeed2,
std::vector< float > &  vclusterDEtaToSeed0,
std::vector< float > &  vclusterDEtaToSeed1,
std::vector< float > &  vclusterDEtaToSeed2,
std::vector< int > &  veleIPhi,
std::vector< int > &  veleIEta,
std::vector< float > &  veleCryPhi,
std::vector< float > &  veleCryEta 
)
inline

Definition at line 145 of file ElectronRegressionValueMapProducer.cc.

References reco::deltaPhi(), reco::deltaR(), edm::detail::isnan(), EcalClusterLocal::localCoordsEB(), EcalClusterLocal::localCoordsEE(), bookConverter::max, cmsHarvester::sep, and mathSSE::sqrt().

171  {
172  LazyTools* tools = static_cast<LazyTools*>(tools_tocast);
173 
174  const auto& the_sc = iEle->superCluster();
175  const auto& theseed = the_sc->seed();
176 
177  std::vector<float> vCov = tools->localCovariances( *theseed );
178 
179  const float eMax = tools->eMax( *theseed );
180  const float e2nd = tools->e2nd( *theseed );
181  const float eTop = tools->eTop( *theseed );
182  const float eLeft = tools->eLeft( *theseed );
183  const float eRight = tools->eRight( *theseed );
184  const float eBottom = tools->eBottom( *theseed );
185 
186  float dummy;
187  int iPhi;
188  int iEta;
189  float cryPhi;
190  float cryEta;
191  EcalClusterLocal _ecalLocal;
192  if (iEle->isEB())
193  _ecalLocal.localCoordsEB(*theseed, iSetup, cryEta, cryPhi, iEta, iPhi, dummy, dummy);
194  else
195  _ecalLocal.localCoordsEE(*theseed, iSetup, cryEta, cryPhi, iEta, iPhi, dummy, dummy);
196 
197  double see = (isnan(vCov[0]) ? 0. : sqrt(vCov[0]));
198  double spp = (isnan(vCov[2]) ? 0. : sqrt(vCov[2]));
199  double sep;
200  if (see*spp > 0)
201  sep = vCov[1] / (see * spp);
202  else if (vCov[1] > 0)
203  sep = 1.0;
204  else
205  sep = -1.0;
206 
207  vsigmaIEtaIPhi.push_back(sep);
208  veMax.push_back(eMax);
209  ve2nd.push_back(e2nd);
210  veTop.push_back(eTop);
211  veBottom.push_back(eBottom);
212  veLeft.push_back(eLeft);
213  veRight.push_back(eRight);
214  veleIPhi.push_back(iPhi);
215  veleIEta.push_back(iEta);
216  veleCryPhi.push_back(cryPhi);
217  veleCryEta.push_back(cryEta);
218 
219  // loop over all clusters that aren't the seed
220  auto clusend = the_sc->clustersEnd();
221  int numberOfClusters = the_sc->clusters().size();
222 
223  std::vector<float> _clusterRawEnergy;
224  _clusterRawEnergy.resize(std::max(3, numberOfClusters), 0);
225  std::vector<float> _clusterDEtaToSeed;
226  _clusterDEtaToSeed.resize(std::max(3, numberOfClusters), 0);
227  std::vector<float> _clusterDPhiToSeed;
228  _clusterDPhiToSeed.resize(std::max(3, numberOfClusters), 0);
229  float _clusterMaxDR = 999.;
230  float _clusterMaxDRDPhi = 999.;
231  float _clusterMaxDRDEta = 999.;
232  float _clusterMaxDRRawEnergy = 0.;
233 
234  size_t iclus = 0;
235  float maxDR = 0;
237  for( auto clus = the_sc->clustersBegin(); clus != clusend; ++clus ) {
238  pclus = *clus;
239 
240  if(theseed == pclus )
241  continue;
242  _clusterRawEnergy.push_back(pclus->energy());
243  _clusterDPhiToSeed.push_back(reco::deltaPhi(pclus->phi(),theseed->phi()));
244  _clusterDEtaToSeed.push_back(pclus->eta() - theseed->eta());
245 
246  // find cluster with max dR
247  if(reco::deltaR(*pclus, *theseed) > maxDR) {
248  maxDR = reco::deltaR(*pclus, *theseed);
249  _clusterMaxDR = maxDR;
250  _clusterMaxDRDPhi = _clusterDPhiToSeed[iclus];
251  _clusterMaxDRDEta = _clusterDEtaToSeed[iclus];
252  _clusterMaxDRRawEnergy = _clusterRawEnergy[iclus];
253  }
254 
255  ++iclus;
256  }
257 
258  vclusterMaxDR.push_back(_clusterMaxDR);
259  vclusterMaxDRDPhi.push_back(_clusterMaxDRDPhi);
260  vclusterMaxDRDEta.push_back(_clusterMaxDRDEta);
261  vclusterMaxDRRawEnergy.push_back(_clusterMaxDRRawEnergy);
262  vclusterRawEnergy0.push_back(_clusterRawEnergy[0]);
263  vclusterRawEnergy1.push_back(_clusterRawEnergy[1]);
264  vclusterRawEnergy2.push_back(_clusterRawEnergy[2]);
265  vclusterDPhiToSeed0.push_back(_clusterDPhiToSeed[0]);
266  vclusterDPhiToSeed1.push_back(_clusterDPhiToSeed[1]);
267  vclusterDPhiToSeed2.push_back(_clusterDPhiToSeed[2]);
268  vclusterDEtaToSeed0.push_back(_clusterDEtaToSeed[0]);
269  vclusterDEtaToSeed1.push_back(_clusterDEtaToSeed[1]);
270  vclusterDEtaToSeed2.push_back(_clusterDEtaToSeed[2]);
271 }
double deltaR(const T1 &t1, const T2 &t2)
Definition: deltaR.h:48
bool isnan(float x)
Definition: math.h:13
T sqrt(T t)
Definition: SSEVec.h:48
double deltaPhi(double phi1, double phi2)
Definition: deltaPhi.h:12
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
static const edm::ParameterSetDescriptionFillerPluginFactory::PMaker<edm::ParameterSetDescriptionFiller< ElectronRegressionValueMapProducer > > s_filler__LINE__ ( "ElectronRegressionValueMapProducer"  )
static
static const edm::MakerPluginFactory ::PMaker< edm::WorkerMaker< ElectronRegressionValueMapProducer > > s_maker__LINE__ ( "ElectronRegressionValueMapProducer"  )
static