FairShip
Loading...
Searching...
No Matches
HNLPythia8Generator.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_HNLPYTHIA8GENERATOR_H_
6#define SHIPGEN_HNLPYTHIA8GENERATOR_H_
7
8#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
9#include "Generator.h"
10#include "Pythia8/Pythia.h"
11#include "TFile.h"
12#include "TROOT.h"
13#include "TRandom1.h"
14#include "TRandom3.h"
15#include "TString.h"
16#include "TTree.h"
17
18class FairPrimaryGenerator;
19
20class PyTr1Rng : public Pythia8::RndmEngine {
21 public:
22 PyTr1Rng() { rng = new TRandom1(gRandom->GetSeed()); };
23 virtual ~PyTr1Rng() {};
24
25 Double_t flat() { return rng->Rndm(); };
26
27 private:
28 TRandom1* rng;
29};
30
31class PyTr3Rng : public Pythia8::RndmEngine {
32 public:
33 PyTr3Rng() { rng = new TRandom3(gRandom->GetSeed()); };
34 virtual ~PyTr3Rng() {};
35
36 Double_t flat() { return rng->Rndm(); };
37
38 private:
39 TRandom3* rng;
40};
41
43 public:
46
48 ~HNLPythia8Generator() override;
49
51 Bool_t ReadEvent(FairPrimaryGenerator*) override;
52 void SetParameters(char*);
53 void Print() { fPythia->settings.listAll(); };
54 void List(int id) { fPythia->particleData.list(id); };
55
57
58 Bool_t Init(const char* inFile) override { return Init(inFile, 0); };
59
60 Bool_t Init(const char* inFile, int startEvent) override {
61 LOG(warning) << "Init with files not implemented for HNLPythia8Generator. "
62 "Using default Init() instead";
63 return Init();
64 };
65 Bool_t Init() override;
66
67 void SetMom(Double_t mom) { fMom = mom; };
68 void SetId(Double_t id) { fId = id; };
69 void SetHNLId(Int_t id) { fHNL = id; };
70 void SetLmin(Double_t z) { fLmin = z * 10; }; // Conversion from cm to mm
71 void SetLmax(Double_t z) { fLmax = z * 10; }; // Conversion from cm to mm
72 void SetSmearBeam(Double_t sb) {
73 fsmearBeam = sb * 10;
74 }; // Conversion from cm to mm
75 void SetPaintRadius(Double_t r) {
76 fPaintBeam = r * 10;
77 }; // Conversion from cm to mm
78 void SetfFDs(Double_t z) { fFDs = z; };
79 void UseRandom1() {
80 fUseRandom1 = kTRUE;
81 fUseRandom3 = kFALSE;
82 };
83 void UseRandom3() {
84 fUseRandom1 = kFALSE;
85 fUseRandom3 = kTRUE;
86 };
87
88 void UseDeepCopy() { fDeepCopy = kTRUE; };
89 Int_t nrOfRetries() { return fnRetries; };
90 Pythia8::Pythia* getPythiaInstance() { return fPythia; };
91 Pythia8::Pythia* fPythia;
92 private:
93 std::shared_ptr<Pythia8::RndmEngine> fRandomEngine;
94
95 protected:
96 Double_t fMom; // proton momentum
97 Int_t fHNL; // HNL 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 Double_t fLmin; // m minimum decay position z
102 Double_t fLmax; // m maximum decay position z
103 Int_t fnRetries; // number of events without any HNL
104 Double_t fctau; // hnl lifetime
105 Double_t fFDs; // correction for Pythia6 to match measured Ds production
106 Double_t fsmearBeam; // finite beam size
107 Double_t fPaintBeam; // beam painting radius
108 TFile* fInputFile;
109 TTree* fTree;
111 Float_t hpx[1], hpy[1], hpz[1], hE[1], hM[1], mpx[1], mpy[1], mpz[1], mE[1],
112 hid[1], mid[1];
113 Bool_t fDeepCopy; // not used
114 FairLogger* fLogger;
115};
116
117#endif // SHIPGEN_HNLPYTHIA8GENERATOR_H_
void SetSmearBeam(Double_t sb)
void SetId(Double_t id)
void SetLmin(Double_t z)
TTree * fTree
pointer to a file
Bool_t Init(const char *inFile, int startEvent) override
void SetLmax(Double_t z)
void SetPaintRadius(Double_t r)
Bool_t ReadEvent(FairPrimaryGenerator *) override
void SetMom(Double_t mom)
Bool_t Init() override
std::shared_ptr< Pythia8::RndmEngine > fRandomEngine
Pythia8::Pythia * getPythiaInstance()
Pythia8::Pythia * fPythia
void SetfFDs(Double_t z)
Bool_t Init(const char *inFile) override
virtual ~PyTr1Rng()
TRandom1 * rng
Double_t flat()
virtual ~PyTr3Rng()
TRandom3 * rng
Double_t flat()
virtual Bool_t Init(const char *, int)=0