Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:19:30

0001 // $Id: $
0002 
0003 /*!
0004  * \file PHHepMCGenHelper.h
0005  * \brief
0006  * \author Jin Huang <jhuang@bnl.gov>
0007  * \version $Revision:   $
0008  * \date $Date: $
0009  */
0010 
0011 #ifndef PHHEPMC_PHHEPMCGENHELPER_H
0012 #define PHHEPMC_PHHEPMCGENHELPER_H
0013 
0014 #include <CLHEP/Vector/ThreeVector.h>
0015 
0016 #include <gsl/gsl_rng.h>
0017 
0018 #include <cmath>
0019 #include <limits>
0020 #include <string>
0021 #include <utility>
0022 #include <vector>
0023 
0024 class PHCompositeNode;
0025 class PHHepMCGenEvent;
0026 class PHHepMCGenEventMap;
0027 
0028 namespace HepMC
0029 {
0030   class GenEvent;
0031 }
0032 
0033 /*!
0034  * \brief PHHepMCGenHelper provides service of DST upload of HepMC subevent, vertex assignment and random generator
0035  */
0036 class PHHepMCGenHelper
0037 {
0038  public:
0039   PHHepMCGenHelper();
0040   virtual ~PHHepMCGenHelper();
0041 
0042   //! supported function distributions
0043   enum VTXFUNC
0044   {
0045     //! uniform distribution with half width set via set_vertex_distribution_width()
0046     Uniform,
0047     //! normal distribution with sigma width set via set_vertex_distribution_width()
0048     Gaus
0049   };
0050 
0051   //! toss a new vertex according to a Uniform or Gaus distribution
0052   void set_vertex_distribution_function(VTXFUNC x, VTXFUNC y, VTXFUNC z, VTXFUNC t);
0053 
0054   //! set the mean value of the vertex distribution, use PHENIX units of cm, ns
0055   void set_vertex_distribution_mean(const double x, const double y, const double z, const double t);
0056 
0057   //! set the width of the vertex distribution function about the mean, use PHENIX units of cm, ns
0058   void set_vertex_distribution_width(const double x, const double y, const double z, const double t);
0059 
0060   //! embedding ID for the event
0061   //! positive ID is the embedded event of interest, e.g. jetty event from pythia
0062   //! negative IDs are backgrounds, .e.g out of time pile up collisions
0063   //! Usually, ID = 0 means the primary Au+Au collision background
0064   int get_embedding_id() const { return _embedding_id; }
0065   //
0066   //! embedding ID for the event
0067   //! positive ID is the embedded event of interest, e.g. jetty event from pythia
0068   //! negative IDs are backgrounds, .e.g out of time pile up collisions
0069   //! Usually, ID = 0 means the primary Au+Au collision background
0070   void set_embedding_id(int id) { _embedding_id = id; }
0071   //
0072   //! reuse vertex from another PHHepMCGenEvent with embedding_id = src_embedding_id Additional smearing and shift possible with set_vertex_distribution_*()
0073   void set_reuse_vertex(int src_embedding_id)
0074   {
0075     _reuse_vertex = true;
0076     _reuse_vertex_embedding_id = src_embedding_id;
0077   }
0078 
0079   //
0080   //! init interface nodes
0081   virtual int create_node_tree(PHCompositeNode *topNode);
0082 
0083   //! choice of reference version of the PHHepMCGenEvent
0084   static const PHHepMCGenEvent *get_PHHepMCGenEvent_template();
0085 
0086   //! send HepMC::GenEvent to DST tree. This function takes ownership of evt
0087   PHHepMCGenEvent *insert_event(HepMC::GenEvent *evt);
0088 
0089   const PHHepMCGenEventMap *get_geneventmap() const
0090   {
0091     return _geneventmap;
0092   }
0093 
0094   PHHepMCGenEventMap *get_geneventmap()
0095   {
0096     return _geneventmap;
0097   }
0098 
0099   gsl_rng *get_random_generator()
0100   {
0101     return RandomGenerator;
0102   }
0103 
0104   void set_geneventmap(PHHepMCGenEventMap *geneventmap)
0105   {
0106     _geneventmap = geneventmap;
0107   }
0108 
0109   //! Beam angle in lab polar coordinate.
0110   //! @param[in] beamA_theta beamA, in pair of Theta-Phi. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0111   //! @param[in] beamA_phi beamA, in pair of Theta-Phi. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0112   //! @param[in] beamB_theta  beamB, in pair of Theta-Phi. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0113   //! @param[in] beamB_phi  beamB, in pair of Theta-Phi. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0114   void set_beam_direction_theta_phi(
0115       const double beamA_theta,
0116       const double beamA_phi,
0117       const double beamB_theta,
0118       const double beamB_phi)
0119   {
0120     m_beam_direction_theta_phi = {{beamA_theta, beamA_phi}, {beamB_theta, beamB_phi}};
0121   }
0122 
0123   //! Beam angle divergence in accelerator beam coordinate.
0124   //! @param[in] beamA_divergence_h beamA, horizontal divergence Gaussian Sigma. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0125   //! @param[in] beamA_divergence_v beamA, vertical divergence Gaussian Sigma. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0126   //! @param[in] beamB_divergence_h beamB, horizontal divergence Gaussian Sigma. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0127   //! @param[in] beamB_divergence_v beamB, vertical divergence Gaussian Sigma. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0128   void set_beam_angular_divergence_hv(
0129       const double beamA_divergence_h,
0130       const double beamA_divergence_v,
0131       const double beamB_divergence_h,
0132       const double beamB_divergence_v)
0133   {
0134     m_beam_angular_divergence_hv = {{beamA_divergence_h, beamA_divergence_v}, {beamB_divergence_h, beamB_divergence_v}};
0135   }
0136 
0137   //! Central beam angle shift as linear function of longitudinal vertex position, d_shift/dz,
0138   //! which is used to represent leading order effect of crab cavity momentum kick on the beam bunch
0139   //! @param[in] beamA_h beamA, horizontal angle dh/dz. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0140   //! @param[in] beamA_v beamA, vertical angle dv/dz. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0141   //! @param[in] beamB_h beamB, horizontal angle dh/dz. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0142   //! @param[in] beamB_v beamB, vertical angle dv/dz. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0143   void set_beam_angular_z_coefficient_hv(
0144       const double beamA_h,
0145       const double beamA_v,
0146       const double beamB_h,
0147       const double beamB_v)
0148   {
0149     m_beam_angular_z_coefficient_hv = {{beamA_h, beamA_v}, {beamB_h, beamB_v}};
0150   }
0151 
0152   //! simulate bunch interaction instead of applying vertex distributions
0153   void use_beam_bunch_sim(bool b) { m_use_beam_bunch_sim = b; }
0154 
0155   //! Beam bunch geometry as 3D Gauss width
0156   //! First element is beamA, in vector of Gaussian Sigma H,V,Longitudinal
0157   //! Second element is beamB, in vector of Gaussian Sigma H,V,Longitudinal
0158   void set_beam_bunch_width(const std::vector<double> &beamA, const std::vector<double> &beamB);
0159 
0160   void CopySettings(PHHepMCGenHelper &helper);
0161 
0162   //! copy setting to helper_dest
0163   void CopySettings(PHHepMCGenHelper *helper_dest);
0164 
0165   //! copy setting from helper_src
0166   void CopyHelperSettings(PHHepMCGenHelper *helper_src);
0167 
0168   void Print(const std::string &what = "ALL") const;
0169 
0170   void PHHepMCGenHelper_Verbosity(int v) { m_verbosity = v; }
0171 
0172   int PHHepMCGenHelper_Verbosity() { return m_verbosity; }
0173 
0174  protected:
0175   //! Record the translation,boost,rotation for HepMC frame to lab frame according to collision settings
0176   void HepMC2Lab_boost_rotation_translation(PHHepMCGenEvent *genevent);
0177 
0178   //! move vertex in translation according to vertex settings
0179   void move_vertex(PHHepMCGenEvent *genevent);
0180 
0181   //! generate vertx with bunch interaction according to
0182   //! https://github.com/eic/documents/blob/d06b5597a0a89dcad215bab50fe3eefa17a097a5/reports/general/Note-Simulations-BeamEffects.pdf
0183   //! \return pair of bunch local z position for beam A and beam B
0184   std::pair<double, double> generate_vertx_with_bunch_interaction(PHHepMCGenEvent *genevent);
0185 
0186  private:
0187   gsl_rng *RandomGenerator{nullptr};
0188 
0189   double smear(const double position, const double width, VTXFUNC dist) const;
0190 
0191   //! function to convert spherical coordinate to Hep3Vector in x-y-z
0192   static CLHEP::Hep3Vector pair2Hep3Vector(const std::pair<double, double> &theta_phi);
0193 
0194   VTXFUNC _vertex_func_x{Gaus};
0195   VTXFUNC _vertex_func_y{Gaus};
0196   VTXFUNC _vertex_func_z{Gaus};
0197   VTXFUNC _vertex_func_t{Gaus};
0198 
0199   double _vertex_x{0.};
0200   double _vertex_y{0.};
0201   double _vertex_z{0.};
0202   double _vertex_t{0.};
0203 
0204   double _vertex_width_x{0.};
0205   double _vertex_width_y{0.};
0206   double _vertex_width_z{0.};
0207   double _vertex_width_t{0.};
0208 
0209   //! Beam angle in lab polar coordinate.
0210   //! First element is beamA, in pair of Theta-Phi. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0211   //! Second element is beamB, in pair of Theta-Phi. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0212   std::pair<std::pair<double, double>, std::pair<double, double>> m_beam_direction_theta_phi = {
0213       {0, 0},    //+z beam
0214       {M_PI, 0}  //-z beam
0215   };
0216 
0217   //! Beam angle divergence in accelerator beam coordinate.
0218   //! First element is beamA, in pair of Gaussian Sigma_H Sigma_V. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0219   //! Second element is beamB, in pair of Gaussian Sigma_H Sigma_V. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0220   std::pair<std::pair<double, double>, std::pair<double, double>> m_beam_angular_divergence_hv = {
0221       {0, 0},  //+z beam
0222       {0, 0}   //-z beam
0223   };
0224 
0225   //! Central beam angle shift as linear function of longitudinal vertex position, d_shift/dz,
0226   //! which is used to represent leading order effect of crab cavity momentum kick on the beam bunch
0227   //! First element is beamA, in pair of dh/dz, dv/dz. BeamA is aimed to +z direction in the HepMC event generator's coordinate
0228   //! Second element is beamB, in pair of dh/dz, dv/dz. BeamA is aimed to -z direction in the HepMC event generator's coordinate
0229   std::pair<std::pair<double, double>, std::pair<double, double>> m_beam_angular_z_coefficient_hv = {
0230       {0, 0},  //+z beam
0231       {0, 0}   //-z beam
0232   };
0233 
0234   //! positive ID is the embedded event of interest, e.g. jetty event from pythia
0235   //! negative IDs are backgrounds, .e.g out of time pile up collisions
0236   //! Usually, ID = 0 means the primary Au+Au collision background
0237   int _embedding_id{0};
0238 
0239   //! whether reuse vertex from another PHHepMCGenEvent
0240   bool _reuse_vertex{false};
0241 
0242   //! if _reuse_vertex, which embedding_id provide the source vertex. Additional smearing and shift possible with set_vertex_distribution_*()
0243   int _reuse_vertex_embedding_id{std::numeric_limits<int>::min()};
0244 
0245   //! pointer to the output container
0246   PHHepMCGenEventMap *_geneventmap{nullptr};
0247 
0248   //! verbosity
0249   int m_verbosity{0};
0250 
0251   //! simulate bunch interaction instead of applying vertex distributions
0252   bool m_use_beam_bunch_sim{false};
0253 
0254   //! Beam bunch geometry as 3D Gauss width
0255   //! First element is beamA, in vector of Gaussian Sigma H,V,Longitudinal
0256   //! Second element is beamB, in vector of Gaussian Sigma H,V,Longitudinal
0257   std::pair<std::vector<double>, std::vector<double>> m_beam_bunch_width = {
0258       {0, 0, 0},  //+z beam
0259       {0, 0, 0}   //-z beam
0260   };
0261 
0262   //! use m_beam_bunch_width to calculate horizontal and vertical collision width
0263   //! \param[in] hv_index 0: horizontal. 1: vertical
0264   double get_collision_width(unsigned int hv_index);
0265 };
0266 
0267 #endif /* PHHEPMC_PHHEPMCGENHELPER_H */