CMS 3D CMS Logo

List of all members | Public Member Functions | Private Attributes
EgammaHLTHcalIsolationDoubleCone Class Reference

#include <RecoEgamma/EgammaHLTAlgos/interface/EgammaHLTHcalIsolationDoubleCone.h>

Public Member Functions

 EgammaHLTHcalIsolationDoubleCone (double egHcalIso_PtMin, double egHcalIso_ConeSize, double egHcalIso_Exclusion)
 
float getConeSize () const
 Get isolation cone size. More...
 
float getExclusion () const
 Get exclusion region. More...
 
float getptMin () const
 Get pt cut for hcal hits. More...
 
float isolPtSum (const reco::RecoCandidate *recocandidate, const HBHERecHitCollection *hbhe, const HFRecHitCollection *hf, const CaloGeometry *geometry) const
 

Private Attributes

const float conesize
 
const float exclusion
 
const float ptMin
 

Detailed Description

Description: sum pt hcal hits in cone around egamma candidate but exlude central cone mostly identical to EgammaHLTHcalIsolation, but with an inner exclusion cone

Usage: <usage>

Definition at line 37 of file EgammaHLTHcalIsolationDoubleCone.h.

Constructor & Destructor Documentation

EgammaHLTHcalIsolationDoubleCone::EgammaHLTHcalIsolationDoubleCone ( double  egHcalIso_PtMin,
double  egHcalIso_ConeSize,
double  egHcalIso_Exclusion 
)
inline

Definition at line 42 of file EgammaHLTHcalIsolationDoubleCone.h.

References photonIsolationHIProducer_cfi::hbhe, photonIsolationHIProducer_cfi::hf, and isolPtSum().

42  :
43  ptMin(egHcalIso_PtMin),conesize(egHcalIso_ConeSize),exclusion(egHcalIso_Exclusion){
44  /*
45  std::cout << "EgammaHLTHcalIsolation instance:"
46  << " ptMin=" << ptMin << "|" << ptMinG
47  << " conesize="<< conesize << "|" << conesizeG
48  << std::endl;
49  */
50  }

Member Function Documentation

float EgammaHLTHcalIsolationDoubleCone::getConeSize ( ) const
inline

Get isolation cone size.

Definition at line 59 of file EgammaHLTHcalIsolationDoubleCone.h.

References conesize.

float EgammaHLTHcalIsolationDoubleCone::getExclusion ( ) const
inline

Get exclusion region.

Definition at line 61 of file EgammaHLTHcalIsolationDoubleCone.h.

References exclusion.

float EgammaHLTHcalIsolationDoubleCone::getptMin ( ) const
inline

Get pt cut for hcal hits.

Definition at line 57 of file EgammaHLTHcalIsolationDoubleCone.h.

References ptMin.

float EgammaHLTHcalIsolationDoubleCone::isolPtSum ( const reco::RecoCandidate recocandidate,
const HBHERecHitCollection hbhe,
const HFRecHitCollection hf,
const CaloGeometry geometry 
) const

Definition at line 24 of file EgammaHLTHcalIsolationDoubleCone.cc.

References edm::SortedCollection< T, SORT >::begin(), conesize, edm::SortedCollection< T, SORT >::end(), PVValHelper::eta, exclusion, JetChargeProducer_cfi::exp, CaloGeometry::getPosition(), phi, PI, ptMin, funct::sin(), reco::RecoCandidate::superCluster(), and TWOPI.

Referenced by EgammaHLTHcalIsolationDoubleCone(), and EgammaHLTHcalIsolationDoubleConeProducers::produce().

24  {
25 
26  float hcalIsol=0.;
27 
28  float candSCphi = recocandidate->superCluster()->phi();
29  float candSCeta = recocandidate->superCluster()->eta();
30  if(candSCphi<0) candSCphi+=TWOPI;
31  float conesizeSquared=conesize*conesize;
32  float exclusionSquared= exclusion*exclusion;
33 
34  for(HBHERecHitCollection::const_iterator hbheItr = hbhe->begin(); hbheItr != hbhe->end(); ++hbheItr){
35  double HcalHit_eta=geometry->getPosition(hbheItr->id()).eta(); //Attention getpos
36  if(fabs(HcalHit_eta-candSCeta)<conesize) {
37  float HcalHit_pth=hbheItr->energy()*sin(2*atan(exp(-HcalHit_eta)));
38  if(HcalHit_pth>ptMin) {
39  double HcalHit_phi=geometry->getPosition(hbheItr->id()).phi();
40  float deltaeta=fabs(HcalHit_eta-candSCeta);
41  if(HcalHit_phi<0) HcalHit_phi+=TWOPI;
42  float deltaphi=fabs(HcalHit_phi-candSCphi);
43  if(deltaphi>TWOPI) deltaphi-=TWOPI;
44  if(deltaphi>PI) deltaphi=TWOPI-deltaphi;
45  float newDelta= (deltaphi*deltaphi+ deltaeta*deltaeta);
46  if(newDelta<conesizeSquared && newDelta>exclusionSquared ) hcalIsol+=HcalHit_pth;
47  }
48  }
49  }
50 
51  for(HFRecHitCollection::const_iterator hfItr = hf->begin(); hfItr != hf->end(); ++hfItr){
52  double HcalHit_eta=geometry->getPosition(hfItr->id()).eta(); //Attention getpos
53  if(fabs(HcalHit_eta-candSCeta)<conesize) {
54  float HcalHit_pth=hfItr->energy()*sin(2*atan(exp(-HcalHit_eta)));
55  if(HcalHit_pth>ptMin) {
56  double HcalHit_phi=geometry->getPosition(hfItr->id()).phi();
57  float deltaeta=fabs(HcalHit_eta-candSCeta);
58  float deltaphi;
59  if(HcalHit_phi<0) HcalHit_phi+=TWOPI;
60  if(candSCphi<0) candSCphi+=TWOPI;
61  deltaphi=fabs(HcalHit_phi-candSCphi);
62  if(deltaphi>TWOPI) deltaphi-=TWOPI;
63  if(deltaphi>PI) deltaphi=TWOPI-deltaphi;
64  float newDelta= (deltaphi*deltaphi+ deltaeta*deltaeta);
65  if(newDelta<conesizeSquared && newDelta>exclusionSquared ) hcalIsol+=HcalHit_pth;
66  }
67  }
68  }
69 
70 
71  return hcalIsol;
72 
73 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
std::vector< HBHERecHit >::const_iterator const_iterator
GlobalPoint getPosition(const DetId &id) const
Get the position of a given detector id.
Definition: CaloGeometry.cc:74
const_iterator end() const
const_iterator begin() const
virtual reco::SuperClusterRef superCluster() const
reference to a SuperCluster

Member Data Documentation

const float EgammaHLTHcalIsolationDoubleCone::conesize
private

Definition at line 68 of file EgammaHLTHcalIsolationDoubleCone.h.

Referenced by getConeSize(), and isolPtSum().

const float EgammaHLTHcalIsolationDoubleCone::exclusion
private

Definition at line 69 of file EgammaHLTHcalIsolationDoubleCone.h.

Referenced by getExclusion(), and isolPtSum().

const float EgammaHLTHcalIsolationDoubleCone::ptMin
private

Definition at line 67 of file EgammaHLTHcalIsolationDoubleCone.h.

Referenced by getptMin(), and isolPtSum().