CMS 3D CMS Logo

Public Member Functions | Public Attributes | Private Attributes

HiGammaJetSignalDef Class Reference

#include <HiPhotonType.h>

List of all members.

Public Member Functions

double getDeltaPhi (const reco::Candidate &track1, const reco::Candidate &track2)
double getDeltaR (const reco::Candidate &track1, const reco::Candidate &track2)
 HiGammaJetSignalDef (const reco::GenParticleCollection *sigPartic)
 HiGammaJetSignalDef ()
bool IsIsolated (const reco::GenParticle &pp)
bool IsIsolatedJP (const reco::GenParticle &pp)
bool IsIsolatedPP (const reco::GenParticle &pp)

Public Attributes

double PI

Private Attributes

const reco::GenParticleCollectionfSigParticles

Detailed Description

Definition at line 23 of file HiPhotonType.h.


Constructor & Destructor Documentation

HiGammaJetSignalDef::HiGammaJetSignalDef ( )

Definition at line 63 of file HiPhotonType.cc.

References PI.

                                          :
  fSigParticles(0)
{   PI = 3.141592653589793238462643383279502884197169399375105820974945;
}
HiGammaJetSignalDef::HiGammaJetSignalDef ( const reco::GenParticleCollection sigPartic)

Definition at line 67 of file HiPhotonType.cc.

References fSigParticles, and PI.

                                                                                        :
  fSigParticles(0)
{
  using namespace std;

  fSigParticles =  sigParticles;
  PI = 3.141592653589793238462643383279502884197169399375105820974945;

}

Member Function Documentation

double HiGammaJetSignalDef::getDeltaPhi ( const reco::Candidate track1,
const reco::Candidate track2 
)

Definition at line 230 of file HiPhotonType.cc.

References dPhi(), reco::Candidate::phi(), and PI.

{
  double dPhi = track1.phi() - track2.phi();

  while(dPhi >= PI)       dPhi -= (2.0*PI);
  while(dPhi < (-1.0*PI)) dPhi += (2.0*PI);

  return dPhi;
}
double HiGammaJetSignalDef::getDeltaR ( const reco::Candidate track1,
const reco::Candidate track2 
)

Definition at line 219 of file HiPhotonType.cc.

References dPhi(), reco::Candidate::eta(), reco::Candidate::phi(), PI, and mathSSE::sqrt().

Referenced by IsIsolated(), IsIsolatedJP(), and IsIsolatedPP().

{
  double dEta = track1.eta() - track2.eta();
  double dPhi = track1.phi() - track2.phi();

  while(dPhi >= PI)       dPhi -= (2.0*PI);
  while(dPhi < (-1.0*PI)) dPhi += (2.0*PI);

  return sqrt(dEta*dEta+dPhi*dPhi);
}
bool HiGammaJetSignalDef::IsIsolated ( const reco::GenParticle pp)

Definition at line 77 of file HiPhotonType.cc.

References abs, reco::GenParticle::collisionId(), reco::LeafCandidate::et(), fSigParticles, getDeltaR(), i, AlCaHLTBitMon_ParallelJobs::p, reco::LeafCandidate::pdgId(), reco::LeafCandidate::pt(), dt_dqm_sourceclient_common_cff::reco, and reco::LeafCandidate::status().

{
  using namespace std;
  using namespace edm;
  using namespace reco;
  // Check if a given particle is isolated.

  double  EtCone = 0;
  double  EtPhoton = 0;
  double  PtMaxHadron = 0;
  double  cone = 0.5;
  const int maxindex = (int)fSigParticles->size();
  for(int i=0; i < maxindex ; ++i) {

    const GenParticle &p=(*fSigParticles)[i];

    if(p.status()!=1)
      continue;
    if (p.collisionId() != pp.collisionId())
      continue;

    int apid= abs(p.pdgId());
    if(apid>11 &&  apid<20)
      continue; //get rid of muons and neutrinos

    if(getDeltaR(p,pp)>cone)
      continue;

    double et=p.et();
    double pt = p.pt();
    EtCone+=et;
    bool isHadron = false;
    if ( fabs(pp.pdgId())>100 && fabs(pp.pdgId()) != 310)
      isHadron = true;

    if(apid>100 && apid!=310 && pt>PtMaxHadron)
      {
        if ( (isHadron == true) && (fabs(pp.pt()-pt)<0.001) && (pp.pdgId()==p.pdgId()) )
          continue;
        else
          PtMaxHadron=pt;
      }

  }
  EtPhoton = pp.et();

  // isolation cuts
  if(EtCone-EtPhoton > 5+EtPhoton/20)
    return kFALSE;
  if(PtMaxHadron > 4.5+EtPhoton/40)
    return kFALSE;

  return kTRUE;

}
bool HiGammaJetSignalDef::IsIsolatedJP ( const reco::GenParticle pp)

Definition at line 179 of file HiPhotonType.cc.

References reco::GenParticle::collisionId(), reco::LeafCandidate::et(), fSigParticles, getDeltaR(), i, AlCaHLTBitMon_ParallelJobs::p, dt_dqm_sourceclient_common_cff::reco, and reco::LeafCandidate::status().

{
  using namespace std;
  using namespace edm;
  using namespace reco;

  // Check if a given particle is isolated. 

    double  EtCone = 0;
    double  EtPhoton = 0;
    double  cone = 0.4;

    const int maxindex = (int)fSigParticles->size();
    for(int i=0; i < maxindex ; ++i) {

      const GenParticle &p=(*fSigParticles)[i];

      if(p.status()!=1)
        continue;
      if (p.collisionId() != pp.collisionId())
        continue;

      if(getDeltaR(p,pp)>cone)
        continue;

      double et=p.et();
      //  double pt = p.pt();           
      EtCone+=et;

    }
    EtPhoton = pp.et();

    // isolation cuts   

      if(EtCone-EtPhoton > 5)
        return kFALSE;
      return kTRUE;

}
bool HiGammaJetSignalDef::IsIsolatedPP ( const reco::GenParticle pp)

Definition at line 133 of file HiPhotonType.cc.

References reco::GenParticle::collisionId(), reco::LeafCandidate::et(), fSigParticles, getDeltaR(), i, AlCaHLTBitMon_ParallelJobs::p, dt_dqm_sourceclient_common_cff::reco, and reco::LeafCandidate::status().

{
  using namespace std;
  using namespace edm;
  using namespace reco;

  // Check if a given particle is isolated. 

    double  EtCone = 0;
    double  EtPhoton = 0;
    double  cone = 0.4;

    const int maxindex = (int)fSigParticles->size();
    for(int i=0; i < maxindex ; ++i) {

      const GenParticle &p=(*fSigParticles)[i];

      if(p.status()!=1)
        continue;
      if (p.collisionId() != pp.collisionId())
        continue;

      //int apid= abs(p.pdgId());
      //  if(apid>11 &&  apid<20)
      //  continue; //get rid of muons and neutrinos 

        if(getDeltaR(p,pp)>cone)
          continue;
        double et=p.et();
        //  double pt = p.pt();
        EtCone+=et;

    }

    EtPhoton = pp.et();

    // isolation cuts 

      if(EtCone-EtPhoton > 5)
        return kFALSE;
      //if(PtMaxHadron > 4.5+EtPhoton/40)
      //  return kFALSE;

      return kTRUE;

}

Member Data Documentation

Definition at line 40 of file HiPhotonType.h.

Referenced by HiGammaJetSignalDef(), IsIsolated(), IsIsolatedJP(), and IsIsolatedPP().

Definition at line 37 of file HiPhotonType.h.

Referenced by getDeltaPhi(), getDeltaR(), and HiGammaJetSignalDef().