CMS 3D CMS Logo

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