evtgen is hosted by Hepforge, IPPP Durham
EvtGen  2.0.0
Monte Carlo generator of particle decays, in particular the weak decays of heavy flavour particles such as B mesons.
EvtPythiaEngine.hh
Go to the documentation of this file.
1 
2 /***********************************************************************
3 * Copyright 1998-2020 CERN for the benefit of the EvtGen authors *
4 * *
5 * This file is part of EvtGen. *
6 * *
7 * EvtGen is free software: you can redistribute it and/or modify *
8 * it under the terms of the GNU General Public License as published by *
9 * the Free Software Foundation, either version 3 of the License, or *
10 * (at your option) any later version. *
11 * *
12 * EvtGen is distributed in the hope that it will be useful, *
13 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
14 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
15 * GNU General Public License for more details. *
16 * *
17 * You should have received a copy of the GNU General Public License *
18 * along with EvtGen. If not, see <https://www.gnu.org/licenses/>. *
19 ***********************************************************************/
20 
21 #ifdef EVTGEN_PYTHIA
22 
23 #ifndef EVTPYTHIAENGINE_HH
24 #define EVTPYTHIAENGINE_HH
25 
27 #include "EvtGenBase/EvtId.hh"
30 
32 
34 
35 #include "Pythia8/ParticleData.h"
36 #include "Pythia8/Pythia.h"
37 
38 #include <map>
39 #include <memory>
40 #include <string>
41 #include <vector>
42 
43 // Description: Interface to the Pytha 8 external generator
44 
45 class EvtPythiaEngine : public EvtAbsExternalGen {
46  public:
47  EvtPythiaEngine( std::string xmlDir = "./xmldoc",
48  bool convertPhysCodes = false, bool useEvtGenRandom = true );
49 
50  virtual ~EvtPythiaEngine();
51 
52  bool doDecay( EvtParticle* theMother ) override;
53 
54  void initialise() override;
55 
56  protected:
57  private:
58  void updateParticleLists();
59  void updatePhysicsParameters();
60 
61  void createPythiaParticle( EvtId& particleId, int PDGCode );
62  bool validPDGCode( int PDGCode );
63  void updatePythiaDecayTable( EvtId& particleId, int aliasInt, int PDGCode );
64  void storeDaughterInfo( EvtParticle* theParticle, int startInt );
65 
66  void clearDaughterVectors();
67  void clearPythiaModeMap();
68 
69  void createDaughterEvtParticles( EvtParticle* theParent );
70 
71  int getModeInt( EvtDecayBase* decayModel );
72 
73  std::unique_ptr<Pythia8::Pythia> _genericPythiaGen;
74  std::unique_ptr<Pythia8::Pythia> _aliasPythiaGen;
75  Pythia8::Pythia* _thePythiaGenerator;
76 
77  std::vector<int> _daugPDGVector;
78  std::vector<EvtVector4R> _daugP4Vector;
79 
80  typedef std::map<int, std::vector<int>> PythiaModeMap;
81  PythiaModeMap _pythiaModeMap;
82 
83  bool _convertPhysCodes, _initialised, _useEvtGenRandom;
84 
85  std::unique_ptr<EvtPythiaRandom> _evtgenRandom;
86 
87  std::map<int, int> _addedPDGCodes;
88 };
89 
90 #endif
91 
92 #endif
Definition: EvtId.hh:27
virtual void initialise()=0
virtual bool doDecay(EvtParticle *theMother)=0