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