SND@LHC Software
Loading...
Searching...
No Matches
Pythia8Generator Class Reference

#include <Pythia8Generator.h>

Inheritance diagram for Pythia8Generator:
Collaboration diagram for Pythia8Generator:

Public Member Functions

 Pythia8Generator ()
 
virtual ~Pythia8Generator ()
 
Bool_t ReadEvent (FairPrimaryGenerator *)
 
void SetParameters (char *)
 
void Print ()
 
virtual Bool_t Init ()
 
void SetMom (Double_t mom)
 
void SetId (Double_t id)
 
void UseRandom1 ()
 
void UseRandom3 ()
 
void UseExternalFile (const char *x, Int_t i)
 
void SetfFDs (Double_t z)
 
void SetTarget (TString s, Double_t x, Double_t y)
 
Int_t nrOfRetries ()
 

Protected Member Functions

 ClassDef (Pythia8Generator, 3)
 

Protected Attributes

Double_t fMom
 
Int_t fId
 
Bool_t fUseRandom1
 
Bool_t fUseRandom3
 
const char * fextFile
 
Float_t hpx [1]
 
Float_t hpy [1]
 
Float_t hpz [1]
 
Float_t hE [1]
 
Float_t hM [1]
 
Float_t mpx [1]
 
Float_t mpy [1]
 
Float_t mpz [1]
 
Float_t mE [1]
 
Float_t hid [1]
 
Float_t mid [1]
 
Float_t ck [1]
 
Float_t ancestors [16]
 
Float_t subprocCodes [16]
 
Int_t fNevents
 
Int_t fn
 
Int_t firstEvent
 
Int_t fShipEventNr
 
TFile * fInputFile
 
FairLogger * fLogger
 pointer to a file
 
TTree * fTree
 don't make it persistent, magic ROOT command
 
Pythia8::Pythia * fPythia
 
Double_t fFDs
 
Int_t fnRetries
 
GenieGeneratorfMaterialInvestigator
 
TString targetName
 
Double_t xOff
 
Double_t yOff
 
Double_t start [3]
 
Double_t end [3]
 
Double_t bparam
 
Double_t mparam [10]
 
Double_t startZ
 
Double_t endZ
 
Double_t maxCrossSection
 

Private Attributes

Pythia8::RndmEngine * fRandomEngine
 

Detailed Description

Definition at line 13 of file Pythia8Generator.h.

Constructor & Destructor Documentation

◆ Pythia8Generator()

Pythia8Generator::Pythia8Generator ( )

default constructor

Definition at line 18 of file Pythia8Generator.cxx.

19{
20 fUseRandom1 = kFALSE;
21 fUseRandom3 = kTRUE;
22 fId = 2212; // proton
23 fMom = 400; // proton
24 fFDs = 7.7/10.4; // correction for Pythia6 to match measured Ds production
25 fextFile = "";
26 fInputFile = NULL;
27 targetName = "";
28 xOff=0; yOff=0;
29 fPythia = new Pythia8::Pythia();
30}
const char * fextFile
Pythia8::Pythia * fPythia

◆ ~Pythia8Generator()

Pythia8Generator::~Pythia8Generator ( )
virtual

destructor

Definition at line 140 of file Pythia8Generator.cxx.

141{
142}

Member Function Documentation

◆ ClassDef()

Pythia8Generator::ClassDef ( Pythia8Generator  ,
 
)
protected

◆ Init()

Bool_t Pythia8Generator::Init ( )
virtual

Definition at line 34 of file Pythia8Generator.cxx.

35{
38 if (fextFile && *fextFile) {
39 if (0 == strncmp("/eos",fextFile,4) ) {
40 TString tmp = gSystem->Getenv("EOSSHIP");
42 fInputFile = TFile::Open(tmp);
43 LOGF(info, "Open external file with charm or beauty hadrons on eos: %s", tmp.Data());
44 if (!fInputFile) {
45 LOG(FATAL) << "Error opening input file. You may have forgotten to provide a krb5 token. Try kinit username@lxplus.cern.ch";
46 return kFALSE; }
47 }else{
48 LOGF(info, "Open external file with charm or beauty hadrons: %s", fextFile);
49 fInputFile = new TFile(fextFile);
50 if (!fInputFile) {
51 LOG(FATAL) << "Error opening input file";
52 return kFALSE; }
53 }
54 if (fInputFile->IsZombie()) {
55 LOG(FATAL) << "File is corrupted";
56 return kFALSE; }
57 fTree = (TTree *)fInputFile->Get("pythia6");
58 fNevents = fTree->GetEntries();
59 fn = firstEvent;
60 fTree->SetBranchAddress("id",&hid); // particle id
61 fTree->SetBranchAddress("px",&hpx); // momentum
62 fTree->SetBranchAddress("py",&hpy);
63 fTree->SetBranchAddress("pz",&hpz);
64 fTree->SetBranchAddress("E",&hE);
65 fTree->SetBranchAddress("M",&hM);
66 fTree->SetBranchAddress("mid",&mid); // mother
67 fTree->SetBranchAddress("mpx",&mpx); // momentum
68 fTree->SetBranchAddress("mpy",&mpy);
69 fTree->SetBranchAddress("mpz",&mpz);
70 fTree->SetBranchAddress("mE",&mE);
71 if (fTree->GetBranch("k")){
72 fTree->SetBranchAddress("k",&ck);
73 if (fTree->GetBranch("a0")){
74 for(Int_t i=0; i<16; i++){
75 TString na = "a";na+=i;
76 fTree->SetBranchAddress(na,&(ancestors[i]));
77 TString ns = "s";ns+=i;
78 fTree->SetBranchAddress(ns,&(subprocCodes[i]));
79 }
80 }
81 }else{
82 ck[0]=1;subprocCodes[0]=88;ancestors[0]=2212;
83 }
84 fPythia->readString("ProcessLevel:all = off");
85// find all long lived particles in pythia
86 Int_t n = 1;
87 while(n!=0){
88 n = fPythia->particleData.nextId(n);
89#ifndef PYTHIA8_V
90 Pythia8::ParticleDataEntry* p = fPythia->particleData.particleDataEntryPtr(n);
91#else
92#if PYTHIA8_V<8309
93 Pythia8::ParticleDataEntry* p = fPythia->particleData.particleDataEntryPtr(n);
94#else
95 Pythia8::ParticleDataEntryPtr p = fPythia->particleData.particleDataEntryPtr(n);
96#endif
97#endif
98 if (p->tau0()>1){
99 std::string particle = std::to_string(n)+":mayDecay = false";
100 fPythia->readString(particle);
101 LOGF(info, "Made %s stable for Pythia, should decay in Geant4", p->name().c_str());
102 }
103 }
104 } else {
105 fPythia->setRndmEnginePtr(fRandomEngine);
106 fPythia->settings.mode("Beams:idA", fId);
107 fPythia->settings.mode("Beams:idB", 2212);
108 fPythia->settings.mode("Beams:frameType", 2);
109 fPythia->settings.parm("Beams:eA",fMom);
110 fPythia->settings.parm("Beams:eB",0.);
111 }
112 fPythia->init();
113 if (targetName!=""){
115 TGeoVolume* top = gGeoManager->GetTopVolume();
116 TGeoNode* target = top->FindNode(targetName);
117 if (!target){
118 LOGF(error, "target not found, %s, program will crash", targetName.Data());
119 }
120 Double_t z_middle = target->GetMatrix()->GetTranslation()[2];
121 TGeoBBox* sha = (TGeoBBox*)target->GetVolume()->GetShape();
122 startZ = z_middle - sha->GetDZ();
123 endZ = z_middle + sha->GetDZ();
124 start[0]=xOff;
125 start[1]=yOff;
126 start[2]=startZ;
127 end[0]=xOff;
128 end[1]=yOff;
129 end[2]=endZ;
130//find maximum interaction length
133 }
134 return kTRUE;
135}
Double_t MeanMaterialBudget(const Double_t *start, const Double_t *end, Double_t *mparam)
Double_t mparam[10]
Pythia8::RndmEngine * fRandomEngine
Float_t subprocCodes[16]
Float_t ancestors[16]
GenieGenerator * fMaterialInvestigator
TTree * fTree
don't make it persistent, magic ROOT command
int i
Definition ShipAna.py:86
int ns
Definition hepunit.py:95

◆ nrOfRetries()

Int_t Pythia8Generator::nrOfRetries ( )
inline

Definition at line 37 of file Pythia8Generator.h.

37{ return fnRetries; };

◆ Print()

void Pythia8Generator::Print ( )

◆ ReadEvent()

Bool_t Pythia8Generator::ReadEvent ( FairPrimaryGenerator *  cpg)

public method ReadEvent

Definition at line 146 of file Pythia8Generator.cxx.

147{
148 Double_t x,y,z,px,py,pz,dl,e,tof;
149 Int_t im,id,key;
150 fnRetries =0;
151// take charm hadrons from external file
152// correct eventually for too much primary Ds produced by pythia6
153 key=0;
154 bool l = true;
155 while(l){
156 if (fn==fNevents) {LOG(WARNING) << "End of input file. Rewind.";}
157 fTree->GetEntry((fn+1)%fNevents);
158// check that this and next entry is charm, otherwise continue reading
159 l = false;
160 if (int(mE[0])== 0){ l = true; }
161 bool isDs = false;
162 if ( int(fabs(hid[0]) ) == 431){ isDs = true; }
163 fTree->GetEntry(fn%fNevents);
164 if (int(mE[0])== 0){ l = true; }
165 fn++;
166 if ( int(fabs(hid[0]) ) == 431 || isDs ){
167 Double_t rnr = gRandom->Uniform(0,1);
168 if( rnr>fFDs ) {
169 l = true;
170 fn++;
171 }
172 }
173 }
174 Double_t zinter=0;
175 if (targetName!=""){
176// calculate primary proton interaction point:
177// loop over trajectory between start and end to pick an interaction point, copied from GenieGenerator and adapted to hadrons
178 Double_t prob2int = -1.;
179 Double_t rndm = 0.;
180 Double_t sigma;
181 Int_t count=0;
182 Double_t zinterStart = start[2];
183// simulate more downstream interaction points for interactions down in the cascade
184 Int_t nInter = ck[0]; if (nInter>16){nInter=16;}
185 for( Int_t nI=0; nI<nInter; nI++){
186 // if (!subprocCodes[nI]<90){continue;} //if process is not inelastic, go to next. Changed by taking now collision length
187 prob2int = -1.;
188 Int_t intLengthFactor = 1; // for nucleons
189 if (TMath::Abs(ancestors[nI]) < 1000){intLengthFactor = 1.16;} // for mesons
190 // Fe: nuclear /\ 16.77 cm pion 20.42 cm f=1.22
191 // W: nuclear /\ 9.946 cm pion 11.33 cm f=1.14
192 // Mo: nuclear /\ 15.25 cm pion 17.98 cm f=1.18
193 while (prob2int<rndm) {
194 //place x,y,z uniform along path
195 zinter = gRandom->Uniform(zinterStart,end[2]);
196 Double_t point[3]={xOff,yOff,zinter};
198 Double_t interLength = mparam[8] * intLengthFactor * 1.7; // 1.7 = interaction length / collision length from PDG Tables
199 TGeoNode *node = gGeoManager->FindNode(point[0],point[1],point[2]);
200 TGeoMaterial *mat = 0;
201 if (node && !gGeoManager->IsOutside()) {
202 mat = node->GetVolume()->GetMaterial();
203 Double_t n = mat->GetDensity()/mat->GetA();
204 sigma = 1./(n*mat->GetIntLen())/mbarn; // no need to multiply with intLengthFactor, will cancel in next equation.
205 prob2int = TMath::Exp(-interLength)*sigma/maxCrossSection;
206 }else{
207 prob2int=0.;
208 }
209 rndm = gRandom->Uniform(0.,1.);
210 count+=1;
211 }
212 zinterStart = zinter;
213 }
214 zinter = zinter*cm;
215 }
216 for(Int_t c=0; c<2; c++){
217 if(c>0){
218 fTree->GetEntry(fn%fNevents);
219 fn++;
220 }
221 fPythia->event.reset();
222 id = (Int_t)hid[0];
223 fPythia->event.append( id, 1, 0, 0, hpx[0], hpy[0], hpz[0], hE[0], hM[0], 0., 9. );
224//simulate displaced vertex, Pythia8 will not do it
225 Double_t tau0 = fPythia->particleData.tau0(id); // ctau in mm
226 dl = gRandom->Exp(tau0) / hM[0]; // mm/GeV
227 fPythia->next();
228 Int_t addedParticles = 0;
229 if(c==0){
230// original particle responsible for interaction, "mother of charm" from external file
231 px=mpx[0];
232 py=mpy[0];
233 pz=mpz[0];
234 x=0.;
235 y=0.;
236 z=zinter;
237 id=mid[0];
238 cpg->AddTrack(id,px,py,pz,x/cm,y/cm,z/cm,-1,false);
239 addedParticles +=1;
240 }
241 for(Int_t ii=1; ii<fPythia->event.size(); ii++){
242 id = fPythia->event[ii].id();
243 Bool_t wanttracking=false;
244 if(fPythia->event[ii].isFinal()){ wanttracking=true; }
245 if (ii>1){
246 z = fPythia->event[ii].zProd()+dl*fPythia->event[1].pz()+zinter;
247 x = fPythia->event[ii].xProd()+dl*fPythia->event[1].px();
248 y = fPythia->event[ii].yProd()+dl*fPythia->event[1].py();
249 tof = fPythia->event[ii].tProd()/ (10*c_light) + dl*fPythia->event[1].e()/cm/c_light;
250 }else{
251 z = fPythia->event[ii].zProd()+zinter;
252 x = fPythia->event[ii].xProd();
253 y = fPythia->event[ii].yProd();
254 tof = fPythia->event[ii].tProd() / (10*c_light) ; // to go from mm to s
255 }
256 pz = fPythia->event[ii].pz();
257 px = fPythia->event[ii].px();
258 py = fPythia->event[ii].py();
259 e = fPythia->event[ii].e();
260 im = fPythia->event[ii].mother1()+key;
261
262 if (ii==1){im = 0;}
263 cpg->AddTrack(id,px,py,pz,x/cm,y/cm,z/cm,im,wanttracking,e,tof,1.);
264 addedParticles+=1;
265 }
266 key+=addedParticles-1; // pythia counts from 1
267 }
268 counter+=1;
269// now the underyling event
270 bool lx = true;
271 while(lx){
272 fTree->GetEntry(fn%fNevents);
273 lx = false;
274 if (mE[0] == 0){
275 lx = true;
276 fn++;
277 cpg->AddTrack((Int_t)hid[0],hpx[0],hpy[0],hpz[0],(mpx[0]+fPythia->event[0].xProd())/cm,
278 (mpy[0]+fPythia->event[0].yProd())/cm,
279 (mpz[0]+fPythia->event[0].zProd())/cm+zinter/cm,-1,true);
280 // mpx,mpy,mpz are the vertex coordinates with respect to charm hadron, first particle in Pythia is (system)
281 }
282 }
283
284 return kTRUE;
285}
const Double_t c_light
const Double_t mbarn
Int_t counter
Double_t cm
c
Definition hnl.py:100

◆ SetfFDs()

void Pythia8Generator::SetfFDs ( Double_t  z)
inline

Definition at line 35 of file Pythia8Generator.h.

35{ fFDs = z; };

◆ SetId()

void Pythia8Generator::SetId ( Double_t  id)
inline

Definition at line 31 of file Pythia8Generator.h.

31{ fId = id; };

◆ SetMom()

void Pythia8Generator::SetMom ( Double_t  mom)
inline

Definition at line 30 of file Pythia8Generator.h.

30{ fMom = mom; };

◆ SetParameters()

void Pythia8Generator::SetParameters ( char *  )

◆ SetTarget()

void Pythia8Generator::SetTarget ( TString  s,
Double_t  x,
Double_t  y 
)
inline

Definition at line 36 of file Pythia8Generator.h.

◆ UseExternalFile()

void Pythia8Generator::UseExternalFile ( const char *  x,
Int_t  i 
)
inline

Definition at line 34 of file Pythia8Generator.h.

34{ fextFile = x; firstEvent=i; };

◆ UseRandom1()

void Pythia8Generator::UseRandom1 ( )
inline

Definition at line 32 of file Pythia8Generator.h.

32{ fUseRandom1 = kTRUE; fUseRandom3 = kFALSE; };

◆ UseRandom3()

void Pythia8Generator::UseRandom3 ( )
inline

Definition at line 33 of file Pythia8Generator.h.

33{ fUseRandom1 = kFALSE; fUseRandom3 = kTRUE; };

Member Data Documentation

◆ ancestors

Float_t Pythia8Generator::ancestors[16]
protected

Definition at line 51 of file Pythia8Generator.h.

◆ bparam

Double_t Pythia8Generator::bparam
protected

Definition at line 66 of file Pythia8Generator.h.

◆ ck

Float_t Pythia8Generator::ck[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ end

Double_t Pythia8Generator::end[3]
protected

Definition at line 65 of file Pythia8Generator.h.

◆ endZ

Double_t Pythia8Generator::endZ
protected

Definition at line 69 of file Pythia8Generator.h.

◆ fextFile

const char* Pythia8Generator::fextFile
protected

Definition at line 49 of file Pythia8Generator.h.

◆ fFDs

Double_t Pythia8Generator::fFDs
protected

Definition at line 57 of file Pythia8Generator.h.

◆ fId

Int_t Pythia8Generator::fId
protected

Definition at line 46 of file Pythia8Generator.h.

◆ fInputFile

TFile* Pythia8Generator::fInputFile
protected

Definition at line 53 of file Pythia8Generator.h.

◆ firstEvent

Int_t Pythia8Generator::firstEvent
protected

Definition at line 52 of file Pythia8Generator.h.

◆ fLogger

FairLogger* Pythia8Generator::fLogger
protected

pointer to a file

Definition at line 54 of file Pythia8Generator.h.

◆ fMaterialInvestigator

GenieGenerator* Pythia8Generator::fMaterialInvestigator
protected

Definition at line 59 of file Pythia8Generator.h.

◆ fMom

Double_t Pythia8Generator::fMom
protected

Definition at line 45 of file Pythia8Generator.h.

◆ fn

Int_t Pythia8Generator::fn
protected

Definition at line 52 of file Pythia8Generator.h.

◆ fNevents

Int_t Pythia8Generator::fNevents
protected

Definition at line 52 of file Pythia8Generator.h.

◆ fnRetries

Int_t Pythia8Generator::fnRetries
protected

Definition at line 58 of file Pythia8Generator.h.

◆ fPythia

Pythia8::Pythia* Pythia8Generator::fPythia
protected

Definition at line 56 of file Pythia8Generator.h.

◆ fRandomEngine

Pythia8::RndmEngine* Pythia8Generator::fRandomEngine
private

Definition at line 41 of file Pythia8Generator.h.

◆ fShipEventNr

Int_t Pythia8Generator::fShipEventNr
protected

Definition at line 52 of file Pythia8Generator.h.

◆ fTree

TTree* Pythia8Generator::fTree
protected

don't make it persistent, magic ROOT command

Definition at line 55 of file Pythia8Generator.h.

◆ fUseRandom1

Bool_t Pythia8Generator::fUseRandom1
protected

Definition at line 47 of file Pythia8Generator.h.

◆ fUseRandom3

Bool_t Pythia8Generator::fUseRandom3
protected

Definition at line 48 of file Pythia8Generator.h.

◆ hE

Float_t Pythia8Generator::hE[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ hid

Float_t Pythia8Generator::hid[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ hM

Float_t Pythia8Generator::hM[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ hpx

Float_t Pythia8Generator::hpx[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ hpy

Float_t Pythia8Generator::hpy[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ hpz

Float_t Pythia8Generator::hpz[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ maxCrossSection

Double_t Pythia8Generator::maxCrossSection
protected

Definition at line 70 of file Pythia8Generator.h.

◆ mE

Float_t Pythia8Generator::mE[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ mid

Float_t Pythia8Generator::mid[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ mparam

Double_t Pythia8Generator::mparam[10]
protected

Definition at line 67 of file Pythia8Generator.h.

◆ mpx

Float_t Pythia8Generator::mpx[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ mpy

Float_t Pythia8Generator::mpy[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ mpz

Float_t Pythia8Generator::mpz[1]
protected

Definition at line 50 of file Pythia8Generator.h.

◆ start

Double_t Pythia8Generator::start[3]
protected

Definition at line 64 of file Pythia8Generator.h.

◆ startZ

Double_t Pythia8Generator::startZ
protected

Definition at line 68 of file Pythia8Generator.h.

◆ subprocCodes

Float_t Pythia8Generator::subprocCodes[16]
protected

Definition at line 51 of file Pythia8Generator.h.

◆ targetName

TString Pythia8Generator::targetName
protected

Definition at line 61 of file Pythia8Generator.h.

◆ xOff

Double_t Pythia8Generator::xOff
protected

Definition at line 62 of file Pythia8Generator.h.

◆ yOff

Double_t Pythia8Generator::yOff
protected

Definition at line 63 of file Pythia8Generator.h.


The documentation for this class was generated from the following files: