CMS 3D CMS Logo

EgammaHLTEcalIsolation Class Reference

Description: sum Et of all island basic clusters in cone around candidate. More...

#include <RecoEgamma/EgammaHLTAlgos/interface/EgammaHLTEcalIsolation.h>

List of all members.

Public Member Functions

 EgammaHLTEcalIsolation (double egEcalIso_EtMin, double egEcalIso_ConeSize, int SC_algo_type)
float getConeSize ()
 Get isolation cone size.
float getetMin ()
 Get Et cut for ecal hits.
float isolPtSum (const reco::RecoCandidate *recocandidate, const std::vector< const reco::SuperCluster * > sclusters, const std::vector< const reco::BasicCluster * > bclusters)

Private Attributes

int algoType_
double conesize
double etMin


Detailed Description

Description: sum Et of all island basic clusters in cone around candidate.

Definition at line 27 of file EgammaHLTEcalIsolation.h.


Constructor & Destructor Documentation

EgammaHLTEcalIsolation::EgammaHLTEcalIsolation ( double  egEcalIso_EtMin,
double  egEcalIso_ConeSize,
int  SC_algo_type 
) [inline]

Definition at line 33 of file EgammaHLTEcalIsolation.h.

00033                                                                                               : 
00034     etMin(egEcalIso_EtMin), conesize(egEcalIso_ConeSize), algoType_(SC_algo_type) {
00035       /*
00036       std::cout << "EgammaHLTEcalIsolation instance:"
00037       << " ptMin=" << etMin
00038       << " conesize=" << conesize
00039       << std::endl;
00040       */
00041 
00042     }


Member Function Documentation

float EgammaHLTEcalIsolation::getConeSize (  )  [inline]

Get isolation cone size.

Definition at line 51 of file EgammaHLTEcalIsolation.h.

References conesize.

00051 { return conesize; }

float EgammaHLTEcalIsolation::getetMin (  )  [inline]

Get Et cut for ecal hits.

Definition at line 49 of file EgammaHLTEcalIsolation.h.

References etMin.

00049 { return etMin; }

float EgammaHLTEcalIsolation::isolPtSum ( const reco::RecoCandidate recocandidate,
const std::vector< const reco::SuperCluster * >  sclusters,
const std::vector< const reco::BasicCluster * >  bclusters 
)

Definition at line 23 of file EgammaHLTEcalIsolation.cc.

References reco::BasicCluster::algo(), algoType_, reco::BasicCluster::chi2(), reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), conesize, reco::CaloCluster::energy(), reco::CaloCluster::eta(), etMin, funct::exp(), reco::CaloCluster::phi(), PI, reco::SuperCluster::seed(), funct::sin(), funct::sqrt(), reco::RecoCandidate::superCluster(), and TWOPI.

Referenced by EgammaHLTEcalIsolationProducersRegional::produce().

00025                                                                                            {
00026 
00027   float ecalIsol=0.;
00028   
00029   float candSCphi = recocandidate->superCluster()->phi();
00030   float candSCeta = recocandidate->superCluster()->eta();
00031 
00032   
00033   // match the photon hybrid supercluster with those with Algo==0 (island)
00034   float delta1=1000.;
00035   float deltacur=1000.;
00036   const reco::SuperCluster *matchedsupercluster=0;
00037   bool MATCHEDSC = false;
00038 
00039 
00040 
00041   for(std::vector<const reco::SuperCluster*>::const_iterator scItr = sclusters.begin(); scItr != sclusters.end(); ++scItr){
00042     
00043     
00044     const reco::SuperCluster *supercluster = *scItr;
00045     
00046     float SCphi = supercluster->phi();
00047     float SCeta = supercluster->eta();
00048    
00049     if(supercluster->seed()->algo() == algoType_){
00050       float deltaphi;
00051       if(candSCphi<0) candSCphi+=TWOPI;
00052       if(SCphi<0) SCphi+=TWOPI;
00053       deltaphi=fabs(candSCphi-SCphi);
00054       if(deltaphi>TWOPI) deltaphi-=TWOPI;
00055       if(deltaphi>PI) deltaphi=TWOPI-deltaphi;
00056       float deltaeta=fabs(SCeta-candSCeta);
00057       deltacur = sqrt(deltaphi*deltaphi+ deltaeta*deltaeta);
00058       
00059       if (deltacur < delta1) {
00060         delta1=deltacur;
00061         matchedsupercluster = supercluster;
00062         MATCHEDSC = true;
00063       }
00064     }
00065   }
00066 
00067   const reco::BasicCluster *cluster= 0;
00068 
00069   //loop over basic clusters
00070   for(std::vector<const reco::BasicCluster*>::const_iterator cItr = bclusters.begin(); cItr != bclusters.end(); ++cItr){
00071  
00072     cluster = *cItr;
00073     float ebc_bcchi2 = cluster->chi2();
00074     int   ebc_bcalgo = cluster->algo();
00075     float ebc_bce    = cluster->energy();
00076     float ebc_bceta  = cluster->eta();
00077     float ebc_bcphi  = cluster->phi();
00078     float ebc_bcet   = ebc_bce*sin(2*atan(exp(ebc_bceta)));
00079     float newDelta;
00080 
00081 
00082     if (ebc_bcet > etMin && ebc_bcalgo == algoType_ ) {
00083       if (ebc_bcchi2 < 30.) {
00084         
00085         if(MATCHEDSC){
00086           bool inSuperCluster = false;
00087 
00088           reco::basicCluster_iterator theEclust = matchedsupercluster->clustersBegin();
00089           // loop over the basic clusters of the matched supercluster
00090           for(;theEclust != matchedsupercluster->clustersEnd();
00091               theEclust++) {
00092             if (&(**theEclust) ==  cluster) inSuperCluster = true;
00093           }
00094           if (!inSuperCluster) {
00095             float deltaphi;
00096             if(ebc_bcphi<0) ebc_bcphi+=TWOPI;
00097             if(candSCphi<0) candSCphi+=TWOPI;
00098             deltaphi=fabs(ebc_bcphi-candSCphi);
00099             if(deltaphi>TWOPI) deltaphi-=TWOPI;
00100             if(deltaphi>PI) deltaphi=TWOPI-deltaphi;
00101             float deltaeta=fabs(ebc_bceta-candSCeta);
00102             newDelta= sqrt(deltaphi*deltaphi+ deltaeta*deltaeta);
00103             if(newDelta < conesize) {
00104               ecalIsol+=ebc_bcet;
00105             }
00106           }
00107         }
00108       } // matches ebc_bcchi2
00109     } // matches ebc_bcet && ebc_bcalgo
00110 
00111   }
00112 
00113 
00114  return ecalIsol;
00115 
00116 }


Member Data Documentation

int EgammaHLTEcalIsolation::algoType_ [private]

Definition at line 61 of file EgammaHLTEcalIsolation.h.

Referenced by isolPtSum().

double EgammaHLTEcalIsolation::conesize [private]

Definition at line 60 of file EgammaHLTEcalIsolation.h.

Referenced by getConeSize(), and isolPtSum().

double EgammaHLTEcalIsolation::etMin [private]

Definition at line 59 of file EgammaHLTEcalIsolation.h.

Referenced by getetMin(), and isolPtSum().


The documentation for this class was generated from the following files:
Generated on Tue Jun 9 18:20:10 2009 for CMSSW by  doxygen 1.5.4