CMS 3D CMS Logo

Public Member Functions

EcalClusterLocal Class Reference

#include <EcalClusterLocal.h>

List of all members.

Public Member Functions

 EcalClusterLocal ()
void localCoordsEB (const reco::BasicCluster &bclus, const edm::EventSetup &es, float &etacry, float &phicry, int &ieta, int &iphi, float &thetatilt, float &phitilt) const
void localCoordsEE (const reco::BasicCluster &bclus, const edm::EventSetup &es, float &xcry, float &ycry, int &ix, int &iy, float &thetatilt, float &phitilt) const
 ~EcalClusterLocal ()

Detailed Description

Function to compute local coordinates of Ecal clusters (adapted from RecoEcal/EgammaCoreTools/plugins/EcalClusterLocal) $Id: EcalClusterLocal.h $Date: $Revision:

Author:
Josh Bendavid, MIT, 2011

Definition at line 19 of file EcalClusterLocal.h.


Constructor & Destructor Documentation

EcalClusterLocal::EcalClusterLocal ( )

Definition at line 16 of file EcalClusterLocal.cc.

{}
EcalClusterLocal::~EcalClusterLocal ( )

Definition at line 19 of file EcalClusterLocal.cc.

{}

Member Function Documentation

void EcalClusterLocal::localCoordsEB ( const reco::BasicCluster bclus,
const edm::EventSetup es,
float &  etacry,
float &  phicry,
int &  ieta,
int &  iphi,
float &  thetatilt,
float &  phitilt 
) const

Definition at line 22 of file EcalClusterLocal.cc.

References deltaR(), DetId::Ecal, EcalBarrel, reco::tau::disc::Eta(), PV3DBase< T, PVType, FrameType >::eta(), first, relativeConstraints::geom, edm::EventSetup::get(), CaloSubdetectorGeometry::getGeometry(), TruncatedPyramid::getPhiAxis(), TruncatedPyramid::getPosition(), TruncatedPyramid::getThetaAxis(), EBDetId::ieta(), EBDetId::iphi(), funct::log(), PV3DBase< T, PVType, FrameType >::phi(), colinearityKinematic::Phi, Phi_mpi_pi(), Pi, PV3DBase< T, PVType, FrameType >::theta(), and X0.

Referenced by EGEnergyCorrector::CorrectedEnergyWithError().

{
  
  assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalBarrel);
  
  edm::ESHandle<CaloGeometry> pG;
  es.get<CaloGeometryRecord>().get(pG); 
  
  const CaloSubdetectorGeometry* geom=pG->getSubdetectorGeometry(DetId::Ecal,EcalBarrel);//EcalBarrel = 1
  
  const math::XYZPoint position_ = bclus.position(); 
  double Theta = -position_.theta()+0.5*TMath::Pi();
  double Eta = position_.eta();
  double Phi = TVector2::Phi_mpi_pi(position_.phi());
  
  //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
  // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
  const float X0 = 0.89; const float T0 = 7.4;
  double depth = X0 * (T0 + log(bclus.energy()));
  
  //find max energy crystal
  std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
  float drmin = 999.;
  EBDetId crystalseed;
  //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
  for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {    
    
    EBDetId crystal(crystals_vector[icry].first);
        
    const CaloCellGeometry* cell=geom->getGeometry(crystal);
    GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
    double EtaCentr = center_pos.eta();
    double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());

    float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
    if (dr<drmin) {
      drmin = dr;
      crystalseed = crystal;
    }
    
  }
  
  ieta = crystalseed.ieta();
  iphi = crystalseed.iphi();
  
  // Get center cell position from shower depth
  const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
  const TruncatedPyramid *cpyr = dynamic_cast<const TruncatedPyramid*>(cell);

  thetatilt = cpyr->getThetaAxis();
  phitilt = cpyr->getPhiAxis();

  GlobalPoint center_pos = cpyr->getPosition(depth);
  
  double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
  double PhiWidth = (TMath::Pi()/180.);
  phicry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
    //Some flips to take into account ECAL barrel symmetries:
  if (ieta<0) phicry *= -1.;  
  
  double ThetaCentr = -center_pos.theta()+0.5*TMath::Pi();
  double ThetaWidth = (TMath::Pi()/180.)*TMath::Cos(ThetaCentr);
  etacry = (Theta-ThetaCentr)/ThetaWidth;    
  //flip to take into account ECAL barrel symmetries:
  if (ieta<0) etacry *= -1.;  
  
  return;

}
void EcalClusterLocal::localCoordsEE ( const reco::BasicCluster bclus,
const edm::EventSetup es,
float &  xcry,
float &  ycry,
int &  ix,
int &  iy,
float &  thetatilt,
float &  phitilt 
) const

Definition at line 92 of file EcalClusterLocal.cc.

References deltaR(), DetId::Ecal, EcalEndcap, reco::tau::disc::Eta(), PV3DBase< T, PVType, FrameType >::eta(), first, relativeConstraints::geom, edm::EventSetup::get(), CaloSubdetectorGeometry::getGeometry(), TruncatedPyramid::getPhiAxis(), TruncatedPyramid::getPosition(), TruncatedPyramid::getThetaAxis(), EEDetId::ix(), EEDetId::iy(), funct::log(), PV3DBase< T, PVType, FrameType >::phi(), colinearityKinematic::Phi, Phi_mpi_pi(), PV3DBase< T, PVType, FrameType >::x(), X, X0, and PV3DBase< T, PVType, FrameType >::y().

{
    
  assert(bclus.hitsAndFractions().at(0).first.subdetId()==EcalEndcap);
  
  edm::ESHandle<CaloGeometry> pG;
  es.get<CaloGeometryRecord>().get(pG); 
  
  const CaloSubdetectorGeometry* geom=pG->getSubdetectorGeometry(DetId::Ecal,EcalEndcap);//EcalBarrel = 1
  
  const math::XYZPoint position_ = bclus.position(); 
  //double Theta = -position_.theta()+0.5*TMath::Pi();
  double Eta = position_.eta();
  double Phi = TVector2::Phi_mpi_pi(position_.phi());
  double X = position_.x();
  double Y = position_.y();
  
  //Calculate expected depth of the maximum shower from energy (like in PositionCalc::Calculate_Location()):
  // The parameters X0 and T0 are hardcoded here because these values were used to calculate the corrections:
  const float X0 = 0.89; float T0 = 1.2;
  //different T0 value if outside of preshower coverage
  if (TMath::Abs(bclus.eta())<1.653) T0 = 3.1;
  
  double depth = X0 * (T0 + log(bclus.energy()));
  
  //find max energy crystal
  std::vector< std::pair<DetId, float> > crystals_vector = bclus.hitsAndFractions();
  float drmin = 999.;
  EEDetId crystalseed;
  //printf("starting loop over crystals, etot = %5f:\n",bclus.energy());
  for (unsigned int icry=0; icry!=crystals_vector.size(); ++icry) {    
    
    EEDetId crystal(crystals_vector[icry].first);
        
    const CaloCellGeometry* cell=geom->getGeometry(crystal);
    GlobalPoint center_pos = (dynamic_cast<const TruncatedPyramid*>(cell))->getPosition(depth);
    double EtaCentr = center_pos.eta();
    double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());

    float dr = reco::deltaR(Eta,Phi,EtaCentr,PhiCentr);
    if (dr<drmin) {
      drmin = dr;
      crystalseed = crystal;
    }
    
  }
  
  ix = crystalseed.ix();
  iy = crystalseed.iy();
  
  // Get center cell position from shower depth
  const CaloCellGeometry* cell=geom->getGeometry(crystalseed);
  const TruncatedPyramid *cpyr = dynamic_cast<const TruncatedPyramid*>(cell);

  thetatilt = cpyr->getThetaAxis();
  phitilt = cpyr->getPhiAxis();

  GlobalPoint center_pos = cpyr->getPosition(depth);
  
  double XCentr = center_pos.x();
  double XWidth = 2.59;
  xcry = (X-XCentr)/XWidth;
 
  
  double YCentr = center_pos.y();
  double YWidth = 2.59;
  ycry = (Y-YCentr)/YWidth;
  
  return;

}