CMS 3D CMS Logo

Public Member Functions | Private Attributes

EgammaHLTEcalIsolation Class Reference

#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.

                                                                                              : 
    etMin(egEcalIso_EtMin), conesize(egEcalIso_ConeSize), algoType_(SC_algo_type) {
      /*
      std::cout << "EgammaHLTEcalIsolation instance:"
      << " ptMin=" << etMin
      << " conesize=" << conesize
      << std::endl;
      */

    }

Member Function Documentation

float EgammaHLTEcalIsolation::getConeSize ( ) [inline]

Get isolation cone size.

Definition at line 51 of file EgammaHLTEcalIsolation.h.

References conesize.

{ return conesize; }
float EgammaHLTEcalIsolation::getetMin ( ) [inline]

Get Et cut for ecal hits.

Definition at line 49 of file EgammaHLTEcalIsolation.h.

References etMin.

{ 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 algoType_, reco::SuperCluster::clustersBegin(), reco::SuperCluster::clustersEnd(), conesize, reco::CaloCluster::eta(), etMin, create_public_lumi_plots::exp, reco::CaloCluster::phi(), PI, reco::SuperCluster::seed(), funct::sin(), mathSSE::sqrt(), reco::RecoCandidate::superCluster(), and TWOPI.

Referenced by EgammaHLTEcalIsolationProducersRegional::produce().

                                                                                           {

  float ecalIsol=0.;
  
  float candSCphi = recocandidate->superCluster()->phi();
  float candSCeta = recocandidate->superCluster()->eta();

  
  // match the photon hybrid supercluster with those with Algo==0 (island)
  float delta1=1000.;
  float deltacur=1000.;
  const reco::SuperCluster *matchedsupercluster=0;
  bool MATCHEDSC = false;



  for(std::vector<const reco::SuperCluster*>::const_iterator scItr = sclusters.begin(); scItr != sclusters.end(); ++scItr){
    
    
    const reco::SuperCluster *supercluster = *scItr;
    
    float SCphi = supercluster->phi();
    float SCeta = supercluster->eta();
   
    if(supercluster->seed()->algo() == algoType_){
      float deltaphi;
      if(candSCphi<0) candSCphi+=TWOPI;
      if(SCphi<0) SCphi+=TWOPI;
      deltaphi=fabs(candSCphi-SCphi);
      if(deltaphi>TWOPI) deltaphi-=TWOPI;
      if(deltaphi>PI) deltaphi=TWOPI-deltaphi;
      float deltaeta=fabs(SCeta-candSCeta);
      deltacur = sqrt(deltaphi*deltaphi+ deltaeta*deltaeta);
      
      if (deltacur < delta1) {
        delta1=deltacur;
        matchedsupercluster = supercluster;
        MATCHEDSC = true;
      }
    }
  }

  const reco::BasicCluster *cluster= 0;

  //loop over basic clusters
  for(std::vector<const reco::BasicCluster*>::const_iterator cItr = bclusters.begin(); cItr != bclusters.end(); ++cItr){
 
    cluster = *cItr;
//    float ebc_bcchi2 = cluster->chi2(); //chi2 for SC was useless and it is removed in 31x
    int   ebc_bcalgo = cluster->algo();
    float ebc_bce    = cluster->energy();
    float ebc_bceta  = cluster->eta();
    float ebc_bcphi  = cluster->phi();
    float ebc_bcet   = ebc_bce*sin(2*atan(exp(ebc_bceta)));
    float newDelta;


    if (ebc_bcet > etMin && ebc_bcalgo == algoType_ ) {
      //  if (ebc_bcchi2 < 30.) {
        
        if(MATCHEDSC){
          bool inSuperCluster = false;

          reco::CaloCluster_iterator theEclust = matchedsupercluster->clustersBegin();
          // loop over the basic clusters of the matched supercluster
          for(;theEclust != matchedsupercluster->clustersEnd();
              theEclust++) {
            if (&(**theEclust) ==  cluster) inSuperCluster = true;
          }
          if (!inSuperCluster) {
            float deltaphi;
            if(ebc_bcphi<0) ebc_bcphi+=TWOPI;
            if(candSCphi<0) candSCphi+=TWOPI;
            deltaphi=fabs(ebc_bcphi-candSCphi);
            if(deltaphi>TWOPI) deltaphi-=TWOPI;
            if(deltaphi>PI) deltaphi=TWOPI-deltaphi;
            float deltaeta=fabs(ebc_bceta-candSCeta);
            newDelta= sqrt(deltaphi*deltaphi+ deltaeta*deltaeta);
            if(newDelta < conesize) {
              ecalIsol+=ebc_bcet;
            }
          }
        }
        //  } // matches ebc_bcchi2
    } // matches ebc_bcet && ebc_bcalgo

  }


 return ecalIsol;

}

Member Data Documentation

Definition at line 61 of file EgammaHLTEcalIsolation.h.

Referenced by isolPtSum().

Definition at line 60 of file EgammaHLTEcalIsolation.h.

Referenced by getConeSize(), and isolPtSum().

Definition at line 59 of file EgammaHLTEcalIsolation.h.

Referenced by getetMin(), and isolPtSum().