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 
15 //#include "DataFormats/BTauReco/interface/JetTagFwd.h"
18 
19 #include <iostream>
20 
22 
23 using namespace edm;
24 using namespace std;
25 using namespace reco;
26 
27 
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 }
45 
46 
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 }
54 
55 
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 }
T getParameter(std::string const &) const
int i
Definition: DBlmapReader.cc:9
virtual double et() const
transverse energy
Base class for all types of Jets.
Definition: Jet.h:21
virtual double eta() const
momentum pseudorapidity
int iEvent
Definition: GenABIO.cc:243
std::vector< IsolatedTauTagInfo > IsolatedTauTagInfoCollection
vector< PseudoJet > jets
bool isValid() const
Definition: HandleBase.h:76
HeavyChHiggsToTauNuSkim(const edm::ParameterSet &)
bool getByLabel(InputTag const &tag, Handle< PROD > &result) const
Definition: Event.h:356
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
virtual bool filter(edm::Event &, const edm::EventSetup &)
T const * product() const
Definition: Handle.h:74
UInt_t nEvents
Definition: hcalCalib.cc:43
#define debug
Definition: MEtoEDMFormat.h:34
virtual double phi() const
momentum azimuthal angle
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects