Class Reference for E1039 Core & Analysis Software
SQPrimaryParticleGen.h
Go to the documentation of this file.
1 /*====================================================================
2 Author: Abinash Pun, Kun Liu
3 Nov, 2019
4 Goal: Import the physics generator of E906 experiment (DPSimPrimaryGeneratorAction)
5 from Kun to E1039 experiment in Fun4All framework
6 =========================================================================*/
7 #ifndef __SQPrimaryParticleGen_H__
8 #define __SQPrimaryParticleGen_H__
9 
10 #ifndef __CINT__
11 #include <Pythia8/Pythia.h>
12 #endif
13 #include <LHAPDF/LHAPDF.h>
14 #ifndef __CINT__
15 #include <gsl/gsl_rng.h>
16 #endif
17 
18 #include <string>
19 
20 #include <TGenPhaseSpace.h>
22 
23 class PHCompositeNode;
25 class PHG4InEvent;
26 class PHG4Particle;
27 
28 class PHGenIntegral;
29 class SQEvent;
30 class SQMCEvent;
31 class SQDimuon;
32 class SQDimuonVector;
33 class SQPrimaryVertexGen;
34 
35 //==========
37 
49 {
50 public:
51  SQPrimaryParticleGen(const std::string& name = "PrimaryGen");
53 
54 
55  int Init(PHCompositeNode* topNode);
56  int InitRun(PHCompositeNode* topNode);
57  int End(PHCompositeNode* topNode);
59  //void GeneratePrimaries(G4Event* anEvent);
60 
62 
63  int generateDrellYan(const TVector3& vtx, const double pARatio, double luminosity);
64  int generateJPsi(const TVector3& vtx, const double pARatio, double luminosity);
65  int generatePsip(const TVector3& vtx, const double pARatio, double luminosity);
66  // void generateDarkPhotonFromEta();
67  int generatePythia(const TVector3& vtx, const double pARatio);
68  //int generateCustomDimuon(PHCompositeNode *topNode, TVector3 vtx, const double pARatio);
69 
71 
73  bool generateDimuon(double mass, double xF);
74 
75  //swith for the generators; Abi
76  //@
77  void enablePythia(){_PythiaGen = true;}
78  void enableCustomDimuon(){_CustomDimuon = true;}
79  void enableDrellYanGen(){_DrellYanGen = true;}
80  void enableJPsiGen(){_JPsiGen = true;}
81  void enablePsipGen(){_PsipGen = true;}
82 
83  void set_pdfset(const std::string name) { _pdfset = name; }
84 
85  void set_pT0DY (const double val) { _pT0DY = val; }
86  void set_pTpowDY (const double val) { _pTpowDY = val; }
87  void set_pT0JPsi (const double val) { _pT0JPsi = val; }
88  void set_pTpowJPsi(const double val) { _pTpowJPsi = val; }
89 
91  void set_config_file(const char *cfg_file)
92  {
93  if (cfg_file) _configFile = cfg_file;
94  }
95 
96  void set_xfRange(const double xmin, const double xmax){
97  xfMin = xmin;
98  xfMax = xmax;
99  }
100  void set_massRange(const double mmin, const double mmax){
101  massMin = mmin;
102  massMax = mmax;
103  }
104  //@
105 
106  double CrossSectionDrellYan(const double mass, const double xF, const double pARatio);
107  double CrossSectionDrellYan(const double mass, const double xF, const double x1, const double x2, const double pARatio);
108  double CrossSectionJPsi(const double xF);
109  double CrossSectionPsip(const double xF);
110 
111  private:
112  bool _PythiaGen;
113  bool _CustomDimuon;
114  bool _DrellYanGen;
115  bool _JPsiGen;
116  bool _PsipGen;
117 
118  SQPrimaryVertexGen* _vertexGen;
119  PHG4InEvent *ineve;
120  SQEvent* _evt; //< An output node
121  SQMCEvent* _mcevt; //< An output node
122  SQDimuonVector* _vec_dim; //< An output node
123  PHGenIntegral *_integral_node; //< An output node
124 
125  SQDimuon* _dim_gen; //< To hold the kinematics of a dimuon generated
126 
127  //Pythia generator
128  Pythia8::Pythia ppGen;
129  Pythia8::Pythia pnGen;
130  Pythia8::Pythia _Pythia;
132  std::string _configFile;
133  int read_config(const char *cfg_file = 0);
134 
136  TGenPhaseSpace phaseGen;
137 
139  std::string _pdfset;
140  LHAPDF::PDF* _pdf;
141 
142  // Parameters (being moved from DPGEN)
143  double _pT0DY;
144  double _pTpowDY;
145  double _pT0JPsi;
146  double _pTpowJPsi;
147 
148  //some initializations
149  double massMin = 0.22;
150  double massMax = 10.;
151  double x1Min = 0.;
152  double x1Max = 1.;
153  double x2Min = 0.;
154  double x2Max = 1.;
155  double xfMin = -1.;
156  double xfMax = 1.;
157  double cosThetaMin = -1.;
158  double cosThetaMax = 1. ;
159  double zOffsetMin = -1.;
160  double zOffsetMax = 1.;
161 
162  void InsertMuonPair(const TVector3& vtx);
163  void InsertEventInfo(const double xsec, const double weight, const TVector3& vtx);
164 
165  double _n_gen_acc_evt; //< N of generator-accepted events
166  double _n_proc_evt; //< N of processed events
167  double _weight_sum; //< Sum of weights
168  double _inte_lumi; //< Integrated luminosity
169 };
170 
171 //========
172 
173 
174 
175 
176 #endif
177 
PHGenIntegral.
Definition: PHGenIntegral.h:21
An SQ interface class to hold a list of SQDimuon objects.
An SQ interface class to hold one true or reconstructed dimuon.
Definition: SQDimuon.h:8
An SQ interface class to hold one event header.
Definition: SQEvent.h:17
An SQ interface class to hold one simulated-event header.
Definition: SQMCEvent.h:12
Physics generator imported from E906 software.
void set_pT0DY(const double val)
int InitRun(PHCompositeNode *topNode)
int process_event(PHCompositeNode *topNode)
double CrossSectionDrellYan(const double mass, const double xF, const double x1, const double x2, const double pARatio)
void set_xfRange(const double xmin, const double xmax)
void set_massRange(const double mmin, const double mmax)
void set_pdfset(const std::string name)
bool generateDimuon(double mass, double xF)
Dimuon phase space generator.
void set_pTpowDY(const double val)
int generateJPsi(const TVector3 &vtx, const double pARatio, double luminosity)
SQPrimaryParticleGen(const std::string &name="PrimaryGen")
void set_pTpowJPsi(const double val)
void set_pT0JPsi(const double val)
int generateDrellYan(const TVector3 &vtx, const double pARatio, double luminosity)
Various generators.
int Init(PHCompositeNode *topNode)
int generatePsip(const TVector3 &vtx, const double pARatio, double luminosity)
void set_config_file(const char *cfg_file)
config file for pythia
int generatePythia(const TVector3 &vtx, const double pARatio)
virtual ~SQPrimaryParticleGen()
double CrossSectionDrellYan(const double mass, const double xF, const double pARatio)
int End(PHCompositeNode *topNode)
Called at the end of all processing.
double CrossSectionJPsi(const double xF)
double CrossSectionPsip(const double xF)
Class to generate the event vertex, based on the beam profile and the target+spectrometer materials g...