CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HFRecoEcalCandidateAlgo Class Reference

#include <HFRecoEcalCandidateAlgo.h>

Public Member Functions

 HFRecoEcalCandidateAlgo (bool correct, double e9e25Cut, double intercept2DCut)
 
void produce (const edm::Handle< reco::SuperClusterCollection > &SuperClusters, const reco::HFEMClusterShapeAssociationCollection &AssocShapes, reco::RecoEcalCandidateCollection &RecoECand)
 

Private Member Functions

reco::RecoEcalCandidate correctEPosition (const reco::SuperCluster &original, const reco::HFEMClusterShape &shape)
 

Private Attributes

bool m_correct
 
double m_e9e25Cut
 
double m_intercept2DCut
 

Detailed Description

Author
K. Klapoetke – Minnesota

Definition at line 22 of file HFRecoEcalCandidateAlgo.h.

Constructor & Destructor Documentation

HFRecoEcalCandidateAlgo::HFRecoEcalCandidateAlgo ( bool  correct,
double  e9e25Cut,
double  intercept2DCut 
)

Definition at line 20 of file HFRecoEcalCandidateAlgo.cc.

20  :
21  m_correct(correct),
22  m_e9e25Cut(e9e25Cut),
23  m_intercept2DCut(intercept2DCut){
24 
25 }

Member Function Documentation

RecoEcalCandidate HFRecoEcalCandidateAlgo::correctEPosition ( const reco::SuperCluster original,
const reco::HFEMClusterShape shape 
)
private

Definition at line 27 of file HFRecoEcalCandidateAlgo.cc.

References reco::HFEMClusterShape::CellEta(), reco::HFEMClusterShape::CellPhi(), funct::cos(), reco::CaloCluster::energy(), reco::CaloCluster::eta(), reco::CaloCluster::phi(), and funct::sin().

Referenced by produce().

27  {
28  double energyCorrect=0.7397;//.7515;
29  double etaCorrect=.00938422+0.00682824*sin(6.28318531*shape.CellEta());//.0144225-.00484597*sin(6.17851*shape.CellEta());//0.01139;
30  double phiAmpCorrect=0.00644091;//-0.006483;
31  double phiFreqCorrect=6.28318531;//6.45377;
32 
33  double corEnergy= original.energy()/energyCorrect;
34  double corEta=original.eta();
35  corEta+=(original.eta()>0)?(etaCorrect):(-etaCorrect);
36  double corPhi=original.phi()+phiAmpCorrect*sin(phiFreqCorrect*shape.CellPhi());
37  double corPx=corEnergy*cos(corPhi)/cosh(corEta);
38  double corPy=corEnergy*sin(corPhi)/cosh(corEta);
39  double corPz=corEnergy*tanh(corEta);
40  RecoEcalCandidate corCand(0,
41  math::XYZTLorentzVector(corPx,corPy,corPz,corEnergy),
42  math::XYZPoint(0,0,0));
43 
44  return corCand;
45 }
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double CellPhi() const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:149
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:109
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:152
double CellEta() const
void HFRecoEcalCandidateAlgo::produce ( const edm::Handle< reco::SuperClusterCollection > &  SuperClusters,
const reco::HFEMClusterShapeAssociationCollection AssocShapes,
reco::RecoEcalCandidateCollection RecoECand 
)

Analyze the hits

Definition at line 47 of file HFRecoEcalCandidateAlgo.cc.

References correctEPosition(), funct::cos(), reco::HFEMClusterShape::eCOREe9(), reco::HFEMClusterShape::eLong3x3(), reco::HFEMClusterShape::eLong5x5(), reco::CaloCluster::energy(), reco::HFEMClusterShape::eSeL(), reco::CaloCluster::eta(), edm::AssociationMap< Tag >::find(), i, m_correct, m_e9e25Cut, m_intercept2DCut, reco::CaloCluster::phi(), funct::sin(), and edm::helpers::KeyVal< K, V >::val.

Referenced by HFRecoEcalCandidateProducer::produce().

49  {
50 
51 
52 
53  //get super's and cluster shapes and associations
54  for (unsigned int i=0; i < SuperClusters->size(); ++i) {
55  const SuperCluster& supClus=(*SuperClusters)[i];
57  const HFEMClusterShapeRef clusShapeRef=AssocShapes.find(theClusRef)->val;
58  const HFEMClusterShape& clusShape=*clusShapeRef;
59 
60  // basic candidate
61  double px=supClus.energy()*cos(supClus.phi())/cosh(supClus.eta());
62  double py=supClus.energy()*sin(supClus.phi())/cosh(supClus.eta());
63  double pz=supClus.energy()*tanh(supClus.eta());
64  RecoEcalCandidate theCand(0,
65  math::XYZTLorentzVector(px,py,pz,supClus.energy()),
66  math::XYZPoint(0,0,0));
67 
68  // correct it?
69  if (m_correct)
70  theCand=correctEPosition(supClus,clusShape);
71 
72  double e9e25=clusShape.eLong3x3()/clusShape.eLong5x5();
73  // EMID cuts...
74  //if((clusShape.e9e25()> m_e9e25Cut)&&((clusShape.eCOREe9()-(clusShape.eSeL()*1.125)) > m_intercept2DCut)){
75  if((e9e25> m_e9e25Cut)&&((clusShape.eCOREe9()-(clusShape.eSeL()*1.125)) > m_intercept2DCut)){
76  theCand.setSuperCluster(theClusRef);
77  RecoECand.push_back(theCand);
78  }
79  }
80 }
reco::RecoEcalCandidate correctEPosition(const reco::SuperCluster &original, const reco::HFEMClusterShape &shape)
int i
Definition: DBlmapReader.cc:9
double eLong5x5() const
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
const_iterator find(const key_type &k) const
find element with specified reference key
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:149
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:30
double eLong3x3() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:109
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:13
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:152

Member Data Documentation

bool HFRecoEcalCandidateAlgo::m_correct
private

Definition at line 35 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_e9e25Cut
private

Definition at line 36 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_intercept2DCut
private

Definition at line 37 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().