CMS 3D CMS Logo

Public Member Functions | Private Attributes

reco::tau::RecoTauEnergyRecoveryPlugin Class Reference

Inheritance diagram for reco::tau::RecoTauEnergyRecoveryPlugin:
reco::tau::RecoTauModifierPlugin reco::tau::RecoTauEventHolderPlugin reco::tau::RecoTauNamedPlugin

List of all members.

Public Member Functions

virtual void beginEvent ()
void operator() (PFTau &) const
 RecoTauEnergyRecoveryPlugin (const edm::ParameterSet &)
virtual ~RecoTauEnergyRecoveryPlugin ()

Private Attributes

unsigned corrLevel_
double lev1EtaWindow_
double lev1PhiWindow_
RecoTauQualityCutsqcuts_
edm::InputTag srcVertices_
const reco::VertextheEventVertex_

Detailed Description

Definition at line 31 of file RecoTauEnergyRecoveryPlugin.cc.


Constructor & Destructor Documentation

reco::tau::RecoTauEnergyRecoveryPlugin::RecoTauEnergyRecoveryPlugin ( const edm::ParameterSet cfg) [explicit]

Definition at line 53 of file RecoTauEnergyRecoveryPlugin.cc.

References edm::ParameterSet::getParameter(), qcuts_, and srcVertices_.

  : RecoTauModifierPlugin(cfg),
    theEventVertex_(0),
    corrLevel_(cfg.getParameter<unsigned>("corrLevel")),
    lev1PhiWindow_(cfg.getParameter<double>("lev1PhiWindow")),
    lev1EtaWindow_(cfg.getParameter<double>("lev1EtaWindow"))
{
  edm::ParameterSet cfgQualityCuts = cfg.getParameter<edm::ParameterSet>("qualityCuts");
  srcVertices_ = cfgQualityCuts.getParameter<edm::InputTag>("primaryVertexSrc");
  qcuts_ = new RecoTauQualityCuts(cfgQualityCuts.getParameter<edm::ParameterSet>("signalQualityCuts"));
}
reco::tau::RecoTauEnergyRecoveryPlugin::~RecoTauEnergyRecoveryPlugin ( ) [virtual]

Definition at line 65 of file RecoTauEnergyRecoveryPlugin.cc.

References qcuts_.

{
  delete qcuts_;
}

Member Function Documentation

void reco::tau::RecoTauEnergyRecoveryPlugin::beginEvent ( ) [virtual]
void reco::tau::RecoTauEnergyRecoveryPlugin::operator() ( PFTau tau) const [virtual]

Implements reco::tau::RecoTauModifierPlugin.

Definition at line 108 of file RecoTauEnergyRecoveryPlugin.cc.

References corrLevel_, reco::LeafCandidate::energy(), reco::LeafCandidate::eta(), reco::tau::RecoTauQualityCuts::filter(), reco::tau::getTrack(), reco::PFCandidate::h0, edm::RefToBase< T >::isNonnull(), reco::tau::isTauSignalPFCandidate(), reco::PFTau::jetRef(), lev1EtaWindow_, lev1PhiWindow_, max(), reco::LeafCandidate::p(), reco::LeafCandidate::phi(), reco::Vertex::position(), reco::LeafCandidate::px(), reco::LeafCandidate::py(), reco::LeafCandidate::pz(), qcuts_, reco::BaseTau::setalternatLorentzVect(), mathSSE::sqrt(), reco::tau::square(), theEventVertex_, and listHistos::trackPt.

{
  double tauEnergyCorr = 0.;

  if ( corrLevel_ & 1 && theEventVertex_ ) {

    bool needsCorrLevel1 = false;

    std::vector<reco::PFCandidatePtr> pfJetConstituents = tau.jetRef()->getPFConstituents();
    for ( std::vector<reco::PFCandidatePtr>::const_iterator pfJetConstituent = pfJetConstituents.begin();
          pfJetConstituent != pfJetConstituents.end(); ++pfJetConstituent ) {
      const reco::TrackBaseRef track = getTrack(**pfJetConstituent);
      if ( track.isNonnull() ) {
        double dIP = fabs(track->dxy(theEventVertex_->position()));
        double dZ = fabs(track->dz(theEventVertex_->position()));
      if ( track->pt() > 2.0 && fabs(tau.eta() - (*pfJetConstituent)->eta()) < lev1EtaWindow_ &&
           !isTauSignalPFCandidate(tau, *pfJetConstituent) && (dZ < 0.2 || dIP > 0.10) ) needsCorrLevel1 = true;
      }
    }

    if ( needsCorrLevel1 ) {
      std::vector<reco::PFCandidatePtr> pfJetConstituents = tau.jetRef()->getPFConstituents();
      for ( std::vector<reco::PFCandidatePtr>::const_iterator pfJetConstituent = pfJetConstituents.begin();
            pfJetConstituent != pfJetConstituents.end(); ++pfJetConstituent ) {
        if ( fabs(tau.eta() - (*pfJetConstituent)->eta()) < lev1EtaWindow_ &&
             fabs(tau.phi() - (*pfJetConstituent)->phi()) < lev1PhiWindow_ ) {
          if ( (*pfJetConstituent)->particleId() == reco::PFCandidate::h0 ) {
            tauEnergyCorr += (*pfJetConstituent)->energy();
          } else {
            if ( !isTauSignalPFCandidate(tau, *pfJetConstituent) ) {
              double caloEn = (*pfJetConstituent)->ecalEnergy() + (*pfJetConstituent)->hcalEnergy();
              tauEnergyCorr += caloEn;
            }
          }
        }
      }
    }
  }

  if ( corrLevel_ & 2 && theEventVertex_ ) {

    double leadTrackMom = 0.;
    double leadTrackMomErr = 0.;
    double jetCaloEn = 0.;

    std::vector<reco::PFCandidatePtr> pfJetConstituents = tau.jetRef()->getPFConstituents();
    for ( std::vector<reco::PFCandidatePtr>::const_iterator pfJetConstituent = pfJetConstituents.begin();
          pfJetConstituent != pfJetConstituents.end(); ++pfJetConstituent ) {
      const reco::TrackBaseRef track = getTrack(**pfJetConstituent);
      if ( track.isNonnull() ) {
        double trackPt = track->pt();
        double trackPtErr = track->ptError();
        if ( qcuts_->filter(**pfJetConstituent) &&
             trackPtErr < (0.20*trackPt) && track->normalizedChi2() < 5.0 && track->hitPattern().numberOfValidPixelHits() >= 1 &&
             (trackPt - 3.*trackPtErr) > (*pfJetConstituent)->pt() && trackPt < (3.*tau.jetRef()->pt()) ) {
          if ( track->p() > leadTrackMom ) {
            leadTrackMom = track->p();
            leadTrackMomErr = leadTrackMom*(trackPtErr/trackPt);
          }
        }
      }

      double caloEn = (*pfJetConstituent)->ecalEnergy() + (*pfJetConstituent)->hcalEnergy();
      jetCaloEn += caloEn;
    }

    if ( leadTrackMom > tau.p() ) {
      const double chargedPionMass = 0.13957; // GeV
      double leadTrackEn = sqrt(square(leadTrackMom) + square(chargedPionMass));
      double jetCaloEnErr = 1.00*sqrt(std::max(jetCaloEn, leadTrackEn));
      double combEn = ((1./square(jetCaloEnErr))*jetCaloEn + (1./square(leadTrackMomErr))*leadTrackEn)/
                      ((1./square(jetCaloEnErr)) + (1./square(leadTrackMomErr)));
      tauEnergyCorr = std::max(tauEnergyCorr, combEn - tau.energy());
    }
  }

  if ( tau.energy() > 0. ) {
    double tauEnergy_corrected = tau.energy() + tauEnergyCorr;
    double scaleFactor = tauEnergy_corrected/tau.energy();
    double tauPx_corrected = scaleFactor*tau.px();
    double tauPy_corrected = scaleFactor*tau.py();
    double tauPz_corrected = scaleFactor*tau.pz();
    tau.setalternatLorentzVect(reco::Candidate::LorentzVector(tauPx_corrected, tauPy_corrected, tauPz_corrected, tauEnergy_corrected));
  }
}

Member Data Documentation

Definition at line 47 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by operator()().

Definition at line 50 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by operator()().

Definition at line 49 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by operator()().

Definition at line 42 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by beginEvent(), and RecoTauEnergyRecoveryPlugin().

Definition at line 43 of file RecoTauEnergyRecoveryPlugin.cc.

Referenced by beginEvent(), and operator()().