Belle II Software  release-05-02-19
KKGenInterface.h
1 /**************************************************************************
2  * BASF2 (Belle Analysis Framework 2) *
3  * Copyright(C) 2010 - Belle II Collaboration *
4  * *
5  * Author: The Belle II Collaboration *
6  * Contributors: Martin Ritter, Susanne Koblitz *
7  * *
8  * This software is provided "as is" without any warranty. *
9  **************************************************************************/
10 
11 #pragma once
12 
13 /* Belle 2 headers. */
14 #include <mdst/dataobjects/MCParticleGraph.h>
15 
16 /* C++ headers. */
17 #include <string>
18 #include <unordered_map>
19 
20 #define nmxhep 4000
21 
25 struct hepevt_type {
26  int nevhep;
27  int nhep;
28  int isthep[nmxhep];
29  int idhep[nmxhep];
30  int jmohep[nmxhep][2];
31  int jdahep[nmxhep][2];
32  double phep[nmxhep][5];
33  double vhep[nmxhep][4];
34 };
35 
39 struct pydat2_type {
40  int KCHG[4][500];
41  double PMAS[4][500];
42  double PARF[2000];
43  double VCKM[4][4];
44 };
45 
46 extern hepevt_type hepevt_;
47 extern pydat2_type pydat2_;
48 
49 extern "C" {
50  void kk_init_(const char*, const char*, const char*, int*, const char*);
51  void kk_begin_run_(double*);
52  void kk_init_seed_();
53  void kk_term_(double*, double*);
54  void kk_event_(int*);
55  void kk_getbeam_(double*, double*, double*, double*,
56  double*, double*, double*, double*);
57  void kk_putbeam_(double*, double*, double*, double*,
58  double*, double*, double*, double*);
59  int pycomp_(int&);
60 
61 }
62 
63 namespace Belle2 {
73 
74  public:
75 
80 
85 
89  KKGenInterface(const KKGenInterface& m) = delete;
90 
94  KKGenInterface& operator= (const KKGenInterface& m) = delete;
95 
99  int setup(const std::string& KKdefaultFileName, const std::string& tauinputFileName, const std::string& taudecaytableFileName,
100  const std::string& KKMCOutputFileName);
101 
105  void set_beam_info(TLorentzVector P4_LER, double Espread_LER, TLorentzVector P4_HER,
106  double Espread_HER);
107 
111  int simulateEvent(MCParticleGraph& graph, TVector3 vertex);
112 
116  void term();
117 
118  private:
119 
123  int addParticles2Graph(MCParticleGraph& graph, TVector3 vertex);
124 
129  TVector3 vertex);
130 
134  std::unordered_map<int, int> m_mapPythiaIDtoPDG;
135 
136  };
137 
139 }
Belle2::KKGenInterface
Interface class for using the KKMC generator.
Definition: KKGenInterface.h:72
Belle2::MCParticleGraph
Class to build, validate and sort a particle decay chain.
Definition: MCParticleGraph.h:48
pydat2_type::PARF
double PARF[2000]
parametrization of dd-uu-ss flavor mixing.
Definition: KKGenInterface.h:42
pydat2_type::PMAS
double PMAS[4][500]
particle information such as mass, width...
Definition: KKGenInterface.h:41
pydat2_type::KCHG
int KCHG[4][500]
particle information such as spin, charge...
Definition: KKGenInterface.h:40
pydat2_type
PYDAT2 common block of PYTHIA6.
Definition: KKGenInterface.h:39
Belle2::KKGenInterface::setup
int setup(const std::string &KKdefaultFileName, const std::string &tauinputFileName, const std::string &taudecaytableFileName, const std::string &KKMCOutputFileName)
Setup for KKMC and TAUOLA.
Definition: KKGenInterface.cc:32
Belle2::KKGenInterface::term
void term()
Terminate the generator.
Definition: KKGenInterface.cc:264
Belle2::KKGenInterface::addParticles2Graph
int addParticles2Graph(MCParticleGraph &graph, TVector3 vertex)
Add particles to the MCParticleGraph.
Definition: KKGenInterface.cc:169
hepevt_type::isthep
int isthep[nmxhep]
status code.
Definition: KKGenInterface.h:28
Belle2::KKGenInterface::m_mapPythiaIDtoPDG
std::unordered_map< int, int > m_mapPythiaIDtoPDG
Map between PYTHIA id and PDG codes.
Definition: KKGenInterface.h:134
hepevt_type::jdahep
int jdahep[nmxhep][2]
childreen particles.
Definition: KKGenInterface.h:31
hepevt_type::phep
double phep[nmxhep][5]
four-momentum, mass [GeV].
Definition: KKGenInterface.h:32
Belle2
Abstract base class for different kinds of events.
Definition: MillepedeAlgorithm.h:19
Belle2::KKGenInterface::~KKGenInterface
~KKGenInterface()
Destructor.
Definition: KKGenInterface.h:84
hepevt_type::jmohep
int jmohep[nmxhep][2]
parent particles.
Definition: KKGenInterface.h:30
Belle2::KKGenInterface::set_beam_info
void set_beam_info(TLorentzVector P4_LER, double Espread_LER, TLorentzVector P4_HER, double Espread_HER)
Setup for beams information.
Definition: KKGenInterface.cc:79
Belle2::KKGenInterface::KKGenInterface
KKGenInterface()
Constructor.
Definition: KKGenInterface.h:79
Belle2::KKGenInterface::operator=
KKGenInterface & operator=(const KKGenInterface &m)=delete
Assignment operator, explicitly deleted.
hepevt_type
HEPEVT common block of PYTHIA6.
Definition: KKGenInterface.h:25
hepevt_type::nevhep
int nevhep
serial number.
Definition: KKGenInterface.h:26
pydat2_type::VCKM
double VCKM[4][4]
squared CKM matrix elements.
Definition: KKGenInterface.h:43
hepevt_type::nhep
int nhep
number of particles.
Definition: KKGenInterface.h:27
Belle2::KKGenInterface::simulateEvent
int simulateEvent(MCParticleGraph &graph, TVector3 vertex)
Simulate the events.
Definition: KKGenInterface.cc:124
hepevt_type::vhep
double vhep[nmxhep][4]
vertex [mm].
Definition: KKGenInterface.h:33
Belle2::MCParticleGraph::GraphParticle
Class to represent Particle data in graph.
Definition: MCParticleGraph.h:86
Belle2::KKGenInterface::updateGraphParticle
void updateGraphParticle(int, MCParticleGraph::GraphParticle *gParticle, TVector3 vertex)
Update the MCParticleGraph.
Definition: KKGenInterface.cc:202
hepevt_type::idhep
int idhep[nmxhep]
particle ident KF.
Definition: KKGenInterface.h:29