CMS 3D CMS Logo

 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Groups Pages
BasicGenParticleValidation.cc
Go to the documentation of this file.
1 /*class BasicGenParticleValidation
2  *
3  * Class to fill dqm monitor elements from existing EDM file
4  *
5  */
6 
8 
9 #include "CLHEP/Units/defs.h"
10 #include "CLHEP/Units/PhysicalConstants.h"
12 
13 using namespace edm;
14 
16  : wmanager_(iPSet, consumesCollector()),
17  hepmcCollection_(iPSet.getParameter<edm::InputTag>("hepmcCollection")),
18  genparticleCollection_(iPSet.getParameter<edm::InputTag>("genparticleCollection")),
19  genjetCollection_(iPSet.getParameter<edm::InputTag>("genjetsCollection")),
20  matchPr_(iPSet.getParameter<double>("matchingPrecision")),
21  verbosity_(iPSet.getUntrackedParameter<unsigned int>("verbosity", 0)),
22  signalParticlesOnly_(iPSet.getParameter<bool>("signalParticlesOnly")) {
23  hepmcCollectionToken_ = consumes<HepMCProduct>(hepmcCollection_);
24  genparticleCollectionToken_ = consumes<reco::GenParticleCollection>(genparticleCollection_);
25  genjetCollectionToken_ = consumes<reco::GenJetCollection>(genjetCollection_);
26 }
27 
29 
32  DQMHelper dqm(&i);
33  i.setCurrentFolder("Generator/GenParticles");
35 
36  // Number of analyzed events
37  nEvt = dqm.book1dHisto("nEvt", "n analyzed Events", 1, 0., 1., "", "Number of Events");
38 
40  genPMultiplicity = dqm.book1dHisto("genPMultiplicty",
41  "Log(No. all GenParticles)",
42  50,
43  -1,
44  5,
45  "log_{10}(N_{All GenParticles})",
46  "Number of Events"); //Log
47  //difference in HepMC and reco multiplicity
48  genMatched = dqm.book1dHisto(
49  "genMatched", "Difference reco - matched", 50, -25, 25, "N_{All GenParticles}-N_{Matched}", "Number of Events");
50  //multiple matching
51  multipleMatching = dqm.book1dHisto("multipleMatching",
52  "multiple reco HepMC matching",
53  50,
54  0,
55  50,
56  "N_{multiple reco HepMC matching}",
57  "Number of Events");
58  //momentum difference of matched particles
59  matchedResolution = dqm.book1dHisto("matchedResolution",
60  "log10(momentum difference of matched particles)",
61  70,
62  -10.,
63  -3.,
64  "log_{10}(#DeltaP_{matched Particles})",
65  "Number of Events");
66 
67  // GenJet general distributions
68  genJetMult = dqm.book1dHisto("genJetMult", "GenJet multiplicity", 50, 0, 50, "N_{gen-jets}", "Number of Events");
70  "genJetEnergy", "Log10(GenJet energy)", 60, -1, 5, "log_{10}(E^{gen-jets}) (log_{10}(GeV))", "Number of Events");
71  genJetPt = dqm.book1dHisto(
72  "genJetPt", "Log10(GenJet pt)", 60, -1, 5, "log_{10}(P_{t}^{gen-jets}) (log_{10}(GeV))", "Number of Events");
73  genJetEta = dqm.book1dHisto("genJetEta", "GenJet eta", 220, -11, 11, "#eta^{gen-jets}", "Number of Events");
74  genJetPhi = dqm.book1dHisto("genJetPhi", "GenJet phi", 360, -180, 180, "#phi^{gen-jets} (rad)", "Number of Events");
76  "genJetDeltaEtaMin", "GenJet minimum rapidity gap", 30, 0, 30, "#delta#eta_{min}^{gen-jets}", "Number of Events");
77 
78  genJetPto1 = dqm.book1dHisto(
79  "genJetPto1", "GenJet multiplicity above 1 GeV", 50, 0, 50, "N_{gen-jets P_{t}>1GeV}", "Number of Events");
81  "genJetPto10", "GenJet multiplicity above 10 GeV", 50, 0, 50, "N_{gen-jets P_{t}>10GeV}", "Number of Events");
83  "genJetPto100", "GenJet multiplicity above 100 GeV", 50, 0, 50, "N_{gen-jets P_{t}>100GeV}", "Number of Events");
85  "genJetCentral", "GenJet multiplicity |eta|.lt.2.5", 50, 0, 50, "N_{gen-jets |#eta|#leq2.5}", "Number of Events");
86 
87  genJetTotPt = dqm.book1dHisto("genJetTotPt",
88  "Log10(GenJet total pt)",
89  100,
90  -5,
91  5,
92  "log_{10}(#SigmaP_{t}^{gen-jets}) (log_{10}(GeV))",
93  "Number of Events");
94 
95  return;
96 }
97 
99  unsigned int initSize = 1000;
100 
103  iEvent.getByToken(hepmcCollectionToken_, evt);
104 
105  //Get HepMC EVENT
106  HepMC::GenEvent* myGenEvent = new HepMC::GenEvent(*(evt->GetEvent()));
107 
108  double weight = wmanager_.weight(iEvent);
109 
110  nEvt->Fill(0.5, weight);
111 
112  std::vector<const HepMC::GenParticle*> hepmcGPCollection;
113  std::vector<int> barcodeList;
114  hepmcGPCollection.reserve(initSize);
115  barcodeList.reserve(initSize);
116 
117  //Looping through HepMC::GenParticle collection to search for status 1 particles
118  for (HepMC::GenEvent::particle_const_iterator iter = myGenEvent->particles_begin();
119  iter != myGenEvent->particles_end();
120  ++iter) {
121  if ((*iter)->status() == 1) {
122  hepmcGPCollection.push_back(*iter);
123  barcodeList.push_back((*iter)->barcode());
124  if (verbosity_ > 0) {
125  std::cout << "HepMC " << std::setw(14) << std::fixed << (*iter)->pdg_id() << std::setw(14) << std::fixed
126  << (*iter)->momentum().px() << std::setw(14) << std::fixed << (*iter)->momentum().py()
127  << std::setw(14) << std::fixed << (*iter)->momentum().pz() << std::endl;
128  }
129  }
130  }
131 
132  // Gather information on the reco::GenParticle collection
134  iEvent.getByToken(genparticleCollectionToken_, genParticles);
135 
136  std::vector<const reco::GenParticle*> particles;
137  particles.reserve(initSize);
138  for (reco::GenParticleCollection::const_iterator iter = genParticles->begin(); iter != genParticles->end(); ++iter) {
139  if ((*iter).status() == 1) {
140  particles.push_back(&*iter);
141  if (verbosity_ > 0) {
142  std::cout << "reco " << std::setw(14) << std::fixed << (*iter).pdgId() << std::setw(14) << std::fixed
143  << (*iter).px() << std::setw(14) << std::fixed << (*iter).py() << std::setw(14) << std::fixed
144  << (*iter).pz() << std::endl;
145  }
146  }
147  }
148 
149  unsigned int nReco = particles.size();
150  unsigned int nHepMC = hepmcGPCollection.size();
151 
152  if (signalParticlesOnly_) {
153  for (unsigned int i = 0; i < particles.size(); ++i) {
154  if (particles[i]->collisionId() != 0)
155  nReco -= 1;
156  }
157  }
158 
159  genPMultiplicity->Fill(std::log10(nReco), weight);
160 
161  // Define vector containing index of hepmc corresponding to the reco::GenParticle
162  std::vector<int> hepmcMatchIndex;
163  hepmcMatchIndex.reserve(initSize);
164 
165  // Matching procedure
166 
167  // Check array size consistency
168 
169  if (nReco != nHepMC) {
170  edm::LogWarning("CollectionSizeInconsistency")
171  << "Collection size inconsistency: HepMC::GenParticle = " << nHepMC << " reco::GenParticle = " << nReco;
172  }
173 
174  // Match each HepMC with a reco
175 
176  for (unsigned int i = 0; i < nReco; ++i) {
177  for (unsigned int j = 0; j < nHepMC; ++j) {
178  if (matchParticles(hepmcGPCollection[j], particles[i])) {
179  hepmcMatchIndex.push_back((int)j);
180  if (hepmcGPCollection[j]->momentum().rho() != 0.) {
181  double reso = 1. - particles[i]->p() / hepmcGPCollection[j]->momentum().rho();
182  if (verbosity_ > 0) {
183  std::cout << "Matching momentum: reco = " << particles[i]->p()
184  << " HepMC = " << hepmcGPCollection[j]->momentum().rho() << " resoultion = " << reso << std::endl;
185  }
186  matchedResolution->Fill(std::log10(std::fabs(reso)), weight);
187  }
188  continue;
189  }
190  }
191  }
192 
193  // Check unicity of matching
194 
195  unsigned int nMatched = hepmcMatchIndex.size();
196 
197  if (nMatched != nReco) {
198  edm::LogWarning("IncorrectMatching") << "Incorrect number of matched indexes: GenParticle = " << nReco
199  << " matched indexes = " << nMatched;
200  }
201  genMatched->Fill(int(nReco - nMatched), weight);
202 
203  unsigned int nWrMatch = 0;
204 
205  for (unsigned int i = 0; i < nMatched; ++i) {
206  for (unsigned int j = i + 1; j < nMatched; ++j) {
207  if (hepmcMatchIndex[i] == hepmcMatchIndex[j]) {
208  int theIndex = hepmcMatchIndex[i];
209  edm::LogWarning("DuplicatedMatching")
210  << "Multiple matching occurencies for GenParticle barcode = " << barcodeList[theIndex];
211  nWrMatch++;
212  }
213  }
214  }
215  multipleMatching->Fill(int(nWrMatch), weight);
216 
217  // Gather information in the GenJet collection
219  iEvent.getByToken(genjetCollectionToken_, genJets);
220 
221  int nJets = 0;
222  int nJetso1 = 0;
223  int nJetso10 = 0;
224  int nJetso100 = 0;
225  int nJetsCentral = 0;
226  double totPt = 0.;
227 
228  std::vector<double> jetEta;
229  jetEta.reserve(initSize);
230 
231  for (reco::GenJetCollection::const_iterator iter = genJets->begin(); iter != genJets->end(); ++iter) {
232  nJets++;
233  double pt = (*iter).pt();
234  totPt += pt;
235  if (pt > 1.)
236  nJetso1++;
237  if (pt > 10.)
238  nJetso10++;
239  if (pt > 100.)
240  nJetso100++;
241  double eta = (*iter).eta();
242  if (std::fabs(eta) < 2.5)
243  nJetsCentral++;
244  jetEta.push_back(eta);
245 
246  genJetEnergy->Fill(std::log10((*iter).energy()), weight);
247  genJetPt->Fill(std::log10(pt), weight);
248  genJetEta->Fill(eta, weight);
249  genJetPhi->Fill((*iter).phi() / CLHEP::degree, weight);
250  }
251 
252  genJetMult->Fill(nJets, weight);
253  genJetPto1->Fill(nJetso1, weight);
254  genJetPto10->Fill(nJetso10, weight);
255  genJetPto100->Fill(nJetso100, weight);
256  genJetCentral->Fill(nJetsCentral, weight);
257 
258  genJetTotPt->Fill(std::log10(totPt), weight);
259 
260  double deltaEta = 999.;
261  if (jetEta.size() > 1) {
262  for (unsigned int i = 0; i < jetEta.size(); i++) {
263  for (unsigned int j = i + 1; j < jetEta.size(); j++) {
264  deltaEta = std::min(deltaEta, std::fabs(jetEta[i] - jetEta[j]));
265  }
266  }
267  }
268 
269  genJetDeltaEtaMin->Fill(deltaEta, weight);
270 
271  delete myGenEvent;
272 } //analyze
273 
275  bool state = false;
276 
277  if (hepmcP->pdg_id() != recoP->pdgId())
278  return state;
279  if (std::fabs(hepmcP->momentum().px() - recoP->px()) < std::fabs(matchPr_ * hepmcP->momentum().px()) &&
280  std::fabs(hepmcP->momentum().py() - recoP->py()) < std::fabs(matchPr_ * hepmcP->momentum().py()) &&
281  std::fabs(hepmcP->momentum().pz() - recoP->pz()) < std::fabs(matchPr_ * hepmcP->momentum().pz())) {
282  state = true;
283  }
284 
285  return state;
286 }
double pz() const final
z coordinate of momentum vector
virtual void setCurrentFolder(std::string const &fullpath)
Definition: DQMStore.cc:32
bool getByToken(EDGetToken token, Handle< PROD > &result) const
Definition: Event.h:539
theIndex(0)
bool matchParticles(const HepMC::GenParticle *&, const reco::GenParticle *&)
static const double deltaEta
Definition: CaloConstants.h:8
int pdgId() const final
PDG identifier.
void Fill(long long x)
double px() const final
x coordinate of momentum vector
int iEvent
Definition: GenABIO.cc:224
edm::EDGetTokenT< reco::GenJetCollection > genjetCollectionToken_
T min(T a, T b)
Definition: MathUtil.h:58
edm::EDGetTokenT< edm::HepMCProduct > hepmcCollectionToken_
double py() const final
y coordinate of momentum vector
MonitorElement * book1dHisto(const std::string &name, const std::string &title, int n, double xmin, double xmax, const std::string &xaxis, const std::string &yaxis)
Definition: DQMHelper.cc:7
void analyze(const edm::Event &, const edm::EventSetup &) override
void bookHistograms(DQMStore::IBooker &i, edm::Run const &, edm::EventSetup const &) override
tuple cout
Definition: gather_cfg.py:144
BasicGenParticleValidation(const edm::ParameterSet &)
Log< level::Warning, false > LogWarning
int weight
Definition: histoStyle.py:51
double weight(const edm::Event &)
edm::EDGetTokenT< reco::GenParticleCollection > genparticleCollectionToken_
Definition: Run.h:45