test
CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
HeavyChHiggsToTauNuSkim.cc
Go to the documentation of this file.
1 
11 
14 //#include "DataFormats/BTauReco/interface/JetTagFwd.h"
16 
17 #include <iostream>
18 
20 
21 using namespace edm;
22 using namespace std;
23 using namespace reco;
24 
25 
27 
28  // Local Debug flag
29  debug = iConfig.getParameter<bool>("DebugHeavyChHiggsToTauNuSkim");
30 
31  hltTauToken = consumes<IsolatedTauTagInfoCollection>(iConfig.getParameter<InputTag>("HLTTauCollection"));
32  jetToken = consumes<CaloJetCollection>(iConfig.getParameter<InputTag>("JetTagCollection"));
33  minNumberOfjets = iConfig.getParameter<int>("minNumberOfJets");
34  jetEtMin = iConfig.getParameter<double>("jetEtMin");
35  jetEtaMin = iConfig.getParameter<double>("jetEtaMin");
36  jetEtaMax = iConfig.getParameter<double>("jetEtaMax");
37  minDRFromTau = iConfig.getParameter<double>("minDRFromTau");
38 
39  nEvents = 0;
40  nSelectedEvents = 0;
41 
42 }
43 
44 
46  edm::LogVerbatim("HeavyChHiggsToTauNuSkim")
47  << " Number_events_read " << nEvents
48  << " Number_events_kept " << nSelectedEvents
49  << " Efficiency " << ((double)nSelectedEvents)/((double) nEvents + 0.01) << std::endl;
50 
51 }
52 
53 
55 
56  nEvents++;
57 
59  iEvent.getByToken(hltTauToken, tauTagL3Handle);
60 
61  if ( !tauTagL3Handle.isValid() ) return false;
62 
63 
64  Jet theTau;
65  double maxEt = 0;
66  if (tauTagL3Handle.isValid() ) {
67  const IsolatedTauTagInfoCollection & L3Taus = *(tauTagL3Handle.product());
68  IsolatedTauTagInfoCollection::const_iterator i;
69  for ( i = L3Taus.begin(); i != L3Taus.end(); i++ ) {
70  if (i->discriminator() == 0) continue;
71  Jet taujet = *(i->jet().get());
72  if (taujet.et() > maxEt) {
73  maxEt = taujet.et();
74  theTau = taujet;
75  }
76  }
77  }
78 
79  if (maxEt == 0) return false;
80 
81  // jets
82 
83  Handle<CaloJetCollection> jetHandle;
84  iEvent.getByToken(jetToken,jetHandle);
85 
86  if ( !jetHandle.isValid() ) return false;
87 
88  bool accepted = false;
89 
90  if (jetHandle.isValid() ) {
91  int nJets = 0;
92  const reco::CaloJetCollection & jets = *(jetHandle.product());
93  CaloJetCollection::const_iterator iJet;
94  for (iJet = jets.begin(); iJet!= jets.end(); iJet++ ) {
95  if (iJet->et() > jetEtMin &&
96  iJet->eta() > jetEtaMin &&
97  iJet->eta() < jetEtaMax ) {
98  double DR = deltaR(theTau.eta(),iJet->eta(),theTau.phi(),iJet->phi());
99  if (DR > minDRFromTau) nJets++;
100  }
101  }
102  if (nJets >= minNumberOfjets) {
103  accepted = true;
104  nSelectedEvents++;
105  }
106  }
107  return accepted;
108 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:462
Base class for all types of Jets.
Definition: Jet.h:20
virtual double phi() const final
momentum azimuthal angle
int iEvent
Definition: GenABIO.cc:230
std::vector< IsolatedTauTagInfo > IsolatedTauTagInfoCollection
vector< PseudoJet > jets
bool isValid() const
Definition: HandleBase.h:75
HeavyChHiggsToTauNuSkim(const edm::ParameterSet &)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
#define debug
Definition: HDRShower.cc:19
T const * product() const
Definition: Handle.h:81
virtual bool filter(edm::Event &, const edm::EventSetup &)
virtual double et() const final
transverse energy
virtual double eta() const final
momentum pseudorapidity
UInt_t nEvents
Definition: hcalCalib.cc:42
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects