photos_tauola_test.cxx
00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011 #include "Pythia.h"
00012 #include "HepMCInterface.h"
00013
00014
00015 #include "Generate.h"
00016 #include "HepMCEvent.H"
00017 #include "Setup.H"
00018
00019
00020 #include "Tauola/Tauola.h"
00021 #include "Tauola/TauolaHepMCEvent.h"
00022
00023
00024 #include "Photos/Photos.h"
00025 #include "Photos/PhotosHepMCParticle.h"
00026 #include "Photos/PhotosHepMCEvent.h"
00027 #include "Photos/Log.h"
00028
00029 using namespace std;
00030 using namespace Pythia8;
00031 using namespace Photospp;
00032 using namespace Tauolapp;
00033
00034 unsigned long NumberOfEvents = 10000;
00035 unsigned int EventsToCheck=20;
00036
00037
00038
00039
00040
00041 void checkMomentumConservationInEvent(HepMC::GenEvent *evt)
00042 {
00043
00044
00045 double px=0.0,py=0.0,pz=0.0,e=0.0;
00046
00047 for ( HepMC::GenEvent::particle_const_iterator p = evt->particles_begin();
00048 p != evt->particles_end(); ++p )
00049 {
00050 if( (*p)->status() == 1 )
00051 {
00052 HepMC::FourVector m = (*p)->momentum();
00053 px+=m.px();
00054 py+=m.py();
00055 pz+=m.pz();
00056 e +=m.e();
00057
00058 }
00059 }
00060 cout.precision(6);
00061 cout.setf(ios_base::floatfield);
00062 cout<<endl<<"Vector Sum: "<<px<<" "<<py<<" "<<pz<<" "<<e<<endl;
00063 }
00064
00065 int main(int argc,char **argv)
00066 {
00067
00068
00069 if(argc<5)
00070 {
00071 cout<<endl<<"Usage: "<<argv[0]<<" <pythia_conf> <pythia_mode> <no_events> <tauola_mode> [ <alpha_order> <ScalarNLO_mode> ]"<<endl;
00072 cout<<endl<<" eg. "<<argv[0]<<" pythia_H.conf 0 10000 4 0 0"<<endl;
00073 cout<<endl;
00074 return -1;
00075 }
00076
00077 HepMC::I_Pythia8 ToHepMC;
00078
00079
00080 Pythia pythia;
00081 Event& event = pythia.event;
00082
00083 pythia.readString("HadronLevel:Hadronize = off");
00084 pythia.readString("SpaceShower:QEDshower = off");
00085 pythia.readString("SpaceShower:QEDshowerByL = off");
00086 pythia.readString("SpaceShower:QEDshowerByQ = off");
00087 pythia.readString("PartonLevel:ISR = off");
00088 pythia.readString("PartonLevel:FSR = off");
00089
00090
00091
00092
00093
00094
00095
00096
00097
00098
00099
00100
00101
00102
00103
00104
00105
00106
00107
00108
00109
00110
00111
00112
00113
00114
00115 if(argc>1) pythia.readFile(argv[1]);
00116
00117
00118 if(atoi(argv[2])==1) pythia.init( 11, -11, 91.187);
00119 else if(atoi(argv[2])==3) pythia.init( 11, -11, 500.);
00120 else
00121 {
00122 cout<<"ERROR: Wrong Pythia mode ("<<atoi(argv[4])<<")"<<endl;
00123 cout<<" Only modes '1' and '3' are used by this program."<<endl;
00124 return -1;
00125 }
00126
00127
00128 if(argc>3) NumberOfEvents=atoi(argv[3]);
00129
00130
00131 if(argc>4)
00132 {
00133
00134
00135 Tauola::setSameParticleDecayMode(atoi(argv[4]));
00136 Tauola::setOppositeParticleDecayMode(atoi(argv[4]));
00137 }
00138
00139 Tauola::initialize();
00140 Photos::initialize();
00141
00142 Photos::setExponentiation(true);
00143 Photos::setInfraredCutOff(1.e-6);
00144 Photos::maxWtInterference(3.0);
00145
00146
00147 if( argc>5 && atoi(argv[5]) )
00148 {
00149 Photos::setDoubleBrem(false);
00150 Photos::setExponentiation(false);
00151
00152
00153 if(atoi(argv[2])==1) Photos::setInfraredCutOff(0.01/91.187);
00154 else Photos::setInfraredCutOff(0.01/500.);
00155
00156 Photos::maxWtInterference(2.0);
00157 }
00158
00159
00160 if( argc>6 )
00161 {
00162 Tauola::setEtaK0sPi(1,1,0);
00163
00164
00165 if(atoi(argv[6])) Photos::setMeCorrectionWtForScalar(true);
00166 }
00167
00168 Log::SummaryAtExit();
00169 cout.setf(ios::fixed);
00170
00171 MC_Initialize();
00172
00173
00174 for(unsigned long iEvent = 0; iEvent < NumberOfEvents; ++iEvent)
00175 {
00176 if(iEvent%1000==0) Log::Info()<<"Event: "<<iEvent<<"\t("<<iEvent*(100./NumberOfEvents)<<"%)"<<endl;
00177 if(!pythia.next()) continue;
00178
00179 HepMC::GenEvent * HepMCEvt = new HepMC::GenEvent();
00180 ToHepMC.fill_next_event(event, HepMCEvt);
00181
00182 if(iEvent<EventsToCheck)
00183 {
00184 cout<<" "<<endl;
00185 cout<<"Momentum conservation chceck BEFORE/AFTER Photos"<<endl;
00186 checkMomentumConservationInEvent(HepMCEvt);
00187 }
00188
00189
00190 TauolaHepMCEvent * t_event = new TauolaHepMCEvent(HepMCEvt);
00191
00192
00193 t_event->undecayTaus();
00194 t_event->decayTaus();
00195 delete t_event;
00196
00197
00198 PhotosHepMCEvent evt(HepMCEvt);
00199 evt.process();
00200
00201 if(iEvent<EventsToCheck)
00202 {
00203 checkMomentumConservationInEvent(HepMCEvt);
00204 }
00205
00206
00207 HepMCEvent temp_event(*HepMCEvt,false);
00208 MC_Analyze(&temp_event);
00209
00210
00211 delete HepMCEvt;
00212 }
00213 pythia.statistics();
00214 MC_Finalize();
00215 }