CMS 3D CMS Logo

/afs/cern.ch/work/a/aaltunda/public/www/CMSSW_6_2_7/src/GeneratorInterface/GenFilters/src/MCZll.cc

Go to the documentation of this file.
00001 
00002 #include "GeneratorInterface/GenFilters/interface/MCZll.h"
00003 
00004 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00005 #include <iostream>
00006 
00007 #include "FWCore/MessageLogger/interface/MessageLogger.h"
00008 
00009 using namespace edm;
00010 using namespace std;
00011 
00012 
00013 MCZll::MCZll(const edm::ParameterSet& iConfig) :
00014   label_(iConfig.getUntrackedParameter("moduleLabel",std::string("generator"))), nEvents_(0), nAccepted_(0)
00015 {
00016   leptonFlavour_ = iConfig.getUntrackedParameter<int>("leptonFlavour",11);
00017   leptonPtMin_ = iConfig.getUntrackedParameter<double>("leptonPtMin",5.);
00018   leptonPtMax_ = iConfig.getUntrackedParameter<double>("leptonPtMax",99999.);
00019   leptonEtaMin_ = iConfig.getUntrackedParameter<double>("leptonEtaMin",0.);
00020   leptonEtaMax_ = iConfig.getUntrackedParameter<double>("leptonEtaMax",2.7);
00021   zMassRange_.first = iConfig.getUntrackedParameter<double>("zMassMin",60.);
00022   zMassRange_.second = iConfig.getUntrackedParameter<double>("zMassMax",120.);
00023   filter_ = iConfig.getUntrackedParameter<bool>("filter",true);
00024   std::ostringstream str;
00025   str << "=========================================================\n" 
00026       << "Filter MCZll being constructed with parameters: " 
00027       << "\nleptonFlavour " << leptonFlavour_ 
00028       << "\nleptonPtMin " << leptonPtMin_
00029       << "\nleptonPtMax " << leptonPtMax_
00030       << "\nleptonEtaMin " << leptonEtaMin_
00031       << "\nleptonEtaMax " << leptonEtaMax_
00032       << "\nzMassMin " << zMassRange_.first
00033       << "\nzMassMax " << zMassRange_.second
00034       << "\n=========================================================" ;
00035   edm::LogVerbatim("MCZllInfo") <<  str.str() ;
00036   if (filter_)
00037     produces< HepMCProduct >();
00038 }
00039 
00040 
00041 MCZll::~MCZll()
00042 {
00043  
00044    // do anything here that needs to be done at desctruction time
00045    // (e.g. close files, deallocate resources etc.)
00046 
00047 }
00048 
00049 void MCZll::endJob() 
00050 {
00051   edm::LogVerbatim("MCZllInfo") << "================MCZll report========================================\n" 
00052             << "Events read " << nEvents_ << " Events accepted " << nAccepted_ << "\nEfficiency " << ((double)nAccepted_)/((double)nEvents_) 
00053             << "\n====================================================================" << std::endl;
00054 }
00055 
00056 // ------------ method called to skim the data  ------------
00057 bool MCZll::filter(edm::Event& iEvent, const edm::EventSetup& iSetup)
00058 {
00059   std::auto_ptr<HepMCProduct> bare_product(new HepMCProduct()); 
00060 
00061   nEvents_++;
00062   using namespace edm;
00063   bool accepted = false;
00064   Handle<HepMCProduct> evt;
00065   iEvent.getByLabel(label_, evt);
00066   HepMC::GenEvent * myGenEvent = new  HepMC::GenEvent(*(evt->GetEvent()));
00067   HepMC::GenEvent * zEvent = new HepMC::GenEvent();
00068 
00069   if (myGenEvent->signal_process_id() != 1) 
00070     {
00071       delete myGenEvent;   
00072       delete zEvent;
00073       return false;
00074     }
00075       
00076   
00077   //found a prompt Z
00078   
00079   for ( HepMC::GenEvent::particle_iterator p = myGenEvent->particles_begin();
00080         p != myGenEvent->particles_end(); ++p ) 
00081     {
00082       if ( !accepted && ( (*p)->pdg_id() == 23 ) && (*p)->status() == 3 ) 
00083         { 
00084           accepted=true;
00085           HepMC::GenVertex* zVertex = new HepMC::GenVertex();
00086           HepMC::GenParticle* myZ= new HepMC::GenParticle(*(*p));
00087           zVertex->add_particle_in(myZ);
00088           //      std::cout << (*p)->momentum().invariantMass() << std::endl;
00089           if ((*p)->momentum().m() < zMassRange_.first || (*p)->momentum().m() > zMassRange_.second)
00090             accepted = false; 
00091           std::vector<HepMC::GenParticle*> children;
00092           HepMC::GenVertex* outVertex=(*p)->end_vertex();
00093           for(HepMC::GenVertex::particles_out_const_iterator iter = outVertex->particles_out_const_begin();
00094               iter != outVertex->particles_out_const_end(); iter++)
00095             children.push_back(*iter);
00096           std::vector<HepMC::GenParticle*>::const_iterator aDaughter;
00097           for (aDaughter = children.begin();aDaughter != children.end();aDaughter++)
00098             {
00099               HepMC::GenParticle* myDa= new HepMC::GenParticle(*(*aDaughter));
00100               zVertex->add_particle_out(myDa);
00101               if ((*aDaughter)->status() == 2)
00102                 continue;
00103               //              (*aDaughter)->print();          
00104 
00105               if (! (abs((*aDaughter)->pdg_id()) == abs(leptonFlavour_)) )
00106                 accepted = false; 
00107               //                std::cout << (*aDaughter)->momentum().perp() << " " << (*aDaughter)->momentum().eta() << std::endl;
00108               if ((*aDaughter)->momentum().perp() < leptonPtMin_)
00109                 accepted = false; 
00110               if ((*aDaughter)->momentum().perp() > leptonPtMax_)
00111                 accepted = false; 
00112               if (fabs((*aDaughter)->momentum().eta()) > leptonEtaMax_)
00113                 accepted = false; 
00114               if (fabs((*aDaughter)->momentum().eta()) < leptonEtaMin_)
00115                 accepted = false; 
00116             }
00117           zEvent->add_vertex( zVertex );
00118           if (accepted)
00119             break;
00120           
00121         }
00122 
00123     } 
00124   
00125 
00126   if (accepted)
00127     {
00128       if(zEvent)  
00129         bare_product->addHepMCData(zEvent);
00130       if (filter_)
00131         iEvent.put(bare_product);
00132       nAccepted_++;
00133       //      std::cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<< std::endl;
00134       LogDebug("MCZll") << "Event " << iEvent.id().event()  << " accepted" << std::endl; 
00135       //      std::cout << "+++++++++++++++++++++++++++++++++++++++++++++++++"<< std::endl;
00136       //       myGenEvent->print(); 
00137       delete myGenEvent;   
00138       return true; 
00139     } 
00140 
00141   delete myGenEvent;   
00142   delete zEvent;
00143   return false;
00144 
00145 
00146 }
00147