CMS 3D CMS Logo

HepMCValidationHelper.cc
Go to the documentation of this file.
4 #include <iostream>
5 #include <limits>
6 #include "TLorentzVector.h"
7 
8 //#define DEBUG_HepMCValidationHelper
9 
10 namespace HepMCValidationHelper {
11  void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
12  const HepMC::GenEvent* all,
13  double deltaRcut,
14  std::vector<const HepMC::GenParticle*>& fsrphotons) {
15  std::vector<const HepMC::GenParticle*> status1;
16  allStatus1(all, status1);
17  findFSRPhotons(leptons, status1, deltaRcut, fsrphotons);
18  }
19 
20  void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
21  const std::vector<const HepMC::GenParticle*>& all,
22  double deltaRcut,
23  std::vector<const HepMC::GenParticle*>& fsrphotons) {
24  //find all status 1 photons
25  std::vector<const HepMC::GenParticle*> allphotons;
26  for (unsigned int i = 0; i < all.size(); ++i) {
27  if (all[i]->status() == 1 && all[i]->pdg_id() == 22)
28  allphotons.push_back(all[i]);
29  }
30 
31  //loop over the photons and check the distance wrt the leptons
32  for (unsigned int ipho = 0; ipho < allphotons.size(); ++ipho) {
33  bool close = false;
34  for (unsigned int ilep = 0; ilep < leptons.size(); ++ilep) {
35  if (deltaR(allphotons[ipho]->momentum(), leptons[ilep]->momentum()) < deltaRcut) {
36  close = true;
37  break;
38  }
39  }
40  if (close)
41  fsrphotons.push_back(allphotons[ipho]);
42  }
43  }
44 
45  //returns true if a status 3 particle is a tau or if a status 1 particle is either an electron or a neutrino
47  int status = part->status();
48  unsigned int pdg_id = abs(part->pdg_id());
49  if (status == 2)
50  return pdg_id == 15;
51  else
52  return status == 1 && (pdg_id == 11 || pdg_id == 13);
53  }
54 
55  //returns true if a status 1 particle is a neutrino
57  int status = part->status();
58  unsigned int pdg_id = abs(part->pdg_id());
59  return status == 1 && (pdg_id == 12 || pdg_id == 14 || pdg_id == 16);
60  }
61 
62  //returns true is status 3 particle is tau
63  bool isTau(const HepMC::GenParticle* part) { return part->status() == 2 && abs(part->pdg_id()) == 15; }
64  /*
65  void getTaus(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& taus){
66  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
67  if (abs((*iter)->pdg_id()) == 15) taus.push_back(*iter);
68  }
69  }
70 */
71  // get all status 1 particles
72  void allStatus1(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1) {
73  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter) {
74  if ((*iter)->status() == 1)
75  status1.push_back(*iter);
76  }
77  }
78 
79  void allStatus2(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1) {
80  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter) {
81  if ((*iter)->status() == 2)
82  status1.push_back(*iter);
83  }
84  }
85 
86  void allStatus3(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1) {
87  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter) {
88  if ((*iter)->status() == 3)
89  status1.push_back(*iter);
90  }
91  }
92 
93  void findDescendents(const HepMC::GenParticle* a, std::vector<const HepMC::GenParticle*>& descendents) {
94  HepMC::GenVertex* decayVertex = a->end_vertex();
95  if (!decayVertex)
96  return;
97  HepMC::GenVertex::particles_out_const_iterator ipart;
98  for (ipart = decayVertex->particles_out_const_begin(); ipart != decayVertex->particles_out_const_end(); ++ipart) {
99  if ((*ipart)->status() == 1)
100  descendents.push_back(*ipart);
101  else
102  findDescendents(*ipart, descendents);
103  }
104  }
105 
107  double deltaRcut,
108  double sumPtCut,
109  std::vector<const HepMC::GenParticle*>& pruned) {
110  //get all status 1 particles
111  std::vector<const HepMC::GenParticle*> status1;
112  allStatus1(all, status1);
113  std::vector<const HepMC::GenParticle*> toRemove;
114  //loop on all particles and find candidates to be isolated
115  for (unsigned int i = 0; i < status1.size(); ++i) {
116  //if it is a neutrino is a charged lepton (not a tau) this is a candidate to be isolated
117  if (isNeutrino(status1[i]) || (isChargedLepton(status1[i]) && !isTau(status1[i]))) {
118  //list of particles not to be considered in the isolation computation.
119  //this includes the particle to be isolated and the fsr photons in case of charged lepton
120  std::vector<const HepMC::GenParticle*> leptons;
121  leptons.push_back(status1[i]);
122  std::vector<const HepMC::GenParticle*> removedForIsolation;
123  removedForIsolation.push_back(status1[i]);
124  if (isChargedLepton(status1[i]))
125  findFSRPhotons(leptons, status1, deltaRcut, removedForIsolation);
126 #ifdef DEBUG_HepMCValidationHelper
127  //std::cout << removedForIsolation.size() << " particles to be removed for isolation calculation " << std::endl;
128 #endif
129  //create vector of particles to compute isolation (removing removedForIsolation);
130  std::vector<const HepMC::GenParticle*> forIsolation;
131  std::vector<const HepMC::GenParticle*>::iterator iiso;
132  for (iiso = status1.begin(); iiso != status1.end(); ++iiso) {
133  std::vector<const HepMC::GenParticle*>::const_iterator iremove;
134  bool marked = false;
135  for (iremove = removedForIsolation.begin(); iremove != removedForIsolation.end(); ++iremove) {
136  if ((*iiso)->barcode() == (*iremove)->barcode()) {
137 #ifdef DEBUG_HepMCValidationHelper
138  //std::cout << "removing particle " << **iiso << " from the list of particles to compute isolation" << std::endl;
139 #endif
140  marked = true;
141  break;
142  }
143  }
144  if (!marked)
145  forIsolation.push_back(*iiso);
146  }
147  //now compute isolation
148  double sumIso = 0;
149  for (iiso = forIsolation.begin(); iiso < forIsolation.end(); ++iiso) {
150  if (deltaR(leptons.front()->momentum(), (*iiso)->momentum()) < deltaRcut) {
151  sumIso += (*iiso)->momentum().perp();
152  }
153  }
154  //if isolated remove from the pruned list
155  if (sumIso < sumPtCut) {
156 #ifdef DEBUG_HepMCValidationHelper
157  std::cout << "particle " << *status1[i] << " is considered isolated, with sumPt " << sumIso << std::endl;
158 #endif
159  toRemove.insert(toRemove.end(), removedForIsolation.begin(), removedForIsolation.end());
160  }
161 #ifdef DEBUG_HepMCValidationHelper
162  else {
163  std::cout << "NOT isolated! " << *status1[i] << " is considered not isolated, with sumPt " << sumIso
164  << std::endl;
165  }
166 #endif
167  }
168  }
169  //at this point we have taken care of the electrons and muons, but pruned could still contain the decay products of isolated taus,
170  //we want to remove these as well
171  std::vector<const HepMC::GenParticle*> status2;
172  allStatus2(all, status2);
173  std::vector<const HepMC::GenParticle*> taus;
174  //getTaus(all, taus);
175  for (unsigned int i = 0; i < status2.size(); ++i) {
176  if (isTau(status2[i])) {
177  //check the list we have already for duplicates
178  //there use to be duplicates in some generators (sherpa)
179  bool duplicate = false;
180  TLorentzVector taumomentum(status2[i]->momentum().x(),
181  status2[i]->momentum().y(),
182  status2[i]->momentum().z(),
183  status2[i]->momentum().t());
184  for (unsigned int j = 0; j < taus.size(); ++j) {
185  //compare momenta
186  TLorentzVector othermomentum(
187  taus[j]->momentum().x(), taus[j]->momentum().y(), taus[j]->momentum().z(), taus[j]->momentum().t());
188  othermomentum -= taumomentum;
189  if (status2[i]->pdg_id() == taus[j]->pdg_id() &&
190  othermomentum.E() < 0.1 && //std::numeric_limits<float>::epsilon() &&
191  othermomentum.P() < 0.1) { //std::numeric_limits<float>::epsilon()){
192  duplicate = true;
193  break;
194  }
195  }
196  if (!duplicate)
197  taus.push_back(status2[i]);
198  }
199  }
200  //loop over the taus, find the descendents, remove all these from the list of particles to compute isolation
201  for (unsigned int i = 0; i < taus.size(); ++i) {
202  std::vector<const HepMC::GenParticle*> taudaughters;
203  findDescendents(taus[i], taudaughters);
204  if (taudaughters.empty()) {
205  edm::LogError("HepMCValidationHelper") << "Tau with no daughters. This is a bug. Fix it";
206  abort();
207  }
208  const HepMC::FourVector& taumom = taus[i]->momentum();
209  //remove the daughters from the list of particles to compute isolation
210  std::vector<const HepMC::GenParticle*> forIsolation;
211  std::vector<const HepMC::GenParticle*>::iterator iiso;
212  for (iiso = status1.begin(); iiso < status1.end(); ++iiso) {
213  bool marked = false;
214  std::vector<const HepMC::GenParticle*>::const_iterator iremove;
215  for (iremove = taudaughters.begin(); iremove != taudaughters.end(); ++iremove) {
216  if ((*iiso)->barcode() == (*iremove)->barcode()) {
217 #ifdef DEBUG_HepMCValidationHelper
218 // std::cout << "removing particle " << **iiso << " from the list of particles to compute isolation because it comes from a tau" << std::endl;
219 #endif
220  marked = true;
221  break;
222  }
223  }
224  if (!marked)
225  forIsolation.push_back(*iiso);
226  }
227  //no compute isolation wrt the status 2 tau direction
228  double sumIso = 0;
229  for (iiso = forIsolation.begin(); iiso < forIsolation.end(); ++iiso) {
230  if (deltaR(taumom, (*iiso)->momentum()) < deltaRcut) {
231  sumIso += (*iiso)->momentum().perp();
232  }
233  }
234  //if isolated remove the tau daughters from the pruned list
235  if (sumIso < sumPtCut) {
236 #ifdef DEBUG_HepMCValidationHelper
237  std::cout << "particle " << *taus[i] << " is considered isolated, with sumPt " << sumIso << std::endl;
238 #endif
239  toRemove.insert(toRemove.end(), taudaughters.begin(), taudaughters.end());
240  }
241  }
242 
243  //now actually remove
244  pruned.clear();
245  for (unsigned int i = 0; i < status1.size(); ++i) {
246  bool marked = false;
247  std::vector<const HepMC::GenParticle*>::const_iterator iremove;
248  for (iremove = toRemove.begin(); iremove != toRemove.end(); ++iremove) {
249  if (status1[i]->barcode() == (*iremove)->barcode()) {
250  marked = true;
251  break;
252  }
253  }
254  if (!marked)
255  pruned.push_back(status1[i]);
256  }
257 
258 #ifdef DEBUG_HepMCValidationHelper
259  std::cout << "list of remaining particles:" << std::endl;
260  for (unsigned int i = 0; i < pruned.size(); ++i) {
261  std::cout << *pruned[i] << std::endl;
262  }
263 #endif
264  }
265 
266  //get all visible status1 particles
267  void allVisibleParticles(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& visible) {
268  std::vector<const HepMC::GenParticle*> status1;
269  visible.clear();
270  allStatus1(all, status1);
271  for (unsigned int i = 0; i < status1.size(); ++i) {
272  if (!isNeutrino(status1[i]))
273  visible.push_back(status1[i]);
274  }
275  }
276 
277  //compute generated met
278  TLorentzVector genMet(const HepMC::GenEvent* all, double etamin, double etamax) {
279  std::vector<const HepMC::GenParticle*> visible;
280  allVisibleParticles(all, visible);
281  TLorentzVector momsum(0., 0., 0., 0.);
282  for (unsigned int i = 0; i < visible.size(); ++i) {
283  if (visible[i]->momentum().eta() > etamin && visible[i]->momentum().eta() < etamax) {
284  TLorentzVector mom(visible[i]->momentum().x(),
285  visible[i]->momentum().y(),
286  visible[i]->momentum().z(),
287  visible[i]->momentum().t());
288  momsum += mom;
289  }
290  }
291  TLorentzVector met(-momsum.Px(), -momsum.Py(), 0., momsum.Pt());
292  return met;
293  }
294 } // namespace HepMCValidationHelper
bool isTau(const HepMC::GenParticle *part)
void allStatus3(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status3)
void findFSRPhotons(const std::vector< const HepMC::GenParticle * > &leptons, const std::vector< const HepMC::GenParticle * > &all, double deltaR, std::vector< const HepMC::GenParticle * > &photons)
bool isChargedLepton(const HepMC::GenParticle *part)
void allStatus1(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status1)
void removeIsolatedLeptons(const HepMC::GenEvent *all, double deltaR, double sumPt, std::vector< const HepMC::GenParticle * > &pruned)
void allStatus2(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status2)
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
void allVisibleParticles(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &visible)
void findDescendents(const HepMC::GenParticle *a, std::vector< const HepMC::GenParticle * > &descendents)
met
===> hadronic RAZOR
part
Definition: HCALResponse.h:20
TLorentzVector genMet(const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
bool isNeutrino(const HepMC::GenParticle *part)
double a
Definition: hdecay.h:121