Go to the documentation of this file.00001 #include "GeneratorInterface/GenFilters/interface/PythiaFilterGammaJetWithOutBg.h"
00002 #include "SimDataFormats/GeneratorProducts/interface/HepMCProduct.h"
00003 #include "SimGeneral/HepPDTRecord/interface/ParticleDataTable.h"
00004 #include "FWCore/Framework/interface/EventSetup.h"
00005 #include "FWCore/Framework/interface/ESHandle.h"
00006 #include <iostream>
00007 #include<list>
00008 #include<vector>
00009 #include<cmath>
00010
00011
00012
00013
00014 namespace{
00015
00016 double deltaR2(double eta0, double phi0, double eta, double phi){
00017 double dphi=phi-phi0;
00018 if(dphi>M_PI) dphi-=2*M_PI;
00019 else if(dphi<=-M_PI) dphi+=2*M_PI;
00020 return dphi*dphi+(eta-eta0)*(eta-eta0);
00021 }
00022
00023 double deltaPhi(double phi0, double phi){
00024 double dphi=phi-phi0;
00025 if(dphi>M_PI) dphi-=2*M_PI;
00026 else if(dphi<=-M_PI) dphi+=2*M_PI;
00027 return dphi;
00028 }
00029
00030 class ParticlePtGreater{
00031 public:
00032 int operator()(const HepMC::GenParticle * p1,
00033 const HepMC::GenParticle * p2) const{
00034 return p1->momentum().perp() > p2->momentum().perp();
00035 }
00036 };
00037 }
00038
00039
00040 PythiaFilterGammaJetWithOutBg::PythiaFilterGammaJetWithOutBg(const edm::ParameterSet& iConfig) :
00041 label_(iConfig.getUntrackedParameter("moduleLabel",std::string("generator"))),
00042 etaMax(iConfig.getUntrackedParameter<double>("MaxPhotonEta", 2.8)),
00043 ptSeed(iConfig.getUntrackedParameter<double>("PhotonSeedPt", 5.)),
00044 ptMin(iConfig.getUntrackedParameter<double>("MinPhotonPt")),
00045 ptMax(iConfig.getUntrackedParameter<double>("MaxPhotonPt")),
00046 dphiMin(iConfig.getUntrackedParameter<double>("MinDeltaPhi", -1)/180*M_PI),
00047 detaMax(iConfig.getUntrackedParameter<double>("MaxDeltaEta", 10.)),
00048 etaPhotonCut2(iConfig.getUntrackedParameter<double>("MinPhotonEtaForwardJet", 1.3)),
00049 cone(0.5),ebEtaMax(1.479),
00050 maxnumberofeventsinrun(iConfig.getUntrackedParameter<int>("MaxEvents",10)){
00051
00052 deltaEB=0.01745/2 *5;
00053 deltaEE=2.93/317/2 *5;
00054 theNumberOfSelected = 0;
00055 }
00056
00057
00058 PythiaFilterGammaJetWithOutBg::~PythiaFilterGammaJetWithOutBg(){}
00059
00060
00061
00062 bool PythiaFilterGammaJetWithOutBg::filter(edm::Event& iEvent, const edm::EventSetup& iSetup){
00063
00064
00065
00066
00067
00068
00069
00070
00071 bool accepted = false;
00072 edm::Handle<edm::HepMCProduct> evt;
00073 iEvent.getByLabel(label_, evt);
00074
00075 std::list<const HepMC::GenParticle *> seeds;
00076 const HepMC::GenEvent * myGenEvent = evt->GetEvent();
00077
00078 for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
00079
00080 if ( (*p)->pdg_id()==22 && (*p)->status()==1
00081 && (*p)->momentum().perp() > ptSeed
00082 && std::abs((*p)->momentum().eta()) < etaMax ) seeds.push_back(*p);
00083
00084 }
00085
00086 seeds.sort(ParticlePtGreater());
00087
00088
00089
00090
00091
00092
00093
00094
00095 for(std::list<const HepMC::GenParticle *>::const_iterator is=
00096 seeds.begin(); is!=seeds.end(); is++){
00097
00098 double etaPhoton=(*is)->momentum().eta();
00099 double phiPhoton=(*is)->momentum().phi();
00100
00101
00102
00103
00104
00105
00106
00107
00108 HepMC::GenEvent::particle_const_iterator ppp = myGenEvent->particles_begin();
00109 for(int i=0;i<6;++i) ppp++;
00110 HepMC::GenParticle* particle7 = (*ppp);
00111 ppp++;
00112 HepMC::GenParticle* particle8 = (*ppp);
00113
00114 double dphi7=std::abs(deltaPhi(phiPhoton,
00115 particle7->momentum().phi()));
00116 double dphi=std::abs(deltaPhi(phiPhoton,
00117 particle8->momentum().phi()));
00118
00119 int jetline=8;
00120 int photonline = 7;
00121 if(dphi7>dphi) {
00122 dphi=dphi7;
00123 jetline=7;
00124 photonline = 8;
00125 }
00126
00127
00128
00129
00130
00131
00132
00133
00134 double etaJet = 0.0;
00135 if(jetline==8) etaJet = particle8->momentum().eta();
00136 else etaJet = particle7->momentum().eta();
00137
00138 double eta1=etaJet-detaMax;
00139 double eta2=etaJet+detaMax;
00140 if (eta1>etaPhotonCut2) eta1=etaPhotonCut2;
00141 if (eta2<-etaPhotonCut2) eta2=-etaPhotonCut2;
00142
00143
00144
00145 if (etaPhoton<eta1 ||etaPhoton>eta2) {
00146
00147 continue;
00148 }
00149 bool inEB(0);
00150 double tgx(0);
00151 double tgy(0);
00152 if( std::abs(etaPhoton)<ebEtaMax) inEB=1;
00153 else{
00154 tgx=(*is)->momentum().px()/(*is)->momentum().pz();
00155 tgy=(*is)->momentum().py()/(*is)->momentum().pz();
00156 }
00157
00158 double etPhoton=0;
00159 double etPhotonCharged=0;
00160 double etCone=0;
00161 double etConeCharged=0;
00162 double ptMaxHadron=0;
00163
00164
00165 for ( HepMC::GenEvent::particle_const_iterator p = myGenEvent->particles_begin(); p != myGenEvent->particles_end(); ++p ) {
00166
00167 if ( (*p)->status()!=1 ) continue;
00168 int pid= (*p)->pdg_id();
00169 int apid= std::abs(pid);
00170 if (apid>11 && apid<21) continue;
00171 double eta=(*p)->momentum().eta();
00172 double phi=(*p)->momentum().phi();
00173 if (deltaR2(etaPhoton, phiPhoton, eta, phi)>cone*cone) continue;
00174 double pt=(*p)->momentum().perp();
00175
00176
00177
00178
00179
00180 edm::ESHandle<ParticleDataTable> pdt;
00181 iSetup.getData( pdt );
00182 int charge3 = ((pdt->particle((*p)->pdg_id()))->ID().threeCharge());
00183
00184
00185 etCone+=pt;
00186 if(charge3 && pt<2) etConeCharged+=pt;
00187
00188
00189 if(inEB) {
00190 if( std::abs(eta-etaPhoton)> deltaEB ||
00191 std::abs(deltaPhi(phi,phiPhoton)) > deltaEB) continue;
00192 }
00193 else if( std::abs((*p)->momentum().px()/(*p)->momentum().pz() - tgx)
00194 > deltaEE ||
00195 std::abs((*p)->momentum().py()/(*p)->momentum().pz() - tgy)
00196 > deltaEE) continue;
00197
00198 etPhoton+=pt;
00199 if(charge3 && pt<2) etPhotonCharged+=pt;
00200 if(apid>100 && apid!=310 && pt>ptMaxHadron) ptMaxHadron=pt;
00201
00202 }
00203
00204 if(etPhoton<ptMin ||etPhoton>ptMax) {
00205
00206 continue;
00207 }
00208
00209
00210
00211
00212
00213
00214
00215
00216
00217
00218
00219
00220
00221
00222
00223
00224
00225
00226
00227
00228
00229
00230 accepted=true;
00231 break;
00232
00233 }
00234
00235
00236 if (accepted) {
00237 theNumberOfSelected++;
00238 std::cout<<" Event preselected "<<theNumberOfSelected<<" Proccess ID "<<myGenEvent->signal_process_id()<<std::endl;
00239 return true;
00240 }
00241 else return false;
00242
00243 }
00244