#include <HeavyChHiggsToTauNuSkim.h>
Public Member Functions | |
virtual bool | filter (edm::Event &, const edm::EventSetup &) |
HeavyChHiggsToTauNuSkim (const edm::ParameterSet &) | |
~HeavyChHiggsToTauNuSkim () | |
Private Member Functions | |
double | deltaPhi (double phi1, double phi2) |
double | deltaR (double eta1, double eta2, double phi1, double phi2) |
Private Attributes | |
bool | debug |
edm::InputTag | hltTauLabel |
double | jetEtaMax |
double | jetEtaMin |
double | jetEtMin |
edm::InputTag | jetLabel |
double | minDRFromTau |
int | minNumberOfjets |
int | nEvents |
int | nSelectedEvents |
Filter to select events passing L1 single tau HLT tau+MET 3 offline jets
This class is an EDFilter for heavy H+->taunu events
Definition at line 30 of file HeavyChHiggsToTauNuSkim.h.
HeavyChHiggsToTauNuSkim::HeavyChHiggsToTauNuSkim | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 28 of file HeavyChHiggsToTauNuSkim.cc.
References debug, edm::ParameterSet::getParameter(), and nEvents.
{ // Local Debug flag debug = iConfig.getParameter<bool>("DebugHeavyChHiggsToTauNuSkim"); hltTauLabel = iConfig.getParameter<InputTag>("HLTTauCollection"); jetLabel = iConfig.getParameter<InputTag>("JetTagCollection"); minNumberOfjets = iConfig.getParameter<int>("minNumberOfJets"); jetEtMin = iConfig.getParameter<double>("jetEtMin"); jetEtaMin = iConfig.getParameter<double>("jetEtaMin"); jetEtaMax = iConfig.getParameter<double>("jetEtaMax"); minDRFromTau = iConfig.getParameter<double>("minDRFromTau"); nEvents = 0; nSelectedEvents = 0; }
HeavyChHiggsToTauNuSkim::~HeavyChHiggsToTauNuSkim | ( | ) |
Definition at line 47 of file HeavyChHiggsToTauNuSkim.cc.
References nEvents.
{ edm::LogVerbatim("HeavyChHiggsToTauNuSkim") << " Number_events_read " << nEvents << " Number_events_kept " << nSelectedEvents << " Efficiency " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl; }
double HeavyChHiggsToTauNuSkim::deltaPhi | ( | double | phi1, |
double | phi2 | ||
) | [inline, private] |
Definition at line 39 of file HeavyChHiggsToTauNuSkim.h.
References PI.
Referenced by deltaR().
{ const double PI = 3.1415926535; // in ORCA phi = [0,2pi], in TLorentzVector phi = [-pi,pi]. // With the conversion below deltaPhi works ok despite the // 2*pi difference in phi definitions. if(phi1 < 0) phi1 += 2*PI; if(phi2 < 0) phi2 += 2*PI; double dphi = fabs(phi1-phi2); if(dphi > PI) dphi = 2*PI - dphi; return dphi; }
double HeavyChHiggsToTauNuSkim::deltaR | ( | double | eta1, |
double | eta2, | ||
double | phi1, | ||
double | phi2 | ||
) | [inline, private] |
Definition at line 53 of file HeavyChHiggsToTauNuSkim.h.
References deltaPhi(), and mathSSE::sqrt().
bool HeavyChHiggsToTauNuSkim::filter | ( | edm::Event & | iEvent, |
const edm::EventSetup & | iSetup | ||
) | [virtual] |
Implements edm::EDFilter.
Definition at line 56 of file HeavyChHiggsToTauNuSkim.cc.
References deltaR(), reco::LeafCandidate::et(), reco::LeafCandidate::eta(), edm::Event::getByLabel(), i, edm::HandleBase::isValid(), analyzePatCleaning_cfg::jets, nEvents, reco::LeafCandidate::phi(), and edm::Handle< T >::product().
{ nEvents++; Handle<IsolatedTauTagInfoCollection> tauTagL3Handle; iEvent.getByLabel(hltTauLabel, tauTagL3Handle); if ( !tauTagL3Handle.isValid() ) return false; Jet theTau; double maxEt = 0; if (tauTagL3Handle.isValid() ) { const IsolatedTauTagInfoCollection & L3Taus = *(tauTagL3Handle.product()); IsolatedTauTagInfoCollection::const_iterator i; for ( i = L3Taus.begin(); i != L3Taus.end(); i++ ) { if (i->discriminator() == 0) continue; Jet taujet = *(i->jet().get()); if (taujet.et() > maxEt) { maxEt = taujet.et(); theTau = taujet; } } } if (maxEt == 0) return false; // jets Handle<CaloJetCollection> jetHandle; iEvent.getByLabel(jetLabel,jetHandle); if ( !jetHandle.isValid() ) return false; bool accepted = false; if (jetHandle.isValid() ) { int nJets = 0; const reco::CaloJetCollection & jets = *(jetHandle.product()); CaloJetCollection::const_iterator iJet; for (iJet = jets.begin(); iJet!= jets.end(); iJet++ ) { if (iJet->et() > jetEtMin && iJet->eta() > jetEtaMin && iJet->eta() < jetEtaMax ) { double DR = deltaR(theTau.eta(),iJet->eta(),theTau.phi(),iJet->phi()); if (DR > minDRFromTau) nJets++; } } if (nJets >= minNumberOfjets) { accepted = true; nSelectedEvents++; } } return accepted; }
bool HeavyChHiggsToTauNuSkim::debug [private] |
Definition at line 59 of file HeavyChHiggsToTauNuSkim.h.
Definition at line 61 of file HeavyChHiggsToTauNuSkim.h.
double HeavyChHiggsToTauNuSkim::jetEtaMax [private] |
Definition at line 66 of file HeavyChHiggsToTauNuSkim.h.
double HeavyChHiggsToTauNuSkim::jetEtaMin [private] |
Definition at line 65 of file HeavyChHiggsToTauNuSkim.h.
double HeavyChHiggsToTauNuSkim::jetEtMin [private] |
Definition at line 64 of file HeavyChHiggsToTauNuSkim.h.
Definition at line 62 of file HeavyChHiggsToTauNuSkim.h.
double HeavyChHiggsToTauNuSkim::minDRFromTau [private] |
Definition at line 67 of file HeavyChHiggsToTauNuSkim.h.
int HeavyChHiggsToTauNuSkim::minNumberOfjets [private] |
Definition at line 63 of file HeavyChHiggsToTauNuSkim.h.
int HeavyChHiggsToTauNuSkim::nEvents [private] |
Definition at line 69 of file HeavyChHiggsToTauNuSkim.h.
int HeavyChHiggsToTauNuSkim::nSelectedEvents [private] |
Definition at line 69 of file HeavyChHiggsToTauNuSkim.h.