SND@LHC Software
Loading...
Searching...
No Matches
FixedTargetGenerator.h
Go to the documentation of this file.
1#ifndef FIXEDTARGETGENERATOR_H
2#define FIXEDTARGETGENERATOR_H 1
3
4#include "TROOT.h"
5#include "FairGenerator.h"
6#include "Pythia8/Pythia.h"
7#include "FairLogger.h" // for FairLogger, MESSAGE_ORIGIN
8#include "TTree.h"
9#include "TNtuple.h"
10#include "GenieGenerator.h"
11
12class FairPrimaryGenerator;
13class EvtGenDecays;
14
15class FixedTargetGenerator : public FairGenerator
16{
17 public:
18
21
23 virtual ~FixedTargetGenerator();
24
26 Bool_t ReadEvent(FairPrimaryGenerator*);
27 void SetParameters(char*);
28 void Print();
29
30 virtual Bool_t Init();
31 Bool_t InitForCharmOrBeauty(TString fInName, Int_t nev, Double_t npots=5E13, Int_t nStart=0);
32
33 void SetMom(Double_t mom) { fMom = mom; };
34 void UseRandom1() { fUseRandom1 = kTRUE; fUseRandom3 = kFALSE; };
35 void UseRandom3() { fUseRandom1 = kFALSE; fUseRandom3 = kTRUE; };
36 void SetCharmTarget(Bool_t charmtarget = true) {fcharmtarget = charmtarget;}; //charm geometry uses a different target, default one is usual
37 void SetTarget(TString s, Double_t x,Double_t y ) { targetName = s; xOff=x; yOff=y; };
38 void SetBoost(Double_t f) { fBoost = f; } // boost factor for rare di-muon decays
39 void SetG4only() { G4only = true; } // only run Geant4, no pythia primary interaction
40 void SetTauOnly() { tauOnly = true; } // only have Ds decay to tau
41 void SetJpsiMainly() { JpsiMainly = true; } // let all Jpsi decay to mumu
42 void SetOnlyMuons() { OnlyMuons = true; } // only transport muons
43 void SetDrellYan() { DrellYan = true; } // only generate prompt Z0* processes
44 void SetPhotonCollision() { PhotonCollision = true; } // only generate prompt photon processes
45 void WithEvtGen() { withEvtGen = true;} // use EvtGen as external decayer to Pythia, experimental phase, only works for one Pythia instance
46 void SetChibb(Double_t x) { chibb = x; } // chibb = bbbar over mbias cross section
47 void SetChicc(Double_t x) { chicc = x; } // chicc = ccbar over mbias cross section
48 inline void SetSeed(Double_t seed){fSeed=seed;}
49 inline void SetHeartBeat(Int_t x){heartbeat=x;}
50 inline void SetEnergyCut(Float_t emax) {EMax=emax;}// min energy to be copied to Geant4
51 inline void SetDebug(Bool_t x){Debug=x;}
52 inline void SetOpt4DP(TNtuple* t){withNtuple=kTRUE; fNtuple = t ; }
53 Double_t GetPotForCharm(){return nrpotspill/wspill;}
54 Pythia8::Pythia* GetPythia() {return fPythiaP;}
55 Pythia8::Pythia* GetPythiaN() {return fPythiaN;}
56 private:
57
58 Pythia8::RndmEngine* fRandomEngine;
59
60 protected:
61
62 Double_t fMom; // proton momentum
63 Bool_t fUseRandom1; // flag to use TRandom1
64 Bool_t fUseRandom3; // flag to use TRandom3 (default)
69 FairLogger* fLogger;
70 Pythia8::Pythia* fPythiaN;
71 Pythia8::Pythia* fPythiaP;
72 EvtGenDecays* evtgenN;
73 EvtGenDecays* evtgenP;
75 Bool_t withNtuple;
76 TNtuple* fNtuple;
78 Double_t xOff;
79 Double_t yOff;
80 Double_t start[3];
81 Double_t end[3];
82 Double_t bparam;
83 Double_t mparam[10];
84 Double_t startZ;
85 Double_t endZ;
87 TFile* fin;
88 TNtuple* nTree;
90 Int_t heartbeat;
91
93};
94#endif /* !FIXEDTARGETGENERATOR_H */
void SetCharmTarget(Bool_t charmtarget=true)
void SetTarget(TString s, Double_t x, Double_t y)
Pythia8::Pythia * GetPythiaN()
GenieGenerator * fMaterialInvestigator
TNtuple * fNtuple
special option for Dark Photon physics studies
Pythia8::RndmEngine * fRandomEngine
Bool_t ReadEvent(FairPrimaryGenerator *)
ClassDef(FixedTargetGenerator, 2)
void SetEnergyCut(Float_t emax)
void SetParameters(char *)
void SetMom(Double_t mom)
Pythia8::Pythia * fPythiaP
Pythia8::Pythia * GetPythia()
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)
void SetSeed(Double_t seed)