CMS 3D CMS Logo

Public Member Functions | Private Member Functions | Private Attributes

PFCandWithSuperClusterExtractor Class Reference

#include <PFCandWithSuperClusterExtractor.h>

Inheritance diagram for PFCandWithSuperClusterExtractor:
reco::isodeposit::IsoDepositExtractor

List of all members.

Public Member Functions

virtual reco::IsoDeposit deposit (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &muon) const
virtual reco::IsoDeposit deposit (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Candidate &cand) const
virtual void fillVetos (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::TrackCollection &cand)
 PFCandWithSuperClusterExtractor (const edm::ParameterSet &par)
 PFCandWithSuperClusterExtractor ()
virtual ~PFCandWithSuperClusterExtractor ()

Private Member Functions

reco::IsoDeposit depositFromObject (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Photon &cand) const
reco::IsoDeposit depositFromObject (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::PFCandidate &cand) const
reco::IsoDeposit depositFromObject (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::Track &cand) const
reco::IsoDeposit depositFromObject (const edm::Event &ev, const edm::EventSetup &evSetup, const reco::GsfElectron &cand) const
reco::IsoDeposit::Veto veto (const reco::IsoDeposit::Direction &dir) const

Private Attributes

std::string theDepositLabel
double theDiff_r
double theDiff_z
double theDR_Max
double theDR_Veto
bool theMissHitVetoSuperClusterMatch
edm::InputTag thePFCandTag
bool theVetoSuperClusterMatch

Detailed Description

Definition at line 20 of file PFCandWithSuperClusterExtractor.h.


Constructor & Destructor Documentation

PFCandWithSuperClusterExtractor::PFCandWithSuperClusterExtractor ( ) [inline]

Definition at line 24 of file PFCandWithSuperClusterExtractor.h.

{};
PFCandWithSuperClusterExtractor::PFCandWithSuperClusterExtractor ( const edm::ParameterSet par)

Definition at line 13 of file PFCandWithSuperClusterExtractor.cc.

                                                                                          :
  thePFCandTag(par.getParameter<edm::InputTag>("inputCandView")),
  theDepositLabel(par.getUntrackedParameter<std::string>("DepositLabel")),
  theVetoSuperClusterMatch(par.getParameter<bool>("SCMatch_Veto")),
  theMissHitVetoSuperClusterMatch(par.getParameter<bool>("MissHitSCMatch_Veto")),
  theDiff_r(par.getParameter<double>("Diff_r")),
  theDiff_z(par.getParameter<double>("Diff_z")),
  theDR_Max(par.getParameter<double>("DR_Max")),
  theDR_Veto(par.getParameter<double>("DR_Veto"))
{
  //  std::cout << " Loading PFCandWithSuperClusterExtractor "  << std::endl;
}
virtual PFCandWithSuperClusterExtractor::~PFCandWithSuperClusterExtractor ( ) [inline, virtual]

Definition at line 27 of file PFCandWithSuperClusterExtractor.h.

{}

Member Function Documentation

virtual reco::IsoDeposit PFCandWithSuperClusterExtractor::deposit ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::Track track 
) const [inline, virtual]

make single IsoDeposit based on track as input purely virtual: have to implement in concrete implementations

Implements reco::isodeposit::IsoDepositExtractor.

Definition at line 33 of file PFCandWithSuperClusterExtractor.h.

References depositFromObject().

Referenced by depositFromObject().

                                                                                                 { 
    return depositFromObject(ev, evSetup, muon);
  }
virtual reco::IsoDeposit PFCandWithSuperClusterExtractor::deposit ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::Candidate track 
) const [inline, virtual]

make single IsoDeposit based on a candidate as input purely virtual: have to implement in concrete implementations

Reimplemented from reco::isodeposit::IsoDepositExtractor.

Definition at line 38 of file PFCandWithSuperClusterExtractor.h.

References depositFromObject().

                                                                                                     { 

    const reco::Photon * myPhoton= dynamic_cast<const reco::Photon*>(&cand);
    if(myPhoton)
      return depositFromObject(ev, evSetup,*myPhoton);
    
    const reco::GsfElectron * myElectron = dynamic_cast<const reco::GsfElectron*>(&cand);
    if(myElectron)
      return depositFromObject(ev,evSetup,*myElectron);

    const reco::PFCandidate * myPFCand = dynamic_cast<const reco::PFCandidate*>(&cand);
    return depositFromObject(ev, evSetup,*myPFCand);
  }
IsoDeposit PFCandWithSuperClusterExtractor::depositFromObject ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::Track cand 
) const [private]

Definition at line 105 of file PFCandWithSuperClusterExtractor.cc.

References abs, deltaR(), deposit(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, eta(), reco::TrackBase::eta(), phi, reco::TrackBase::phi(), reco::TrackBase::pt(), theDiff_r, theDiff_z, theDR_Max, theDR_Veto, thePFCandTag, reco::TrackBase::vertex(), veto(), and reco::TrackBase::vz().

{
    reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
    IsoDeposit deposit(candDir );
    deposit.setVeto( veto(candDir) );
    deposit.addCandEnergy(cand.pt());
    Handle< PFCandidateCollection > PFCandH;
    event.getByLabel(thePFCandTag, PFCandH);

    double eta = cand.eta(), phi = cand.phi();
    reco::Particle::Point vtx = cand.vertex();
    for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
        double dR = deltaR(it->eta(), it->phi(), eta, phi);

        if ( (dR < theDR_Max) && (dR > theDR_Veto) &&
                (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
                ((it->vertex() - vtx).Rho() < theDiff_r)) {
            // ok
            reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
            deposit.addDeposit(dirTrk, it->pt());
        }
    }
    
    return deposit;
}
IsoDeposit PFCandWithSuperClusterExtractor::depositFromObject ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::Photon cand 
) const [private]

Definition at line 43 of file PFCandWithSuperClusterExtractor.cc.

References abs, deltaR(), deposit(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, eta(), reco::LeafCandidate::eta(), edm::Ref< C, T, F >::isNonnull(), reco::LeafCandidate::phi(), phi, reco::LeafCandidate::pt(), reco::Photon::superCluster(), theDiff_r, theDiff_z, theDR_Max, theDR_Veto, thePFCandTag, theVetoSuperClusterMatch, reco::LeafCandidate::vertex(), veto(), and reco::LeafCandidate::vz().

Referenced by deposit().

{
    reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
    IsoDeposit deposit(candDir );
    deposit.setVeto( veto(candDir) );
    deposit.addCandEnergy(cand.pt());

    Handle< PFCandidateCollection > PFCandH;
    event.getByLabel(thePFCandTag, PFCandH);

    double eta = cand.eta(), phi = cand.phi();
    reco::Particle::Point vtx = cand.vertex();
    for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
      double dR = deltaR(it->eta(), it->phi(), eta, phi);
      // veto SC 
      if (theVetoSuperClusterMatch && cand.superCluster().isNonnull() && it->superClusterRef().isNonnull() && cand.superCluster() == it->superClusterRef()) continue;
      if ( (dR < theDR_Max) && (dR > theDR_Veto) &&
           (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
           ((it->vertex() - vtx).Rho() < theDiff_r)) {
        // ok
        reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
        deposit.addDeposit(dirTrk, it->pt());
      }
    }
    
    return deposit;
}
IsoDeposit PFCandWithSuperClusterExtractor::depositFromObject ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::PFCandidate cand 
) const [private]

Definition at line 132 of file PFCandWithSuperClusterExtractor.cc.

References abs, deltaR(), deposit(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, eta(), reco::LeafCandidate::eta(), edm::Ref< C, T, F >::isNonnull(), reco::LeafCandidate::phi(), phi, reco::LeafCandidate::pt(), reco::PFCandidate::superClusterRef(), theDiff_r, theDiff_z, theDR_Max, theDR_Veto, thePFCandTag, theVetoSuperClusterMatch, reco::PFCandidate::vertex(), veto(), and reco::PFCandidate::vz().

{
    reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
    IsoDeposit deposit(candDir );
    deposit.setVeto( veto(candDir) );
    deposit.addCandEnergy(cand.pt());
    Handle< PFCandidateCollection > PFCandH;
    event.getByLabel(thePFCandTag, PFCandH);

    double eta = cand.eta(), phi = cand.phi();
    reco::Particle::Point vtx = cand.vertex();
    for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
      // veto SC 
      if (theVetoSuperClusterMatch && cand.superClusterRef().isNonnull() && it->superClusterRef().isNonnull() && cand.superClusterRef() == it->superClusterRef()) continue;
      double dR = deltaR(it->eta(), it->phi(), eta, phi);
      
      if ( (dR < theDR_Max) && (dR > theDR_Veto) &&
           (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
           ((it->vertex() - vtx).Rho() < theDiff_r)) {
        // ok
        reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
        deposit.addDeposit(dirTrk, it->pt());
      }
    }
    
    return deposit;
}
IsoDeposit PFCandWithSuperClusterExtractor::depositFromObject ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::GsfElectron cand 
) const [private]

Definition at line 72 of file PFCandWithSuperClusterExtractor.cc.

References abs, deltaR(), deposit(), PFRecoTauDiscriminationAgainstElectronDeadECAL_cfi::dR, eta(), reco::LeafCandidate::eta(), reco::GsfElectron::gsfTrack(), edm::Ref< C, T, F >::isNonnull(), reco::LeafCandidate::phi(), phi, reco::LeafCandidate::pt(), reco::GsfElectron::superCluster(), theDiff_r, theDiff_z, theDR_Max, theDR_Veto, theMissHitVetoSuperClusterMatch, thePFCandTag, reco::LeafCandidate::vertex(), veto(), and reco::LeafCandidate::vz().

{
    reco::isodeposit::Direction candDir(cand.eta(), cand.phi());
    IsoDeposit deposit(candDir );
    deposit.setVeto( veto(candDir) );
    deposit.addCandEnergy(cand.pt());

    Handle< PFCandidateCollection > PFCandH;
    event.getByLabel(thePFCandTag, PFCandH);

    double eta = cand.eta(), phi = cand.phi();
    reco::Particle::Point vtx = cand.vertex();
    for (PFCandidateCollection::const_iterator it = PFCandH->begin(), ed = PFCandH->end(); it != ed; ++it) {
      double dR = deltaR(it->eta(), it->phi(), eta, phi);
      // If MissHits>0 (possibly reconstructed as a photon in the PF in this case, kill the the photon if sharing the same SC)
      if (cand.gsfTrack()->trackerExpectedHitsInner().numberOfHits()>0 && theMissHitVetoSuperClusterMatch && 
          it->mva_nothing_gamma() > 0.99 && cand.superCluster().isNonnull() 
          && it->superClusterRef().isNonnull() 
          && cand.superCluster() == it->superClusterRef())   continue;
      if ( (dR < theDR_Max) && (dR > theDR_Veto) &&
           (std::abs(it->vz() - cand.vz()) < theDiff_z) &&
           ((it->vertex() - vtx).Rho() < theDiff_r)) {
        // ok
        reco::isodeposit::Direction dirTrk(it->eta(), it->phi());
        deposit.addDeposit(dirTrk, it->pt());
      }
    }
    
    return deposit;
}
virtual void PFCandWithSuperClusterExtractor::fillVetos ( const edm::Event ev,
const edm::EventSetup evSetup,
const reco::TrackCollection tracks 
) [inline, virtual]

fill vetoes: to exclude deposits at IsoDeposit creation stage check concrete extractors if it's no-op !

Implements reco::isodeposit::IsoDepositExtractor.

Definition at line 29 of file PFCandWithSuperClusterExtractor.h.

                                                                       { }
reco::IsoDeposit::Veto PFCandWithSuperClusterExtractor::veto ( const reco::IsoDeposit::Direction dir) const [private]

Member Data Documentation

Definition at line 70 of file PFCandWithSuperClusterExtractor.h.

Definition at line 73 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject().

Definition at line 74 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject().

Definition at line 75 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject().

Definition at line 76 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject(), and veto().

Definition at line 72 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject().

Definition at line 69 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject().

Definition at line 71 of file PFCandWithSuperClusterExtractor.h.

Referenced by depositFromObject().