CMS 3D CMS Logo

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 
26  // Local Debug flag
27  debug = iConfig.getParameter<bool>("DebugHeavyChHiggsToTauNuSkim");
28 
29  hltTauToken = consumes<IsolatedTauTagInfoCollection>(iConfig.getParameter<InputTag>("HLTTauCollection"));
30  jetToken = consumes<CaloJetCollection>(iConfig.getParameter<InputTag>("JetTagCollection"));
31  minNumberOfjets = iConfig.getParameter<int>("minNumberOfJets");
32  jetEtMin = iConfig.getParameter<double>("jetEtMin");
33  jetEtaMin = iConfig.getParameter<double>("jetEtaMin");
34  jetEtaMax = iConfig.getParameter<double>("jetEtaMax");
35  minDRFromTau = iConfig.getParameter<double>("minDRFromTau");
36 
37  nEvents = 0;
38  nSelectedEvents = 0;
39 }
40 
42  edm::LogVerbatim("HeavyChHiggsToTauNuSkim")
43  << " Number_events_read " << nEvents << " Number_events_kept " << nSelectedEvents << " Efficiency "
44  << ((double)nSelectedEvents) / ((double)nEvents + 0.01) << std::endl;
45 }
46 
48  nEvents++;
49 
51  iEvent.getByToken(hltTauToken, tauTagL3Handle);
52 
53  if (!tauTagL3Handle.isValid())
54  return false;
55 
56  Jet theTau;
57  double maxEt = 0;
58  if (tauTagL3Handle.isValid()) {
59  const IsolatedTauTagInfoCollection& L3Taus = *(tauTagL3Handle.product());
60  IsolatedTauTagInfoCollection::const_iterator i;
61  for (i = L3Taus.begin(); i != L3Taus.end(); i++) {
62  if (i->discriminator() == 0)
63  continue;
64  Jet taujet = *(i->jet().get());
65  if (taujet.et() > maxEt) {
66  maxEt = taujet.et();
67  theTau = taujet;
68  }
69  }
70  }
71 
72  if (maxEt == 0)
73  return false;
74 
75  // jets
76 
77  Handle<CaloJetCollection> jetHandle;
78  iEvent.getByToken(jetToken, jetHandle);
79 
80  if (!jetHandle.isValid())
81  return false;
82 
83  bool accepted = false;
84 
85  if (jetHandle.isValid()) {
86  int nJets = 0;
87  const reco::CaloJetCollection& jets = *(jetHandle.product());
88  CaloJetCollection::const_iterator iJet;
89  for (iJet = jets.begin(); iJet != jets.end(); iJet++) {
90  if (iJet->et() > jetEtMin && iJet->eta() > jetEtaMin && iJet->eta() < jetEtaMax) {
91  double DR = deltaR(theTau.eta(), iJet->eta(), theTau.phi(), iJet->phi());
92  if (DR > minDRFromTau)
93  nJets++;
94  }
95  }
96  if (nJets >= minNumberOfjets) {
97  accepted = true;
98  nSelectedEvents++;
99  }
100  }
101  return accepted;
102 }
JetTag.h
mps_fire.i
i
Definition: mps_fire.py:355
MessageLogger.h
edm::Handle::product
T const * product() const
Definition: Handle.h:70
reco::IsolatedTauTagInfoCollection
std::vector< IsolatedTauTagInfo > IsolatedTauTagInfoCollection
Definition: IsolatedTauTagInfo.h:95
edm
HLT enums.
Definition: AlignableModifier.h:19
MssmHbbBtagTriggerMonitor_cfi.jetEtaMax
jetEtaMax
Definition: MssmHbbBtagTriggerMonitor_cfi.py:24
Jet.h
singleTopDQM_cfi.jets
jets
Definition: singleTopDQM_cfi.py:42
reco
fixed size matrix
Definition: AlignmentAlgorithmBase.h:45
edm::Handle
Definition: AssociativeIterator.h:50
debug
#define debug
Definition: HDRShower.cc:19
Jet
Definition: Jet.py:1
cms::dd::accepted
bool accepted(std::vector< std::regex > const &, std::string_view)
PbPb_ZMuSkimMuonDPG_cff.deltaR
deltaR
Definition: PbPb_ZMuSkimMuonDPG_cff.py:63
heavyChHiggsToTauNu_Filter_cfi.jetEtaMin
jetEtaMin
Definition: heavyChHiggsToTauNu_Filter_cfi.py:8
edm::ParameterSet
Definition: ParameterSet.h:36
reco::LeafCandidate::eta
double eta() const final
momentum pseudorapidity
Definition: LeafCandidate.h:152
iEvent
int iEvent
Definition: GenABIO.cc:224
reco::CaloJetCollection
std::vector< CaloJet > CaloJetCollection
collection of CaloJet objects
Definition: CaloJetCollection.h:15
HeavyChHiggsToTauNuSkim::HeavyChHiggsToTauNuSkim
HeavyChHiggsToTauNuSkim(const edm::ParameterSet &)
Definition: HeavyChHiggsToTauNuSkim.cc:25
edm::LogVerbatim
Definition: MessageLogger.h:297
edm::EventSetup
Definition: EventSetup.h:57
edm::ParameterSet::getParameter
T getParameter(std::string const &) const
HeavyChHiggsToTauNuSkim.h
reco::LeafCandidate::et
double et() const final
transverse energy
Definition: LeafCandidate.h:127
HeavyChHiggsToTauNuSkim::filter
bool filter(edm::Event &, const edm::EventSetup &) override
Definition: HeavyChHiggsToTauNuSkim.cc:47
std
Definition: JetResolutionObject.h:76
HeavyChHiggsToTauNuSkim::~HeavyChHiggsToTauNuSkim
~HeavyChHiggsToTauNuSkim() override
Definition: HeavyChHiggsToTauNuSkim.cc:41
reco::LeafCandidate::phi
double phi() const final
momentum azimuthal angle
Definition: LeafCandidate.h:148
heavyionUCCDQM_cfi.maxEt
maxEt
Definition: heavyionUCCDQM_cfi.py:14
heavyChHiggsToTauNu_Filter_cfi.minDRFromTau
minDRFromTau
Definition: heavyChHiggsToTauNu_Filter_cfi.py:6
heavyChHiggsToTauNu_Filter_cfi.jetEtMin
jetEtMin
Definition: heavyChHiggsToTauNu_Filter_cfi.py:11
HepMCProduct.h
edm::HandleBase::isValid
bool isValid() const
Definition: HandleBase.h:70
nEvents
UInt_t nEvents
Definition: hcalCalib.cc:41
edm::Event
Definition: Event.h:73
edm::InputTag
Definition: InputTag.h:15
unpackData-CaloStage1.jetToken
jetToken
Definition: unpackData-CaloStage1.py:163