CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_4_1_8_patch12/src/RecoEgamma/EgammaHFProducers/src/HFRecoEcalCandidateAlgo.cc

Go to the documentation of this file.
00001 
00009 //#includes
00010 #include "DataFormats/Math/interface/LorentzVector.h"
00011 #include "DataFormats/Math/interface/deltaR.h"
00012 #include "DataFormats/Math/interface/deltaPhi.h"
00013 #include "CLHEP/Vector/LorentzVector.h"
00014 #include <iostream>
00015 #include "DataFormats/RecoCandidate/interface/RecoEcalCandidate.h"
00016 #include "RecoEgamma/EgammaHFProducers/interface/HFRecoEcalCandidateAlgo.h"
00017 using namespace std;
00018 using namespace reco;
00019 
00020 HFRecoEcalCandidateAlgo::HFRecoEcalCandidateAlgo(bool correct,double e9e25Cut,double intercept2DCut) :
00021   m_correct(correct), 
00022   m_e9e25Cut(e9e25Cut),
00023   m_intercept2DCut(intercept2DCut){
00024   
00025 }
00026 
00027 RecoEcalCandidate HFRecoEcalCandidateAlgo::correctEPosition(const SuperCluster& original , const HFEMClusterShape& shape) {
00028   double energyCorrect=0.7397;//.7515;
00029   double etaCorrect=.00938422+0.00682824*sin(6.28318531*shape.CellEta());//.0144225-.00484597*sin(6.17851*shape.CellEta());//0.01139;
00030   double phiAmpCorrect=0.00644091;//-0.006483;
00031   double phiFreqCorrect=6.28318531;//6.45377;
00032 
00033   double corEnergy= original.energy()/energyCorrect;
00034   double corEta=original.eta();
00035   corEta+=(original.eta()>0)?(etaCorrect):(-etaCorrect);
00036   double corPhi=original.phi()+phiAmpCorrect*sin(phiFreqCorrect*shape.CellPhi());
00037   double corPx=corEnergy*cos(corPhi)/cosh(corEta);
00038   double corPy=corEnergy*sin(corPhi)/cosh(corEta);
00039   double corPz=corEnergy*tanh(corEta);
00040     RecoEcalCandidate corCand(0,
00041                               math::XYZTLorentzVector(corPx,corPy,corPz,corEnergy),
00042                               math::XYZPoint(0,0,0));
00043 
00044  return corCand;
00045 }
00046 
00047 void HFRecoEcalCandidateAlgo::produce(const edm::Handle<SuperClusterCollection>& SuperClusters,
00048                                       const HFEMClusterShapeAssociationCollection& AssocShapes,
00049                                       RecoEcalCandidateCollection& RecoECand) {
00050   
00051   
00052   
00053   //get super's and cluster shapes and associations 
00054   for (unsigned int i=0; i < SuperClusters->size(); ++i) {
00055     const SuperCluster& supClus=(*SuperClusters)[i];    
00056     reco::SuperClusterRef theClusRef=edm::Ref<SuperClusterCollection>(SuperClusters,i);
00057   const HFEMClusterShapeRef clusShapeRef=AssocShapes.find(theClusRef)->val;
00058     const HFEMClusterShape& clusShape=*clusShapeRef;
00059 
00060     // basic candidate
00061     double px=supClus.energy()*cos(supClus.phi())/cosh(supClus.eta());
00062     double py=supClus.energy()*sin(supClus.phi())/cosh(supClus.eta());
00063     double pz=supClus.energy()*tanh(supClus.eta());
00064     RecoEcalCandidate theCand(0,
00065                               math::XYZTLorentzVector(px,py,pz,supClus.energy()),
00066                               math::XYZPoint(0,0,0));
00067 
00068     // correct it?
00069     if (m_correct)
00070       theCand=correctEPosition(supClus,clusShape);
00071 
00072     double e9e25=clusShape.eLong3x3()/clusShape.eLong5x5();
00073     // EMID cuts...  
00074     //if((clusShape.e9e25()> m_e9e25Cut)&&((clusShape.eCOREe9()-(clusShape.eSeL()*1.125)) > m_intercept2DCut)){
00075     if((e9e25> m_e9e25Cut)&&((clusShape.eCOREe9()-(clusShape.eSeL()*1.125)) > m_intercept2DCut)){
00076       theCand.setSuperCluster(theClusRef);
00077       RecoECand.push_back(theCand);
00078     }
00079   }
00080 }