Class Reference for E1039 Core & Analysis Software
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
PHPythia8.h
Go to the documentation of this file.
1 #ifndef __PHPYTHIA8_H__
2 #define __PHPYTHIA8_H__
3 
4 #include <fun4all/SubsysReco.h>
6 
7 #ifndef __CINT__
8 #include <Pythia8/Pythia.h>
9 #endif
10 
11 #ifndef __CINT__
12 #include <gsl/gsl_rng.h>
13 #endif
14 
15 #include <iostream>
16 #include <string>
17 
18 class PHCompositeNode;
19 class PHHepMCGenEvent;
20 class PHHepMCFilter;
21 class PHGenIntegral;
22 class PHPy8GenTrigger;
23 
24 namespace HepMC
25 {
26 class GenEvent;
27 class Pythia8ToHepMC;
28 };
29 
30 namespace Pythia8
31 {
32 class Pythia;
33 };
34 
35 class PHPythia8 : public SubsysReco
36 {
37  public:
38  PHPythia8(const std::string &name = "PHPythia8");
39  virtual ~PHPythia8();
40 
41  int Init(PHCompositeNode *topNode);
42  int process_event(PHCompositeNode *topNode);
43  int ResetEvent(PHCompositeNode *topNode);
44  int End(PHCompositeNode *topNode);
45 
46  void set_config_file(const char *cfg_file)
47  {
48  if (cfg_file) _configFile = cfg_file;
49  }
50 
51  void print_config() const;
52 
54  void register_trigger(PHPy8GenTrigger *theTrigger);
56  {
57  _triggersOR = true;
58  _triggersAND = false;
59  } // default true
61  {
62  _triggersAND = true;
63  _triggersOR = false;
64  }
65 
67  void process_string(std::string s) { _commands.push_back(s); }
68  void beam_vertex_parameters(double beamX,
69  double beamY,
70  double beamZ,
71  double beamXsigma,
72  double beamYsigma,
73  double beamZsigma)
74  {
75  set_vertex_distribution_mean(beamX, beamY, beamZ, 0);
76  set_vertex_distribution_width(beamXsigma, beamYsigma, beamZsigma, 0);
77  }
78 
81  {
82  hepmc_helper.set_vertex_distribution_function(x, y, z, t);
83  }
84 
86  void set_vertex_distribution_mean(const double x, const double y, const double z, const double t)
87  {
88  hepmc_helper.set_vertex_distribution_mean(x, y, z, t);
89  }
90 
92  void set_vertex_distribution_width(const double x, const double y, const double z, const double t)
93  {
94  hepmc_helper.set_vertex_distribution_width(x, y, z, t);
95  }
96  //
98  void set_reuse_vertex(int src_embedding_id)
99  {
100  hepmc_helper.set_reuse_vertex(src_embedding_id);
101  }
102 
107  int get_embedding_id() const { return hepmc_helper.get_embedding_id(); }
108  //
113  void set_embedding_id(int id) { hepmc_helper.set_embedding_id(id); }
115  void save_integrated_luminosity(const bool b) { _save_integrated_luminosity = b; }
116  private:
117  int read_config(const char *cfg_file = 0);
118  int create_node_tree(PHCompositeNode *topNode);
119  double percent_diff(const double a, const double b) { return abs((a - b) / a); }
120  int _eventcount;
121 
122  // event selection
123  std::vector<PHPy8GenTrigger *> _registeredTriggers;
124  bool _triggersOR;
125  bool _triggersAND;
126 
127 // PYTHIA
128 #ifndef __CINT__
129  Pythia8::Pythia *_pythia;
130 #endif
131 
132  std::string _configFile;
133  std::vector<std::string> _commands;
134 
135  // HepMC
136  HepMC::Pythia8ToHepMC *_pythiaToHepMC;
137 
139  PHHepMCGenHelper hepmc_helper;
140 
142  bool _save_integrated_luminosity;
143 
145  PHGenIntegral *_integral_node;
146 };
147 
148 #endif /* __PHPYTHIA8_H__ */
int get_embedding_id() const
Definition: PHPythia8.h:107
void set_trigger_AND()
Definition: PHPythia8.h:60
void set_trigger_OR()
Definition: PHPythia8.h:55
void set_reuse_vertex(int src_embedding_id)
reuse vertex from another PHHepMCGenEvent with embedding_id = src_embedding_id Additional smearing an...
Definition: PHPythia8.h:98
void set_reuse_vertex(int src_embedding_id)
reuse vertex from another PHHepMCGenEvent with embedding_id = src_embedding_id Additional smearing an...
const double s
void process_string(std::string s)
pass commands directly to PYTHIA8
Definition: PHPythia8.h:67
int get_embedding_id() const
PHHepMCGenHelper provides service of DST upload of HepMC subevent, vertex assignment and random gener...
int End(PHCompositeNode *topNode)
Called at the end of all processing.
Definition: PHPythia8.C:109
void set_vertex_distribution_function(PHHepMCGenHelper::VTXFUNC x, PHHepMCGenHelper::VTXFUNC y, PHHepMCGenHelper::VTXFUNC z, PHHepMCGenHelper::VTXFUNC t)
toss a new vertex according to a Uniform or Gaus distribution
Definition: PHPythia8.h:80
void set_config_file(const char *cfg_file)
Definition: PHPythia8.h:46
void set_vertex_distribution_mean(const double x, const double y, const double z, const double t)
set the mean value of the vertex distribution, use PHENIX units of cm, ns
void set_vertex_distribution_function(VTXFUNC x, VTXFUNC y, VTXFUNC z, VTXFUNC t)
toss a new vertex according to a Uniform or Gaus distribution
void set_vertex_distribution_mean(const double x, const double y, const double z, const double t)
set the mean value of the vertex distribution, use PHENIX units of cm, ns
Definition: PHPythia8.h:86
void print_config() const
Definition: PHPythia8.C:162
void save_integrated_luminosity(const bool b)
whether to store the integrated luminosity and other event statistics to the TOP/RUN/PHGenIntegral no...
Definition: PHPythia8.h:115
int Init(PHCompositeNode *topNode)
Definition: PHPythia8.C:69
void set_vertex_distribution_width(const double x, const double y, const double z, const double t)
set the width of the vertex distribution function about the mean, use PHENIX units of cm...
Definition: PHPythia8.h:92
int process_event(PHCompositeNode *topNode)
Definition: PHPythia8.C:167
void register_trigger(PHPy8GenTrigger *theTrigger)
set event selection criteria
Definition: PHPythia8.C:303
void beam_vertex_parameters(double beamX, double beamY, double beamZ, double beamXsigma, double beamYsigma, double beamZsigma)
Definition: PHPythia8.h:68
int ResetEvent(PHCompositeNode *topNode)
Clean up after each event.
Definition: PHPythia8.C:298
void set_vertex_distribution_width(const double x, const double y, const double z, const double t)
set the width of the vertex distribution function about the mean, use PHENIX units of cm...
virtual ~PHPythia8()
Definition: PHPythia8.C:64
VTXFUNC
supported function distributions
PHPythia8(const std::string &name="PHPythia8")
Definition: PHPythia8.C:32
void set_embedding_id(int id)
Definition: PHPythia8.h:113
void set_embedding_id(int id)
PHGenIntegral.
Definition: PHGenIntegral.h:20