single_photos_gun_example.cxx
00001
00002
00003
00004
00005
00006
00007
00008
00009
00010 #include "Pythia.h"
00011 #include "HepMCInterface.h"
00012
00013
00014 #include "Photos/Photos.h"
00015 #include "Photos/PhotosHepMCParticle.h"
00016 #include "Photos/Log.h"
00017
00018 using namespace std;
00019 using namespace Pythia8;
00020 using namespace Photospp;
00021
00022 int EventsToCheck=20;
00023
00024
00025
00026
00027
00028 void checkMomentumConservationInEvent(HepMC::GenEvent *evt)
00029 {
00030
00031
00032 double px=0.0,py=0.0,pz=0.0,e=0.0;
00033
00034 for ( HepMC::GenEvent::particle_const_iterator p = evt->particles_begin();
00035 p != evt->particles_end(); ++p )
00036 {
00037 if( (*p)->status() == 1 )
00038 {
00039 HepMC::FourVector m = (*p)->momentum();
00040 px+=m.px();
00041 py+=m.py();
00042 pz+=m.pz();
00043 e +=m.e();
00044
00045 }
00046 }
00047 cout.precision(6);
00048 cout.setf(ios_base::floatfield);
00049 cout<<endl<<"Vector Sum: "<<px<<" "<<py<<" "<<pz<<" "<<e<<endl;
00050 }
00051
00052 int main(int argc,char **argv)
00053 {
00054
00055 HepMC::I_Pythia8 ToHepMC;
00056 Pythia pythia;
00057 Event& event = pythia.event;
00058 pythia.readString("WeakSingleBoson:ffbar2gmZ = on");
00059 pythia.readString("23:onMode = off");
00060 pythia.readString("23:onIfAny = 15");
00061 pythia.readString("HadronLevel:Hadronize = off");
00062 pythia.readString("SpaceShower:QEDshower = off");
00063 pythia.readString("SpaceShower:QEDshowerByL = off");
00064 pythia.readString("SpaceShower:QEDshowerByQ = off");
00065 pythia.readString("PartonLevel:ISR = off");
00066 pythia.readString("PartonLevel:FSR = off");
00067 pythia.init( 11, -11, 92.);
00068
00069 Photos::initialize();
00070 Photos::setInfraredCutOff(0.001/200);
00071
00072 int NumberOfEvents = 10000;
00073 if(argc>1) NumberOfEvents=atoi(argv[1]);
00074
00075 int photonAdded=0,twoAdded=0,moreAdded=0,tauCount=0;
00076
00077
00078 for (int iEvent = 0; iEvent < NumberOfEvents; ++iEvent)
00079 {
00080 if(iEvent%(NumberOfEvents/10)==0) Log::Info()<<iEvent<<endl;
00081 if(!pythia.next()) continue;
00082
00083 HepMC::GenEvent * HepMCEvt = new HepMC::GenEvent();
00084 ToHepMC.fill_next_event(event, HepMCEvt);
00085
00086 if(iEvent<EventsToCheck)
00087 {
00088 cout<<" "<<endl;
00089 cout<<"Momentum conservation chceck BEFORE/AFTER Photos"<<endl;
00090 checkMomentumConservationInEvent(HepMCEvt);
00091 }
00092
00093
00094 HepMC::GenParticle *tau=0;
00095 for(HepMC::GenEvent::vertex_const_iterator i = HepMCEvt->vertices_begin();i!=HepMCEvt->vertices_end();i++)
00096 {
00097 for(HepMC::GenVertex::particles_in_const_iterator p=(*i)->particles_in_const_begin();p!=(*i)->particles_in_const_end(); p++)
00098 {
00099 if((*p)->pdg_id()==15) tau=*p;
00100 break;
00101 }
00102 if(tau) break;
00103 }
00104 if(tau)
00105 {
00106 tauCount++;
00107 int buf = -HepMCEvt->particles_size();
00108
00109
00110 Photos::processParticle( new PhotosHepMCParticle(tau) );
00111
00112 buf+=HepMCEvt->particles_size();
00113 if(buf==1) photonAdded++;
00114 else if(buf==2) twoAdded++;
00115 else if(buf>2) moreAdded++;
00116 }
00117
00118 if(iEvent<EventsToCheck)
00119 {
00120 checkMomentumConservationInEvent(HepMCEvt);
00121 }
00122
00123
00124 delete HepMCEvt;
00125 }
00126 pythia.statistics();
00127
00128
00129 cout.precision(2);
00130 cout.setf(ios::fixed);
00131 cout<<endl;
00132 if(tauCount==0)
00133 {
00134 cout<<"Something went wrong with pythia generation."<<endl;
00135 cout<<"No taus were processed."<<endl<<endl;
00136 return 0;
00137 }
00138 cout<<"Summary (single tau decay processing):"<<endl;
00139 cout<<tauCount <<"\ttaus processed"<<endl;
00140 cout<<photonAdded<<"\ttimes one photon added to the decay \t("<<(photonAdded*100./tauCount)<<"%)"<<endl;
00141 cout<<twoAdded <<"\ttimes two photons added to the decay \t("<<(twoAdded*100./tauCount)<<"%)"<<endl;
00142 cout<<moreAdded <<"\ttimes more than two photons added to the decay\t("<<(moreAdded*100./tauCount)<<"%)"<<endl<<endl;
00143 cout<<"(Contrary to results from MC-Tester, these values are technical and infrared unstable)"<<endl<<endl;
00144 cout<<"To proccess different number of events use:"<<endl<<" ./single_photos_gun_example <number_of_events>"<<endl<<endl;
00145 }
00146