Belle II Software  release-06-01-15
KKGenInterface.h
1 /**************************************************************************
2  * basf2 (Belle II Analysis Software Framework) *
3  * Author: The Belle II Collaboration *
4  * *
5  * See git log for contributors and copyright holders. *
6  * This file is licensed under LGPL-3.0, see LICENSE.md. *
7  **************************************************************************/
8 
9 #pragma once
10 
11 /* Belle 2 headers. */
12 #include <mdst/dataobjects/MCParticleGraph.h>
13 
14 /* C++ headers. */
15 #include <string>
16 #include <unordered_map>
17 
18 #define nmxhep 4000
19 
23 struct hepevt_type {
24  int nevhep;
25  int nhep;
26  int isthep[nmxhep];
27  int idhep[nmxhep];
28  int jmohep[nmxhep][2];
29  int jdahep[nmxhep][2];
30  double phep[nmxhep][5];
31  double vhep[nmxhep][4];
32 };
33 
37 struct pydat2_type {
38  int KCHG[4][500];
39  double PMAS[4][500];
40  double PARF[2000];
41  double VCKM[4][4];
42 };
43 
44 extern hepevt_type hepevt_;
45 extern pydat2_type pydat2_;
46 
47 extern "C" {
48  void kk_init_(const char*, const char*, const char*, int*, const char*);
49  void kk_begin_run_(double*);
50  void kk_init_seed_();
51  void kk_term_(double*, double*);
52  void kk_event_(int*);
53  void kk_getbeam_(double*, double*, double*, double*,
54  double*, double*, double*, double*);
55  void kk_putbeam_(double*, double*, double*, double*,
56  double*, double*, double*, double*);
57  int pycomp_(int&);
58 
59 }
60 
61 namespace Belle2 {
71 
72  public:
73 
78 
83 
87  KKGenInterface(const KKGenInterface& m) = delete;
88 
93 
97  int setup(const std::string& KKdefaultFileName, const std::string& tauinputFileName, const std::string& taudecaytableFileName,
98  const std::string& KKMCOutputFileName);
99 
103  void set_beam_info(TLorentzVector P4_LER, double Espread_LER, TLorentzVector P4_HER,
104  double Espread_HER);
105 
109  int simulateEvent(MCParticleGraph& graph, TVector3 vertex);
110 
114  void term();
115 
116  private:
117 
121  int addParticles2Graph(MCParticleGraph& graph, TVector3 vertex);
122 
127  TVector3 vertex);
128 
132  std::unordered_map<int, int> m_mapPythiaIDtoPDG;
133 
134  };
135 
137 }
Interface class for using the KKMC generator.
int simulateEvent(MCParticleGraph &graph, TVector3 vertex)
Simulate the events.
KKGenInterface()
Constructor.
std::unordered_map< int, int > m_mapPythiaIDtoPDG
Map between PYTHIA id and PDG codes.
void term()
Terminate the generator.
KKGenInterface(const KKGenInterface &m)=delete
Copy constructor, explicitly deleted.
int addParticles2Graph(MCParticleGraph &graph, TVector3 vertex)
Add particles to the MCParticleGraph.
void set_beam_info(TLorentzVector P4_LER, double Espread_LER, TLorentzVector P4_HER, double Espread_HER)
Setup for beams information.
~KKGenInterface()
Destructor.
int setup(const std::string &KKdefaultFileName, const std::string &tauinputFileName, const std::string &taudecaytableFileName, const std::string &KKMCOutputFileName)
Setup for KKMC and TAUOLA.
void updateGraphParticle(int, MCParticleGraph::GraphParticle *gParticle, TVector3 vertex)
Update the MCParticleGraph.
KKGenInterface & operator=(const KKGenInterface &m)=delete
Assignment operator, explicitly deleted.
Class to represent Particle data in graph.
Class to build, validate and sort a particle decay chain.
Abstract base class for different kinds of events.
HEPEVT common block of PYTHIA6.
int jmohep[nmxhep][2]
parent particles.
int isthep[nmxhep]
status code.
double vhep[nmxhep][4]
vertex [mm].
int jdahep[nmxhep][2]
childreen particles.
double phep[nmxhep][5]
four-momentum, mass [GeV].
int nevhep
serial number.
int idhep[nmxhep]
particle ident KF.
int nhep
number of particles.
PYDAT2 common block of PYTHIA6.
int KCHG[4][500]
particle information such as spin, charge...
double VCKM[4][4]
squared CKM matrix elements.
double PARF[2000]
parametrization of dd-uu-ss flavor mixing.
double PMAS[4][500]
particle information such as mass, width...