FairShip
Loading...
Searching...
No Matches
FixedTargetGenerator.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_FIXEDTARGETGENERATOR_H_
6#define SHIPGEN_FIXEDTARGETGENERATOR_H_
7
8#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
9#include "Generator.h"
10#include "GenieGenerator.h"
11#include "Pythia8/Pythia.h"
12#include "TNtuple.h"
13#include "TROOT.h"
14#include "TTree.h"
15
16namespace Pythia8 {
17class EvtGenDecays;
18}
19
20class FairPrimaryGenerator;
21
23 public:
26
28 ~FixedTargetGenerator() override;
29
31 Bool_t ReadEvent(FairPrimaryGenerator*) override;
32 void SetParameters(char*);
33 void Print();
34 Bool_t Init(const char* inFile) override { return Init(inFile, 0); };
35
36 Bool_t Init(const char* inFile, int startEvent) override {
37 LOG(warning) << "Init with files not implemented for FixedTargetGenerator. "
38 "Using default Init() instead";
39 return Init();
40 };
42 Bool_t Init() override;
43
44 Bool_t InitForCharmOrBeauty(TString fInName, Int_t nev, Double_t npots = 5E13,
45 Int_t nStart = 0);
46
47 void SetMom(Double_t mom) { fMom = mom; };
48 void UseRandom1() {
49 fUseRandom1 = kTRUE;
50 fUseRandom3 = kFALSE;
51 };
52 void UseRandom3() {
53 fUseRandom1 = kFALSE;
54 fUseRandom3 = kTRUE;
55 };
56 void SetTarget(TString s, Double_t x, Double_t y) {
57 targetName = s;
58 xOff = x;
59 yOff = y;
60 };
61 void SetZoffset(Double_t z) { zOff = z; };
62 void SetXoffset(Double_t x) { xOff = x; };
63 void SetYoffset(Double_t y) { yOff = y; };
64 void SetSmearBeam(Double_t sb) { fsmearBeam = sb; };
65 void SetPaintRadius(Double_t r) { fPaintBeam = r; };
66 void SetTargetCoordinates(Double_t start_z, Double_t end_z) {
67 startZ = start_z;
68 endZ = end_z;
69 targetFromGeometry = true;
70 };
71 void SetBoost(Double_t f) {
72 fBoost = f;
73 } // boost factor for rare di-muon decays
74 void SetG4only() {
75 G4only = true;
76 } // only run Geant4, no pythia primary interaction
77 void SetTauOnly() { tauOnly = true; } // only have Ds decay to tau
78 void SetJpsiMainly() { JpsiMainly = true; } // let all Jpsi decay to mumu
79 void SetOnlyMuons() { OnlyMuons = true; } // only transport muons
80 void SetDrellYan() { DrellYan = true; } // only generate prompt Z0* processes
82 PhotonCollision = true;
83 } // only generate prompt photon processes
84 void WithEvtGen() {
85 withEvtGen = true;
86 } // use EvtGen as external decayer to Pythia, experimental phase, only works
87 // for one Pythia instance
88 void SetChibb(Double_t x) {
89 chibb = x;
90 } // chibb = bbbar over mbias cross section
91 void SetChicc(Double_t x) {
92 chicc = x;
93 } // chicc = ccbar over mbias cross section
94 inline void SetSeed(Double_t seed) { fSeed = seed; }
95 inline void SetHeartBeat(Int_t x) { heartbeat = x; }
96 inline void SetEnergyCut(Float_t emax) {
97 EMax = emax;
98 } // min energy to be copied to Geant4
99 inline void SetDebug(Bool_t x) { Debug = x; }
100 inline void SetOpt4DP(TNtuple* t) {
101 withNtuple = kTRUE;
102 fNtuple = t;
103 }
104 Double_t GetPotForCharm() { return nrpotspill / wspill; }
105 Pythia8::Pythia* GetPythia() { return fPythiaP; }
106 Pythia8::Pythia* GetPythiaN() { return fPythiaN; }
107
108 private:
109 std::shared_ptr<Pythia8::RndmEngine> fRandomEngine;
110
111 protected:
112 Double_t fMom; // proton momentum
113 Bool_t fUseRandom1; // flag to use TRandom1
114 Bool_t fUseRandom3; // flag to use TRandom3 (default)
119 FairLogger* fLogger;
120 Pythia8::Pythia* fPythiaN;
121 Pythia8::Pythia* fPythiaP;
122 Pythia8::EvtGenDecays* evtgenN;
123 Pythia8::EvtGenDecays* evtgenP;
125 Bool_t withNtuple;
126 TNtuple* fNtuple;
128 Double_t xOff;
129 Double_t yOff;
130 Double_t zOff;
131 Double_t fsmearBeam; // finite beam size
132 Double_t fPaintBeam; // beam painting radius
133 Double_t start[3];
134 Double_t end[3];
135 Double_t bparam;
136 Double_t mparam[10];
137 Double_t startZ;
138 Double_t endZ;
139 Bool_t targetFromGeometry; // flag to indicate coordinates set from geometry
141 TFile* fin;
142 TNtuple* nTree;
146};
147#endif // SHIPGEN_FIXEDTARGETGENERATOR_H_
Bool_t ReadEvent(FairPrimaryGenerator *) override
void SetTargetCoordinates(Double_t start_z, Double_t end_z)
void SetTarget(TString s, Double_t x, Double_t y)
Pythia8::Pythia * GetPythiaN()
GenieGenerator * fMaterialInvestigator
TNtuple * fNtuple
special option for Dark Photon physics studies
void SetXoffset(Double_t x)
void SetChibb(Double_t x)
Bool_t Init(const char *inFile, int startEvent) override
void SetEnergyCut(Float_t emax)
void SetSmearBeam(Double_t sb)
void SetYoffset(Double_t y)
void SetParameters(char *)
Bool_t Init(const char *inFile) override
void SetMom(Double_t mom)
void SetZoffset(Double_t z)
Pythia8::Pythia * fPythiaP
Pythia8::Pythia * GetPythia()
void SetOpt4DP(TNtuple *t)
Pythia8::EvtGenDecays * evtgenN
void SetBoost(Double_t f)
void SetPaintRadius(Double_t r)
Pythia8::Pythia * fPythiaN
don't make it persistent, magic ROOT command
Bool_t InitForCharmOrBeauty(TString fInName, Int_t nev, Double_t npots=5E13, Int_t nStart=0)
Pythia8::EvtGenDecays * evtgenP
void SetSeed(Double_t seed)
void SetChicc(Double_t x)
std::shared_ptr< Pythia8::RndmEngine > fRandomEngine
virtual Bool_t Init(const char *, int)=0