StRoot  1
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
taumain_stand_alone_example.c
1 
10 #include "Tauola/Tauola.h"
11 #include "Tauola/TauolaHepMCEvent.h"
12 
13 //HepMC headers
14 #include "HepMC/GenEvent.h"
15 #include "HepMC/GenParticle.h"
16 
17 #include "tauola_print_parameters.h"
18 
19 //#include "tauola_set_seed.h" // uncomment for old releases
20 // and use setSeed( 47238, 985439, 0 );
21 // In new libraries Tauola::setSeed( 47238, 985439, 0 );
22 // is available.
23 
24 using namespace std;
25 using namespace Tauolapp;
26 
30 HepMC::GenEvent * make_simple_tau_event(){
31 
32  HepMC::GenEvent * event = new HepMC::GenEvent();
33 
34  //Create some four vectors for the electrons
35  double ele_mass_sqr = parmas_.amell*parmas_.amell;
36  HepMC::FourVector momentum_e1(0,0,0,0);
37  HepMC::FourVector momentum_e2(0,0,0,0);
38 
39  momentum_e1.setPz(-2); //change these
40  momentum_e2.setPz(3.5); //as needed
41 
42  momentum_e1.setE(sqrt(momentum_e1.pz()*momentum_e1.pz()+ele_mass_sqr));
43  momentum_e2.setE(sqrt(momentum_e2.pz()*momentum_e2.pz()+ele_mass_sqr));
44 
45  //make taus
46  double tau_mass = parmas_.amtau;
47 
48  //make tau's four vector
49  HepMC::FourVector momentum_tau1(0,0,0,0);
50  HepMC::FourVector momentum_tau2(0,0,0,0);
51 
52  //make particles
53  HepMC::GenParticle * e1 = new HepMC::GenParticle(momentum_e1,-11,3);
54  HepMC::GenParticle * e2 = new HepMC::GenParticle(momentum_e2,11,3);
55  HepMC::GenParticle * tau1 = new HepMC::GenParticle(momentum_tau1,-15,1);
56  HepMC::GenParticle * tau2 = new HepMC::GenParticle(momentum_tau2,15,1);
57 
58  //set the masses
59  e1->set_generated_mass(parmas_.amell);
60  e2->set_generated_mass(parmas_.amell);
61  tau1->set_generated_mass(tau_mass);
62  tau2->set_generated_mass(tau_mass);
63 
64  //make the vertex
66  vertex->add_particle_in(e1);
67  vertex->add_particle_in(e2);
68  vertex->add_particle_out(tau1);
69  vertex->add_particle_out(tau2);
70 
71  event->add_vertex(vertex);
72 
73  //calculate center of mass frame
74  HepMC::FourVector cms(0,0,(momentum_e1.pz()+momentum_e2.pz()),
75  momentum_e1.e()+momentum_e2.e());
76 
77  HepMC::GenParticle * cms_particle = new HepMC::GenParticle(cms,0,0);
78 
79  //Make TauolaParticles for doing boosting:
80  TauolaHepMCParticle * cms_boost = new TauolaHepMCParticle(cms_particle);
81 
82  TauolaHepMCParticle first_tau(tau1);
83  TauolaHepMCParticle second_tau(tau2);
84 
85  double tau_energy = 0.5*sqrt(cms.e()*cms.e() - (cms.px()*cms.px()
86  + cms.py()*cms.py()+cms.pz()*cms.pz()));
87 
88  first_tau.setE(tau_energy);
89  first_tau.setPx((1.0/sqrt(2.0))*sqrt(tau_energy*tau_energy-tau_mass*tau_mass));
90  first_tau.setPy((1.0/sqrt(2.0))*sqrt(tau_energy*tau_energy-tau_mass*tau_mass));
91 
92  second_tau.setE(tau_energy);
93  second_tau.setPx(-1*(1.0/sqrt(2.0))*sqrt(tau_energy*tau_energy-tau_mass*tau_mass));
94  second_tau.setPy(-1*(1.0/sqrt(2.0))*sqrt(tau_energy*tau_energy-tau_mass*tau_mass));
95 
96  first_tau.boostFromRestFrame(cms_boost);
97  second_tau.boostFromRestFrame(cms_boost);
98 
99  //clean up
100  delete cms_boost;
101  delete cms_particle;
102 
103  return event;
104 }
105 
107 int main(void){
108 
109  int NumberOfEvents = 10;
110 
111  //These three lines are not really necessary since they are the default
112  Tauola::setDecayingParticle(15);
113  Tauola::setSameParticleDecayMode(0);
114  Tauola::setOppositeParticleDecayMode(0);
115 
116  // Change TAUOLA-FORTRAN random generator seed:
117  // Tauola::setSeed( 47238, 985439, 0 );
118 
119  // for older releases use setSeed( 47238, 985439, 0 ); // and tauola_set_seed.h
120 
121  Tauola::initialize();
122 
123  tauola_print_parameters(); // Prints TAUOLA parameters (residing inside its library): e.g. to test user interface
124 
125  // Begin event loop. Generate event.
126  for (int iEvent = 0; iEvent < NumberOfEvents; ++iEvent) {
127 
128  // Convert event record to Tau format
129 
130  HepMC::GenEvent * event = make_simple_tau_event();
131 
132  cout << "BEFORE:"<<endl;
133  event->print();
134  TauolaHepMCEvent * t_event = new TauolaHepMCEvent(event);
135  t_event->decayTaus();
136 
137  cout << "AFTER:"<<endl;
138  event->print();
139 
140  //clean up
141  delete event;
142  delete t_event;
143  }
144 
145  // This is an access to old FORTRAN info on generated tau sample.
146  // That is why it refers to old version number (eg. 2.7) for TAUOLA.
147  //Tauola::summary();
148 }
149 
void set_generated_mass(const double &m)
define the actual generated mass
Definition: GenParticle.cc:240
void add_particle_in(GenParticle *inparticle)
add incoming particle
Definition: GenVertex.cc:273
GenVertex contains information about decay vertices.
Definition: GenVertex.h:52
The GenEvent class is the core of HepMC.
Definition: GenEvent.h:155
void add_particle_out(GenParticle *outparticle)
add outgoing particle
Definition: GenVertex.cc:284
FourVector is a simple representation of a physics 4 vector.
Definition: SimpleVector.h:42
The GenParticle class contains information about generated particles.
Definition: GenParticle.h:60