FairShip
Loading...
Searching...
No Matches
Pythia8Generator.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_PYTHIA8GENERATOR_H_
6#define SHIPGEN_PYTHIA8GENERATOR_H_
7
8#include <string>
9#include <vector>
10
11#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
12#include "Generator.h"
13#include "GenieGenerator.h"
14#include "Pythia8/Pythia.h"
15#include "TROOT.h"
16#include "TTree.h"
17
18class FairPrimaryGenerator;
19
21 public:
24
26 ~Pythia8Generator() override;
27
29 Bool_t ReadEvent(FairPrimaryGenerator*) override;
30 void SetParameters(char*);
31 void Print();
32
34 Bool_t Init() override;
35 Bool_t Init(const char* inFile) override { return Init(inFile, 0); };
36
37 Bool_t Init(const char* inFile, int startEvent) override {
38 LOG(warning) << "Init with files not implemented for Pythia8Generator. "
39 "Using default Init() instead";
40 return Init();
41 };
42
43 void SetMom(Double_t mom) { fMom = mom; };
44 void SetId(Double_t id) { fId = id; };
45 void UseRandom1() {
46 fUseRandom1 = kTRUE;
47 fUseRandom3 = kFALSE;
48 };
49 void UseRandom3() {
50 fUseRandom1 = kFALSE;
51 fUseRandom3 = kTRUE;
52 };
53
54 void SetfFDs(Double_t z) { fFDs = z; };
55 void SetTarget(TString s, Double_t x, Double_t y) {
56 targetName = s;
57 xOff = x;
58 yOff = y;
59 };
60 void SetTargetCoordinates(Double_t start_z, Double_t end_z) {
61 startZ = start_z;
62 endZ = end_z;
63 targetFromGeometry = true;
64 };
65 Int_t nrOfRetries() { return fnRetries; };
66
67 private:
68 std::shared_ptr<Pythia8::RndmEngine> fRandomEngine;
69
70 protected:
71 Double_t fMom; // proton momentum
72 Int_t fId; // target type
73 Bool_t fUseRandom1; // flag to use TRandom1
74 Bool_t fUseRandom3; // flag to use TRandom3 (default)
75 Float_t hpx[1], hpy[1], hpz[1], hE[1], hM[1], mpx[1], mpy[1], mpz[1], mE[1],
76 hid[1], mid[1], ck[1];
77 Float_t ancestors[16], subprocCodes[16];
79 TFile* fInputFile;
80 FairLogger* fLogger;
81 TTree* fTree;
82 Pythia8::Pythia* fPythia;
83 Double_t fFDs; // correction for Pythia6 to match measured Ds production
84 Int_t fnRetries; //
86 TString targetName;
87 Double_t xOff;
88 Double_t yOff;
89 Double_t start[3];
90 Double_t end[3];
91 Double_t bparam;
92 Double_t mparam[10];
93 Double_t startZ;
94 Double_t endZ;
95 Bool_t targetFromGeometry; // flag to indicate coordinates set from geometry
97};
98
99#endif // SHIPGEN_PYTHIA8GENERATOR_H_ /* !PNDP8GENERATOR_H */
void SetMom(Double_t mom)
~Pythia8Generator() override
Double_t mparam[10]
Double_t maxCrossSection
void SetTarget(TString s, Double_t x, Double_t y)
void SetParameters(char *)
Bool_t Init() override
void SetTargetCoordinates(Double_t start_z, Double_t end_z)
void SetId(Double_t id)
Float_t subprocCodes[16]
Bool_t Init(const char *inFile) override
FairLogger * fLogger
pointer to a file
Bool_t ReadEvent(FairPrimaryGenerator *) override
Float_t ancestors[16]
Bool_t Init(const char *inFile, int startEvent) override
Pythia8::Pythia * fPythia
void SetfFDs(Double_t z)
GenieGenerator * fMaterialInvestigator
TTree * fTree
don't make it persistent, magic ROOT command
std::shared_ptr< Pythia8::RndmEngine > fRandomEngine
virtual Bool_t Init(const char *, int)=0