Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-08-05 08:16:00

0001 #ifndef PHPYTHIA8_PHPYTHIA8_H
0002 #define PHPYTHIA8_PHPYTHIA8_H
0003 
0004 #include <fun4all/SubsysReco.h>
0005 
0006 #include <phhepmc/PHHepMCGenHelper.h>
0007 
0008 #include <algorithm>  // for max
0009 #include <cmath>
0010 #include <memory>
0011 #include <string>
0012 #include <vector>
0013 
0014 class PHCompositeNode;
0015 class PHGenIntegral;
0016 class PHPy8GenTrigger;
0017 
0018 namespace HepMC
0019 {
0020   class Pythia8ToHepMC;
0021 }  // namespace HepMC
0022 
0023 namespace Pythia8
0024 {
0025   class Pythia;
0026 }
0027 
0028 class PHPythia8 : public SubsysReco, public PHHepMCGenHelper
0029 {
0030  public:
0031 
0032   //! constructor
0033   explicit PHPythia8(const std::string &name = "PHPythia8");
0034 
0035   //! destructor
0036   ~PHPythia8() override = default;
0037 
0038   int Init(PHCompositeNode *topNode) override;
0039   int process_event(PHCompositeNode *topNode) override;
0040   int End(PHCompositeNode *topNode) override;
0041 
0042   void set_config_file(const std::string &cfg_file)
0043   {
0044     m_ConfigFileName = cfg_file;
0045   }
0046 
0047   void print_config() const;
0048 
0049   /// set event selection criteria
0050   void register_trigger(PHPy8GenTrigger *theTrigger);
0051   void set_trigger_OR()
0052   {
0053     m_TriggersOR = true;
0054     m_TriggersAND = false;
0055   }  // default true
0056   void set_trigger_AND()
0057   {
0058     m_TriggersAND = true;
0059     m_TriggersOR = false;
0060   }
0061 
0062   /// pass commands directly to PYTHIA8
0063   void process_string(const std::string &s) { m_Commands.push_back(s); }
0064   void beam_vertex_parameters(double beamX,
0065                               double beamY,
0066                               double beamZ,
0067                               double beamXsigma,
0068                               double beamYsigma,
0069                               double beamZsigma)
0070   {
0071     set_vertex_distribution_mean(beamX, beamY, beamZ, 0);
0072     set_vertex_distribution_width(beamXsigma, beamYsigma, beamZsigma, 0);
0073   }
0074 
0075   void save_event_weight(const bool b) { m_SaveEventWeightFlag = b; }
0076   void save_integrated_luminosity(const bool b) { m_SaveIntegratedLuminosityFlag = b; }
0077 
0078  private:
0079   int read_config(const std::string &cfg_file);
0080   int create_node_tree(PHCompositeNode *topNode) final;
0081   double percent_diff(const double a, const double b) { return fabs((a - b) / a); }
0082   int m_EventCount = 0;
0083 
0084   // event selection
0085   std::vector<PHPy8GenTrigger *> m_RegisteredTriggers;
0086   bool m_TriggersOR{true};
0087   bool m_TriggersAND{false};
0088 
0089   // PYTHIA
0090   std::unique_ptr<Pythia8::Pythia> m_Pythia8;
0091 
0092   std::string m_ConfigFileName{"phpythia8.cfg"};
0093   std::vector<std::string> m_Commands;
0094 
0095   // HepMC
0096   std::unique_ptr<HepMC::Pythia8ToHepMC> m_Pythia8ToHepMC;
0097 
0098   //! whether to store the overall event weight into the HepMC weights
0099   bool m_SaveEventWeightFlag{true};
0100 
0101   //! whether to store the integrated luminosity and other event statistics to the TOP/RUN/PHGenIntegral node
0102   bool m_SaveIntegratedLuminosityFlag{true};
0103 
0104   //! pointer to data node saving the integrated luminosity
0105   PHGenIntegral *m_IntegralNode{};
0106 };
0107 
0108 #endif /* PHPYTHIA8_PHPYTHIA8_H */