CMS 3D CMS Logo

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, double intercept2DSlope, const std::vector< double > &e1e9Cut, const std::vector< double > &eCOREe9Cut, const std::vector< double > &eSeLCut, const reco::HFValueStruct hfvv)
 
void produce (const edm::Handle< reco::SuperClusterCollection > &SuperClusters, const reco::HFEMClusterShapeAssociationCollection &AssocShapes, reco::RecoEcalCandidateCollection &RecoECand, int nvtx)
 

Private Member Functions

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

Private Attributes

bool m_correct
 
bool m_correctForPileup
 
double m_e1e9Cuthi
 
double m_e1e9Cutlo
 
double m_e9e25Cut
 
double m_eCOREe9Cuthi
 
double m_eCOREe9Cutlo
 
int m_era
 
double m_eSeLCuthi
 
double m_eSeLCutlo
 
reco::HFValueStruct m_hfvv
 
double m_intercept2DCut
 
double m_intercept2DSlope
 

Detailed Description

Author
K. Klapoetke – Minnesota

Definition at line 23 of file HFRecoEcalCandidateAlgo.h.

Constructor & Destructor Documentation

HFRecoEcalCandidateAlgo::HFRecoEcalCandidateAlgo ( bool  correct,
double  e9e25Cut,
double  intercept2DCut,
double  intercept2DSlope,
const std::vector< double > &  e1e9Cut,
const std::vector< double > &  eCOREe9Cut,
const std::vector< double > &  eSeLCut,
const reco::HFValueStruct  hfvv 
)

Definition at line 22 of file HFRecoEcalCandidateAlgo.cc.

27  :
28 
29  m_correct(correct),
33  m_e1e9Cuthi(e1e9Cut[1]),
35  m_eSeLCuthi(eSeLCut[1]),
36  m_e1e9Cutlo(e1e9Cut[0]),
38  m_eSeLCutlo(eSeLCut[0]),
39  m_era(4),
40  m_hfvv(hfvv)
41 {
42 
43 }

Member Function Documentation

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

Definition at line 45 of file HFRecoEcalCandidateAlgo.cc.

References reco::HFEMClusterShape::CellEta(), reco::HFEMClusterShape::CellPhi(), funct::cos(), reco::dp, reco::HFEMClusterShape::eLong3x3(), reco::HFValueStruct::EnCor(), reco::CaloCluster::energy(), reco::CaloCluster::eta(), GetRecoTauVFromDQM_MC_cff::kk, cmsBatch::log, m_hfvv, p1, p2, p3, reco::CaloCluster::phi(), reco::HFValueStruct::PUIntercept(), reco::HFValueStruct::PUSlope(), FWPFMaths::sgn(), and funct::sin().

Referenced by produce().

45  {
46  double corEta=original.eta();
47  //piece-wise log energy correction to eta
48  double logel=log(shape.eLong3x3()/100.0);
49  double lowC= 0.00911;
50  double hiC = 0.0193;
51  double logSlope = 0.0146;
52  double logIntercept=-0.00988;
53  double sgn=(original.eta()>0)?(1.0):(-1.0);
54  if(logel<=1.25)corEta+=(lowC)*sgn;
55  else if((logel>1.25)&&(logel<=2))corEta+=(logIntercept+logSlope*logel)*sgn;
56  else if(logel>2)corEta+=hiC*sgn;
57  //poly-3 cell eta correction to eta (de)
58  double p0= 0.0125;
59  double p1=-0.0475;
60  double p2= 0.0732;
61  double p3=-0.0506;
62  double de= sgn*(p0+(p1*(shape.CellEta()))+(p2*(shape.CellEta()*shape.CellEta()))+(p3*(shape.CellEta()*shape.CellEta()*shape.CellEta())));
63  corEta+=de;
64 
65  double corPhi=original.phi();
66  //poly-3 cell phi correction to phi (dp)
67  p0= 0.0115;
68  p1=-0.0394;
69  p2= 0.0486;
70  p3=-0.0335;
71  double dp =p0+(p1*(shape.CellPhi()))+(p2*(shape.CellPhi()*shape.CellPhi()))+(p3*(shape.CellPhi()*shape.CellPhi()*shape.CellPhi()));
72  corPhi+=dp;
73 
74 
75  double corEnergy= original.energy();
76  //energy corrections
77  //flat correction for using only long fibers
78  double energyCorrect=0.7397;
79  corEnergy= corEnergy/energyCorrect;
80  //corection based on ieta for pileup and general
81  //find ieta
82  int ieta=0;
83  double etabounds[30]={
84  2.846,2.957,3.132,3.307,3.482,
85  3.657,3.833,4.006,4.184,4.357,
86  4.532,4.709,4.882,5.184};
87  for (int kk=0;kk<=12;kk++){
88  if((fabs(corEta) < etabounds[kk+1])&&(fabs(corEta) > etabounds[kk])){
89  ieta = (corEta > 0)?(kk+29):(-kk-29);
90  }
91  }
92  //check if ieta is in bounds, should be [-41, -29] U [29, 41]
93  if (ieta != 0) corEnergy=(m_hfvv.PUSlope(ieta)*1.0*(nvtx-1)+m_hfvv.PUIntercept(ieta)*1.0)*corEnergy*1.0*m_hfvv.EnCor(ieta);
94 
95  //re-calculate the energy vector
96  double corPx=corEnergy*cos(corPhi)/cosh(corEta);
97  double corPy=corEnergy*sin(corPhi)/cosh(corEta);
98  double corPz=corEnergy*tanh(corEta);
99 
100  //make the candidate
101  RecoEcalCandidate corCand(0,
102  math::XYZTLorentzVector(corPx,corPy,corPz,corEnergy),
103  math::XYZPoint(0,0,0));
104 
105  return corCand;
106 }
float sgn(float val)
Definition: FWPFMaths.cc:9
Sin< T >::type sin(const T &t)
Definition: Sin.h:22
double CellPhi() const
double PUIntercept(int ieta) const
double eta() const
pseudorapidity of cluster centroid
Definition: CaloCluster.h:166
double EnCor(int ieta) const
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double eLong3x3() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:124
double p2[4]
Definition: TauolaWrapper.h:90
double PUSlope(int ieta) const
auto dp
Definition: deltaR.h:22
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double p1[4]
Definition: TauolaWrapper.h:89
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:169
double CellEta() const
double p3[4]
Definition: TauolaWrapper.h:91
void HFRecoEcalCandidateAlgo::produce ( const edm::Handle< reco::SuperClusterCollection > &  SuperClusters,
const reco::HFEMClusterShapeAssociationCollection AssocShapes,
reco::RecoEcalCandidateCollection RecoECand,
int  nvtx 
)

Analyze the hits

Definition at line 108 of file HFRecoEcalCandidateAlgo.cc.

References correctEPosition(), funct::cos(), reco::HFEMClusterShape::eCOREe9(), reco::HFEMClusterShape::eLong1x1(), reco::HFEMClusterShape::eLong3x3(), reco::HFEMClusterShape::eLong5x5(), reco::CaloCluster::energy(), reco::HFEMClusterShape::eSeL(), hf_egamma::eSeLCorrected(), reco::HFEMClusterShape::eShort3x3(), reco::CaloCluster::eta(), edm::AssociationMap< Tag >::find(), mps_fire::i, m_correct, m_e1e9Cutlo, m_e9e25Cut, m_eCOREe9Cuthi, m_eCOREe9Cutlo, m_eSeLCuthi, m_eSeLCutlo, m_intercept2DCut, m_intercept2DSlope, reco::CaloCluster::phi(), and funct::sin().

Referenced by HFRecoEcalCandidateProducer::produce(), and HLTHFRecoEcalCandidateProducer::produce().

111  {
112 
113 
114 
115  //get super's and cluster shapes and associations
116  for (unsigned int i=0; i < SuperClusters->size(); ++i) {
117  const SuperCluster& supClus=(*SuperClusters)[i];
118  reco::SuperClusterRef theClusRef=edm::Ref<SuperClusterCollection>(SuperClusters,i);
119  const HFEMClusterShapeRef clusShapeRef=AssocShapes.find(theClusRef)->val;
120  const HFEMClusterShape& clusShape=*clusShapeRef;
121 
122  // basic candidate
123  double px=supClus.energy()*cos(supClus.phi())/cosh(supClus.eta());
124  double py=supClus.energy()*sin(supClus.phi())/cosh(supClus.eta());
125  double pz=supClus.energy()*tanh(supClus.eta());
126  RecoEcalCandidate theCand(0,
127  math::XYZTLorentzVector(px,py,pz,supClus.energy()),
128  math::XYZPoint(0,0,0));
129 
130  // correct it?
131  if (m_correct)
132  theCand=correctEPosition(supClus,clusShape,nvtx);
133 
134  double e9e25=clusShape.eLong3x3()/clusShape.eLong5x5();
135  double e1e9=clusShape.eLong1x1()/clusShape.eLong3x3();
136  double eSeL=hf_egamma::eSeLCorrected(clusShape.eShort3x3(),clusShape.eLong3x3(),4);
137  double var2d=(clusShape.eCOREe9()-(eSeL*m_intercept2DSlope));
138 
139  bool isAcceptable=true;
140  isAcceptable=isAcceptable && (e9e25> m_e9e25Cut);
141  isAcceptable=isAcceptable && (var2d > m_intercept2DCut);
142  isAcceptable=isAcceptable && ((e1e9< m_e1e9Cuthi)&&(e1e9> m_e1e9Cutlo));
143  isAcceptable=isAcceptable && ((clusShape.eCOREe9()< m_eCOREe9Cuthi)&&(clusShape.eCOREe9()> m_eCOREe9Cutlo));
144  isAcceptable=isAcceptable && ((clusShape.eSeL()<m_eSeLCuthi)&&(clusShape.eSeL()> m_eSeLCutlo));
145 
146 
147  if(isAcceptable){
148 
149  theCand.setSuperCluster(theClusRef);
150  RecoECand.push_back(theCand);
151  }
152  }
153 }
double eLong5x5() const
double eLong1x1() 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:166
XYZTLorentzVectorD XYZTLorentzVector
Lorentz vector with cylindrical internal representation using pseudorapidity.
Definition: LorentzVector.h:29
double eShort3x3() const
double eLong3x3() const
Cos< T >::type cos(const T &t)
Definition: Cos.h:22
double energy() const
cluster energy
Definition: CaloCluster.h:124
XYZPointD XYZPoint
point in space with cartesian internal representation
Definition: Point3D.h:12
double eSeLCorrected(double es, double el, double pc, double px, double py)
reco::RecoEcalCandidate correctEPosition(const reco::SuperCluster &original, const reco::HFEMClusterShape &shape, int nvtx)
double phi() const
azimuthal angle of cluster centroid
Definition: CaloCluster.h:169

Member Data Documentation

bool HFRecoEcalCandidateAlgo::m_correct
private

Definition at line 43 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

bool HFRecoEcalCandidateAlgo::m_correctForPileup
private

Definition at line 54 of file HFRecoEcalCandidateAlgo.h.

double HFRecoEcalCandidateAlgo::m_e1e9Cuthi
private

Definition at line 47 of file HFRecoEcalCandidateAlgo.h.

double HFRecoEcalCandidateAlgo::m_e1e9Cutlo
private

Definition at line 50 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_e9e25Cut
private

Definition at line 44 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_eCOREe9Cuthi
private

Definition at line 48 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_eCOREe9Cutlo
private

Definition at line 51 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

int HFRecoEcalCandidateAlgo::m_era
private

Definition at line 53 of file HFRecoEcalCandidateAlgo.h.

double HFRecoEcalCandidateAlgo::m_eSeLCuthi
private

Definition at line 49 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_eSeLCutlo
private

Definition at line 52 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

reco::HFValueStruct HFRecoEcalCandidateAlgo::m_hfvv
private

Definition at line 55 of file HFRecoEcalCandidateAlgo.h.

Referenced by correctEPosition().

double HFRecoEcalCandidateAlgo::m_intercept2DCut
private

Definition at line 45 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

double HFRecoEcalCandidateAlgo::m_intercept2DSlope
private

Definition at line 46 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().