CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

HFRecoEcalCandidateAlgo Class Reference

#include <HFRecoEcalCandidateAlgo.h>

List of all members.

Public Member Functions

 HFRecoEcalCandidateAlgo (bool correct, double e9e25Cut, double intercept2DCut, const std::vector< double > &e1e9Cut, const std::vector< double > &eCOREe9Cut, const std::vector< double > &eSeLCut)
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_e1e9Cuthi
double m_e1e9Cutlo
double m_e9e25Cut
double m_eCOREe9Cuthi
double m_eCOREe9Cutlo
double m_eSeLCuthi
double m_eSeLCutlo
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,
const std::vector< double > &  e1e9Cut,
const std::vector< double > &  eCOREe9Cut,
const std::vector< double > &  eSeLCut 
)

Definition at line 20 of file HFRecoEcalCandidateAlgo.cc.

                                                                                   :

  m_correct(correct), 
  m_e9e25Cut(e9e25Cut),
  m_intercept2DCut(intercept2DCut),
  m_e1e9Cuthi(e1e9Cut[1]),
  m_eCOREe9Cuthi(eCOREe9Cut[1]),
  m_eSeLCuthi(eSeLCut[1]),
  m_e1e9Cutlo(e1e9Cut[0]),
  m_eCOREe9Cutlo(eCOREe9Cut[0]),
  m_eSeLCutlo(eSeLCut[0]){

}

Member Function Documentation

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

Definition at line 37 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().

                                                                                                                        {
  double energyCorrect=0.7397;//.7515;
  double etaCorrect=.00938422+0.00682824*sin(6.28318531*shape.CellEta());//.0144225-.00484597*sin(6.17851*shape.CellEta());//0.01139;
  double phiAmpCorrect=0.00644091;//-0.006483;
  double phiFreqCorrect=6.28318531;//6.45377;

  double corEnergy= original.energy()/energyCorrect;
  double corEta=original.eta();
  corEta+=(original.eta()>0)?(etaCorrect):(-etaCorrect);
  double corPhi=original.phi()+phiAmpCorrect*sin(phiFreqCorrect*shape.CellPhi());
  double corPx=corEnergy*cos(corPhi)/cosh(corEta);
  double corPy=corEnergy*sin(corPhi)/cosh(corEta);
  double corPz=corEnergy*tanh(corEta);
    RecoEcalCandidate corCand(0,
                              math::XYZTLorentzVector(corPx,corPy,corPz,corEnergy),
                              math::XYZPoint(0,0,0));

 return corCand;
}
void HFRecoEcalCandidateAlgo::produce ( const edm::Handle< reco::SuperClusterCollection > &  SuperClusters,
const reco::HFEMClusterShapeAssociationCollection AssocShapes,
reco::RecoEcalCandidateCollection RecoECand 
)

Analyze the hits

Definition at line 57 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(), reco::CaloCluster::eta(), edm::AssociationMap< Tag >::find(), i, m_correct, m_e1e9Cutlo, m_e9e25Cut, m_eCOREe9Cuthi, m_eCOREe9Cutlo, m_eSeLCuthi, m_eSeLCutlo, m_intercept2DCut, reco::CaloCluster::phi(), and funct::sin().

Referenced by HFRecoEcalCandidateProducer::produce().

                                                                              {
  
  
  
  //get super's and cluster shapes and associations 
  for (unsigned int i=0; i < SuperClusters->size(); ++i) {
    const SuperCluster& supClus=(*SuperClusters)[i];    
    reco::SuperClusterRef theClusRef=edm::Ref<SuperClusterCollection>(SuperClusters,i);
  const HFEMClusterShapeRef clusShapeRef=AssocShapes.find(theClusRef)->val;
    const HFEMClusterShape& clusShape=*clusShapeRef;

    // basic candidate
    double px=supClus.energy()*cos(supClus.phi())/cosh(supClus.eta());
    double py=supClus.energy()*sin(supClus.phi())/cosh(supClus.eta());
    double pz=supClus.energy()*tanh(supClus.eta());
    RecoEcalCandidate theCand(0,
                              math::XYZTLorentzVector(px,py,pz,supClus.energy()),
                              math::XYZPoint(0,0,0));

    // correct it?
    if (m_correct)
      theCand=correctEPosition(supClus,clusShape);

    double e9e25=clusShape.eLong3x3()/clusShape.eLong5x5();
    double e1e9=clusShape.eLong1x1()/clusShape.eLong3x3();
    // EMID cuts...  
    //if((clusShape.e9e25()> m_e9e25Cut)&&((clusShape.eCOREe9()-(clusShape.eSeL()*1.125)) > m_intercept2DCut)){
    //  if((e9e25> m_e9e25Cut)&&((clusShape.eCOREe9()-(clusShape.eSeL()*1.125)) > m_intercept2DCut)){
    double var2d=(clusShape.eCOREe9()-(clusShape.eSeL()*1.125));
 
    bool isAcceptable=true;
    isAcceptable=isAcceptable && (e9e25> m_e9e25Cut);
    isAcceptable=isAcceptable && (var2d > m_intercept2DCut);
    isAcceptable=isAcceptable && ((e1e9< m_e1e9Cuthi)&&(e1e9> m_e1e9Cutlo));
    isAcceptable=isAcceptable && ((clusShape.eCOREe9()< m_eCOREe9Cuthi)&&(clusShape.eCOREe9()>  m_eCOREe9Cutlo));
    isAcceptable=isAcceptable && ((clusShape.eSeL()<m_eSeLCuthi)&&(clusShape.eSeL()>  m_eSeLCutlo));
   
    
    if(isAcceptable){

      theCand.setSuperCluster(theClusRef);
      RecoECand.push_back(theCand);
    }
  }
}

Member Data Documentation

Definition at line 39 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 42 of file HFRecoEcalCandidateAlgo.h.

Definition at line 45 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 40 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 43 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 46 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 44 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 47 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().

Definition at line 41 of file HFRecoEcalCandidateAlgo.h.

Referenced by produce().