Go to the documentation of this file.00001 #include "GeneratorInterface/GenFilters/interface/PythiaFilterEMJet.h"
00002 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00003 #include <iostream>
00004 #include<list>
00005 #include<map>
00006 #include<vector>
00007 #include<cmath>
00008
00009 using namespace edm;
00010 using namespace std;
00011
00012 namespace{
00013
00014 double deltaR2(double eta0, double phi0, double eta, double phi){
00015 double dphi=phi-phi0;
00016 if(dphi>M_PI) dphi-=2*M_PI;
00017 else if(dphi<=-M_PI) dphi+=2*M_PI;
00018 return dphi*dphi+(eta-eta0)*(eta-eta0);
00019 }
00020
00021 double deltaPhi(double phi0, double phi){
00022 double dphi=phi-phi0;
00023 if(dphi>M_PI) dphi-=2*M_PI;
00024 else if(dphi<=-M_PI) dphi+=2*M_PI;
00025 return dphi;
00026 }
00027
00028 class ParticlePtGreater{
00029 public:
00030 int operator()(const HepMC::GenParticle * p1,
00031 const HepMC::GenParticle * p2) const{
00032 return p1->momentum().perp() > p2->momentum().perp();
00033 }
00034 };
00035 }
00036
00037
00038 PythiaFilterEMJet::PythiaFilterEMJet(const edm::ParameterSet& iConfig) :
00039 label_(iConfig.getUntrackedParameter("moduleLabel",std::string("generator"))),
00040 etaMin(iConfig.getUntrackedParameter<double>("MinEMEta", 0)),
00041 eTSumMin(iConfig.getUntrackedParameter<double>("ETSumMin", 50.)),
00042 pTMin(iConfig.getUntrackedParameter<double>("MinEMpT", 5.)),
00043 etaMax(iConfig.getUntrackedParameter<double>("MaxEMEta", 2.7)),
00044 eTSumMax(iConfig.getUntrackedParameter<double>("ETSumMax", 100.)),
00045 pTMax(iConfig.getUntrackedParameter<double>("MaxEMpT", 999999.)),
00046 ebEtaMax(1.479),
00047 maxnumberofeventsinrun(iConfig.getUntrackedParameter<int>("MaxEvents",3000000)){
00048
00049 deltaEB=0.01745/2 *5;
00050 deltaEE=2.93/317/2 *5;
00051 theNumberOfTestedEvt = 0;
00052 theNumberOfSelected = 0;
00053
00054 cout << " Max Events : " << maxnumberofeventsinrun << endl;
00055 cout << " Cut Definition: " << endl;
00056 cout << " MinEMEta = " << etaMin << endl;
00057 cout << " ETSumMin = " << eTSumMin << endl;
00058 cout << " MinEMpT = " << pTMin << endl;
00059 cout << " MaxEMEta = " << etaMax << endl;
00060 cout << " ETSumMax = " << eTSumMax << endl;
00061 cout << " MaxEMpT = " << pTMax << endl;
00062 cout << " 5x5 crystal cone around EM axis in ECAL Barrel = " << deltaEB << endl;
00063 cout << " 5x5 crystal cone around EM axis in ECAL Endcap = " << deltaEE << endl;
00064
00065 }
00066
00067 PythiaFilterEMJet::~PythiaFilterEMJet()
00068 {
00069 std::cout << "Total number of tested events = " << theNumberOfTestedEvt << std::endl;
00070 std::cout << "Total number of accepted events = " << theNumberOfSelected << std::endl;
00071 }
00072
00073
00074
00075 bool PythiaFilterEMJet::filter(edm::Event& iEvent, const edm::EventSetup& iSetup){
00076
00077 if(theNumberOfSelected>=maxnumberofeventsinrun) {
00078 throw cms::Exception("endJob") << "we have reached the maximum number of events ";
00079 }
00080
00081 theNumberOfTestedEvt++;
00082 if(theNumberOfTestedEvt%1000 == 0) cout << "Number of tested events = " << theNumberOfTestedEvt << endl;
00083
00084 bool accepted = false;
00085 Handle<edm::HepMCProduct> evt;
00086 iEvent.getByLabel(label_, evt);
00087
00088 list<const HepMC::GenParticle *> EM_seeds;
00089 const HepMC::GenEvent * myGenEvent = evt->GetEvent();
00090
00091 int particle_id = 1;
00092 int EM_id = -1;
00093
00094
00095 for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
00096 p != myGenEvent->particles_end();
00097 ++p ) {
00098
00099
00100 if ( (abs((*p)->pdg_id())==11 || (*p)->pdg_id()==22)
00101 && (*p)->status()==1
00102 && (*p)->momentum().perp() > pTMin
00103 && (*p)->momentum().perp() < pTMax
00104 && fabs((*p)->momentum().eta()) < etaMax
00105 && fabs((*p)->momentum().eta()) > etaMin ) {
00106 EM_id = particle_id;
00107 EM_seeds.push_back(*p);
00108 }
00109 particle_id++;
00110 }
00111
00112 EM_seeds.sort(ParticlePtGreater());
00113
00114 double etaEMClus=0;
00115 double phiEMClus=0;
00116 double ptEMClus=0;
00117 for(std::list<const HepMC::GenParticle *>::const_iterator is=EM_seeds.begin();
00118 is!= EM_seeds.end();
00119 ++is){
00120 double etaEM=(*is)->momentum().eta();
00121 double phiEM=(*is)->momentum().phi();
00122 double ptEM=(*is)->momentum().perp();
00123
00124 etaEMClus = etaEM;
00125 phiEMClus = phiEM;
00126 ptEMClus = ptEM;
00127
00128
00129 bool inEB(0);
00130 double tgx(0);
00131 double tgy(0);
00132 if( std::abs(etaEM)<ebEtaMax ) inEB=1;
00133 else{
00134 tgx=(*is)->momentum().px()/(*is)->momentum().pz();
00135 tgy=(*is)->momentum().py()/(*is)->momentum().pz();
00136 }
00137
00138
00139
00140
00141 for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin();
00142 p != myGenEvent->particles_end();
00143 ++p ) {
00144
00145 if (((*p)->status()==1 && (*p)->pdg_id() == 22) ||
00146 ((*p)->status()==1 && abs((*p)->pdg_id()) == 11) )
00147
00148
00149
00150
00151
00152 {
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165
00166
00167
00168
00169
00170
00171
00172 double pt=(*p)->momentum().perp();
00173 if (pt == ptEM) continue ;
00174 double eta=(*p)->momentum().eta();
00175 double phi=(*p)->momentum().phi();
00176
00177 if(inEB) {
00178 if( std::abs(eta-etaEM)> deltaEB ||
00179 std::abs(deltaPhi(phi,phiEM)) > deltaEB) continue;
00180 }
00181 else if( std::abs((*p)->momentum().px()/(*p)->momentum().pz() - tgx)
00182 > deltaEE ||
00183 std::abs((*p)->momentum().py()/(*p)->momentum().pz() - tgy)
00184 > deltaEE) continue;
00185 ptEMClus += pt ;
00186 if(inEB)
00187 {
00188 etaEMClus += (eta-etaEMClus)*pt/ptEMClus ;
00189 phiEMClus += deltaPhi(phi,phiEM)*pt/ptEMClus;
00190 }
00191 else
00192 {
00193 etaEMClus += ((*p)->momentum().px()/(*p)->momentum().pz() - tgx)*pt/ptEMClus ;
00194 phiEMClus += ((*p)->momentum().py()/(*p)->momentum().pz() - tgy)*pt/ptEMClus;
00195 }
00196 }
00197 }
00198 if( ptEMClus > eTSumMin)
00199 accepted = true ;
00200 }
00201
00202 if (accepted) {
00203 theNumberOfSelected++;
00204 cout << "========> Event preselected " << theNumberOfSelected
00205 << " Proccess ID " << myGenEvent->signal_process_id() << endl;
00206 return true;
00207 }
00208 else return false;
00209 }
00210