CMS 3D CMS Logo

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