CMS 3D CMS Logo

TagAndProbeBtagTriggerMonitor.cc
Go to the documentation of this file.
1 #include <cmath>
2 
4 
6 
8 
10 
11 //#include "TLorentzVector.h"
12 
13 
14 // -----------------------------
15 // constructors and destructor
16 // -----------------------------
17 
19 {
20  folderName_ = iConfig.getParameter<std::string>("dirname");
21  processname_ = iConfig.getParameter<std::string>("processname");
22  triggerobjbtag_ = iConfig.getParameter<std::string>("triggerobjbtag");
23  jetPtmin_ = iConfig.getParameter<double>("jetPtMin");
24  jetEtamax_ = iConfig.getParameter<double>("jetEtaMax");
25  tagBtagmin_ = iConfig.getParameter<double>("tagBtagMin");
26  probeBtagmin_ = iConfig.getParameter<double>("probeBtagMin");
27  triggerSummaryLabel_ = iConfig.getParameter<edm::InputTag>("triggerSummary");
28  triggerSummaryToken_ = consumes <trigger::TriggerEvent> (triggerSummaryLabel_);
29  offlineBtagToken_ = consumes<reco::JetTagCollection> (iConfig.getParameter<edm::InputTag>("offlineBtag"));
30  genTriggerEventFlag_ = new GenericTriggerEventFlag(iConfig.getParameter<edm::ParameterSet>("genericTriggerEventPSet"),consumesCollector(), *this);
31 
32  jetPtbins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetPt");
33  jetEtabins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetEta");
34  jetPhibins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetPhi");
35  jetBtagbins_ = iConfig.getParameter<edm::ParameterSet>("histoPSet").getParameter<std::vector<double> >("jetBtag");
36 
37 }
38 
40 {
42 
43 }
44 
46  edm::Run const & iRun,
47  edm::EventSetup const & iSetup)
48 {
49 
50  std::string currentFolder = folderName_ ;
51  ibooker.setCurrentFolder(currentFolder);
52 
53  int ptnbins = jetPtbins_.size()-1;
54  int etanbins = jetEtabins_.size()-1;
55  int phinbins = jetPhibins_.size()-1;
56  int btagnbins = jetBtagbins_.size()-1;
57 
58  std::vector<float> fptbins(jetPtbins_.begin(),jetPtbins_.end());
59  std::vector<float> fetabins(jetEtabins_.begin(),jetEtabins_.end());
60  std::vector<float> fphibins(jetPhibins_.begin(),jetPhibins_.end());
61  std::vector<float> fbtagbins(jetBtagbins_.begin(),jetBtagbins_.end());
62 
63  float * ptbins = &fptbins[0];
64  float * etabins = &fetabins[0];
65  float * phibins = &fphibins[0];
66  float * btagbins = &fbtagbins[0];
67 
68 
69  pt_jet1_ = ibooker.book1D("pt_jet1","pt_jet1",ptnbins, ptbins);
70  pt_jet2_ = ibooker.book1D("pt_jet2","pt_jet2",ptnbins, ptbins);
71  eta_jet1_ = ibooker.book1D("eta_jet1","eta_jet1",etanbins, etabins);
72  eta_jet2_ = ibooker.book1D("eta_jet2","eta_jet2",etanbins, etabins);
73  phi_jet1_ = ibooker.book1D("phi_jet1","phi_jet1",phinbins, phibins);
74  phi_jet2_ = ibooker.book1D("phi_jet2","phi_jet2",phinbins, phibins);
75  eta_phi_jet1_ = ibooker.book2D("eta_phi_jet1","eta_phi_jet1",etanbins, etabins,phinbins, phibins);
76  eta_phi_jet2_ = ibooker.book2D("eta_phi_jet2","eta_phi_jet2",etanbins, etabins,phinbins, phibins);
77 
78  pt_probe_ = ibooker.book1D("pt_probe","pt_probe",ptnbins, ptbins);
79  pt_probe_match_ = ibooker.book1D("pt_probe_match","pt_probe_match",ptnbins, ptbins);
80  eta_probe_ = ibooker.book1D("eta_probe","eta_probe",etanbins, etabins);
81  eta_probe_match_ = ibooker.book1D("eta_probe_match","eta_probe_match",etanbins, etabins);
82  phi_probe_ = ibooker.book1D("phi_probe","phi_probe",phinbins, phibins);
83  phi_probe_match_ = ibooker.book1D("phi_probe_match","phi_probe_match",phinbins, phibins);
84  eta_phi_probe_ = ibooker.book2D("eta_phi_probe","eta_phi_probe",etanbins, etabins,phinbins, phibins);
85  eta_phi_probe_match_ = ibooker.book2D("eta_phi_probe_match","eta_phi_match",etanbins, etabins,phinbins, phibins);
86 
87  discr_offline_btag_jet1_ = ibooker.book1D("discr_offline_btag_jet1","discr_offline_btag_jet1",btagnbins, btagbins);
88  discr_offline_btag_jet2_ = ibooker.book1D("discr_offline_btag_jet2","discr_offline_btag_jet2",btagnbins, btagbins);
89 
90  // Initialize the GenericTriggerEventFlag
92 }
93 
95 {
96 
97 // bool accept = false;
98  bool match1 = false;
99  bool match2 = false;
100 
101  edm::Handle<reco::JetTagCollection> offlineJetTagPFHandler;
102  iEvent.getByToken(offlineBtagToken_, offlineJetTagPFHandler);
103 
104  if( ! offlineJetTagPFHandler.isValid()) return;
105 
106  // applying selection for event; tag & probe -> selection for all events
107  if (genTriggerEventFlag_->on() && ! genTriggerEventFlag_->accept( iEvent, iSetup) )
108  {
109 
110  reco::JetTagCollection jettags = *(offlineJetTagPFHandler.product());
111  if ( jettags.size() > 1 )
112  {
113 
114  const reco::Jet jet1 = *(jettags.key(0).get());
115  const reco::Jet jet2 = *(jettags.key(1).get());
116 
117  float btag1 = jettags.value(0);
118  float btag2 = jettags.value(1);
119 
120  if ( jet1.pt() > jetPtmin_ && jet2.pt() > jetPtmin_ && fabs(jet1.eta()) < jetEtamax_ && fabs(jet2.eta()) < jetEtamax_ )
121  {
122  if ( btag1 > tagBtagmin_ && btag2 > probeBtagmin_ )
123  {
124  pt_jet1_ -> Fill(jet1.pt());
125  pt_jet2_ -> Fill(jet2.pt());
126  eta_jet1_ -> Fill(jet1.eta());
127  eta_jet2_ -> Fill(jet2.eta());
128  phi_jet1_ -> Fill(jet1.phi());
129  phi_jet2_ -> Fill(jet2.phi());
130  eta_phi_jet1_ -> Fill(jet1.eta(),jet1.phi());
131  eta_phi_jet2_ -> Fill(jet2.eta(),jet2.phi());
132  if ( btag1 < 0 ) btag1 = -0.5;
133  if ( btag2 < 0 ) btag2 = -0.5;
134  discr_offline_btag_jet1_ -> Fill(btag1);
135  discr_offline_btag_jet2_ -> Fill(btag2);
136 
137 // trigger objects matching
138  std::vector<trigger::TriggerObject> onlinebtags;
139  edm::Handle<trigger::TriggerEvent> triggerEventHandler;
140  iEvent.getByToken(triggerSummaryToken_, triggerEventHandler);
141  const unsigned int filterIndex(triggerEventHandler->filterIndex(edm::InputTag(triggerobjbtag_,"",processname_)));
142  if ( filterIndex < triggerEventHandler->sizeFilters() )
143  {
144  const trigger::Keys& keys(triggerEventHandler->filterKeys(filterIndex));
145  const trigger::TriggerObjectCollection & triggerObjects = triggerEventHandler->getObjects();
146  for ( auto & key : keys )
147  {
148  onlinebtags.reserve(onlinebtags.size()+keys.size());
149  onlinebtags.push_back(triggerObjects[key]);
150  }
151  }
152  for ( auto const & to : onlinebtags )
153  {
154  if ( reco::deltaR(jet1,to) ) match1 = true;
155  if ( reco::deltaR(jet2,to) ) match2 = true;
156  }
157 
158  if ( match1 ) // jet1 is the tag
159  {
160  pt_probe_ -> Fill(jet2.pt());
161  eta_probe_ -> Fill(jet2.eta());
162  phi_probe_ -> Fill(jet2.phi());
163  eta_phi_probe_ -> Fill(jet2.eta(),jet2.phi());
164  if ( match2 ) // jet2 is the probe
165  {
166  pt_probe_match_ -> Fill(jet2.pt());
167  eta_probe_match_ -> Fill(jet2.eta());
168  phi_probe_match_ -> Fill(jet2.phi());
169  eta_phi_probe_match_ -> Fill(jet2.eta(),jet2.phi());
170  }
171  }
172  } // offline jets btag
173  } // offline jets kinematic selection
174  } // at least two offline jets
175  } // accept trigger
176 
177 }
179 {
180 }
181 
182 
183 // Define this as a plug-in
T getParameter(std::string const &) const
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
Base class for all types of Jets.
Definition: Jet.h:20
const Keys & filterKeys(trigger::size_type index) const
Definition: TriggerEvent.h:111
trigger::size_type filterIndex(const edm::InputTag &filterTag) const
find index of filter in data-member vector from filter tag
Definition: TriggerEvent.h:123
double pt() const final
transverse momentum
edm::EDGetTokenT< reco::JetTagCollection > offlineBtagToken_
Provides a code based selection for trigger and DCS information in order to have no failing filters i...
GenericTriggerEventFlag * genTriggerEventFlag_
int iEvent
Definition: GenABIO.cc:224
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:268
edm::EDGetTokenT< trigger::TriggerEvent > triggerSummaryToken_
const TriggerObjectCollection & getObjects() const
Definition: TriggerEvent.h:98
void Fill(HcalDetId &id, double val, std::vector< TH2F > &depth)
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
bool accept(const edm::Event &event, const edm::EventSetup &setup)
To be called from analyze/filter() methods.
void dqmBeginRun(edm::Run const &run, edm::EventSetup const &iSetup) override
MonitorElement * book1D(Args &&...args)
Definition: DQMStore.h:106
CVal::value_type const value(size_type i) const
bool isValid() const
Definition: HandleBase.h:74
KeyRef key(size_type i) const
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
std::vector< TriggerObject > TriggerObjectCollection
collection of trigger physics objects (e.g., all isolated muons)
Definition: TriggerObject.h:81
T const * product() const
Definition: Handle.h:74
MonitorElement * book2D(Args &&...args)
Definition: DQMStore.h:109
std::vector< size_type > Keys
TagAndProbeBtagTriggerMonitor(const edm::ParameterSet &)
void initRun(const edm::Run &run, const edm::EventSetup &setup)
To be called from beginRun() methods.
void analyze(edm::Event const &iEvent, edm::EventSetup const &iSetup) override
double phi() const final
momentum azimuthal angle
Definition: Run.h:45
size_type size() const