CMS 3D CMS Logo

Functions
HepMCValidationHelper Namespace Reference

Functions

void allStatus1 (const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status1)
 
void allStatus2 (const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status2)
 
void allStatus3 (const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status3)
 
void allVisibleParticles (const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &visible)
 
void findDescendents (const HepMC::GenParticle *a, std::vector< const HepMC::GenParticle * > &descendents)
 
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 findFSRPhotons (const std::vector< const HepMC::GenParticle * > &leptons, const HepMC::GenEvent *all, double deltaR, std::vector< const HepMC::GenParticle * > &photons)
 
TLorentzVector genMet (const HepMC::GenEvent *all, double etamin=-9999., double etamax=9999.)
 
template<class T >
bool GreaterByE (const T &a1, const T &a2)
 
bool isChargedLepton (const HepMC::GenParticle *part)
 
bool isNeutrino (const HepMC::GenParticle *part)
 
bool isTau (const HepMC::GenParticle *part)
 
void removeIsolatedLeptons (const HepMC::GenEvent *all, double deltaR, double sumPt, std::vector< const HepMC::GenParticle * > &pruned)
 
bool sortByE (const HepMC::GenParticle *a, const HepMC::GenParticle *b)
 
bool sortByPseudoRapidity (const HepMC::GenParticle *a, const HepMC::GenParticle *b)
 
bool sortByPt (const HepMC::GenParticle *a, const HepMC::GenParticle *b)
 
bool sortByRapidity (const HepMC::GenParticle *a, const HepMC::GenParticle *b)
 

Function Documentation

void HepMCValidationHelper::allStatus1 ( const HepMC::GenEvent *  all,
std::vector< const HepMC::GenParticle * > &  status1 
)

Definition at line 72 of file HepMCValidationHelper.cc.

Referenced by allVisibleParticles(), MBUEandQCDValidation::analyze(), findFSRPhotons(), removeIsolatedLeptons(), and sortByPseudoRapidity().

72  {
73  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
74  if ((*iter)->status() == 1) status1.push_back(*iter);
75  }
76  }
void HepMCValidationHelper::allStatus2 ( const HepMC::GenEvent *  all,
std::vector< const HepMC::GenParticle * > &  status2 
)

Definition at line 78 of file HepMCValidationHelper.cc.

Referenced by removeIsolatedLeptons(), and sortByPseudoRapidity().

78  {
79  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
80  if ((*iter)->status() == 2) status1.push_back(*iter);
81  }
82  }
void HepMCValidationHelper::allStatus3 ( const HepMC::GenEvent *  all,
std::vector< const HepMC::GenParticle * > &  status3 
)

Definition at line 84 of file HepMCValidationHelper.cc.

Referenced by sortByPseudoRapidity().

84  {
85  for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
86  if ((*iter)->status() == 3) status1.push_back(*iter);
87  }
88  }
void HepMCValidationHelper::allVisibleParticles ( const HepMC::GenEvent *  all,
std::vector< const HepMC::GenParticle * > &  visible 
)

Definition at line 248 of file HepMCValidationHelper.cc.

References allStatus1(), mps_fire::i, and isNeutrino().

Referenced by genMet(), and sortByPseudoRapidity().

248  {
249  std::vector<const HepMC::GenParticle*> status1;
250  visible.clear();
251  allStatus1(all, status1);
252  for (unsigned int i = 0; i < status1.size(); ++i){
253  if (!isNeutrino(status1[i])) visible.push_back(status1[i]);
254  }
255  }
void allStatus1(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status1)
bool isNeutrino(const Candidate &part)
Definition: pdgIdUtils.h:25
void HepMCValidationHelper::findDescendents ( const HepMC::GenParticle *  a,
std::vector< const HepMC::GenParticle * > &  descendents 
)

Definition at line 90 of file HepMCValidationHelper.cc.

Referenced by removeIsolatedLeptons(), and sortByPseudoRapidity().

90  {
91  HepMC::GenVertex* decayVertex = a->end_vertex();
92  if (!decayVertex) return;
93  HepMC::GenVertex::particles_out_const_iterator ipart;
94  for (ipart = decayVertex->particles_out_const_begin(); ipart != decayVertex->particles_out_const_end(); ++ipart ){
95  if ((*ipart)->status() == 1) descendents.push_back(*ipart);
96  else findDescendents(*ipart, descendents);
97  }
98  }
double a
Definition: hdecay.h:121
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)
void HepMCValidationHelper::findFSRPhotons ( const std::vector< const HepMC::GenParticle * > &  leptons,
const std::vector< const HepMC::GenParticle * > &  all,
double  deltaR,
std::vector< const HepMC::GenParticle * > &  photons 
)

Definition at line 21 of file HepMCValidationHelper.cc.

References deltaR(), mps_fire::i, and mps_update::status.

Referenced by WValidation::analyze(), DrellYanValidation::analyze(), findFSRPhotons(), removeIsolatedLeptons(), and sortByPseudoRapidity().

24  {
25 
26  //find all status 1 photons
27  std::vector<const HepMC::GenParticle*> allphotons;
28  for (unsigned int i = 0; i < all.size(); ++i){
29  if (all[i]->status()==1 && all[i]->pdg_id()==22) 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) fsrphotons.push_back(allphotons[ipho]);
42  }
43  }
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void HepMCValidationHelper::findFSRPhotons ( const std::vector< const HepMC::GenParticle * > &  leptons,
const HepMC::GenEvent *  all,
double  deltaR,
std::vector< const HepMC::GenParticle * > &  photons 
)

Definition at line 11 of file HepMCValidationHelper.cc.

References allStatus1(), and findFSRPhotons().

14  {
15 
16  std::vector<const HepMC::GenParticle*> status1;
17  allStatus1(all, status1);
18  findFSRPhotons(leptons, status1, deltaRcut, fsrphotons);
19  }
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 allStatus1(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status1)
TLorentzVector HepMCValidationHelper::genMet ( const HepMC::GenEvent *  all,
double  etamin = -9999.,
double  etamax = 9999. 
)

Definition at line 258 of file HepMCValidationHelper.cc.

References allVisibleParticles(), mps_fire::i, RazorAnalyzer::met, and lumiQTWidget::t.

Referenced by GenSpecificAlgo::addInfo(), WValidation::analyze(), HLTJetMETValidation::analyze(), pat::GenMETExtractor::produce(), HLTTauMCProducer::produce(), l1t::GenToInputProducer::produce(), and sortByPseudoRapidity().

258  {
259  std::vector<const HepMC::GenParticle*> visible;
260  allVisibleParticles(all, visible);
261  TLorentzVector momsum(0., 0., 0., 0.);
262  for (unsigned int i = 0; i < visible.size(); ++i){
263  if (visible[i]->momentum().eta() > etamin && visible[i]->momentum().eta() < etamax) {
264  TLorentzVector mom(visible[i]->momentum().x(), visible[i]->momentum().y(), visible[i]->momentum().z(), visible[i]->momentum().t());
265  momsum += mom;
266  }
267  }
268  TLorentzVector met(-momsum.Px(), -momsum.Py(), 0., momsum.Pt());
269  return met;
270  }
float float float z
void allVisibleParticles(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &visible)
met
===> hadronic RAZOR
template<class T >
bool HepMCValidationHelper::GreaterByE ( const T a1,
const T a2 
)
inline

Definition at line 9 of file HepMCValidationHelper.h.

9 {return a1.E() > a2.E();}
bool HepMCValidationHelper::isChargedLepton ( const HepMC::GenParticle *  part)

Definition at line 46 of file HepMCValidationHelper.cc.

References funct::abs(), HiggsValidation_cfi::pdg_id, and mps_update::status.

Referenced by removeIsolatedLeptons(), and sortByPseudoRapidity().

46  {
47  int status = part->status();
48  unsigned int pdg_id = abs(part->pdg_id());
49  if (status == 2) return pdg_id == 15;
50  else return status == 1 && (pdg_id == 11 || pdg_id == 13 ) ;
51  }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
part
Definition: HCALResponse.h:20
bool HepMCValidationHelper::isNeutrino ( const HepMC::GenParticle *  part)

Definition at line 54 of file HepMCValidationHelper.cc.

References funct::abs(), HiggsValidation_cfi::pdg_id, and mps_update::status.

Referenced by allVisibleParticles(), removeIsolatedLeptons(), and sortByPseudoRapidity().

54  {
55  int status = part->status();
56  unsigned int pdg_id = abs(part->pdg_id());
57  return status == 1 && (pdg_id == 12 || pdg_id == 14 || pdg_id == 16) ;
58  }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
part
Definition: HCALResponse.h:20
bool HepMCValidationHelper::isTau ( const HepMC::GenParticle *  part)

Definition at line 61 of file HepMCValidationHelper.cc.

References funct::abs().

Referenced by removeIsolatedLeptons(), and sortByPseudoRapidity().

61  {
62  return part->status() == 2 && abs(part->pdg_id()) == 15;
63  }
Abs< T >::type abs(const T &t)
Definition: Abs.h:22
part
Definition: HCALResponse.h:20
void HepMCValidationHelper::removeIsolatedLeptons ( const HepMC::GenEvent *  all,
double  deltaR,
double  sumPt,
std::vector< const HepMC::GenParticle * > &  pruned 
)

Definition at line 100 of file HepMCValidationHelper.cc.

References allStatus1(), allStatus2(), gather_cfg::cout, deltaR(), findDescendents(), findFSRPhotons(), mps_fire::i, isChargedLepton(), isNeutrino(), isTau(), AK4GenJetFlavourInfos_cfi::leptons, HiggsValidation_cfi::pdg_id, lumiQTWidget::t, and nano_cff::taus.

Referenced by MBUEandQCDValidation::analyze(), and sortByPseudoRapidity().

100  {
101  //get all status 1 particles
102  std::vector<const HepMC::GenParticle*> status1;
103  allStatus1(all, status1);
104  std::vector<const HepMC::GenParticle*> toRemove;
105  //loop on all particles and find candidates to be isolated
106  for (unsigned int i = 0; i < status1.size(); ++i){
107  //if it is a neutrino is a charged lepton (not a tau) this is a candidate to be isolated
108  if (isNeutrino(status1[i]) || (isChargedLepton(status1[i]) && !isTau(status1[i]))){
109  //list of particles not to be considered in the isolation computation.
110  //this includes the particle to be isolated and the fsr photons in case of charged lepton
111  std::vector<const HepMC::GenParticle*> leptons;
112  leptons.push_back(status1[i]);
113  std::vector<const HepMC::GenParticle*> removedForIsolation;
114  removedForIsolation.push_back(status1[i]);
115  if (isChargedLepton(status1[i])) findFSRPhotons(leptons, status1, deltaRcut, removedForIsolation);
116 #ifdef DEBUG_HepMCValidationHelper
117  //std::cout << removedForIsolation.size() << " particles to be removed for isolation calculation " << std::endl;
118 #endif
119  //create vector of particles to compute isolation (removing removedForIsolation);
120  std::vector<const HepMC::GenParticle*> forIsolation;
121  std::vector<const HepMC::GenParticle*>::iterator iiso;
122  for (iiso = status1.begin(); iiso != status1.end(); ++iiso){
123  std::vector<const HepMC::GenParticle*>::const_iterator iremove;
124  bool marked = false;
125  for ( iremove = removedForIsolation.begin(); iremove != removedForIsolation.end(); ++iremove){
126  if ((*iiso)->barcode() == (*iremove)->barcode()) {
127 #ifdef DEBUG_HepMCValidationHelper
128  //std::cout << "removing particle " << **iiso << " from the list of particles to compute isolation" << std::endl;
129 #endif
130  marked = true;
131  break;
132  }
133  }
134  if (!marked) forIsolation.push_back(*iiso);
135  }
136  //now compute isolation
137  double sumIso = 0;
138  for (iiso = forIsolation.begin(); iiso < forIsolation.end(); ++iiso){
139  if (deltaR(leptons.front()->momentum(), (*iiso)->momentum()) < deltaRcut){
140  sumIso += (*iiso)->momentum().perp();
141  }
142  }
143  //if isolated remove from the pruned list
144  if (sumIso < sumPtCut) {
145 #ifdef DEBUG_HepMCValidationHelper
146  std::cout << "particle " << *status1[i] << " is considered isolated, with sumPt " << sumIso << std::endl;
147 #endif
148  toRemove.insert(toRemove.end(), removedForIsolation.begin(), removedForIsolation.end());
149  }
150 #ifdef DEBUG_HepMCValidationHelper
151  else {
152  std::cout << "NOT isolated! " << *status1[i] << " is considered not isolated, with sumPt " << sumIso << std::endl;
153  }
154 #endif
155  }
156  }
157  //at this point we have taken care of the electrons and muons, but pruned could still contain the decay products of isolated taus,
158  //we want to remove these as well
159  std::vector<const HepMC::GenParticle*> status2;
160  allStatus2(all, status2);
161  std::vector<const HepMC::GenParticle*> taus;
162  //getTaus(all, taus);
163  for (unsigned int i = 0; i < status2.size(); ++i){
164  if (isTau(status2[i])) {
165  //check the list we have already for duplicates
166  //there use to be duplicates in some generators (sherpa)
167  bool duplicate = false;
168  TLorentzVector taumomentum(status2[i]->momentum().x(), status2[i]->momentum().y(), status2[i]->momentum().z(), status2[i]->momentum().t());
169  for (unsigned int j = 0; j < taus.size(); ++j){
170  //compare momenta
171  TLorentzVector othermomentum(taus[j]->momentum().x(), taus[j]->momentum().y(), taus[j]->momentum().z(), taus[j]->momentum().t());
172  othermomentum-=taumomentum;
173  if (status2[i]->pdg_id() == taus[j]->pdg_id() &&
174  othermomentum.E() < 0.1 &&//std::numeric_limits<float>::epsilon() &&
175  othermomentum.P() < 0.1 ){ //std::numeric_limits<float>::epsilon()){
176  duplicate = true;
177  break;
178  }
179  }
180  if (!duplicate) taus.push_back(status2[i]);
181  }
182  }
183  //loop over the taus, find the descendents, remove all these from the list of particles to compute isolation
184  for (unsigned int i = 0; i < taus.size(); ++i){
185  std::vector<const HepMC::GenParticle*> taudaughters;
186  findDescendents(taus[i], taudaughters);
187  if ( taudaughters.empty() ) {
188  edm::LogError("HepMCValidationHelper") << "Tau with no daughters. This is a bug. Fix it";
189  abort();
190  }
191  const HepMC::FourVector& taumom = taus[i]->momentum();
192  //remove the daughters from the list of particles to compute isolation
193  std::vector<const HepMC::GenParticle*> forIsolation;
194  std::vector<const HepMC::GenParticle*>::iterator iiso;
195  for (iiso = status1.begin(); iiso < status1.end(); ++iiso){
196  bool marked = false;
197  std::vector<const HepMC::GenParticle*>::const_iterator iremove;
198  for ( iremove = taudaughters.begin(); iremove != taudaughters.end(); ++iremove){
199  if ((*iiso)->barcode() == (*iremove)->barcode()) {
200 #ifdef DEBUG_HepMCValidationHelper
201 // std::cout << "removing particle " << **iiso << " from the list of particles to compute isolation because it comes from a tau" << std::endl;
202 #endif
203  marked = true;
204  break;
205  }
206  if (!marked) forIsolation.push_back(*iiso);
207  }
208  }
209  //no compute isolation wrt the status 2 tau direction
210  double sumIso = 0;
211  for (iiso = forIsolation.begin(); iiso < forIsolation.end(); ++iiso){
212  if (deltaR(taumom, (*iiso)->momentum()) < deltaRcut){
213  sumIso += (*iiso)->momentum().perp();
214  }
215  }
216  //if isolated remove the tau daughters from the pruned list
217  if (sumIso < sumPtCut) {
218 #ifdef DEBUG_HepMCValidationHelper
219  std::cout << "particle " << *taus[i] << " is considered isolated, with sumPt " << sumIso << std::endl;
220 #endif
221  toRemove.insert(toRemove.end(), taudaughters.begin(), taudaughters.end());
222  }
223  }
224 
225  //now actually remove
226  pruned.clear();
227  for (unsigned int i = 0; i < status1.size(); ++i){
228  bool marked = false;
229  std::vector<const HepMC::GenParticle*>::const_iterator iremove;
230  for ( iremove = toRemove.begin(); iremove != toRemove.end(); ++iremove){
231  if (status1[i]->barcode() == (*iremove)->barcode()){
232  marked = true;
233  break;
234  }
235  }
236  if (!marked) pruned.push_back(status1[i]);
237  }
238 
239 #ifdef DEBUG_HepMCValidationHelper
240  std::cout<< "list of remaining particles:" <<std::endl;
241  for (unsigned int i = 0; i < pruned.size(); ++i){
242  std::cout << *pruned[i] << std::endl;
243  }
244 #endif
245  }
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)
float float float z
bool isNeutrino(const Candidate &part)
Definition: pdgIdUtils.h:25
void allStatus2(const HepMC::GenEvent *all, std::vector< const HepMC::GenParticle * > &status2)
double deltaR(double eta1, double eta2, double phi1, double phi2)
Definition: TreeUtility.cc:17
void findDescendents(const reco::GenParticleRef &base, reco::GenParticleRefVector &descendents, int status, int pdgId=0)
find all descendents of a given status and pdgId (recursive)
bool isTau(const Candidate &part)
Definition: pdgIdUtils.h:15
bool HepMCValidationHelper::sortByE ( const HepMC::GenParticle *  a,
const HepMC::GenParticle *  b 
)
inline

Definition at line 15 of file HepMCValidationHelper.h.

15 {return a->momentum().e() > b->momentum().e(); }
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
bool HepMCValidationHelper::sortByPseudoRapidity ( const HepMC::GenParticle *  a,
const HepMC::GenParticle *  b 
)
inline
bool HepMCValidationHelper::sortByPt ( const HepMC::GenParticle *  a,
const HepMC::GenParticle *  b 
)
inline

Definition at line 12 of file HepMCValidationHelper.h.

Referenced by WValidation::analyze(), and DrellYanValidation::analyze().

12 {return a->momentum().perp() > b->momentum().perp(); }
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121
bool HepMCValidationHelper::sortByRapidity ( const HepMC::GenParticle *  a,
const HepMC::GenParticle *  b 
)
inline

Definition at line 18 of file HepMCValidationHelper.h.

References cmsBatch::log.

18  {
19  const HepMC::FourVector& amom = a->momentum();
20  const HepMC::FourVector& bmom = b->momentum();
21  double rapa = 0.5 * std::log( (amom.e() + amom.z()) / (amom.e() - amom.z()) );
22  double rapb = 0.5 * std::log( (bmom.e() + bmom.z()) / (bmom.e() - bmom.z()) );
23  return rapa > rapb;
24  }
double b
Definition: hdecay.h:120
double a
Definition: hdecay.h:121