SND@LHC Software
Loading...
Searching...
No Matches
HNLPythia8Generator.h
Go to the documentation of this file.
1#ifndef PNDH8GENERATOR_H
2#define PNDH8GENERATOR_H 1
3
4#include "TROOT.h"
5#include "TString.h"
6#include "TFile.h"
7#include "TTree.h"
8#include "FairGenerator.h"
9#include "Pythia8/Pythia.h"
10#include "TRandom1.h"
11#include "TRandom3.h"
12#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
13
14class FairPrimaryGenerator;
15
16class PyTr1Rng : public Pythia8::RndmEngine
17{
18 public:
19 PyTr1Rng() { rng = new TRandom1(gRandom->GetSeed()); };
20 virtual ~PyTr1Rng() {};
21
22 Double_t flat() { return rng->Rndm(); };
23
24 private:
25 TRandom1 *rng;
26};
27
28class PyTr3Rng : public Pythia8::RndmEngine
29{
30 public:
31 PyTr3Rng() { rng = new TRandom3(gRandom->GetSeed()); };
32 virtual ~PyTr3Rng() {};
33
34 Double_t flat() { return rng->Rndm(); };
35
36 private:
37 TRandom3 *rng;
38};
39
40
41
42class HNLPythia8Generator : public FairGenerator
43{
44 public:
45
48
50 virtual ~HNLPythia8Generator();
51
53 Bool_t ReadEvent(FairPrimaryGenerator*);
54 void SetParameters(char*);
55 void Print(){fPythia->settings.listAll(); };
56 void List(int id){fPythia->particleData.list(id);};
57
58 virtual Bool_t Init();
59
60 void SetMom(Double_t mom) { fMom = mom; };
61 void SetId(Double_t id) { fId = id; };
62 void SetHNLId(Int_t id) { fHNL = id; };
63 void SetLmin(Double_t z) { fLmin = z*10; };
64 void SetLmax(Double_t z) { fLmax = z*10; };
65 void SetSmearBeam(Double_t sb) { fsmearBeam = sb; };
66 void SetfFDs(Double_t z) { fFDs = z; };
67 void UseRandom1() { fUseRandom1 = kTRUE; fUseRandom3 = kFALSE; };
68 void UseRandom3() { fUseRandom1 = kFALSE; fUseRandom3 = kTRUE; };
69 void UseExternalFile(const char* x, Int_t i){ fextFile = x; firstEvent=i; };
70 void UseDeepCopy(){ fDeepCopy = kTRUE; };
71 Int_t nrOfRetries(){ return fnRetries; };
72 Pythia8::Pythia* getPythiaInstance(){return fPythia;};
73 Pythia8::Pythia* fPythia;
74 private:
75
76 Pythia8::RndmEngine* fRandomEngine;
77
78 protected:
79
80 Double_t fMom; // proton momentum
81 Int_t fHNL; // HNL ID
82 Int_t fId; // target type
83 Bool_t fUseRandom1; // flag to use TRandom1
84 Bool_t fUseRandom3; // flag to use TRandom3 (default)
85 Double_t fLmin; // m minimum decay position z
86 Double_t fLmax; // m maximum decay position z
87 Int_t fnRetries; // number of events without any HNL
88 Double_t fctau; // hnl lifetime
89 Double_t fFDs; // correction for Pythia6 to match measured Ds production
90 Double_t fsmearBeam; // finite beam size
91 const char* fextFile; // read charm and beauty hadrons from external file, decay with Pythia
92 TFile* fInputFile;
93 TTree* fTree;
95 Float_t hpx[1], hpy[1], hpz[1], hE[1],hM[1],mpx[1], mpy[1], mpz[1], mE[1],hid[1], mid[1];
96 Bool_t fDeepCopy; // not used
97 FairLogger* fLogger;
98
100};
101
102#endif /* !PNDH8GENERATOR_H */
ClassDef(HNLPythia8Generator, 6)
don't make it persistent, magic ROOT command
Pythia8::RndmEngine * fRandomEngine
void SetSmearBeam(Double_t sb)
void SetId(Double_t id)
void SetLmin(Double_t z)
TTree * fTree
pointer to a file
void UseExternalFile(const char *x, Int_t i)
void SetLmax(Double_t z)
void SetMom(Double_t mom)
Bool_t ReadEvent(FairPrimaryGenerator *)
Pythia8::Pythia * getPythiaInstance()
Pythia8::Pythia * fPythia
void SetfFDs(Double_t z)
virtual ~PyTr1Rng()
TRandom1 * rng
Double_t flat()
virtual ~PyTr3Rng()
TRandom3 * rng
Double_t flat()