![]() |
![]() |
00001 #include "GeneratorInterface/GenFilters/interface/MCDijetResonance.h" 00002 00003 #include "SimDataFormats/HepMCProduct/interface/HepMCProduct.h" 00004 #include <iostream> 00005 00006 #include "FWCore/MessageLogger/interface/MessageLogger.h" 00007 00008 using namespace edm; 00009 using namespace std; 00010 00011 00012 MCDijetResonance::MCDijetResonance(const edm::ParameterSet& iConfig) : 00013 label_(iConfig.getUntrackedParameter("moduleLabel",std::string("source"))) 00014 { 00015 //here do whatever other initialization is needed 00016 00017 //Get dijet process type we will select 00018 dijetProcess = iConfig.getUntrackedParameter<string>( "dijetProcess" ); 00019 cout << "Dijet Resonance Filter Selecting Process = " << dijetProcess << endl; 00020 00021 maxQuarkID = 4; //Maximum |ID| of Light Quark = charm quark gives u, d, s, and c decays of Zprime. 00022 bosonID = 21; //Gluon 00023 00024 //Number of events and Number Accepted 00025 nEvents = 0; 00026 nAccepted = 0; 00027 } 00028 00029 00030 MCDijetResonance::~MCDijetResonance() 00031 { 00032 00033 // do anything here that needs to be done at desctruction time 00034 // (e.g. close files, deallocate resources etc.) 00035 00036 } 00037 00038 void MCDijetResonance::endJob() 00039 { 00040 edm::LogVerbatim("MCDijetResonanceInfo") << "================MCDijetResonance report========================================\n" 00041 << "Events read " << nEvents << " Events accepted " << nAccepted << "\nEfficiency " << ((double)nAccepted)/((double)nEvents) 00042 << "\n====================================================================" << std::endl; 00043 } 00044 00045 // ------------ method called to skim the data ------------ 00046 bool MCDijetResonance::filter(edm::Event& iEvent, const edm::EventSetup& iSetup) 00047 { 00048 nEvents++; 00049 cout << endl << "Event=" << nEvents << endl; 00050 using namespace edm; 00051 Handle<HepMCProduct> evt; 00052 iEvent.getByLabel(label_, evt); 00053 HepMC::GenEvent * myGenEvent = new HepMC::GenEvent(*(evt->GetEvent())); 00054 00055 //If process is not the desired primary process, cleanup and reject the event. 00056 if (dijetProcess == "ZprimeLightQuarks" && myGenEvent->signal_process_id() != 141){ 00057 // Wanted a Z' but didn't find it, so reject event. 00058 delete myGenEvent; 00059 return false; 00060 } 00061 if (dijetProcess == "QstarQuarkGluon" && myGenEvent->signal_process_id() != 147 && myGenEvent->signal_process_id() != 148){ 00062 // Wanted a q* but didn't find it, so reject event. 00063 delete myGenEvent; 00064 return false; 00065 } 00066 00067 //found a dijet resonance 00068 00069 //debug 00070 // int count = 0; 00071 //for ( HepMC::GenEvent::particle_iterator p = myGenEvent->particles_begin(); 00072 // p != myGenEvent->particles_end() & count<13; ++p ) 00073 // { 00074 //std::cout << count << ": ID=" << (*p)->pdg_id() << ", status=" << (*p)->status() << ", mass=" << (*p)->momentum().invariantMass() << ", pt=" <<(*p)->momentum().perp() << std::endl; 00075 //count++; 00076 //} 00077 00078 for ( HepMC::GenEvent::particle_iterator p = myGenEvent->particles_begin(); 00079 p != myGenEvent->particles_end(); ++p ) 00080 { 00081 //Find resonance particle and check that it is part of hard collision 00082 if ( (*p)->status() == 3 && ( (dijetProcess == "ZprimeLightQuarks" && (*p)->pdg_id() == 32) || 00083 (dijetProcess == "QstarQuarkGluon" && abs((*p)->pdg_id()) == 4000001) || 00084 (dijetProcess == "QstarQuarkGluon" && abs((*p)->pdg_id()) == 4000002) ) ) 00085 { 00086 // The next two particles are the outgoing particles from the resonance decay 00087 p++; 00088 int ID1 = (*p)->pdg_id(); 00089 p++; 00090 int ID2 = (*p)->pdg_id(); 00091 00092 //Check for the process we want 00093 if( (dijetProcess == "ZprimeLightQuarks" && abs(ID1) <= maxQuarkID && abs(ID2) <= maxQuarkID ) || 00094 (dijetProcess == "QstarQuarkGluon" && ( ID1 == bosonID || ID2 == bosonID ) ) ) 00095 { 00096 //cout << "dijet resonance " << dijetProcess << " found " << endl; 00097 nAccepted++; 00098 delete myGenEvent; 00099 return true; 00100 } 00101 else 00102 { 00103 delete myGenEvent; 00104 return false; 00105 } 00106 } 00107 00108 } 00109 00110 delete myGenEvent; 00111 return false; 00112 }