CMS 3D CMS Logo

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 
68 
70 }
71 
73  //iBooker.setCurrentFolder(_dirName);
74  _Hists.Book(iBooker, _dirName);
75 };
76 
78  using namespace edm;
79  using namespace std;
80  using namespace l1extra;
81  using namespace reco;
82 
88  Handle<GenEventInfoProduct> genEvtInfoProduct;
90 
91  bool GotEverything=true;
92 
93  if(!iEvent.getByToken(_GenSource, GenParticles)) GotEverything=false;
94  if(!iEvent.getByToken(_L1MuonBXSource, MuonsBX)) GotEverything=false;
95  if(!iEvent.getByToken(_L1EGammaBXSource, EGammasBX)) GotEverything=false;
96  if(!iEvent.getByToken(_L1TauBXSource, TausBX)) GotEverything=false;
97  if(!iEvent.getByToken(_L1JetBXSource, JetsBX)) GotEverything=false;
98  if(!iEvent.getByToken(_srcToken, genEvtInfoProduct)) GotEverything=false;
99  if(!iEvent.getByToken(_L1GenJetSource, GenJets)) GotEverything=false;
100 
101  if(!GotEverything) return;
102 
103  /*
104  std::string moduleName = "";
105  if( genEvtInfoProduct.isValid() ) {
106  const edm::Provenance& prov = iEvent.getProvenance(genEvtInfoProduct.id());
107  moduleName = edm::moduleName(prov);
108  //cout<<" generator name: "<<moduleName<<endl;
109  }
110  */
111 
112  _Hists.NEvents++;
113 
114  int nL1Muons = 0, nL1EGammas = 0, nL1Taus = 0, nL1Jets = 0;
115  if(MuonsBX->getFirstBX()>=0) nL1Muons = MuonsBX->size(0);
116  if(EGammasBX->getFirstBX()>=0) nL1EGammas = EGammasBX->size(0);
117  if(TausBX->getFirstBX()>=0) nL1Taus = TausBX->size(0);
118  if(JetsBX->getFirstBX()>=0) nL1Jets = JetsBX->size(0);
119 
121  _Hists.FillNumber(L1ValidatorHists::Type::Egamma, nL1EGammas);
124 
125  //For gen jet
126 
127  for(auto &Genjet : *GenJets ){
128 
129  // eta within calorimeter acceptance 4.7
130  if(fabs((&Genjet)->eta())>4.7) continue;
131 
132  // only consider the gen jet with pt greater than 10 GeV
133  if((&Genjet)->pt()<10.0) continue;
134 
135  double minDR = 999.0;
136 
137  // match L1T object
138  const l1t::Jet *L1Part=nullptr;
139  for(int iBx = JetsBX->getFirstBX(); iBx<=JetsBX->getLastBX(); ++iBx){
140  if(iBx>0) continue;
141  for(std::vector<l1t::Jet>::const_iterator jet = JetsBX->begin(iBx); jet != JetsBX->end(iBx); ++jet){
142  double idR = reco::deltaR((&Genjet)->eta(), (&Genjet)->phi(), jet->eta(), jet->phi());
143  if( idR < minDR ){
144  minDR = idR;
145  L1Part = &(*jet);
146  }
147  }
148  }
149  _Hists.Fill(L1ValidatorHists::Type::Jet, &Genjet, L1Part);
150  }
151 
152 
153  for(uint i=0; i < GenParticles->size(); i++){
154  const GenParticle *GenPart = &GenParticles->at(i);
155 
156  int pdg = GenPart->pdgId(), status = GenPart->status();
157 
158  double minDR = 999.0;
159 
160  // only consider the gen particle with pt greater than 10 GeV
161  if(GenPart->pt()<10.0) continue;
162 
164  if(status==1 && abs(pdg)==13){ //Muon
165 
166  // eta within tracker acceptance 2.4
167  if(fabs(GenPart->eta())>2.4) continue;
168 
169  // match L1T object
170  const l1t::Muon *L1Part=nullptr;
171  for(int iBx = MuonsBX->getFirstBX(); iBx<=MuonsBX->getLastBX(); ++iBx){
172  if(iBx>0) continue;
173  for(std::vector<l1t::Muon>::const_iterator mu = MuonsBX->begin(iBx); mu != MuonsBX->end(iBx); ++mu){
174  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), mu->eta(), mu->phi());
175  if(idR < minDR ){
176  minDR = idR;
177  L1Part = &(*mu);
178  }
179 
180  }
181  _Hists.Fill(L1ValidatorHists::Type::Muon, GenPart, L1Part);
182  }
183 
184 
186  } else if(status==1 && (abs(pdg)==11 || pdg==22)){ //Egamma
187 
188  // eta within EM calorimeter acceptance 2.5
189  if(fabs(GenPart->eta())>2.5) continue;
190 
191  // exclude the calorimeter barrel and endcap overlap region
192  if(fabs(GenPart->eta())>1.4442 && fabs(GenPart->eta())<1.5660) continue;
193 
194  // match L1T object
195  const l1t::EGamma *L1Part=nullptr;
196  for(int iBx = EGammasBX->getFirstBX(); iBx<=EGammasBX->getLastBX(); ++iBx){
197  if(iBx>0) continue;
198  for(std::vector<l1t::EGamma>::const_iterator eg = EGammasBX->begin(iBx); eg != EGammasBX->end(iBx); ++eg){
199  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), eg->eta(), eg->phi());
200  if(idR < minDR ){
201  minDR = idR;
202  L1Part = &(*eg);
203  }
204  }
205  }
206  _Hists.Fill(L1ValidatorHists::Type::Egamma, GenPart, L1Part);
207 
208 
210  } else if(status==2 && abs(pdg)==15){ //Tau
211 
212  // eta within tracker acceptance 2.4
213  if(fabs(GenPart->eta())>2.4) continue;
214 
215  // match L1T object
216  const l1t::Tau *L1Part=nullptr;
217  for(int iBx = TausBX->getFirstBX(); iBx<=TausBX->getLastBX(); ++iBx){
218  if(iBx>0) continue;
219  for(std::vector<l1t::Tau>::const_iterator tau = TausBX->begin(iBx); tau != TausBX->end(iBx); ++tau){
220  double idR = reco::deltaR(GenPart->eta(), GenPart->phi(), tau->eta(), tau->phi());
221  if(idR < minDR ){
222  minDR = idR;
223  L1Part = &(*tau);
224  }
225  }
226  }
227  _Hists.Fill(L1ValidatorHists::Type::Tau, GenPart, L1Part);
228 
229  }
230 
231  }
232 
233 }
234 
235 //The next three are exactly the same, but apparently inheritance doesn't work like I thought it did.
236 const reco::LeafCandidate *L1Validator::FindBest(const reco::GenParticle *GenPart, const std::vector<l1extra::L1EmParticle> *Collection1, const std::vector<l1extra::L1EmParticle> *Collection2=nullptr){
237  const reco::LeafCandidate *BestPart=nullptr;
238  double BestDR=999.;
239 
240  for(uint i=0; i < Collection1->size(); i++){
241  const reco::LeafCandidate *ThisPart = &Collection1->at(i);
242  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
243  if(ThisDR < BestDR){
244  BestDR = ThisDR;
245  BestPart = ThisPart;
246  }
247  }
248 
249  if(Collection2==nullptr) return BestPart;
250 
251  for(uint i=0; i < Collection2->size(); i++){
252  const reco::LeafCandidate *ThisPart = &Collection2->at(i);
253  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
254  if(ThisDR < BestDR){
255  BestDR = ThisDR;
256  BestPart = ThisPart;
257  }
258  }
259 
260  return BestPart;
261 }
262 
263 const reco::LeafCandidate *L1Validator::FindBest(const reco::GenParticle *GenPart, const std::vector<l1extra::L1JetParticle> *Collection1, const std::vector<l1extra::L1JetParticle> *Collection2=nullptr){
264  const reco::LeafCandidate *BestPart=nullptr;
265  double BestDR=999.;
266 
267  for(uint i=0; i < Collection1->size(); i++){
268  const reco::LeafCandidate *ThisPart = &Collection1->at(i);
269  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
270  if(ThisDR < BestDR){
271  BestDR = ThisDR;
272  BestPart = ThisPart;
273  }
274  }
275 
276  if(Collection2==nullptr) return BestPart;
277 
278  for(uint i=0; i < Collection2->size(); i++){
279  const reco::LeafCandidate *ThisPart = &Collection2->at(i);
280  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
281  if(ThisDR < BestDR){
282  BestDR = ThisDR;
283  BestPart = ThisPart;
284  }
285  }
286 
287  return BestPart;
288 }
289 
290 const reco::LeafCandidate *L1Validator::FindBest(const reco::GenParticle *GenPart, const std::vector<l1extra::L1MuonParticle> *Collection1){
291  const reco::LeafCandidate *BestPart=nullptr;
292  double BestDR=999.;
293 
294  for(uint i=0; i < Collection1->size(); i++){
295  const reco::LeafCandidate *ThisPart = &Collection1->at(i);
296  double ThisDR = reco::deltaR(GenPart->eta(), GenPart->phi(), ThisPart->eta(), ThisPart->phi());
297  if(ThisDR < BestDR){
298  BestDR = ThisDR;
299  BestPart = ThisPart;
300  }
301  }
302 
303  return BestPart;
304 }
305 
306 
307 // ------------ method fills 'descriptions' with the allowed parameters for the module ------------
309  //The following says we do not know what parameters are allowed so do no validation
310  // Please change this to state exactly what you do use, even if it is no parameters
312  desc.setUnknown();
313  descriptions.addDefault(desc);
314 }
315 
316 //define this as a plug-in
void Fill(int, const reco::LeafCandidate *, const reco::LeafCandidate *)
const_iterator end(int bx) const
T getParameter(std::string const &) const
edm::EDGetTokenT< reco::GenJetCollection > _L1GenJetSource
Definition: L1Validator.h:90
int pdgId() const final
PDG identifier.
unsigned size(int bx) const
double eta() const final
momentum pseudorapidity
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:519
#define DEFINE_FWK_MODULE(type)
Definition: MakerMacros.h:17
edm::EDGetTokenT< reco::GenParticleCollection > _GenSource
Definition: L1Validator.h:84
Definition: Tau.h:21
edm::EDGetTokenT< GenEventInfoProduct > _srcToken
Definition: L1Validator.h:89
static void fillDescriptions(edm::ConfigurationDescriptions &descriptions)
Definition: L1Validator.cc:308
double pt() const final
transverse momentum
L1Validator(const edm::ParameterSet &)
Definition: L1Validator.cc:54
~L1Validator() override
Definition: L1Validator.cc:69
void FillNumber(int, int)
Definition: Jet.h:21
int iEvent
Definition: GenABIO.cc:230
void Book(DQMStore::IBooker &, std::string dirname)
void addDefault(ParameterSetDescription const &psetDescription)
edm::EDGetTokenT< l1t::TauBxCollection > _L1TauBXSource
Definition: L1Validator.h:87
edm::EDGetTokenT< l1t::JetBxCollection > _L1JetBXSource
Definition: L1Validator.h:88
std::string _dirName
Definition: L1Validator.h:81
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
auto deltaR(const T1 &t1, const T2 &t2) -> decltype(t1.eta())
Definition: deltaR.h:28
const int mu
Definition: Constants.h:22
edm::EDGetTokenT< l1t::MuonBxCollection > _L1MuonBXSource
Definition: L1Validator.h:85
Definition: Muon.h:21
void analyze(const edm::Event &, const edm::EventSetup &) override
Definition: L1Validator.cc:77
void bookHistograms(DQMStore::IBooker &, edm::Run const &, edm::EventSetup const &) override
Definition: L1Validator.cc:72
int getFirstBX() const
def uint(string)
L1ValidatorHists _Hists
Definition: L1Validator.h:92
const reco::LeafCandidate * FindBest(const reco::GenParticle *, const std::vector< l1extra::L1EmParticle > *, const std::vector< l1extra::L1EmParticle > *)
Definition: L1Validator.cc:236
fixed size matrix
HLT enums.
int getLastBX() const
int status() const final
status word
edm::EDGetTokenT< l1t::EGammaBxCollection > _L1EGammaBXSource
Definition: L1Validator.h:86
const_iterator begin(int bx) const
double phi() const final
momentum azimuthal angle
Definition: Run.h:43