CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
List of all members | Public Member Functions | Private Member Functions | Private Attributes
HeavyChHiggsToTauNuSkim Class Reference

#include <HeavyChHiggsToTauNuSkim.h>

Inheritance diagram for HeavyChHiggsToTauNuSkim:
edm::EDFilter edm::ProducerBase edm::EDConsumerBase edm::ProductRegistryHelper

Public Member Functions

virtual bool filter (edm::Event &, const edm::EventSetup &)
 
 HeavyChHiggsToTauNuSkim (const edm::ParameterSet &)
 
 ~HeavyChHiggsToTauNuSkim ()
 
- Public Member Functions inherited from edm::EDFilter
 EDFilter ()
 
ModuleDescription const & moduleDescription () const
 
virtual ~EDFilter ()
 
- Public Member Functions inherited from edm::ProducerBase
 ProducerBase ()
 
void registerProducts (ProducerBase *, ProductRegistry *, ModuleDescription const &)
 
std::function< void(BranchDescription
const &)> 
registrationCallback () const
 used by the fwk to register list of products More...
 
virtual ~ProducerBase ()
 
- Public Member Functions inherited from edm::EDConsumerBase
 EDConsumerBase ()
 
ProductHolderIndexAndSkipBit indexFrom (EDGetToken, BranchType, TypeID const &) const
 
void itemsMayGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
void itemsToGet (BranchType, std::vector< ProductHolderIndexAndSkipBit > &) const
 
std::vector
< ProductHolderIndexAndSkipBit >
const & 
itemsToGetFromEvent () const
 
void labelsForToken (EDGetToken iToken, Labels &oLabels) const
 
bool registeredToConsume (ProductHolderIndex, bool, BranchType) const
 
bool registeredToConsumeMany (TypeID const &, BranchType) const
 
void updateLookup (BranchType iBranchType, ProductHolderIndexHelper const &)
 
virtual ~EDConsumerBase ()
 

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
 

Additional Inherited Members

- Public Types inherited from edm::EDFilter
typedef EDFilter ModuleType
 
- Public Types inherited from edm::ProducerBase
typedef
ProductRegistryHelper::TypeLabelList 
TypeLabelList
 
- Static Public Member Functions inherited from edm::EDFilter
static const std::string & baseType ()
 
static void fillDescriptions (ConfigurationDescriptions &descriptions)
 
static void prevalidate (ConfigurationDescriptions &)
 
- Protected Member Functions inherited from edm::ProducerBase
void callWhenNewProductsRegistered (std::function< void(BranchDescription const &)> const &func)
 
- Protected Member Functions inherited from edm::EDConsumerBase
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > consumes (edm::InputTag const &tag)
 
EDGetToken consumes (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken consumes (TypeToGet const &id, edm::InputTag const &tag)
 
ConsumesCollector consumesCollector ()
 Use a ConsumesCollector to gather consumes information from helper functions. More...
 
template<typename ProductType , BranchType B = InEvent>
void consumesMany ()
 
void consumesMany (const TypeToGet &id)
 
template<BranchType B>
void consumesMany (const TypeToGet &id)
 
template<typename ProductType , BranchType B = InEvent>
EDGetTokenT< ProductType > mayConsume (edm::InputTag const &tag)
 
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 
template<BranchType B>
EDGetToken mayConsume (const TypeToGet &id, edm::InputTag const &tag)
 

Detailed Description

Filter to select events passing L1 single tau HLT tau+MET 3 offline jets

Author
Sami Lehti - HIP Helsinki

This class is an EDFilter for heavy H+->taunu events

Author
Sami Lehti - HIP Helsinki

Definition at line 30 of file HeavyChHiggsToTauNuSkim.h.

Constructor & Destructor Documentation

HeavyChHiggsToTauNuSkim::HeavyChHiggsToTauNuSkim ( const edm::ParameterSet iConfig)
explicit

Definition at line 28 of file HeavyChHiggsToTauNuSkim.cc.

References debug, edm::ParameterSet::getParameter(), and nEvents.

28  {
29 
30  // Local Debug flag
31  debug = iConfig.getParameter<bool>("DebugHeavyChHiggsToTauNuSkim");
32 
33  hltTauLabel = iConfig.getParameter<InputTag>("HLTTauCollection");
34  jetLabel = iConfig.getParameter<InputTag>("JetTagCollection");
35  minNumberOfjets = iConfig.getParameter<int>("minNumberOfJets");
36  jetEtMin = iConfig.getParameter<double>("jetEtMin");
37  jetEtaMin = iConfig.getParameter<double>("jetEtaMin");
38  jetEtaMax = iConfig.getParameter<double>("jetEtaMax");
39  minDRFromTau = iConfig.getParameter<double>("minDRFromTau");
40 
41  nEvents = 0;
42  nSelectedEvents = 0;
43 
44 }
T getParameter(std::string const &) const
HeavyChHiggsToTauNuSkim::~HeavyChHiggsToTauNuSkim ( )

Definition at line 47 of file HeavyChHiggsToTauNuSkim.cc.

References nEvents.

47  {
48  edm::LogVerbatim("HeavyChHiggsToTauNuSkim")
49  << " Number_events_read " << nEvents
50  << " Number_events_kept " << nSelectedEvents
51  << " Efficiency " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl;
52 
53 }

Member Function Documentation

double HeavyChHiggsToTauNuSkim::deltaPhi ( double  phi1,
double  phi2 
)
inlineprivate

Definition at line 39 of file HeavyChHiggsToTauNuSkim.h.

References PI.

Referenced by deltaR().

39  {
40  const double PI = 3.1415926535;
41  // in ORCA phi = [0,2pi], in TLorentzVector phi = [-pi,pi].
42  // With the conversion below deltaPhi works ok despite the
43  // 2*pi difference in phi definitions.
44  if(phi1 < 0) phi1 += 2*PI;
45  if(phi2 < 0) phi2 += 2*PI;
46 
47  double dphi = fabs(phi1-phi2);
48 
49  if(dphi > PI) dphi = 2*PI - dphi;
50  return dphi;
51  }
#define PI
double HeavyChHiggsToTauNuSkim::deltaR ( double  eta1,
double  eta2,
double  phi1,
double  phi2 
)
inlineprivate

Definition at line 53 of file HeavyChHiggsToTauNuSkim.h.

References deltaPhi(), and mathSSE::sqrt().

53  {
54  double dphi = deltaPhi(phi1,phi2);
55  double deta = fabs(eta1-eta2);
56  return sqrt(dphi*dphi + deta*deta);
57  }
double deltaPhi(double phi1, double phi2)
T sqrt(T t)
Definition: SSEVec.h:48
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(), fwrapper::jets, nEvents, reco::LeafCandidate::phi(), and edm::Handle< T >::product().

56  {
57 
58  nEvents++;
59 
61  iEvent.getByLabel(hltTauLabel, tauTagL3Handle);
62 
63  if ( !tauTagL3Handle.isValid() ) return false;
64 
65 
66  Jet theTau;
67  double maxEt = 0;
68  if (tauTagL3Handle.isValid() ) {
69  const IsolatedTauTagInfoCollection & L3Taus = *(tauTagL3Handle.product());
70  IsolatedTauTagInfoCollection::const_iterator i;
71  for ( i = L3Taus.begin(); i != L3Taus.end(); i++ ) {
72  if (i->discriminator() == 0) continue;
73  Jet taujet = *(i->jet().get());
74  if (taujet.et() > maxEt) {
75  maxEt = taujet.et();
76  theTau = taujet;
77  }
78  }
79  }
80 
81  if (maxEt == 0) return false;
82 
83  // jets
84 
85  Handle<CaloJetCollection> jetHandle;
86  iEvent.getByLabel(jetLabel,jetHandle);
87 
88  if ( !jetHandle.isValid() ) return false;
89 
90  bool accepted = false;
91 
92  if (jetHandle.isValid() ) {
93  int nJets = 0;
94  const reco::CaloJetCollection & jets = *(jetHandle.product());
95  CaloJetCollection::const_iterator iJet;
96  for (iJet = jets.begin(); iJet!= jets.end(); iJet++ ) {
97  if (iJet->et() > jetEtMin &&
98  iJet->eta() > jetEtaMin &&
99  iJet->eta() < jetEtaMax ) {
100  double DR = deltaR(theTau.eta(),iJet->eta(),theTau.phi(),iJet->phi());
101  if (DR > minDRFromTau) nJets++;
102  }
103  }
104  if (nJets >= minNumberOfjets) {
105  accepted = true;
106  nSelectedEvents++;
107  }
108  }
109  return accepted;
110 }
virtual double et() const GCC11_FINAL
transverse energy
int i
Definition: DBlmapReader.cc:9
double deltaR(double eta1, double eta2, double phi1, double phi2)
Base class for all types of Jets.
Definition: Jet.h:20
virtual float phi() const GCC11_FINAL
momentum azimuthal angle
vector< PseudoJet > jets
bool isValid() const
Definition: HandleBase.h:76
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:390
virtual float eta() const GCC11_FINAL
momentum pseudorapidity
T const * product() const
Definition: Handle.h:81
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects

Member Data Documentation

bool HeavyChHiggsToTauNuSkim::debug
private

Definition at line 59 of file HeavyChHiggsToTauNuSkim.h.

edm::InputTag HeavyChHiggsToTauNuSkim::hltTauLabel
private

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.

edm::InputTag HeavyChHiggsToTauNuSkim::jetLabel
private

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.