Belle II Software development
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/* Basf2 headers. */
12#include <mdst/dataobjects/MCParticleGraph.h>
13#include <framework/utilities/ConditionalGaussGenerator.h>
14
15/* C++ headers. */
16#include <string>
17#include <unordered_map>
18
19#define nmxhep 4000
20
25 int nevhep;
26 int nhep;
27 int isthep[nmxhep];
28 int idhep[nmxhep];
29 int jmohep[nmxhep][2];
30 int jdahep[nmxhep][2];
31 double phep[nmxhep][5];
32 double vhep[nmxhep][4];
33};
34
39 int KCHG[4][500];
40 double PMAS[4][500];
41 double PARF[2000];
42 double VCKM[4][4];
43};
44
45extern hepevt_type hepevt_;
46extern pydat2_type pydat2_;
47
48extern "C" {
49 void kk_init_(const char*, const char*, const char*, int*, const char*);
50 void kk_begin_run_(double*, double*);
51 void kk_init_seed_();
52 void kk_term_(double*, double*);
53 void kk_event_(int*);
54 void kk_shifttaudecayvtx_();
55 int pycomp_(int&);
56
57}
58
59namespace Belle2 {
69
70 public:
71
76
81
85 KKGenInterface(const KKGenInterface& m) = delete;
86
91
95 int setup(const std::string& KKdefaultFileName, const std::string& tauinputFileName, const std::string& taudecaytableFileName,
96 const std::string& KKMCOutputFileName);
97
101 void set_beam_info(double Ecms0, double Ecms0Spread);
102
106 int simulateEvent(MCParticleGraph& graph, const ConditionalGaussGenerator& lorentzGenerator, ROOT::Math::XYZVector vertex);
107
111 void term();
112
113 private:
114
118 int addParticles2Graph(MCParticleGraph& graph, ROOT::Math::XYZVector vertex);
119
124 ROOT::Math::XYZVector vertex);
125
129 std::unordered_map<int, int> m_mapPythiaIDtoPDG;
130
131 };
132
134}
Class implementing n-dimensional random number generator from Gaussian distribution where the first c...
Interface class for using the KKMC generator.
KKGenInterface & operator=(const KKGenInterface &m)=delete
Assignment operator, explicitly deleted.
int simulateEvent(MCParticleGraph &graph, const ConditionalGaussGenerator &lorentzGenerator, ROOT::Math::XYZVector 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, ROOT::Math::XYZVector vertex)
Add particles to the MCParticleGraph.
~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, ROOT::Math::XYZVector vertex)
Update the MCParticleGraph.
void set_beam_info(double Ecms0, double Ecms0Spread)
Setup for beams information.
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...