CMS 3D CMS Logo

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