CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_5/src/Validation/EventGenerator/src/HepMCValidationHelper.cc

Go to the documentation of this file.
00001 #include "Validation/EventGenerator/interface/HepMCValidationHelper.h"
00002 #include "DataFormats/Math/interface/deltaR.h"
00003 #include <iostream>
00004 #include <cassert>
00005 #include <limits>
00006 #include "TLorentzVector.h"
00007 
00008 //#define DEBUG_HepMCValidationHelper
00009 
00010 namespace HepMCValidationHelper {
00011   void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons, 
00012                       const HepMC::GenEvent* all,
00013                       double deltaRcut,
00014                       std::vector<const HepMC::GenParticle*>& fsrphotons){
00015 
00016     std::vector<const HepMC::GenParticle*> status1;
00017     allStatus1(all, status1);
00018     findFSRPhotons(leptons, status1, deltaRcut, fsrphotons);
00019   }
00020   
00021   void findFSRPhotons(const std::vector<const HepMC::GenParticle*>& leptons,
00022                       const std::vector<const HepMC::GenParticle*>& all,
00023                       double deltaRcut,
00024                       std::vector<const HepMC::GenParticle*>& fsrphotons){
00025                         
00026     //find all status 1 photons
00027     std::vector<const HepMC::GenParticle*> allphotons;
00028     for (unsigned int i = 0; i < all.size(); ++i){
00029       if  (all[i]->status()==1 && all[i]->pdg_id()==22) allphotons.push_back(all[i]);
00030     }
00031 
00032     //loop over the photons and check the distance wrt the leptons
00033     for (unsigned int ipho = 0; ipho < allphotons.size(); ++ipho){
00034       bool close = false;
00035       for (unsigned int ilep = 0; ilep < leptons.size(); ++ilep){
00036         if (deltaR(allphotons[ipho]->momentum(), leptons[ilep]->momentum()) < deltaRcut){
00037           close = true;
00038           break;
00039         }  
00040       }
00041       if (close) fsrphotons.push_back(allphotons[ipho]);
00042     }
00043   }
00044 
00045   //returns true if a status 3 particle is a tau or if a status 1 particle is either an electron or a neutrino
00046   bool isChargedLepton(const HepMC::GenParticle* part){
00047     int status = part->status();
00048     unsigned int pdg_id = abs(part->pdg_id());
00049     if (status == 2) return pdg_id == 15;
00050     else return status == 1 &&  (pdg_id == 11 || pdg_id == 13 ) ;
00051   }
00052 
00053   //returns true if a status 1 particle is a neutrino
00054   bool isNeutrino(const HepMC::GenParticle* part ) {
00055       int status = part->status();
00056       unsigned int pdg_id = abs(part->pdg_id());
00057       return status == 1 &&  (pdg_id == 12 || pdg_id == 14 || pdg_id == 16) ;
00058   } 
00059 
00060   //returns true is status 3 particle is tau
00061   bool isTau(const HepMC::GenParticle* part){
00062     return part->status() == 2 && abs(part->pdg_id()) == 15;
00063   }
00064 /* 
00065   void getTaus(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& taus){
00066     for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
00067       if (abs((*iter)->pdg_id()) == 15) taus.push_back(*iter);
00068     }
00069   }
00070 */
00071   // get all status 1 particles
00072   void allStatus1(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1){
00073     for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
00074       if ((*iter)->status() == 1) status1.push_back(*iter);
00075     }
00076   }
00077 
00078   void allStatus2(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1){
00079     for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
00080       if ((*iter)->status() == 2) status1.push_back(*iter);
00081     }
00082   }
00083 
00084   void allStatus3(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& status1){
00085     for (HepMC::GenEvent::particle_const_iterator iter = all->particles_begin(); iter != all->particles_end(); ++iter){
00086       if ((*iter)->status() == 3) status1.push_back(*iter);
00087     }
00088   }
00089 
00090   void findDescendents(const HepMC::GenParticle* a, std::vector<const HepMC::GenParticle*>& descendents){
00091     HepMC::GenVertex* decayVertex = a->end_vertex();
00092     if (!decayVertex) return;
00093     HepMC::GenVertex::particles_out_const_iterator ipart;
00094     for (ipart = decayVertex->particles_out_const_begin(); ipart != decayVertex->particles_out_const_end(); ++ipart ){
00095       if ((*ipart)->status() == 1) descendents.push_back(*ipart);
00096       else findDescendents(*ipart, descendents);
00097     }
00098   }
00099 
00100   void removeIsolatedLeptons(const HepMC::GenEvent* all, double deltaRcut, double sumPtCut, std::vector<const HepMC::GenParticle*>& pruned) {
00101     //get all status 1 particles
00102     std::vector<const HepMC::GenParticle*> status1;
00103     allStatus1(all, status1);
00104     std::vector<const HepMC::GenParticle*> toRemove;
00105     //loop on all particles and find candidates to be isolated
00106     for (unsigned int i = 0; i < status1.size(); ++i){
00107       //if it is a neutrino is a charged lepton (not a tau) this is a candidate to be isolated
00108       if (isNeutrino(status1[i]) || (isChargedLepton(status1[i]) && !isTau(status1[i]))){
00109         //list of particles not to be considered in the isolation computation.
00110         //this includes the particle to be isolated and the fsr photons in case of charged lepton
00111         std::vector<const HepMC::GenParticle*> leptons;
00112         leptons.push_back(status1[i]);
00113         std::vector<const HepMC::GenParticle*> removedForIsolation;
00114         removedForIsolation.push_back(status1[i]);
00115         if (isChargedLepton(status1[i])) findFSRPhotons(leptons, status1, deltaRcut, removedForIsolation);
00116 #ifdef  DEBUG_HepMCValidationHelper
00117         //std::cout << removedForIsolation.size() << " particles to be removed for isolation calculation " << std::endl;
00118 #endif
00119         //create vector of particles to compute isolation (removing removedForIsolation);
00120         std::vector<const HepMC::GenParticle*> forIsolation;
00121         std::vector<const HepMC::GenParticle*>::iterator iiso;
00122         for (iiso = status1.begin(); iiso != status1.end(); ++iiso){
00123           std::vector<const HepMC::GenParticle*>::const_iterator iremove;
00124           bool marked = false;
00125           for ( iremove = removedForIsolation.begin(); iremove != removedForIsolation.end(); ++iremove){
00126             if ((*iiso)->barcode() == (*iremove)->barcode()) {
00127 #ifdef  DEBUG_HepMCValidationHelper
00128               //std::cout << "removing particle " << **iiso << " from the list of particles to compute isolation" << std::endl;
00129 #endif
00130               marked = true;  
00131               break;
00132             }
00133           }
00134           if (!marked) forIsolation.push_back(*iiso);
00135         }
00136         //now compute isolation
00137         double sumIso = 0;
00138         for (iiso = forIsolation.begin(); iiso < forIsolation.end(); ++iiso){
00139            if (deltaR(leptons.front()->momentum(), (*iiso)->momentum()) < deltaRcut){
00140               sumIso += (*iiso)->momentum().perp();
00141            }
00142         }
00143         //if isolated remove from the pruned list
00144         if (sumIso < sumPtCut) {
00145 #ifdef  DEBUG_HepMCValidationHelper
00146            std::cout << "particle " <<  *status1[i] <<  " is considered isolated, with sumPt " << sumIso << std::endl; 
00147 #endif
00148            toRemove.insert(toRemove.end(), removedForIsolation.begin(), removedForIsolation.end()); 
00149         }
00150 #ifdef  DEBUG_HepMCValidationHelper        
00151         else {
00152           std::cout << "NOT isolated! " <<  *status1[i] <<  " is considered not isolated, with sumPt " << sumIso << std::endl;
00153         }
00154 #endif        
00155       } 
00156     }
00157     //at this point we have taken care of the electrons and muons, but pruned could  still contain the decay products of isolated taus, 
00158     //we want to remove these as well
00159     std::vector<const HepMC::GenParticle*> status2;
00160     allStatus2(all, status2);
00161     std::vector<const HepMC::GenParticle*> taus;
00162     //getTaus(all, taus);
00163     for (unsigned int i = 0; i < status2.size(); ++i){
00164       if (isTau(status2[i])) {
00165         //check the list we have already for duplicates
00166         //there use to be duplicates in some generators (sherpa) 
00167         bool duplicate = false; 
00168         TLorentzVector taumomentum(status2[i]->momentum().x(), status2[i]->momentum().y(), status2[i]->momentum().z(), status2[i]->momentum().t());
00169         for (unsigned int j = 0; j < taus.size(); ++j){
00170           //compare momenta
00171           TLorentzVector othermomentum(taus[j]->momentum().x(), taus[j]->momentum().y(), taus[j]->momentum().z(), taus[j]->momentum().t());
00172           othermomentum-=taumomentum;
00173           if (status2[i]->pdg_id() == taus[j]->pdg_id() &&
00174               othermomentum.E() < 0.1 &&//std::numeric_limits<float>::epsilon() &&
00175               othermomentum.P() < 0.1 ){ //std::numeric_limits<float>::epsilon()){
00176             duplicate = true;
00177             break;
00178           }    
00179         }
00180         if (!duplicate) taus.push_back(status2[i]);
00181       }  
00182     }
00183     //loop over the taus, find the descendents, remove all these from the list of particles to compute isolation   
00184     for (unsigned int i = 0; i < taus.size(); ++i){
00185       std::vector<const HepMC::GenParticle*> taudaughters;
00186       findDescendents(taus[i], taudaughters);
00187       assert(taudaughters.size()>0);
00188       const HepMC::FourVector& taumom = taus[i]->momentum();
00189       //remove the daughters from the list of particles to compute isolation
00190       std::vector<const HepMC::GenParticle*> forIsolation;
00191       std::vector<const HepMC::GenParticle*>::iterator iiso;
00192       for (iiso = status1.begin(); iiso < status1.end(); ++iiso){
00193         bool marked = false;
00194         std::vector<const HepMC::GenParticle*>::const_iterator iremove;
00195         for ( iremove = taudaughters.begin(); iremove != taudaughters.end(); ++iremove){
00196           if ((*iiso)->barcode() == (*iremove)->barcode()) {
00197 #ifdef  DEBUG_HepMCValidationHelper
00198 //            std::cout << "removing particle " << **iiso << " from the list of particles to compute isolation because it comes from a tau" << std::endl;
00199 #endif
00200             marked = true;
00201             break;
00202           }
00203           if (!marked) forIsolation.push_back(*iiso);
00204         }
00205       }
00206       //no compute isolation wrt the status 2 tau direction
00207       double sumIso = 0;
00208       for (iiso = forIsolation.begin(); iiso < forIsolation.end(); ++iiso){
00209          if (deltaR(taumom, (*iiso)->momentum()) < deltaRcut){
00210             sumIso += (*iiso)->momentum().perp();
00211          }
00212       }
00213       //if isolated remove the tau daughters from the pruned list
00214       if (sumIso < sumPtCut) {
00215 #ifdef  DEBUG_HepMCValidationHelper
00216         std::cout << "particle " <<  *taus[i] <<  " is considered isolated, with sumPt " << sumIso << std::endl;
00217 #endif
00218         toRemove.insert(toRemove.end(), taudaughters.begin(), taudaughters.end());
00219       }
00220     }
00221 
00222     //now actually remove
00223     pruned.clear();
00224     for (unsigned int i = 0; i < status1.size(); ++i){
00225       bool marked = false; 
00226       std::vector<const HepMC::GenParticle*>::const_iterator iremove;
00227       for ( iremove = toRemove.begin(); iremove != toRemove.end(); ++iremove){
00228         if (status1[i]->barcode() == (*iremove)->barcode()){
00229           marked = true;
00230           break;
00231         }
00232       }
00233       if (!marked) pruned.push_back(status1[i]);
00234     }
00235 
00236 #ifdef  DEBUG_HepMCValidationHelper    
00237     std::cout<< "list of remaining particles:" <<std::endl;   
00238     for (unsigned int i = 0; i < pruned.size(); ++i){
00239       std::cout << *pruned[i] << std::endl;
00240     }
00241 #endif    
00242   } 
00243 
00244   //get all visible status1 particles
00245   void allVisibleParticles(const HepMC::GenEvent* all, std::vector<const HepMC::GenParticle*>& visible)  {
00246     std::vector<const HepMC::GenParticle*> status1;
00247     visible.clear();    
00248     allStatus1(all, status1);
00249     for (unsigned int i = 0; i < status1.size(); ++i){
00250       if (!isNeutrino(status1[i])) visible.push_back(status1[i]);
00251     }    
00252   }
00253 
00254   //compute generated met
00255   TLorentzVector genMet(const HepMC::GenEvent* all, double etamin, double etamax){
00256     std::vector<const HepMC::GenParticle*> visible;
00257     allVisibleParticles(all, visible);
00258     TLorentzVector momsum(0., 0., 0., 0.);      
00259     for (unsigned int i = 0; i < visible.size(); ++i){
00260       if (visible[i]->momentum().eta() > etamin && visible[i]->momentum().eta() < etamax) {
00261         TLorentzVector mom(visible[i]->momentum().x(), visible[i]->momentum().y(), visible[i]->momentum().z(), visible[i]->momentum().t());     
00262         momsum += mom;  
00263       }                         
00264     }
00265     TLorentzVector met(-momsum.Px(), -momsum.Py(), 0., momsum.Pt());    
00266     return met; 
00267   }              
00268 }