CMS 3D CMS Logo

RazorVarAnalyzer.cc
Go to the documentation of this file.
2 
5 
8 
11 
12 #include <cmath>
13 
14 // A simple constructor which takes as inoput only the name of the PF jet
15 // collection
17  : ScoutingAnalyzerBase(conf),
18  m_eleCollectionTag(conf.getUntrackedParameter<edm::InputTag>("eleCollectionName",
19  edm::InputTag("hltPixelMatchElectronsActivity"))),
20  m_jetCollectionTag(
21  conf.getUntrackedParameter<edm::InputTag>("jetCollectionName", edm::InputTag("hltCaloJetIDPassed"))),
22  m_muCollectionTag(
23  conf.getUntrackedParameter<edm::InputTag>("muCollectionName", edm::InputTag("hltL3MuonCandidates"))),
24  m_razorVarCollectionTag(conf.getUntrackedParameter<edm::InputTag>("razorVarCollectionName")) {
25  // set Token(-s)
26  m_jetCollectionTagToken_ = consumes<reco::CaloJetCollection>(
27  conf.getUntrackedParameter<edm::InputTag>("jetCollectionName", edm::InputTag("hltCaloJetIDPassed")));
28  m_muCollectionTagToken_ = consumes<std::vector<reco::RecoChargedCandidate>>(
29  conf.getUntrackedParameter<edm::InputTag>("muCollectionName", edm::InputTag("hltL3MuonCandidates")));
30  m_eleCollectionTagToken_ = consumes<reco::ElectronCollection>(
31  conf.getUntrackedParameter<edm::InputTag>("eleCollectionName", edm::InputTag("hltPixelMatchElectronsActivity")));
33  consumes<std::vector<double>>(conf.getUntrackedParameter<edm::InputTag>("razorVarCollectionName"));
34 }
35 
37 
38 // Function to book the Monitoring Elements.
41  // the full inclusive histograms
43  iBooker, "rsqMRFullyInc", "M_{R} vs R^{2} (All Events)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
45  iBooker, "rsqMRInc4J", "M_{R} vs R^{2} (>= 4j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
47  iBooker, "rsqMRInc6J", "M_{R} vs R^{2} (>= 6j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
49  iBooker, "rsqMRInc8J", "M_{R} vs R^{2} (>= 8j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
51  iBooker, "rsqMRInc10J", "M_{R} vs R^{2} (>= 10j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
53  iBooker, "rsqMRInc12J", "M_{R} vs R^{2} (>= 12j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
55  iBooker, "rsqMRInc14J", "M_{R} vs R^{2} (>= 14j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
56 
57  // the by box histograms
59  iBooker, "rsqMREleMu", "M_{R} vs R^{2} (EleMu box)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
61  iBooker, "rsqMRMuMu", "M_{R} vs R^{2} (MuMu box)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
63  iBooker, "rsqMREleEle", "M_{R} vs R^{2} (EleEle box)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
65  iBooker, "rsqMRMu", "M_{R} vs R^{2} (Mu box)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
67  iBooker, "rsqMREle", "M_{R} vs R^{2} (Ele box)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
69  iBooker, "rsqMRHad", "M_{R} vs R^{2} (Had box)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
70 
71  // the by box histograms
73  iBooker, "rsqMRMuMJ", "M_{R} vs R^{2} (Mu box >= 4j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
75  iBooker, "rsqMREleMJ", "M_{R} vs R^{2} (Ele box >= 5j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
77  iBooker, "rsqMRHadMJ", "M_{R} vs R^{2} (Had box >= 6j)", 400, 0., 4000., 50, 0., 1., "M_{R} [GeV]", "R^{2}");
78 }
79 
80 // Usual analyze method
82  // count the number of jets with a minimal selection
84  iEvent.getByToken(m_jetCollectionTagToken_, calojets_handle);
85 
86  unsigned int njets = 0;
87  for (reco::CaloJetCollection::const_iterator it = calojets_handle->begin(); it != calojets_handle->end(); ++it) {
88  if (it->pt() >= 30. && fabs(it->eta()) <= 3.0) {
89  njets++;
90  }
91  }
92 
93  // count the number of muons
95  iEvent.getByToken(m_muCollectionTagToken_, muon_handle);
96 
97  unsigned int nmu_loose = 0;
98  unsigned int nmu_tight = 0;
99  if (muon_handle.isValid()) {
100  for (std::vector<reco::RecoChargedCandidate>::const_iterator it = muon_handle->begin(); it != muon_handle->end();
101  ++it) {
102  if (it->pt() >= 15 && fabs(it->eta()) <= 2.1)
103  nmu_tight++;
104  if (it->pt() >= 10 && fabs(it->eta()) <= 2.4)
105  nmu_loose++;
106  }
107  }
108 
109  // count the number of electrons
111  iEvent.getByToken(m_eleCollectionTagToken_, ele_handle);
112 
113  unsigned int nele_loose = 0;
114  unsigned int nele_tight = 0;
115  if (ele_handle.isValid()) {
116  for (reco::ElectronCollection::const_iterator it = ele_handle->begin(); it != ele_handle->end(); ++it) {
117  if (it->pt() >= 20 && fabs(it->eta()) <= 2.5)
118  nele_tight++;
119  if (it->pt() >= 10 && fabs(it->eta()) <= 2.5)
120  nele_loose++;
121  }
122  }
123 
124  // now get the box number:
125  // {'MuEle':0,'MuMu':1,'EleEle':2,'Mu':3,'Ele':4,'Had':5}
126  unsigned int box_num = 5;
127  if (nmu_tight > 0 && nele_tight > 0) {
128  box_num = 0;
129  } else if (nmu_tight > 0 && nmu_loose > 1) {
130  box_num = 1;
131  } else if (nele_tight > 0 && nele_loose > 1) {
132  box_num = 2;
133  } else if (nmu_tight > 0) {
134  box_num = 3;
135  } else if (nele_tight > 0) {
136  box_num = 4;
137  }
138 
139  edm::Handle<std::vector<double>> razorvar_handle;
140  iEvent.getByToken(m_razorVarCollectionTagToken_, razorvar_handle);
141  if (razorvar_handle->size() > 1) {
142  const double MR = razorvar_handle->at(0);
143  const double R = razorvar_handle->at(1);
144  m_rsqMRFullyInc->Fill(MR, R * R);
145  if (njets >= 4)
146  m_rsqMRInc4J->Fill(MR, R * R);
147  if (njets >= 6)
148  m_rsqMRInc6J->Fill(MR, R * R);
149  if (njets >= 8)
150  m_rsqMRInc8J->Fill(MR, R * R);
151  if (njets >= 10)
152  m_rsqMRInc10J->Fill(MR, R * R);
153  if (njets >= 12)
154  m_rsqMRInc12J->Fill(MR, R * R);
155  if (njets >= 14)
156  m_rsqMRInc14J->Fill(MR, R * R);
157 
158  // now fill the boxes
159  if (box_num == 0)
160  m_rsqMREleMu->Fill(MR, R * R);
161  if (box_num == 1)
162  m_rsqMRMuMu->Fill(MR, R * R);
163  if (box_num == 2)
164  m_rsqMREleEle->Fill(MR, R * R);
165  if (box_num == 3)
166  m_rsqMRMu->Fill(MR, R * R);
167  if (box_num == 4)
168  m_rsqMREle->Fill(MR, R * R);
169  if (box_num == 5)
170  m_rsqMRHad->Fill(MR, R * R);
171 
172  // finally the multijet boxes - think ttbar
173  // muon boxes: muons are not in jets
174  if (box_num == 3 && njets >= 4)
175  m_rsqMRMuMJ->Fill(MR, R * R);
176  // ele boxes: electrons are in jets
177  else if (box_num == 4 && njets >= 5)
178  m_rsqMREleMJ->Fill(MR, R * R);
179  // fill the Had box
180  else if (box_num == 5 && njets >= 6)
181  m_rsqMRHadMJ->Fill(MR, R * R);
182  }
183 }
MonitorElement * m_rsqMRInc12J
MonitorElement * m_rsqMRHadMJ
T getUntrackedParameter(std::string const &, T const &) const
edm::EDGetTokenT< std::vector< double > > m_razorVarCollectionTagToken_
edm::EDGetTokenT< reco::ElectronCollection > m_eleCollectionTagToken_
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
MonitorElement * m_rsqMREle
MonitorElement * m_rsqMRInc6J
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:517
void prepareBooking(DQMStore::IBooker &)
void Fill(long long x)
MonitorElement * m_rsqMREleMJ
MonitorElement * m_rsqMRInc8J
MonitorElement * m_rsqMRMu
int iEvent
Definition: GenABIO.cc:224
MonitorElement * m_rsqMRMuMJ
edm::EDGetTokenT< std::vector< reco::RecoChargedCandidate > > m_muCollectionTagToken_
RazorVarAnalyzer(const edm::ParameterSet &)
void analyze(const edm::Event &, const edm::EventSetup &) override
MonitorElement * m_rsqMREleMu
MonitorElement * m_rsqMRMuMu
bool isValid() const
Definition: HandleBase.h:74
MonitorElement * m_rsqMRInc14J
MonitorElement * m_rsqMRInc4J
MonitorElement * bookH2withSumw2(DQMStore::IBooker &, const std::string &name, const std::string &title, int nchX, double lowX, double highX, int nchY, double lowY, double highY, const std::string &titleX="", const std::string &titleY="", Option_t *option="COLZ")
HLT enums.
MonitorElement * m_rsqMREleEle
~RazorVarAnalyzer() override
edm::EDGetTokenT< reco::CaloJetCollection > m_jetCollectionTagToken_
MonitorElement * m_rsqMRHad
MonitorElement * m_rsqMRInc10J
MonitorElement * m_rsqMRFullyInc
Definition: Run.h:45