@@ -168,6 +168,149 @@ class O2_GeneratorParamPsiMidY : public GeneratorTGenerator
168168 GeneratorParam * paramPsi = nullptr ;
169169};
170170
171+ class O2_GeneratorParamJpsiMidY_5TeV : public GeneratorTGenerator
172+ {
173+
174+ public :
175+ O2_GeneratorParamJpsiMidY_5TeV () : GeneratorTGenerator ("paramJpsi" )
176+ {
177+ paramJpsi = new GeneratorParam (1 , -1 , PtJPsipp5TeV , YJPsipp5TeV , V2JPsipp5TeV , IpJPsipp5TeV );
178+ paramJpsi -> SetMomentumRange (0. , 1.e6 );
179+ paramJpsi -> SetPtRange (0. , 1000. );
180+ paramJpsi -> SetYRange (-1.0 , 1.0 );
181+ paramJpsi -> SetPhiRange (0. , 360. );
182+ paramJpsi -> SetDecayer (new TPythia6Decayer ()); // Pythia
183+ paramJpsi -> SetForceDecay (kNoDecay ); // particle left undecayed
184+ setTGenerator (paramJpsi );
185+ };
186+
187+ ~O2_GeneratorParamJpsiMidY_5TeV ()
188+ {
189+ delete paramJpsi ;
190+ };
191+
192+ Bool_t Init () override
193+ {
194+ GeneratorTGenerator ::Init ();
195+ paramJpsi -> Init ();
196+ return true;
197+ }
198+
199+ void SetNSignalPerEvent (Int_t nsig ) { paramJpsi -> SetNumberParticles (nsig ); }
200+
201+ //-------------------------------------------------------------------------//
202+ static Double_t PtJPsipp5TeV (const Double_t * px , const Double_t * /*dummy*/ )
203+ {
204+ // JPSi pt at 5.02 TeV: https: // www.hepdata.net/record/ins1735351
205+ //
206+ const Double_t kC = 1774.9 ;
207+ const Double_t kpt0 = 3.38452 ;
208+ const Double_t kn = 2.77889 ;
209+ Double_t pt = px [0 ];
210+
211+ return kC * pt / TMath ::Power ((1. + (pt / kpt0 ) * (pt / kpt0 )), kn );
212+ }
213+
214+ //-------------------------------------------------------------------------//
215+ static Double_t YJPsipp5TeV (const Double_t * py , const Double_t * /*dummy*/ )
216+ {
217+ // Taken the same as: jpsi y in pp at 13 TeV, tuned on data, prompt jpsi ALICE+LHCb, 13 TeV
218+ Double_t y = * py ;
219+ Float_t p0 , p1 , p2 ;
220+ p0 = 7.79382e+00 ;
221+ p1 = 2.87827e-06 ;
222+ p2 = 4.41847e+00 ;
223+ return p0 * TMath ::Exp (- (1. / 2. ) * TMath ::Power (((y - p1 ) / p2 ), 2 ));
224+ }
225+
226+ //-------------------------------------------------------------------------//
227+ static Double_t V2JPsipp5TeV (const Double_t * /*dummy*/ , const Double_t * /*dummy*/ )
228+ {
229+ // jpsi v2
230+ return 0. ;
231+ }
232+
233+ //-------------------------------------------------------------------------//
234+ static Int_t IpJPsipp5TeV (TRandom * )
235+ {
236+ return 443 ;
237+ }
238+
239+ private :
240+ GeneratorParam * paramJpsi = nullptr ;
241+ };
242+
243+ class O2_GeneratorParamPsiMidY_5TeV : public GeneratorTGenerator
244+ {
245+
246+ public :
247+ O2_GeneratorParamPsiMidY_5TeV () : GeneratorTGenerator ("ParamPsi" )
248+ {
249+ paramPsi = new GeneratorParam (1 , -1 , PtPsipp5TeV , YPsipp5TeV , V2Psipp5TeV , IpPsipp5TeV );
250+ paramPsi -> SetMomentumRange (0. , 1.e6 ); // Momentum range added from me
251+ paramPsi -> SetPtRange (0. , 1000. ); // transverse of momentum range
252+ paramPsi -> SetYRange (-1.0 , 1.0 ); // rapidity range
253+ paramPsi -> SetPhiRange (0. , 360. ); // phi range
254+ paramPsi -> SetDecayer (new TPythia6Decayer ()); // Pythia decayer
255+ paramPsi -> SetForceDecay (kNoDecay ); // particle left undecayed
256+ setTGenerator (paramPsi ); // Setting parameters to ParamPsi for Psi(2S)
257+ };
258+
259+ ~O2_GeneratorParamPsiMidY_5TeV ()
260+ {
261+ delete paramPsi ;
262+ };
263+
264+ Bool_t Init () override
265+ {
266+ GeneratorTGenerator ::Init ();
267+ paramPsi -> Init ();
268+ return true;
269+ }
270+ void SetNSignalPerEvent (Int_t nsig ) { paramPsi -> SetNumberParticles (nsig ); }
271+
272+ //-------------------------------------------------------------------------//
273+ static Double_t PtPsipp5TeV (const Double_t * px , const Double_t * /*dummy*/ )
274+ {
275+ // Same as JPsi at 5.02 TeV since ratio is almost flat in pT: https: // www.hepdata.net/record/ins1735351
276+ //
277+ const Double_t kC = 1774.9 ;
278+ const Double_t kpt0 = 3.38452 ;
279+ const Double_t kn = 2.77889 ;
280+ Double_t pt = px [0 ];
281+
282+ return kC * pt / TMath ::Power ((1. + (pt / kpt0 ) * (pt / kpt0 )), kn );
283+ }
284+
285+ //-------------------------------------------------------------------------//
286+ static Double_t YPsipp5TeV (const Double_t * py , const Double_t * /*dummy*/ )
287+ {
288+ // Taken same as jpsi y in pp at 13 TeV, tuned on data, prompt jpsi ALICE+LHCb, 13 TeV
289+ Double_t y = * py ;
290+ Float_t p0 , p1 , p2 ;
291+ p0 = 7.79382e+00 ;
292+ p1 = 2.87827e-06 ;
293+ p2 = 4.41847e+00 ;
294+ return p0 * TMath ::Exp (- (1. / 2. ) * TMath ::Power (((y - p1 ) / p2 ), 2 ));
295+ }
296+
297+ //-------------------------------------------------------------------------//
298+ static Double_t V2Psipp5TeV (const Double_t * /*dummy*/ , const Double_t * /*dummy*/ )
299+ {
300+ // jpsi v2
301+ return 0. ;
302+ }
303+
304+ //-------------------------------------------------------------------------//
305+ static Int_t IpPsipp5TeV (TRandom * )
306+ {
307+ return 100443 ;
308+ }
309+
310+ private :
311+ GeneratorParam * paramPsi = nullptr ;
312+ };
313+
171314class O2_GeneratorParamJpsiFwdY : public GeneratorTGenerator
172315{
173316
@@ -982,6 +1125,31 @@ FairGenerator* GeneratorCocktailPromptCharmoniaToMuonEvtGen_pp5TeV()
9821125 return genCocktailEvtGen ;
9831126}
9841127
1128+ FairGenerator * GeneratorCocktailPromptCharmoniaToElectronEvtGen_pp5TeV ()
1129+ {
1130+
1131+ auto genCocktailEvtGen = new o2 ::eventgen ::GeneratorEvtGen < GeneratorCocktail > ();
1132+
1133+ auto genJpsi = new o2 ::eventgen ::O2_GeneratorParamJpsiMidY_5TeV ;
1134+ genJpsi -> SetNSignalPerEvent (1 ); // 1 J/psi generated per event by GeneratorParam
1135+ auto genPsi = new o2 ::eventgen ::O2_GeneratorParamPsiMidY_5TeV ;
1136+ genPsi -> SetNSignalPerEvent (1 ); // 1 Psi(2S) generated per event by GeneratorParam
1137+ genCocktailEvtGen -> AddGenerator (genJpsi , 1 ); // add J/psi generator
1138+ genCocktailEvtGen -> AddGenerator (genPsi , 1 ); // add Psi(2S) generator
1139+
1140+ TString pdgs = "443;100443" ;
1141+ std ::string spdg ;
1142+ TObjArray * obj = pdgs .Tokenize (";" );
1143+ genCocktailEvtGen -> SetSizePdg (obj -> GetEntriesFast ());
1144+ for (int i = 0 ; i < obj -> GetEntriesFast (); i ++ ) {
1145+ spdg = obj -> At (i )-> GetName ();
1146+ genCocktailEvtGen -> AddPdg (std ::stoi (spdg ), i );
1147+ printf ("PDG %d \n" , std ::stoi (spdg ));
1148+ }
1149+ genCocktailEvtGen -> SetForceDecay (kEvtDiElectron );
1150+
1151+ return genCocktailEvtGen ;
1152+ }
9851153
9861154FairGenerator *
9871155 GeneratorParamPromptPsiToJpsiPiPiEvtGen_pp13TeV (TString pdgs = "100443" )
0 commit comments