StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
tauola_photos_pythia_example.cxx
1 
9 //pythia header files
10 #include "Pythia8/Pythia.h"
11 #include "Pythia8Plugins/HepMC2.h"
12 
13 //MC-TESTER header files
14 #include "Generate.h"
15 #include "HepMCEvent.H"
16 #include "Setup.H"
17 
18 //PHOTOS header files
19 #include "Photos/Photos.h"
20 #include "Photos/PhotosHepMCEvent.h"
21 #include "Photos/Log.h"
22 
23 //TAUOLA header files
24 #include "Tauola/Tauola.h"
25 #include "Tauola/TauolaHepMCEvent.h"
26 
27 using namespace std;
28 using namespace Pythia8;
29 using namespace Photospp;
30 using namespace Tauolapp;
31 
32 unsigned long NumberOfEvents = 10000;
33 unsigned int EventsToCheck=20;
34 
35 // elementary test of HepMC typically executed before
36 // detector simulation based on http://home.fnal.gov/~mrenna/HCPSS/HCPSShepmc.html
37 // similar test was performed in Fortran
38 // we perform it before and after Photos (for the first several events)
39 void checkMomentumConservationInEvent(HepMC::GenEvent *evt)
40 {
41  //cout<<"List of stable particles: "<<endl;
42 
43  double px=0.0,py=0.0,pz=0.0,e=0.0;
44 
46  p != evt->particles_end(); ++p )
47  {
48  if( (*p)->status() == 1 )
49  {
50  HepMC::FourVector m = (*p)->momentum();
51  px+=m.px();
52  py+=m.py();
53  pz+=m.pz();
54  e +=m.e();
55  //(*p)->print();
56  }
57  }
58  cout.precision(6);
59  cout.setf(ios_base::floatfield);
60  cout<<endl<<"Vector Sum: "<<px<<" "<<py<<" "<<pz<<" "<<e<<endl;
61 }
62 
63 int main(int argc,char **argv)
64 {
65  HepMC::Pythia8ToHepMC ToHepMC;
66 
67  // Initialization of pythia
68  Pythia pythia;
69  Event& event = pythia.event;
70 
71  pythia.readFile("tauola_photos_pythia_example.conf");
72  pythia.init();
73 
74  // TAUOLA and PHOTOS initialization
75  Tauola::initialize();
76  Photos::initialize();
77 
78  Photos::setInfraredCutOff(0.01/91.187); // 10MeV for scale to M_Z=91.187
79  //Photos::setDoubleBrem(false);
80  //Photos::setExponentiation(false);
81 
82  Log::SummaryAtExit();
83  cout.setf(ios::fixed);
84  //Log::LogInfo(false) //To turn printing of last five events and pythia statistics off
85 
86  // Example setup - suppress processing of whole Z0 decay,
87  // leaving only the Z0 -> tau+ tau- decay and whole branch starting
88  // from tau- to be processed
89  //Photos::suppressBremForBranch(0,23);
90  //Photos::forceBremForDecay (2,23,15,-15);
91  //Photos::forceBremForBranch(0,15);
92 
93  // Force mass of electron and positron to be 0.000511
94  //Photos::forceMass(11,0.000511);
95 
96  // Force mass of electron and positron to be taken
97  // from event record instead of being calculated from 4-vector
98  //Photos::forceMassFromEventRecord(11);
99 
100  // Exclude particles with given status code from being processed
101  // or taken into account during momentum conservation calculation
102  //Photos::ignoreParticlesWithStatus(3);
103 
104  // Remove status code from the list of ignored status codes
105  //Photos::DeIgnoreParticlesWithStatus(3);
106 
107  // Force writing history decay products for vertices
108  // modified i.e. with added photons. These particles will
109  // have the provided status code. Photos will ignore
110  // all particles with this status code.
111  //Photos::createHistoryEntries(true,3);
112 
113  MC_Initialize();
114 
115  // Begin event loop
116  for (unsigned long iEvent = 0; iEvent < NumberOfEvents; ++iEvent)
117  {
118  if(iEvent%1000==0) Log::Info()<<"Event: "<<iEvent<<"\t("<<iEvent*(100./NumberOfEvents)<<"%)"<<endl;
119  if(!pythia.next()) continue;
120 
121  // Convert event record to HepMC
122  HepMC::GenEvent * HepMCEvt = new HepMC::GenEvent();
123  ToHepMC.fill_next_event(event, HepMCEvt);
124 
125  // Run TAUOLA on the event
126  TauolaHepMCEvent * t_event = new TauolaHepMCEvent(HepMCEvt);
127 
128  // We may want to undecay previously decayed taus.
129  //t_event->undecayTaus();
130  t_event->decayTaus();
131  delete t_event;
132 
133  //Log::LogPhlupa(2,4);
134 
135  if(iEvent<EventsToCheck)
136  {
137  cout<<" "<<endl;
138  cout<<"Momentum conservation chceck BEFORE/AFTER Photos"<<endl;
139  checkMomentumConservationInEvent(HepMCEvt);
140  }
141 
142  // Run PHOTOS on the event
143  PhotosHepMCEvent evt(HepMCEvt);
144  evt.process();
145 
146  if(iEvent<EventsToCheck)
147  {
148  checkMomentumConservationInEvent(HepMCEvt);
149  }
150 
151  // Run MC-TESTER on the event
152  HepMCEvent temp_event(*HepMCEvt,false);
153  MC_Analyze(&temp_event);
154 
155  // Print out last 5 events
156  if(iEvent>=NumberOfEvents-5)
157  {
158  Log::RedirectOutput(Log::Info());
159  //pythia.event.list();
160  HepMCEvt->print();
161  Log::RevertOutput();
162  }
163 
164  // Clean up
165  delete HepMCEvt;
166  }
167 
168  Log::RedirectOutput(Log::Info());
169  pythia.stat();
170  Log::RevertOutput();
171 
172  MC_Finalize();
173 }
const particle iterator
Definition: GenEvent.h:464
particle_const_iterator particles_begin() const
begin particle iteration
Definition: GenEvent.h:507
double e() const
return E
Definition: SimpleVector.h:73
double px() const
return px
Definition: SimpleVector.h:70
The GenEvent class is the core of HepMC.
Definition: GenEvent.h:155
double pz() const
return pz
Definition: SimpleVector.h:72
void print(std::ostream &ostr=std::cout) const
dumps to ostr
Definition: GenEvent.cc:277
particle_const_iterator particles_end() const
end particle iteration
Definition: GenEvent.h:511
FourVector is a simple representation of a physics 4 vector.
Definition: SimpleVector.h:42
double py() const
return py
Definition: SimpleVector.h:71