CMS 3D CMS Logo

/data/refman/pasoursint/CMSSW_5_3_0/src/RecoParticleFlow/PFClusterTools/src/PFPhotonClusters.cc

Go to the documentation of this file.
00001 #include "DataFormats/EcalDetId/interface/EBDetId.h"
00002 #include "DataFormats/EcalDetId/interface/EEDetId.h"
00003 #include "DataFormats/ParticleFlowReco/interface/PFRecHit.h"
00004 #include "RecoParticleFlow/PFClusterTools/interface/PFPhotonClusters.h"
00005 
00006 #include <TMath.h>
00007 #include <TVector2.h>
00008 PFPhotonClusters::PFPhotonClusters(PFClusterRef PFClusterRef):
00009   PFClusterRef_(PFClusterRef)
00010 {
00011   if(PFClusterRef_->layer()==PFLayer:: ECAL_BARREL )isEB_=true;
00012   else isEB_=false;
00013   SetSeed();
00014   PFCrystalCoor();
00015   for(int i=0; i<5; ++i)
00016     for(int j=0; j<5; ++j)e5x5_[i][j]=0;
00017   FillClusterShape();
00018   FillClusterWidth();
00019 }
00020 
00021 void PFPhotonClusters::SetSeed(){
00022   double PFSeedE=0;
00023   math::XYZVector axis;
00024   math::XYZVector position;
00025   DetId idseed;
00026   const std::vector< reco::PFRecHitFraction >& PFRecHits=
00027     PFClusterRef_->recHitFractions();
00028   
00029   for(std::vector< reco::PFRecHitFraction >::const_iterator it = PFRecHits.begin();
00030       it != PFRecHits.end(); ++it){
00031     const PFRecHitRef& RefPFRecHit = it->recHitRef();
00032     double frac=it->fraction();
00033     float E= RefPFRecHit->energy()* frac;
00034     if(E>PFSeedE){
00035       PFSeedE=E;  
00036       axis=RefPFRecHit->getAxisXYZ();
00037       position=RefPFRecHit->position();
00038       idseed = RefPFRecHit->detId();
00039     }
00040   }
00041   idseed_=idseed;
00042   seedPosition_=position;
00043   seedAxis_=axis;
00044 }
00045 
00046 void PFPhotonClusters::PFCrystalCoor(){
00047   if(PFClusterRef_->layer()==PFLayer:: ECAL_BARREL ){//is Barrel
00048     isEB_=true;
00049     EBDetId EBidSeed=EBDetId(idseed_.rawId());
00050     CrysIEta_=EBidSeed.ieta();
00051     CrysIPhi_=EBidSeed.iphi();
00052     double depth = PFClusterRef_->getDepthCorrection(PFClusterRef_->energy(), false, false);
00053     math::XYZVector center_pos = seedPosition_+depth*seedAxis_;
00054     //Crystal Coordinates:
00055     double Pi=TMath::Pi();
00056     float Phi=PFClusterRef_->position().phi(); 
00057     double Theta = -(PFClusterRef_->position().theta())+0.5* Pi;
00058     double PhiCentr = TVector2::Phi_mpi_pi(center_pos.phi());
00059     double PhiWidth = (Pi/180.);
00060     double PhiCry = (TVector2::Phi_mpi_pi(Phi-PhiCentr))/PhiWidth;
00061     double ThetaCentr = -center_pos.theta()+0.5*Pi;
00062     double ThetaWidth = (Pi/180.)*cos(ThetaCentr);
00063     
00064     double EtaCry = (Theta-ThetaCentr)/ThetaWidth; 
00065     CrysEta_=EtaCry;
00066     CrysPhi_=PhiCry;
00067     
00068     if(abs(CrysIEta_)==1 || abs(CrysIEta_)==2 )
00069       CrysIEtaCrack_=abs(CrysIEta_);
00070     if(abs(CrysIEta_)>2 && abs(CrysIEta_)<24)
00071       CrysIEtaCrack_=3;
00072     if(abs(CrysIEta_)==24)
00073       CrysIEtaCrack_=4;
00074     if(abs(CrysIEta_)==25)
00075       CrysIEtaCrack_=5;
00076     if(abs(CrysIEta_)==26)
00077       CrysIEtaCrack_=6;
00078     if(abs(CrysIEta_)==27)
00079       CrysIEtaCrack_=7;
00080     if(abs(CrysIEta_)>27 &&  abs(CrysIEta_)<44)
00081       CrysIEtaCrack_=8;
00082     if(abs(CrysIEta_)==44)
00083       CrysIEtaCrack_=9;
00084     if(abs(CrysIEta_)==45)
00085       CrysIEtaCrack_=10;
00086     if(abs(CrysIEta_)==46)
00087       CrysIEtaCrack_=11;
00088     if(abs(CrysIEta_)==47)
00089       CrysIEtaCrack_=12;
00090     if(abs(CrysIEta_)>47 &&  abs(CrysIEta_)<64)
00091       CrysIEtaCrack_=13;
00092     if(abs(CrysIEta_)==64)
00093       CrysIEtaCrack_=14;
00094     if(abs(CrysIEta_)==65)
00095         CrysIEtaCrack_=15;
00096     if(abs(CrysIEta_)==66)
00097       CrysIEtaCrack_=16;
00098     if(abs(CrysIEta_)==67)
00099       CrysIEtaCrack_=17;
00100     if(abs(CrysIEta_)>67 &&  abs(CrysIEta_)<84)
00101       CrysIEtaCrack_=18;
00102     if(abs(CrysIEta_)==84)
00103       CrysIEtaCrack_=19;
00104     if(abs(CrysIEta_)==85)
00105       CrysIEtaCrack_=20;
00106   }
00107   else{
00108     isEB_=false;
00109     EEDetId EEidSeed=EEDetId(idseed_.rawId());
00110     CrysIX_=EEidSeed.ix();
00111     CrysIY_=EEidSeed.iy();
00112     float X0 = 0.89; float T0 = 1.2;
00113     if(fabs(PFClusterRef_->eta())>1.653)T0=3.1;
00114     double depth = X0 * (T0 + log(PFClusterRef_->energy()));
00115     math::XYZVector center_pos=(seedPosition_)+depth*seedAxis_;
00116     double XCentr = center_pos.x();
00117     double YCentr = center_pos.y();
00118     double XWidth = 2.59;
00119     double YWidth = 2.59;
00120     
00121     CrysX_=(PFClusterRef_->x()-XCentr)/XWidth;
00122     CrysY_=(PFClusterRef_->y()-YCentr)/YWidth;
00123   }
00124 }
00125 
00126 void PFPhotonClusters::FillClusterShape(){
00127   const std::vector< reco::PFRecHitFraction >& PFRecHits=PFClusterRef_->recHitFractions();  
00128   for(std::vector< reco::PFRecHitFraction >::const_iterator it = PFRecHits.begin(); it != PFRecHits.end(); ++it){
00129     const PFRecHitRef& RefPFRecHit = it->recHitRef();
00130     DetId id=RefPFRecHit->detId();
00131     double frac=it->fraction();
00132     float E=RefPFRecHit->energy()*frac;
00133     if(isEB_){  
00134       int deta=EBDetId::distanceEta(id,idseed_);
00135       int dphi=EBDetId::distancePhi(id,idseed_);
00136       if(abs(deta)>2 ||abs(dphi)>2)continue;
00137       //f(abs(dphi)<=2 && abs(deta)<=2){
00138       EBDetId EBidSeed=EBDetId(idseed_.rawId());
00139       EBDetId EBid=EBDetId(id.rawId());
00140       int ind1=EBidSeed.ieta()-EBid.ieta();
00141       int ind2=EBid.iphi()-EBidSeed.iphi();
00142       if(EBidSeed.ieta() * EBid.ieta() > 0){
00143         ind1=EBid.ieta()-EBidSeed.ieta();
00144       }
00145       else{ //near EB+ EB-
00146         ind1=(1-(EBidSeed.ieta()-EBid.ieta())); 
00147       }
00148       int iEta=ind1+2;
00149       int iPhi=ind2+2;
00150       //std::cout<<"IEta, IPhi "<<iEta<<", "<<iPhi<<std::endl;
00151       e5x5_[iEta][iPhi]=E;
00152     }
00153     else{
00154       int dx=EEDetId::distanceX(id,idseed_);
00155       int dy=EEDetId::distanceY(id,idseed_);
00156       if(abs(dx)>2 ||abs(dy>2))continue;
00157       EEDetId EEidSeed=EEDetId(idseed_.rawId());
00158       EEDetId EEid=EEDetId(id.rawId());
00159       int ind1=EEid.ix()-EEidSeed.ix();
00160       int ind2=EEid.iy()-EEidSeed.iy();
00161       int ix=ind1+2;
00162       int iy=ind2+2;
00163       //std::cout<<"IX, IY "<<ix<<", "<<iy<<std::endl;      
00164       e5x5_[ix][iy]=E;
00165     }
00166   }
00167 }
00168 
00169 void PFPhotonClusters::FillClusterWidth(){
00170   double numeratorEtaWidth = 0.;
00171   double numeratorPhiWidth = 0.;
00172   double numeratorEtaPhiWidth = 0.;
00173   double ClustEta=PFClusterRef_->eta();
00174   double ClustPhi=PFClusterRef_->phi();
00175   const std::vector< reco::PFRecHitFraction >& PFRecHits=PFClusterRef_->recHitFractions();  
00176   for(std::vector< reco::PFRecHitFraction >::const_iterator it = PFRecHits.begin(); it != PFRecHits.end(); ++it){
00177     const PFRecHitRef& RefPFRecHit = it->recHitRef();  
00178     float E=RefPFRecHit->energy() * it->fraction();
00179     double dEta = RefPFRecHit->position().eta() - ClustEta;     
00180     double dPhi = RefPFRecHit->position().phi() - ClustPhi;
00181     if (dPhi > + TMath::Pi()) { dPhi = TMath::TwoPi() - dPhi; }
00182     if (dPhi < - TMath::Pi()) { dPhi = TMath::TwoPi() + dPhi; }
00183     numeratorEtaWidth += E * dEta * dEta;
00184     numeratorPhiWidth += E * dPhi * dPhi;
00185     numeratorEtaPhiWidth += E * fabs(dPhi) * fabs(dEta);
00186   }
00187   double denominator=PFClusterRef_->energy();
00188   sigetaeta_ = sqrt(numeratorEtaWidth / denominator);
00189   sigphiphi_ = sqrt(numeratorPhiWidth / denominator);
00190   sigetaphi_ = sqrt(numeratorEtaPhiWidth / denominator);
00191 }