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

#include <DPPythia8Generator.h>

Inheritance diagram for DPPythia8Generator:
Collaboration diagram for DPPythia8Generator:

Public Member Functions

 DPPythia8Generator ()
 
virtual ~DPPythia8Generator ()
 
Bool_t ReadEvent (FairPrimaryGenerator *)
 
void SetParameters (char *)
 
void Print ()
 
void List (int id)
 
virtual Bool_t Init ()
 
void SetMom (Double_t mom)
 
Double_t GetMom ()
 
void SetId (Double_t id)
 
void SetDPId (Int_t id)
 
Int_t GetDPId ()
 
void SetLmin (Double_t z)
 
void SetLmax (Double_t z)
 
void SetSmearBeam (Double_t sb)
 
void SetfFDs (Double_t z)
 
void UseRandom1 ()
 
void UseRandom3 ()
 
void UseExternalFile (const char *x, Int_t i)
 
void SetPbrem (TH2F *pdf)
 
Bool_t IsPbrem ()
 
void SetDY ()
 
Double_t MinDPMass ()
 
void SetMinDPMass (Double_t m)
 
void UseDeepCopy ()
 
Int_t nrOfRetries ()
 
Int_t nrOfDP ()
 
Pythia8::Pythia * getPythiaInstance ()
 

Public Attributes

Pythia8::Pythia * fPythia
 

Protected Member Functions

 ClassDef (DPPythia8Generator, 2)
 don't make it persistent, magic ROOT command
 

Protected Attributes

Double_t fMom
 
Int_t fDP
 
Int_t fId
 
Bool_t fUseRandom1
 
Bool_t fUseRandom3
 
Bool_t fpbrem
 
TH2F * fpbremPDF
 
Bool_t fdy
 
Double_t fDPminM
 
Double_t fLmin
 
Double_t fLmax
 
Int_t fnRetries
 
Int_t fnDPtot
 
Double_t fctau
 
Double_t fFDs
 
Double_t fsmearBeam
 
const char * fextFile
 
TFile * fInputFile
 
TTree * fTree
 pointer to a file
 
Int_t fNevents
 
Int_t fn
 
Int_t firstEvent
 
Int_t fShipEventNr
 
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]
 
Bool_t fDeepCopy
 
FairLogger * fLogger
 

Private Attributes

Pythia8::RndmEngine * fRandomEngine
 

Detailed Description

Definition at line 23 of file DPPythia8Generator.h.

Constructor & Destructor Documentation

◆ DPPythia8Generator()

DPPythia8Generator::DPPythia8Generator ( )

default constructor

Definition at line 14 of file DPPythia8Generator.cxx.

15{
16 //fHadDecay = false;
17 fId = 2212; // proton
18 fMom = 400; // proton
19 fDP = 9900015; // DP pdg code
20 fLmin = 5000.*cm; // mm minimum decay position z ROOT units !
21 fLmax = 12000.*cm; // mm maximum decay position z
22 fFDs = 7.7/10.4; // correction for Pythia6 to match measured Ds production
23 fpbrem = kFALSE;
24 fpbremPDF = 0;
25 fdy = kFALSE;
26 fDPminM = 0.5;
27 fextFile = "";
28 fInputFile = NULL;
29 fnRetries = 0;
30 fnDPtot = 0;
31 fShipEventNr = 0;
32 fPythia = new Pythia8::Pythia();
33 //fPythiaHadDecay = new Pythia8::Pythia();
34}
Double_t cm
Pythia8::Pythia * fPythia

◆ ~DPPythia8Generator()

DPPythia8Generator::~DPPythia8Generator ( )
virtual

destructor

Definition at line 140 of file DPPythia8Generator.cxx.

141{
142}

Member Function Documentation

◆ ClassDef()

DPPythia8Generator::ClassDef ( DPPythia8Generator  ,
 
)
protected

don't make it persistent, magic ROOT command

◆ GetDPId()

Int_t DPPythia8Generator::GetDPId ( )
inline

Definition at line 50 of file DPPythia8Generator.h.

50{ return fDP; };

◆ GetMom()

Double_t DPPythia8Generator::GetMom ( )
inline

Definition at line 47 of file DPPythia8Generator.h.

47{ return fMom; };

◆ getPythiaInstance()

Pythia8::Pythia * DPPythia8Generator::getPythiaInstance ( )
inline

Definition at line 75 of file DPPythia8Generator.h.

75{return fPythia;};

◆ Init()

Bool_t DPPythia8Generator::Init ( )
virtual

Definition at line 38 of file DPPythia8Generator.cxx.

39{
42 fPythia->setRndmEnginePtr(fRandomEngine);
43 //fPythiaHadDecay->setRndmEnginePtr(fRandomEngine);
44 fn = 0;
45 if (fextFile && *fextFile){
46 /*if (0 == strncmp("/eos",fextFile,4) ) {
47 if (0 == strncmp("/eos",fextFile,4) ) {
48 TString tmp = gSystem->Getenv("EOSSHIP");
49 tmp+=fextFile;
50 fInputFile = TFile::Open(tmp);
51 fLogger->Info(MESSAGE_ORIGIN,""Open external file with charm or beauty hadrons on eos: %s",tmp->Data());
52 if (!fInputFile) {
53 fLogger->Fatal(MESSAGE_ORIGIN, "Error opening input file. You may have forgotten to provide a krb5 token. Try kinit username@lxplus.cern.ch");
54 return kFALSE; }
55 }else{
56 fLogger->Info(MESSAGE_ORIGIN,"Open external file with charm or beauty hadrons: %s",fextFile);
57 fInputFile = new TFile(fextFile);
58 if (!fInputFile) {
59 fLogger->Fatal(MESSAGE_ORIGIN, "Error opening input file");
60 return kFALSE; }
61 }
62 if (fInputFile->IsZombie()) {
63 fLogger->Fatal(MESSAGE_ORIGIN, "File is corrupted");
64 return kFALSE; }
65 fTree = (TTree *)fInputFile->Get("pythia6");
66 fNevents = fTree->GetEntries();
67 fn = firstEvent;
68 fTree->SetBranchAddress("id",&hid); // particle id
69 fTree->SetBranchAddress("px",&hpx); // momentum
70 fTree->SetBranchAddress("py",&hpy);
71 fTree->SetBranchAddress("pz",&hpz);
72 fTree->SetBranchAddress("E",&hE);
73 fTree->SetBranchAddress("M",&hM);
74 fTree->SetBranchAddress("mid",&mid); // mother
75 fTree->SetBranchAddress("mpx",&mpx); // momentum
76 fTree->SetBranchAddress("mpy",&mpy);
77 fTree->SetBranchAddress("mpz",&mpz);
78 fTree->SetBranchAddress("mE",&mE);
79 */
80 }
81 else if (!fpbrem){
82 if ( debug ){std::cout<<"Beam Momentum "<<fMom<<std::endl;}
83 fPythia->settings.mode("Beams:idA", fId);
84 fPythia->settings.mode("Beams:idB", 2212);
85 fPythia->settings.mode("Beams:frameType", 2);
86 fPythia->settings.parm("Beams:eA",fMom);
87 fPythia->settings.parm("Beams:eB",0.);
88
89 if (fdy) fPythia->settings.parm("PhaseSpace:mHatMin",fDPminM);
90
91 }
92 else {
93 if (!fpbremPDF) {
94 //std::cout << " Failed in retrieving dark photon PDF for production by proton bremstrahlung! Exiting..." << std::endl;
95 LOG(FATAL) << "Failed in retrieving dark photon PDF for production by proton bremstrahlung!";
96 return kFALSE;
97 }
98 }
99 /*if (fHadDecay) {
100 std::cout << " ******************************** " << std::endl
101 << " ** Initialise Pythia for e+e-->hadrons " << std::endl
102 << " ******************************** " << std::endl
103 << " Mass of the A: " << fPythia->particleData.m0(fDP) << " GeV" << std::endl;
104 fPythiaHadDecay->settings.mode("Beams:idA", 11);
105 fPythiaHadDecay->settings.mode("Beams:idB", -11);
106 fPythiaHadDecay->settings.mode("Beams:frameType", 1);
107 fPythiaHadDecay->settings.parm("Beams:eCM",10);
108 fPythiaHadDecay->readString("WeakSingleBoson:ffbar2ffbar(s:gm) = on");
109 //fPythiaHadDecay->readString("23:isResonance = true")
110 //force to hadrons
111 fPythiaHadDecay->readString("23:onMode = off");
112 fPythiaHadDecay->readString("23:onIfAny = 1 2 3 4 5");
113 }*/
114 TDatabasePDG* pdgBase = TDatabasePDG::Instance();
115 Double_t root_ctau = pdgBase->GetParticle(fDP)->Lifetime();
116 //fPythia->particleData.readString("4900023:useBreitWigner = false");
117 if ( debug ){
118 std::cout<<"Final particle parameters for PDGID " << fDP << ":" << std::endl;
119 List(fDP);
120 }
121 if ( debug ){std::cout<<"tau root PDG database "<<root_ctau<< "[s] ctau root = " << root_ctau*3e10 << "[cm]"<<std::endl;}
122 fctau = fPythia->particleData.tau0(fDP); //* 3.3333e-12
123 if ( debug ){std::cout<<"ctau pythia "<<fctau<<"[mm]"<<std::endl;}
124 int initPass = fPythia->init();
125 if ( debug ){std::cout<<"Pythia initialisation bool: " << initPass << std::endl;}
126
127 if (!initPass) {
128 LOG(FATAL) << "Pythia initialisation failed";
129 return kFALSE;
130 }
131
132 return kTRUE;
133 //if (fHadDecay) fPythiaHadDecay->init();
134 //return kTRUE;
135}
Pythia8::RndmEngine * fRandomEngine

◆ IsPbrem()

Bool_t DPPythia8Generator::IsPbrem ( )
inline

Definition at line 62 of file DPPythia8Generator.h.

62{ return fpbrem; };

◆ List()

void DPPythia8Generator::List ( int  id)
inline

Definition at line 37 of file DPPythia8Generator.h.

37{fPythia->particleData.list(id);};

◆ MinDPMass()

Double_t DPPythia8Generator::MinDPMass ( )
inline

Definition at line 67 of file DPPythia8Generator.h.

67{ return fDPminM; };

◆ nrOfDP()

Int_t DPPythia8Generator::nrOfDP ( )
inline

Definition at line 74 of file DPPythia8Generator.h.

74{ return fnDPtot; };

◆ nrOfRetries()

Int_t DPPythia8Generator::nrOfRetries ( )
inline

Definition at line 73 of file DPPythia8Generator.h.

73{ return fnRetries; };

◆ Print()

void DPPythia8Generator::Print ( )
inline

Definition at line 36 of file DPPythia8Generator.h.

36{fPythia->settings.listAll(); };

◆ ReadEvent()

Bool_t DPPythia8Generator::ReadEvent ( FairPrimaryGenerator *  cpg)

public method ReadEvent

Definition at line 146 of file DPPythia8Generator.cxx.

147{
148 Double_t tp,tS,zp,xp,yp,zS,xS,yS,pz,px,py,e,w;
149 Double_t tm,zm,xm,ym,pmz,pmx,pmy,em;
150 //Double_t td,zd,xd,yd;
151 Int_t im;
152// take DP decay of Pythia, move it to the SHiP decay region
153// recalculate decay time
154// weight event with exp(-t_ship/tau_DP) / exp(-t_pythia/tau_DP)
155
156 int iDP = 0; // index of the chosen DP, also ensures that at least 1 DP is produced
157 std::vector<int> dec_chain; // pythia indices of the particles to be stored on the stack
158 std::vector<int> dpvec; // pythia indices of DP particles
159 bool hadDecay = false;
160 do {
161
162 if (fextFile && *fextFile){
163 // take charm or beauty hadron from external file
164 // correct for too much Ds produced by pythia6
165 /*bool x = true;
166 while(x){
167 if (fn==fNevents) {fLogger->Warning(MESSAGE_ORIGIN, "End of input file. Rewind.");}
168 fTree->GetEntry(fn%fNevents);
169 fn++;
170 if ( int(fabs(hid[0]) ) != 431){ x = false; }
171 else {
172 Double_t rnr = gRandom->Uniform(0,1);
173 if( rnr<fFDs ) { x = false; };
174 //std::cout<<"what is x "<<x<<" id "<<int(fabs(hid[0]))<<" rnr " << rnr <<" "<< fFDs <<std::endl ;
175 }
176 }
177 fPythia->event.reset();
178 fPythia->event.append( (Int_t)hid[0], 1, 0, 0, hpx[0], hpy[0], hpz[0], hE[0], hM[0], 0., 9. );
179 */
180 }
181 //bit for proton brem
182 if (fpbrem){
183 fPythia->event.reset();
184 double dpmom = 0;
185 double thetain = 0;
186 fpbremPDF->GetRandom2(dpmom, thetain);
187 double dpm = fPythia->particleData.m0(fDP);
188 double dpe = sqrt(dpmom*dpmom+dpm*dpm);
189 double phiin = 2. * M_PI * gRandom->Rndm();
190
191 if ( debug > 1){std::cout << " Adding DP gun with p "
192 << dpmom
193 << " m " << dpm
194 << " e " << dpe
195 << " theta,phi " << thetain << "," << phiin << std::endl << std::flush;}
196 fPythia->event.append( fDP, 1, 0, 0, dpmom * sin(thetain) * cos(phiin), dpmom * sin(thetain) * sin(phiin), dpmom * cos(thetain), dpe, dpm);
197 }
198
199 if (!fPythia->next()) LOG(FATAL) << "fPythia->next() failed";
200
201 //fPythia->event.list();
202 for(int i=0; i<fPythia->event.size(); i++){
203 // find all DP
204 if (abs(fPythia->event[i].id())==fDP){
205 dpvec.push_back( i );
206 }
207 }
208 iDP = dpvec.size();
209 fnDPtot += iDP;
210 if ( iDP == 0 ){
211 //fLogger->Info(MESSAGE_ORIGIN,"Event without DP. Retry.");
212 //fPythia->event.list();
213 fnRetries+=1; // can happen if phasespace does not allow meson to decay to DP
214 }else{
215 //for mesons, could have more than one ... but for DY prod, need to take the last one...
216 //int r = int( gRandom->Uniform(0,iDP) );
217 int r = iDP-1;
218 // std::cout << " ----> debug 2 " << r << std::endl;
219 int i = dpvec[r];
220 // production vertex
221 zp =fPythia->event[i].zProd();
222 xp =fPythia->event[i].xProd();
223 yp =fPythia->event[i].yProd();
224 tp =fPythia->event[i].tProd();
225 // momentum
226 pz =fPythia->event[i].pz();
227 px =fPythia->event[i].px();
228 py =fPythia->event[i].py();
229 e =fPythia->event[i].e();
230 // old decay vertex
231 //Int_t ida =fPythia->event[i].daughter1();
232 std::cout << " Debug: decay product of A: "
233 << fPythia->event[fPythia->event[i].daughter1()].id() << " "
234 << fPythia->event[fPythia->event[i].daughter2()].id()
235 << std::endl;
236 //hadDecay = abs(fPythia->event[ida].id())!=11 && abs(fPythia->event[ida].id())!=13 && abs(fPythia->event[ida].id())!=15;
237
238 //zd =fPythia->event[ida].zProd();
239 //xd =fPythia->event[ida].xProd();
240 //yd =fPythia->event[ida].yProd();
241 //td =fPythia->event[ida].tProd();
242 // new decay vertex
243 Double_t LS = gRandom->Uniform(fLmin,fLmax); // mm, G4 and Pythia8 units
244 Double_t p = TMath::Sqrt(px*px+py*py+pz*pz);
245 Double_t lam = LS/p;
246 xS = xp + lam * px;
247 yS = yp + lam * py;
248 zS = zp + lam * pz;
249 Double_t gam = e/TMath::Sqrt(e*e-p*p);
250 Double_t beta = p/e;
251 tS = tp + LS/beta; // units ? [mm/c] + [mm/beta] (beta is dimensionless speed, and c=1 here)
252 // if one would use [s], then tS = tp/(cm*c_light) + (LS/cm)/(beta*c_light) = tS/(cm*c_light) i.e. units look consisent
253 w = TMath::Exp(-LS/(beta*gam*fctau))*( (fLmax-fLmin)/(beta*gam*fctau) );
254 im = (Int_t)fPythia->event[i].mother1();
255 zm =fPythia->event[im].zProd();
256 xm =fPythia->event[im].xProd();
257 ym =fPythia->event[im].yProd();
258 pmz =fPythia->event[im].pz();
259 pmx =fPythia->event[im].px();
260 pmy =fPythia->event[im].py();
261 em =fPythia->event[im].e();
262 tm =fPythia->event[im].tProd();
263 // foresee finite beam size
264 Double_t dx=0;
265 Double_t dy=0;
266 if (fsmearBeam>0){
267 Double_t test = fsmearBeam*fsmearBeam;
268 Double_t Rsq = test+1.;
269 while(Rsq>test){
270 dx = gRandom->Uniform(-1.,1.) * fsmearBeam;
271 dy = gRandom->Uniform(-1.,1.) * fsmearBeam;
272 Rsq = dx*dx+dy*dy;
273 }
274 }
275 if (fextFile && *fextFile){
276 // take grand mother particle from input file, to know if primary or secondary production
277 //cpg->AddTrack((Int_t)mid[0],mpx[0],mpy[0],mpz[0],xm/cm+dx,ym/cm+dy,zm/cm,-1,false,mE[0],0.,1.);
278 //cpg->AddTrack((Int_t)fPythia->event[im].id(),pmx,pmy,pmz,xm/cm+dx,ym/cm+dy,zm/cm,0,false,em,tm/cm/c_light,w); // convert pythia's (x,y,z[mm], t[mm/c]) to ([cm], [s])
279 //cpg->AddTrack(fDP, px, py, pz, xp/cm+dx,yp/cm+dy,zp/cm, 1,false,e,tp/cm/c_light,w);
280 }else{
281 cpg->AddTrack((Int_t)fPythia->event[im].id(),pmx,pmy,pmz,xm/cm+dx,ym/cm+dy,zm/cm,-1,false,em,tm/cm/c_light,w); // convert pythia's (x,y,z[mm], t[mm/c]) to ([cm], [s])
282 cpg->AddTrack(fDP, px, py, pz, xp/cm+dx,yp/cm+dy,zp/cm, 0,false,e,tp/cm/c_light,w);
283 }
284 // bookkeep the indices of stored particles
285 dec_chain.push_back( im );
286 dec_chain.push_back( i );
287 if (debug>1) std::cout << std::endl << " insert mother id " << im << " pdg=" <<fPythia->event[im].id() << " pmz = " << pmz << " [GeV], zm = " << zm << " [mm] tm = " << tm << " [mm/c]" << std::endl;
288 if (debug>1) std::cout << " ----> insert DP id " << i << " pdg=" << fDP << " pz = " << pz << " [GeV] zp = " << zp << " [mm] tp = " << tp << " [mm/c]" << std::endl;
289 iDP = i;
290 }
291 } while ( iDP == 0 ); // ----------- avoid rare empty events w/o any DP's produced
292
293 if (fShipEventNr%100==0) {
294 LOGF(info, "ship event %i / pythia event-nr %i", fShipEventNr, fn);
295 }
296 fShipEventNr += 1;
297 // fill a container with pythia indices of the DP decay chain
298 if (debug>1) std::cout << "Filling daughter particles" << std::endl;
299 //if (!hadDecay){
300 for(int k=0; k<fPythia->event.size(); k++){
301 // if daughter of DP, copy
302 if (debug>1) std::cout <<k<< " pdg =" <<fPythia->event[k].id() << " mum " << fPythia->event[k].mother1() << std::endl;
303 im =fPythia->event[k].mother1();
304 while (im>0){
305 if ( im == iDP ){break;} // pick the decay products of only 1 chosen DP
306 // if ( abs(fPythia->event[im].id())==fDP && im == iDP ){break;}
307 else {im =fPythia->event[im].mother1();}
308 }
309 if (im < 1) {
310 if (debug>1) std::cout << "reject" << std::endl;
311 continue;
312 }
313 if (debug>1) std::cout << "accept" << std::endl;
314 dec_chain.push_back( k );
315 }
316 //}
317 /*else {
318 //get decay from e+e- collision.....
319 fPythiaHadDecay->settings.parm("Beams:eCM",20);
320 fPythiaHadDecay->next();
321 for (int k=0; k<fPythiaHadDecay->event.size(); k++){
322 fPythia->event.append( fPythiaHadDecay->event[k].id(),fPythiaHadDecay->event[k].status() ,fPythiaHadDecay->event[k].mother1() , fPythiaHadDecay->event[k].mother2(), fPythiaHadDecay->event[k].daughter1(), fPythiaHadDecay->event[k].daughter2(), fPythiaHadDecay->event[k].col(), fPythiaHadDecay->event[k].acol(), fPythiaHadDecay->event[k].px(), fPythiaHadDecay->event[k].py(), fPythiaHadDecay->event[k].pz(), fPythiaHadDecay->event[k].e(), fPythiaHadDecay->event[k].m(), 0., 9. );
323 dec_chain.push_back( fPythia->event.size()-1 );
324 std::cout << " Adding decay product: " << k << " "
325 << fPythiaHadDecay->event[k].id() << " "
326 << fPythiaHadDecay->event[k].status() << " "
327 << fPythiaHadDecay->event[k].mother1() << " "
328 << fPythiaHadDecay->event[k].mother2() << " "
329 << std::endl;
330 }
331 }*/
332
333 // go over daughters and store them on the stack, starting from 2 to account for DP and its mother
334 for(std::vector<int>::iterator it = dec_chain.begin() + 2; it != dec_chain.end(); ++it){
335 // pythia index of the paricle to store
336 int k = *it;
337 // find mother position on the output stack: impy -> im
338 int impy =fPythia->event[k].mother1();
339 std::vector<int>::iterator itm = std::find( dec_chain.begin(), dec_chain.end(), impy);
340 im =-1; // safety
341 if ( itm != dec_chain.end() )
342 im = itm - dec_chain.begin(); // convert iterator into sequence number
343
344 Bool_t wanttracking=false;
345 if(fPythia->event[k].isFinal()){ wanttracking=true;}
346 pz =fPythia->event[k].pz();
347 px =fPythia->event[k].px();
348 py =fPythia->event[k].py();
349 e =fPythia->event[k].e();
350 if (fextFile && *fextFile){im+=1;};
351 cpg->AddTrack((Int_t)fPythia->event[k].id(),px,py,pz,xS/cm,yS/cm,zS/cm,im,wanttracking,e,tS/cm/c_light,w);
352 // std::cout <<k<< " insert pdg =" <<fPythia->event[k].id() << " pz = " << pz << " [GeV] zS = " << zS << " [mm] tS = " << tS << "[mm/c]" << std::endl;
353 }
354 return kTRUE;
355}
const Double_t c_light
int i
Definition ShipAna.py:86

◆ SetDPId()

void DPPythia8Generator::SetDPId ( Int_t  id)
inline

Definition at line 49 of file DPPythia8Generator.h.

◆ SetDY()

void DPPythia8Generator::SetDY ( )
inline

Definition at line 63 of file DPPythia8Generator.h.

63 {
64 fdy = kTRUE;
65 };

◆ SetfFDs()

void DPPythia8Generator::SetfFDs ( Double_t  z)
inline

Definition at line 54 of file DPPythia8Generator.h.

◆ SetId()

void DPPythia8Generator::SetId ( Double_t  id)
inline

Definition at line 48 of file DPPythia8Generator.h.

48{ fId = id; };

◆ SetLmax()

void DPPythia8Generator::SetLmax ( Double_t  z)
inline

Definition at line 52 of file DPPythia8Generator.h.

52{ fLmax = z*10; };

◆ SetLmin()

void DPPythia8Generator::SetLmin ( Double_t  z)
inline

Definition at line 51 of file DPPythia8Generator.h.

51{ fLmin = z*10; };

◆ SetMinDPMass()

void DPPythia8Generator::SetMinDPMass ( Double_t  m)
inline

Definition at line 68 of file DPPythia8Generator.h.

68 {
69 fDPminM = m;
70 };
Double_t m

◆ SetMom()

void DPPythia8Generator::SetMom ( Double_t  mom)
inline

Definition at line 46 of file DPPythia8Generator.h.

46{ fMom = mom; };

◆ SetParameters()

void DPPythia8Generator::SetParameters ( char *  par)

Definition at line 357 of file DPPythia8Generator.cxx.

358{
359 // Set Parameters
360 fPythia->readString(par);
361 if ( debug ){std::cout<<"fPythia->readString(\""<<par<<"\")"<<std::endl;}
362}

◆ SetPbrem()

void DPPythia8Generator::SetPbrem ( TH2F *  pdf)
inline

Definition at line 58 of file DPPythia8Generator.h.

58 {
59 fpbrem = kTRUE;
60 fpbremPDF = pdf;
61 };

◆ SetSmearBeam()

void DPPythia8Generator::SetSmearBeam ( Double_t  sb)
inline

Definition at line 53 of file DPPythia8Generator.h.

53{ fsmearBeam = sb; };

◆ UseDeepCopy()

void DPPythia8Generator::UseDeepCopy ( )
inline

Definition at line 72 of file DPPythia8Generator.h.

72{ fDeepCopy = kTRUE; };

◆ UseExternalFile()

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

◆ UseRandom1()

void DPPythia8Generator::UseRandom1 ( )
inline

Definition at line 55 of file DPPythia8Generator.h.

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

◆ UseRandom3()

void DPPythia8Generator::UseRandom3 ( )
inline

Definition at line 56 of file DPPythia8Generator.h.

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

Member Data Documentation

◆ fctau

Double_t DPPythia8Generator::fctau
protected

Definition at line 98 of file DPPythia8Generator.h.

◆ fDeepCopy

Bool_t DPPythia8Generator::fDeepCopy
protected

Definition at line 106 of file DPPythia8Generator.h.

◆ fDP

Int_t DPPythia8Generator::fDP
protected

Definition at line 86 of file DPPythia8Generator.h.

◆ fDPminM

Double_t DPPythia8Generator::fDPminM
protected

Definition at line 93 of file DPPythia8Generator.h.

◆ fdy

Bool_t DPPythia8Generator::fdy
protected

Definition at line 92 of file DPPythia8Generator.h.

◆ fextFile

const char* DPPythia8Generator::fextFile
protected

Definition at line 101 of file DPPythia8Generator.h.

◆ fFDs

Double_t DPPythia8Generator::fFDs
protected

Definition at line 99 of file DPPythia8Generator.h.

◆ fId

Int_t DPPythia8Generator::fId
protected

Definition at line 87 of file DPPythia8Generator.h.

◆ fInputFile

TFile* DPPythia8Generator::fInputFile
protected

Definition at line 102 of file DPPythia8Generator.h.

◆ firstEvent

Int_t DPPythia8Generator::firstEvent
protected

Definition at line 104 of file DPPythia8Generator.h.

◆ fLmax

Double_t DPPythia8Generator::fLmax
protected

Definition at line 95 of file DPPythia8Generator.h.

◆ fLmin

Double_t DPPythia8Generator::fLmin
protected

Definition at line 94 of file DPPythia8Generator.h.

◆ fLogger

FairLogger* DPPythia8Generator::fLogger
protected

Definition at line 107 of file DPPythia8Generator.h.

◆ fMom

Double_t DPPythia8Generator::fMom
protected

Definition at line 85 of file DPPythia8Generator.h.

◆ fn

Int_t DPPythia8Generator::fn
protected

Definition at line 104 of file DPPythia8Generator.h.

◆ fnDPtot

Int_t DPPythia8Generator::fnDPtot
protected

Definition at line 97 of file DPPythia8Generator.h.

◆ fNevents

Int_t DPPythia8Generator::fNevents
protected

Definition at line 104 of file DPPythia8Generator.h.

◆ fnRetries

Int_t DPPythia8Generator::fnRetries
protected

Definition at line 96 of file DPPythia8Generator.h.

◆ fpbrem

Bool_t DPPythia8Generator::fpbrem
protected

Definition at line 90 of file DPPythia8Generator.h.

◆ fpbremPDF

TH2F* DPPythia8Generator::fpbremPDF
protected

Definition at line 91 of file DPPythia8Generator.h.

◆ fPythia

Pythia8::Pythia* DPPythia8Generator::fPythia

Definition at line 76 of file DPPythia8Generator.h.

◆ fRandomEngine

Pythia8::RndmEngine* DPPythia8Generator::fRandomEngine
private

Definition at line 80 of file DPPythia8Generator.h.

◆ fShipEventNr

Int_t DPPythia8Generator::fShipEventNr
protected

Definition at line 104 of file DPPythia8Generator.h.

◆ fsmearBeam

Double_t DPPythia8Generator::fsmearBeam
protected

Definition at line 100 of file DPPythia8Generator.h.

◆ fTree

TTree* DPPythia8Generator::fTree
protected

pointer to a file

Definition at line 103 of file DPPythia8Generator.h.

◆ fUseRandom1

Bool_t DPPythia8Generator::fUseRandom1
protected

Definition at line 88 of file DPPythia8Generator.h.

◆ fUseRandom3

Bool_t DPPythia8Generator::fUseRandom3
protected

Definition at line 89 of file DPPythia8Generator.h.

◆ hE

Float_t DPPythia8Generator::hE[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ hid

Float_t DPPythia8Generator::hid[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ hM

Float_t DPPythia8Generator::hM[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ hpx

Float_t DPPythia8Generator::hpx[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ hpy

Float_t DPPythia8Generator::hpy[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ hpz

Float_t DPPythia8Generator::hpz[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ mE

Float_t DPPythia8Generator::mE[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ mid

Float_t DPPythia8Generator::mid[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ mpx

Float_t DPPythia8Generator::mpx[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ mpy

Float_t DPPythia8Generator::mpy[1]
protected

Definition at line 105 of file DPPythia8Generator.h.

◆ mpz

Float_t DPPythia8Generator::mpz[1]
protected

Definition at line 105 of file DPPythia8Generator.h.


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