CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
L1Validator.cc
Go to the documentation of this file.
1 // -*- C++ -*-
2 //
3 // Package: L1T
4 // Class: L1Validator
5 //
14 //
15 // Original Author: Scott Wilbur
16 // Created: Wed, 28 Aug 2013 09:42:55 GMT
17 // $Id$
18 //
19 //
20 
21 #include <string>
22 
24 
26 
28 
29 #include "TFile.h"
30 
31 // defining as a macro instead of a function because inheritance doesn't work:
32 #define FINDRECOPART(TYPE, COLLECTION1, COLLECTION2) \
33  const TYPE *RecoPart = NULL; \
34  double BestDist = 999.; \
35  for (uint i = 0; i < COLLECTION1->size(); i++) { \
36  const TYPE *ThisPart = &COLLECTION1->at(i); \
37  double ThisDist = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi()); \
38  if (ThisDist < 1.0 && ThisDist < BestDist) { \
39  BestDist = ThisDist; \
40  RecoPart = ThisPart; \
41  } \
42  } \
43  if (COLLECTION1.product() != COLLECTION2.product()) { \
44  for (uint i = 0; i < COLLECTION2->size(); i++) { \
45  const TYPE *ThisPart = &COLLECTION2->at(i); \
46  double ThisDist = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi()); \
47  if (ThisDist < 1.0 && ThisDist < BestDist) { \
48  BestDist = ThisDist; \
49  RecoPart = ThisPart; \
50  } \
51  } \
52  }
53 
55  _dirName = iConfig.getParameter<std::string>("dirName");
56  _GenSource = consumes<reco::GenParticleCollection>(iConfig.getParameter<edm::InputTag>("GenSource"));
57 
58  _L1MuonBXSource = consumes<l1t::MuonBxCollection>(iConfig.getParameter<edm::InputTag>("L1MuonBXSource"));
59  _L1EGammaBXSource = consumes<l1t::EGammaBxCollection>(iConfig.getParameter<edm::InputTag>("L1EGammaBXSource"));
60  _L1TauBXSource = consumes<l1t::TauBxCollection>(iConfig.getParameter<edm::InputTag>("L1TauBXSource"));
61  _L1JetBXSource = consumes<l1t::JetBxCollection>(iConfig.getParameter<edm::InputTag>("L1JetBXSource"));
62  _srcToken = mayConsume<GenEventInfoProduct>(iConfig.getParameter<edm::InputTag>("srcToken"));
63  _L1GenJetSource = consumes<reco::GenJetCollection>(iConfig.getParameter<edm::InputTag>("L1GenJetSource"));
64 
65  //_fileName = iConfig.getParameter<std::string>("fileName");
66 }
67 
69 
71  // iBooker.setCurrentFolder(_dirName);
72  _Hists.Book(iBooker, _dirName);
73 };
74 
76  using namespace edm;
77  using namespace std;
78  using namespace l1extra;
79  using namespace reco;
80 
81  Handle<GenParticleCollection> GenParticles;
86  Handle<GenEventInfoProduct> genEvtInfoProduct;
88 
89  bool GotEverything = true;
90 
91  if (!iEvent.getByToken(_GenSource, GenParticles))
92  GotEverything = false;
93  if (!iEvent.getByToken(_L1MuonBXSource, MuonsBX))
94  GotEverything = false;
95  if (!iEvent.getByToken(_L1EGammaBXSource, EGammasBX))
96  GotEverything = false;
97  if (!iEvent.getByToken(_L1TauBXSource, TausBX))
98  GotEverything = false;
99  if (!iEvent.getByToken(_L1JetBXSource, JetsBX))
100  GotEverything = false;
101  if (!iEvent.getByToken(_srcToken, genEvtInfoProduct))
102  GotEverything = false;
103  if (!iEvent.getByToken(_L1GenJetSource, GenJets))
104  GotEverything = false;
105 
106  if (!GotEverything)
107  return;
108 
109  /*
110  std::string moduleName = "";
111  if( genEvtInfoProduct.isValid() ) {
112  const edm::Provenance& prov =
113  iEvent.getProvenance(genEvtInfoProduct.id()); moduleName =
114  edm::moduleName(prov);
115  //cout<<" generator name: "<<moduleName<<endl;
116  }
117  */
118 
119  _Hists.NEvents++;
120 
121  int nL1Muons = 0, nL1EGammas = 0, nL1Taus = 0, nL1Jets = 0;
122  if (MuonsBX->getFirstBX() >= 0)
123  nL1Muons = MuonsBX->size(0);
124  if (EGammasBX->getFirstBX() >= 0)
125  nL1EGammas = EGammasBX->size(0);
126  if (TausBX->getFirstBX() >= 0)
127  nL1Taus = TausBX->size(0);
128  if (JetsBX->getFirstBX() >= 0)
129  nL1Jets = JetsBX->size(0);
130 
132  _Hists.FillNumber(L1ValidatorHists::Type::Egamma, nL1EGammas);
135 
136  // For gen jet
137 
138  for (auto &Genjet : *GenJets) {
139  // eta within calorimeter acceptance 4.7
140  if (fabs((&Genjet)->eta()) > 4.7)
141  continue;
142 
143  // only consider the gen jet with pt greater than 10 GeV
144  if ((&Genjet)->pt() < 10.0)
145  continue;
146 
147  double minDR = 999.0;
148 
149  // match L1T object
150  const l1t::Jet *L1Part = nullptr;
151  for (int iBx = JetsBX->getFirstBX(); iBx <= JetsBX->getLastBX(); ++iBx) {
152  if (iBx > 0)
153  continue;
154  for (std::vector<l1t::Jet>::const_iterator jet = JetsBX->begin(iBx); jet != JetsBX->end(iBx); ++jet) {
155  double idR = reco::deltaR((&Genjet)->eta(), (&Genjet)->phi(), jet->eta(), jet->phi());
156  if (idR < minDR) {
157  minDR = idR;
158  L1Part = &(*jet);
159  }
160  }
161  }
162  _Hists.Fill(L1ValidatorHists::Type::Jet, &Genjet, L1Part);
163  }
164 
165  for (uint i = 0; i < GenParticles->size(); i++) {
166  const GenParticle *GenPart = &GenParticles->at(i);
167 
168  int pdg = GenPart->pdgId(), status = GenPart->status();
169 
170  double minDR = 999.0;
171 
172  // only consider the gen particle with pt greater than 10 GeV
173  if (GenPart->pt() < 10.0)
174  continue;
175 
177  if (status == 1 && abs(pdg) == 13) { // Muon
178 
179  // eta within tracker acceptance 2.4
180  if (fabs(GenPart->eta()) > 2.4)
181  continue;
182 
183  // match L1T object
184  const l1t::Muon *L1Part = nullptr;
185  for (int iBx = MuonsBX->getFirstBX(); iBx <= MuonsBX->getLastBX(); ++iBx) {
186  if (iBx > 0)
187  continue;
188  for (std::vector<l1t::Muon>::const_iterator mu = MuonsBX->begin(iBx); mu != MuonsBX->end(iBx); ++mu) {
189  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), mu->eta(), mu->phi());
190  if (idR < minDR) {
191  minDR = idR;
192  L1Part = &(*mu);
193  }
194  }
195  _Hists.Fill(L1ValidatorHists::Type::Muon, GenPart, L1Part);
196  }
197 
200  } else if (status == 1 && (abs(pdg) == 11 || pdg == 22)) { // Egamma
201 
202  // eta within EM calorimeter acceptance 2.5
203  if (fabs(GenPart->eta()) > 2.5)
204  continue;
205 
206  // exclude the calorimeter barrel and endcap overlap region
207  if (fabs(GenPart->eta()) > 1.4442 && fabs(GenPart->eta()) < 1.5660)
208  continue;
209 
210  // match L1T object
211  const l1t::EGamma *L1Part = nullptr;
212  for (int iBx = EGammasBX->getFirstBX(); iBx <= EGammasBX->getLastBX(); ++iBx) {
213  if (iBx > 0)
214  continue;
215  for (std::vector<l1t::EGamma>::const_iterator eg = EGammasBX->begin(iBx); eg != EGammasBX->end(iBx); ++eg) {
216  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), eg->eta(), eg->phi());
217  if (idR < minDR) {
218  minDR = idR;
219  L1Part = &(*eg);
220  }
221  }
222  }
223  _Hists.Fill(L1ValidatorHists::Type::Egamma, GenPart, L1Part);
224 
227  } else if (status == 2 && abs(pdg) == 15) { // Tau
228 
229  // eta within tracker acceptance 2.4
230  if (fabs(GenPart->eta()) > 2.4)
231  continue;
232 
233  // match L1T object
234  const l1t::Tau *L1Part = nullptr;
235  for (int iBx = TausBX->getFirstBX(); iBx <= TausBX->getLastBX(); ++iBx) {
236  if (iBx > 0)
237  continue;
238  for (std::vector<l1t::Tau>::const_iterator tau = TausBX->begin(iBx); tau != TausBX->end(iBx); ++tau) {
239  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), tau->eta(), tau->phi());
240  if (idR < minDR) {
241  minDR = idR;
242  L1Part = &(*tau);
243  }
244  }
245  }
246  _Hists.Fill(L1ValidatorHists::Type::Tau, GenPart, L1Part);
247  }
248  }
249 }
250 
251 // The next three are exactly the same, but apparently inheritance doesn't work
252 // like I thought it did.
254  const std::vector<l1extra::L1EmParticle> *Collection1,
255  const std::vector<l1extra::L1EmParticle> *Collection2 = nullptr) {
256  const reco::LeafCandidate *BestPart = nullptr;
257  double BestDR = 999.;
258 
259  for (uint i = 0; i < Collection1->size(); i++) {
260  const reco::LeafCandidate *ThisPart = &Collection1->at(i);
261  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
262  if (ThisDR < BestDR) {
263  BestDR = ThisDR;
264  BestPart = ThisPart;
265  }
266  }
267 
268  if (Collection2 == nullptr)
269  return BestPart;
270 
271  for (uint i = 0; i < Collection2->size(); i++) {
272  const reco::LeafCandidate *ThisPart = &Collection2->at(i);
273  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
274  if (ThisDR < BestDR) {
275  BestDR = ThisDR;
276  BestPart = ThisPart;
277  }
278  }
279 
280  return BestPart;
281 }
282 
284  const std::vector<l1extra::L1JetParticle> *Collection1,
285  const std::vector<l1extra::L1JetParticle> *Collection2 = nullptr) {
286  const reco::LeafCandidate *BestPart = nullptr;
287  double BestDR = 999.;
288 
289  for (uint i = 0; i < Collection1->size(); i++) {
290  const reco::LeafCandidate *ThisPart = &Collection1->at(i);
291  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
292  if (ThisDR < BestDR) {
293  BestDR = ThisDR;
294  BestPart = ThisPart;
295  }
296  }
297 
298  if (Collection2 == nullptr)
299  return BestPart;
300 
301  for (uint i = 0; i < Collection2->size(); i++) {
302  const reco::LeafCandidate *ThisPart = &Collection2->at(i);
303  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
304  if (ThisDR < BestDR) {
305  BestDR = ThisDR;
306  BestPart = ThisPart;
307  }
308  }
309 
310  return BestPart;
311 }
312 
314  const std::vector<l1extra::L1MuonParticle> *Collection1) {
315  const reco::LeafCandidate *BestPart = nullptr;
316  double BestDR = 999.;
317 
318  for (uint i = 0; i < Collection1->size(); i++) {
319  const reco::LeafCandidate *ThisPart = &Collection1->at(i);
320  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
321  if (ThisDR < BestDR) {
322  BestDR = ThisDR;
323  BestPart = ThisPart;
324  }
325  }
326 
327  return BestPart;
328 }
329 
330 // ------------ method fills 'descriptions' with the allowed parameters for the
331 // module ------------
333  // The following says we do not know what parameters are allowed so do no
334  // validation
335  // Please change this to state exactly what you do use, even if it is no
336  // parameters
338  desc.setUnknown();
339  descriptions.addDefault(desc);
340 }
341 
342 // define this as a plug-in
void Fill(int, const reco::LeafCandidate *, const reco::LeafCandidate *)
edm::EDGetTokenT< reco::GenJetCollection > _L1GenJetSource
Definition: L1Validator.h:89
double pt() const final
transverse momentum
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:16
edm::EDGetTokenT< reco::GenParticleCollection > _GenSource
Definition: L1Validator.h:83
Definition: Tau.h:20
edm::EDGetTokenT< GenEventInfoProduct > _srcToken
Definition: L1Validator.h:88
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1Validator.cc:332
list status
Definition: mps_update.py:107
int status() const final
status word
L1Validator(const edm::ParameterSet &)
Definition: L1Validator.cc:54
~L1Validator() override
Definition: L1Validator.cc:68
int pdgId() const final
PDG identifier.
void FillNumber(int, int)
Definition: Jet.h:20
int iEvent
Definition: GenABIO.cc:224
void Book(DQMStore::IBooker &, std::string dirname)
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< l1t::TauBxCollection > _L1TauBXSource
Definition: L1Validator.h:86
edm::EDGetTokenT< l1t::JetBxCollection > _L1JetBXSource
Definition: L1Validator.h:87
std::string _dirName
Definition: L1Validator.h:80
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
const int mu
Definition: Constants.h:22
edm::EDGetTokenT< l1t::MuonBxCollection > _L1MuonBXSource
Definition: L1Validator.h:84
Definition: Muon.h:21
constexpr auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:30
tuple GenJets
Definition: jetmet_cfg.py:40
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: L1Validator.cc:75
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: L1Validator.cc:70
T getParameter(std::string const &) const
Definition: ParameterSet.h:303
L1ValidatorHists _Hists
Definition: L1Validator.h:91
constexpr char Jet[]
Definition: modules.cc:9
const reco::LeafCandidate * FindBest(const reco::GenParticle *, const std::vector< l1extra::L1EmParticle > *, const std::vector< l1extra::L1EmParticle > *)
Definition: L1Validator.cc:253
edm::EDGetTokenT< l1t::EGammaBxCollection > _L1EGammaBXSource
Definition: L1Validator.h:85
double phi() const final
momentum azimuthal angle
Definition: Run.h:45
double eta() const final
momentum pseudorapidity