CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_6_1_2_SLHC4_patch1/src/RecoEgamma/EgammaTools/src/EcalClusterLocal.cc

Go to the documentation of this file.
00001 #include "../interface/EcalClusterLocal.h"
00002 
00003 #include "FWCore/Framework/interface/EventSetup.h"
00004 #include "FWCore/Framework/interface/ESHandle.h"
00005 #include "Geometry/Records/interface/CaloGeometryRecord.h"
00006 #include "Geometry/CaloGeometry/interface/CaloGeometry.h"
00007 #include "Geometry/CaloGeometry/interface/TruncatedPyramid.h"
00008 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00009 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00010 #include "Geometry/CaloGeometry/interface/CaloSubdetectorGeometry.h"
00011 #include "DataFormats/Math/interface/deltaR.h"
00012 #include "DataFormats/EgammaReco/interface/BasicCluster.h"
00013 #include "TMath.h"
00014 #include "TVector2.h"
00015 
00016 EcalClusterLocal::EcalClusterLocal()
00017 {}
00018 
00019 EcalClusterLocal::~EcalClusterLocal()
00020 {}
00021 
00022 void EcalClusterLocal::localCoordsEB( const reco::BasicCluster &bclus, const edm::EventSetup &es, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const
00023 {
00024   
00025   assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalBarrel);
00026   
00027   edm::ESHandle<CaloGeometry> pG;
00028   es.get<CaloGeometryRecord>().get(pG); 
00029   
00030   const CaloSubdetectorGeometry* geom=pG->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);//EcalBarrel = 1
00031   
00032   const math::XYZPoint position_ = bclus.position(); 
00033   double Theta = -position_.theta()+0.5*TMath::Pi();
00034   double Eta = position_.eta();
00035   double Phi = TVector2::Phi_mpi_pi(position_.phi());
00036   
00037   //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
00038   // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
00039   const float X0 = 0.89; const float T0 = 7.4;
00040   double depth = X0 * (T0 + log(bclus.energy()));
00041   
00042   //find max energy crystal
00043   std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
00044   float drmin = 999.;
00045   EBDetId crystalseed;
00046   //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
00047   for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {    
00048     
00049     EBDetId crystal(crystals_vector[icry].first);
00050         
00051     const CaloCellGeometry* cell=geom->getGeometry(crystal);
00052     GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
00053     double EtaCentr = center_pos.eta();
00054     double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
00055 
00056     float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
00057     if (dr<drmin) {
00058       drmin = dr;
00059       crystalseed = crystal;
00060     }
00061     
00062   }
00063   
00064   ieta = crystalseed.ieta();
00065   iphi = crystalseed.iphi();
00066   
00067   // Get center cell position from shower depth
00068   const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
00069   const TruncatedPyramid *cpyr = dynamic_cast<const TruncatedPyramid*>(cell);
00070 
00071   thetatilt = cpyr->getThetaAxis();
00072   phitilt = cpyr->getPhiAxis();
00073 
00074   GlobalPoint center_pos = cpyr->getPosition(depth);
00075   
00076   double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
00077   double PhiWidth = (TMath::Pi()/180.);
00078   phicry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
00079     //Some flips to take into account ECAL barrel symmetries:
00080   if (ieta<0) phicry *= -1.;  
00081   
00082   double ThetaCentr = -center_pos.theta()+0.5*TMath::Pi();
00083   double ThetaWidth = (TMath::Pi()/180.)*TMath::Cos(ThetaCentr);
00084   etacry = (Theta-ThetaCentr)/ThetaWidth;    
00085   //flip to take into account ECAL barrel symmetries:
00086   if (ieta<0) etacry *= -1.;  
00087   
00088   return;
00089 
00090 }
00091 
00092 void EcalClusterLocal::localCoordsEE( const reco::BasicCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
00093 {
00094     
00095   assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalEndcap);
00096   
00097   edm::ESHandle<CaloGeometry> pG;
00098   es.get<CaloGeometryRecord>().get(pG); 
00099   
00100   const CaloSubdetectorGeometry* geom=pG->getSubdetectorGeometry(DetId::Ecal,EcalEndcap);//EcalBarrel = 1
00101   
00102   const math::XYZPoint position_ = bclus.position(); 
00103   //double Theta = -position_.theta()+0.5*TMath::Pi();
00104   double Eta = position_.eta();
00105   double Phi = TVector2::Phi_mpi_pi(position_.phi());
00106   double X = position_.x();
00107   double Y = position_.y();
00108   
00109   //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
00110   // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
00111   const float X0 = 0.89; float T0 = 1.2;
00112   //different T0 value if outside of preshower coverage
00113   if (TMath::Abs(bclus.eta())<1.653) T0 = 3.1;
00114   
00115   double depth = X0 * (T0 + log(bclus.energy()));
00116   
00117   //find max energy crystal
00118   std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
00119   float drmin = 999.;
00120   EEDetId crystalseed;
00121   //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
00122   for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {    
00123     
00124     EEDetId crystal(crystals_vector[icry].first);
00125         
00126     const CaloCellGeometry* cell=geom->getGeometry(crystal);
00127     GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
00128     double EtaCentr = center_pos.eta();
00129     double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
00130 
00131     float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
00132     if (dr<drmin) {
00133       drmin = dr;
00134       crystalseed = crystal;
00135     }
00136     
00137   }
00138   
00139   ix = crystalseed.ix();
00140   iy = crystalseed.iy();
00141   
00142   // Get center cell position from shower depth
00143   const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
00144   const TruncatedPyramid *cpyr = dynamic_cast<const TruncatedPyramid*>(cell);
00145 
00146   thetatilt = cpyr->getThetaAxis();
00147   phitilt = cpyr->getPhiAxis();
00148 
00149   GlobalPoint center_pos = cpyr->getPosition(depth);
00150   
00151   double XCentr = center_pos.x();
00152   double XWidth = 2.59;
00153   xcry = (X-XCentr)/XWidth;
00154  
00155   
00156   double YCentr = center_pos.y();
00157   double YWidth = 2.59;
00158   ycry = (Y-YCentr)/YWidth;
00159   
00160   return;
00161 
00162 }
00163 
00164