FairShip
Loading...
Searching...
No Matches
DPPythia8Generator.h
Go to the documentation of this file.
1// SPDX-License-Identifier: LGPL-3.0-or-later
2// SPDX-FileCopyrightText: Copyright CERN for the benefit of the SHiP
3// Collaboration
4
5#ifndef SHIPGEN_DPPYTHIA8GENERATOR_H_
6#define SHIPGEN_DPPYTHIA8GENERATOR_H_
7#define DPP8GENERATOR_H
8// Avoid the inclusion of dlfcn.h by Pyhtia.h that CINT is not able to process
9// #ifdef __CINT__
10// #define _DLFCN_H_
11// #define _DLFCN_H
12// #endif
13
14#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
15#include "Generator.h"
16#include "HNLPythia8Generator.h"
17#include "Pythia8/Pythia.h"
18#include "TFile.h"
19#include "TH2F.h"
20#include "TROOT.h"
21#include "TRandom1.h"
22#include "TRandom3.h"
23#include "TString.h"
24#include "TTree.h"
25
26class FairPrimaryGenerator;
27
29 public:
32
34 ~DPPythia8Generator() override;
35
37 Bool_t ReadEvent(FairPrimaryGenerator*) override;
38 void SetParameters(char*);
39 void Print() { fPythia->settings.listAll(); };
40 void List(int id) { fPythia->particleData.list(id); };
41
42 // void SetDecayToHadrons(){
43 // std::cout << " INFO: Adding decay to hadrons." << std::endl;
44 // fHadDecay = true;
45 // };
46
48 Bool_t Init() override;
49 Bool_t Init(const char* inFile) override { return Init(inFile, 0); };
50
51 Bool_t Init(const char* inFile, int startEvent) override {
52 LOG(warning) << "Init with files not implemented for DPPythia8Generator. "
53 "Using default Init() instead";
54 return Init();
55 };
56
57 void SetMom(Double_t mom) { fMom = mom; };
58 Double_t GetMom() { return fMom; };
59 void SetId(Double_t id) { fId = id; };
60 void SetDPId(Int_t id) { fDP = id; };
61 Int_t GetDPId() { return fDP; };
62 void SetLmin(Double_t z) { fLmin = z * 10; };
63 void SetLmax(Double_t z) { fLmax = z * 10; };
64 void SetSmearBeam(Double_t sb) { fsmearBeam = sb; };
65 void SetfFDs(Double_t z) { fFDs = z; };
66 void UseRandom1() {
67 fUseRandom1 = kTRUE;
68 fUseRandom3 = kFALSE;
69 };
70 void UseRandom3() {
71 fUseRandom1 = kFALSE;
72 fUseRandom3 = kTRUE;
73 };
74
75 void SetPbrem(TH2F* pdf) {
76 fpbrem = kTRUE;
77 fpbremPDF = pdf;
78 };
79 Bool_t IsPbrem() { return fpbrem; };
80 void SetDY() { fdy = kTRUE; };
81
82 Double_t MinDPMass() { return fDPminM; };
83 void SetMinDPMass(Double_t m) { fDPminM = m; };
84
85 void UseDeepCopy() { fDeepCopy = kTRUE; };
86 Int_t nrOfRetries() { return fnRetries; };
87 Int_t nrOfDP() { return fnDPtot; };
88 Pythia8::Pythia* getPythiaInstance() { return fPythia; };
89 Pythia8::Pythia* fPythia;
90 // Pythia8::Pythia* fPythiaHadDecay; //!
91 private:
92 std::shared_ptr<Pythia8::RndmEngine> fRandomEngine;
93
94 protected:
95 // Bool_t fHadDecay; //select hadronic decay
96 Double_t fMom; // proton energy
97 Int_t fDP; // DP ID
98 Int_t fId; // target type
99 Bool_t fUseRandom1; // flag to use TRandom1
100 Bool_t fUseRandom3; // flag to use TRandom3 (default)
101 Bool_t
102 fpbrem; // flag to do proton bremstrahlung production (default is false)
103 TH2F* fpbremPDF; // pointer to TH2 containing PDF(p,theta) to have a dark
104 // photon with momentum p and angle theta to be produced by
105 // pbrem.
106 Bool_t fdy; // flag to do Drell-Yan QCD production
107 Double_t fDPminM; // Minimum mass, in GeV, for the DP produced in ffbar to DP
108 // QCD production.
109 Double_t fLmin; // m minimum decay position z
110 Double_t fLmax; // m maximum decay position z
111 Int_t fnRetries; // number of events without any DP
112 Int_t fnDPtot; // total number of DP from multiple mesons in single collision
113 Double_t fctau; // dark photon lifetime
114 Double_t fFDs; // correction for Pythia6 to match measured Ds production
115 Double_t fsmearBeam; // finite beam size
116 TFile* fInputFile;
117 TTree* fTree;
119 Float_t hpx[1], hpy[1], hpz[1], hE[1], hM[1], mpx[1], mpy[1], mpz[1], mE[1],
120 hid[1], mid[1];
121 Bool_t fDeepCopy; // not used
122 FairLogger* fLogger;
123};
124
125#endif // SHIPGEN_DPPYTHIA8GENERATOR_H_
Double_t m
void SetPbrem(TH2F *pdf)
Bool_t Init(const char *inFile) override
void SetSmearBeam(Double_t sb)
void SetMom(Double_t mom)
void SetDPId(Int_t id)
Pythia8::Pythia * fPythia
Bool_t ReadEvent(FairPrimaryGenerator *) override
Bool_t Init(const char *inFile, int startEvent) override
void SetLmax(Double_t z)
Pythia8::Pythia * getPythiaInstance()
std::shared_ptr< Pythia8::RndmEngine > fRandomEngine
void SetMinDPMass(Double_t m)
void SetId(Double_t id)
void SetLmin(Double_t z)
Bool_t Init() override
TTree * fTree
pointer to a file
void SetfFDs(Double_t z)
virtual Bool_t Init(const char *, int)=0