#include <HiggsAnalysis/Skimming/interface/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 |
InputTag | hltTauLabel |
double | jetEtaMax |
double | jetEtaMin |
double | jetEtMin |
InputTag | jetLabel |
double | minDRFromTau |
int | minNumberOfjets |
int | nEvents |
int | nSelectedEvents |
This class is an EDFilter for heavy H+->taunu events.
Definition at line 33 of file HeavyChHiggsToTauNuSkim.h.
HeavyChHiggsToTauNuSkim::HeavyChHiggsToTauNuSkim | ( | const edm::ParameterSet & | iConfig | ) | [explicit] |
Definition at line 28 of file HeavyChHiggsToTauNuSkim.cc.
References debug, edm::ParameterSet::getParameter(), hltTauLabel, jetEtaMax, jetEtaMin, jetEtMin, jetLabel, minDRFromTau, minNumberOfjets, nEvents, and nSelectedEvents.
00028 { 00029 00030 // Local Debug flag 00031 debug = iConfig.getParameter<bool>("DebugHeavyChHiggsToTauNuSkim"); 00032 00033 hltTauLabel = iConfig.getParameter<InputTag>("HLTTauCollection"); 00034 jetLabel = iConfig.getParameter<InputTag>("JetTagCollection"); 00035 minNumberOfjets = iConfig.getParameter<int>("minNumberOfJets"); 00036 jetEtMin = iConfig.getParameter<double>("jetEtMin"); 00037 jetEtaMin = iConfig.getParameter<double>("jetEtaMin"); 00038 jetEtaMax = iConfig.getParameter<double>("jetEtaMax"); 00039 minDRFromTau = iConfig.getParameter<double>("minDRFromTau"); 00040 00041 nEvents = 0; 00042 nSelectedEvents = 0; 00043 00044 }
HeavyChHiggsToTauNuSkim::~HeavyChHiggsToTauNuSkim | ( | ) |
Definition at line 47 of file HeavyChHiggsToTauNuSkim.cc.
References lat::endl(), nEvents, and nSelectedEvents.
00047 { 00048 edm::LogVerbatim("HeavyChHiggsToTauNuSkim") 00049 << " Number_events_read " << nEvents 00050 << " Number_events_kept " << nSelectedEvents 00051 << " Efficiency " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl; 00052 00053 }
double HeavyChHiggsToTauNuSkim::deltaPhi | ( | double | phi1, | |
double | phi2 | |||
) | [inline, private] |
Definition at line 42 of file HeavyChHiggsToTauNuSkim.h.
References PI.
Referenced by deltaR().
00042 { 00043 const double PI = 3.1415926535; 00044 // in ORCA phi = [0,2pi], in TLorentzVector phi = [-pi,pi]. 00045 // With the conversion below deltaPhi works ok despite the 00046 // 2*pi difference in phi definitions. 00047 if(phi1 < 0) phi1 += 2*PI; 00048 if(phi2 < 0) phi2 += 2*PI; 00049 00050 double dphi = fabs(phi1-phi2); 00051 00052 if(dphi > PI) dphi = 2*PI - dphi; 00053 return dphi; 00054 }
double HeavyChHiggsToTauNuSkim::deltaR | ( | double | eta1, | |
double | eta2, | |||
double | phi1, | |||
double | phi2 | |||
) | [inline, private] |
Definition at line 56 of file HeavyChHiggsToTauNuSkim.h.
References deltaPhi(), and funct::sqrt().
Referenced by filter().
00056 { 00057 double dphi = deltaPhi(phi1,phi2); 00058 double deta = fabs(eta1-eta2); 00059 return sqrt(dphi*dphi + deta*deta); 00060 }
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::Particle::et(), reco::Particle::eta(), edm::Event::getByLabel(), hltTauLabel, i, edm::Handle< T >::isValid(), jetEtaMax, jetEtaMin, jetEtMin, jetLabel, pfTauBenchmarkGeneric_cfi::jets, minDRFromTau, minNumberOfjets, nEvents, nSelectedEvents, reco::Particle::phi(), and edm::Handle< T >::product().
00056 { 00057 00058 nEvents++; 00059 00060 Handle<IsolatedTauTagInfoCollection> tauTagL3Handle; 00061 iEvent.getByLabel(hltTauLabel, tauTagL3Handle); 00062 00063 if ( !tauTagL3Handle.isValid() ) return false; 00064 00065 00066 Jet theTau; 00067 double maxEt = 0; 00068 if (tauTagL3Handle.isValid() ) { 00069 const IsolatedTauTagInfoCollection & L3Taus = *(tauTagL3Handle.product()); 00070 IsolatedTauTagInfoCollection::const_iterator i; 00071 for ( i = L3Taus.begin(); i != L3Taus.end(); i++ ) { 00072 if (i->discriminator() == 0) continue; 00073 Jet taujet = *(i->jet().get()); 00074 if (taujet.et() > maxEt) { 00075 maxEt = taujet.et(); 00076 theTau = taujet; 00077 } 00078 } 00079 } 00080 00081 if (maxEt == 0) return false; 00082 00083 // jets 00084 00085 Handle<CaloJetCollection> jetHandle; 00086 iEvent.getByLabel(jetLabel,jetHandle); 00087 00088 if ( !jetHandle.isValid() ) return false; 00089 00090 bool accepted = false; 00091 00092 if (jetHandle.isValid() ) { 00093 int nJets = 0; 00094 const reco::CaloJetCollection & jets = *(jetHandle.product()); 00095 CaloJetCollection::const_iterator iJet; 00096 for (iJet = jets.begin(); iJet!= jets.end(); iJet++ ) { 00097 if (iJet->et() > jetEtMin && 00098 iJet->eta() > jetEtaMin && 00099 iJet->eta() < jetEtaMax ) { 00100 double DR = deltaR(theTau.eta(),iJet->eta(),theTau.phi(),iJet->phi()); 00101 if (DR > minDRFromTau) nJets++; 00102 } 00103 } 00104 if (nJets >= minNumberOfjets) { 00105 accepted = true; 00106 nSelectedEvents++; 00107 } 00108 } 00109 return accepted; 00110 }
bool HeavyChHiggsToTauNuSkim::debug [private] |
InputTag HeavyChHiggsToTauNuSkim::hltTauLabel [private] |
Definition at line 64 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
double HeavyChHiggsToTauNuSkim::jetEtaMax [private] |
Definition at line 69 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
double HeavyChHiggsToTauNuSkim::jetEtaMin [private] |
Definition at line 68 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
double HeavyChHiggsToTauNuSkim::jetEtMin [private] |
Definition at line 67 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
InputTag HeavyChHiggsToTauNuSkim::jetLabel [private] |
Definition at line 65 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
double HeavyChHiggsToTauNuSkim::minDRFromTau [private] |
Definition at line 70 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
int HeavyChHiggsToTauNuSkim::minNumberOfjets [private] |
Definition at line 66 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), and HeavyChHiggsToTauNuSkim().
int HeavyChHiggsToTauNuSkim::nEvents [private] |
Definition at line 72 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), HeavyChHiggsToTauNuSkim(), and ~HeavyChHiggsToTauNuSkim().
int HeavyChHiggsToTauNuSkim::nSelectedEvents [private] |
Definition at line 72 of file HeavyChHiggsToTauNuSkim.h.
Referenced by filter(), HeavyChHiggsToTauNuSkim(), and ~HeavyChHiggsToTauNuSkim().