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

#include <FixedTargetGenerator.h>

Inheritance diagram for FixedTargetGenerator:
Collaboration diagram for FixedTargetGenerator:

Public Member Functions

 FixedTargetGenerator ()
 
virtual ~FixedTargetGenerator ()
 
Bool_t ReadEvent (FairPrimaryGenerator *)
 
void SetParameters (char *)
 
void Print ()
 
virtual Bool_t Init ()
 
Bool_t InitForCharmOrBeauty (TString fInName, Int_t nev, Double_t npots=5E13, Int_t nStart=0)
 
void SetMom (Double_t mom)
 
void UseRandom1 ()
 
void UseRandom3 ()
 
void SetCharmTarget (Bool_t charmtarget=true)
 
void SetTarget (TString s, Double_t x, Double_t y)
 
void SetBoost (Double_t f)
 
void SetG4only ()
 
void SetTauOnly ()
 
void SetJpsiMainly ()
 
void SetOnlyMuons ()
 
void SetDrellYan ()
 
void SetPhotonCollision ()
 
void WithEvtGen ()
 
void SetChibb (Double_t x)
 
void SetChicc (Double_t x)
 
void SetSeed (Double_t seed)
 
void SetHeartBeat (Int_t x)
 
void SetEnergyCut (Float_t emax)
 
void SetDebug (Bool_t x)
 
void SetOpt4DP (TNtuple *t)
 
Double_t GetPotForCharm ()
 
Pythia8::Pythia * GetPythia ()
 
Pythia8::Pythia * GetPythiaN ()
 

Protected Member Functions

 ClassDef (FixedTargetGenerator, 2)
 

Protected Attributes

Double_t fMom
 
Bool_t fUseRandom1
 
Bool_t fUseRandom3
 
Double_t fSeed
 
Double_t EMax
 
Double_t fBoost
 
Double_t chicc
 
Double_t chibb
 
Double_t wspill
 
Double_t nrpotspill
 
Int_t nEvents
 
Int_t nEntry
 
Int_t pot
 
Int_t nDsprim
 
Int_t ntotprim
 
Bool_t tauOnly
 
Bool_t JpsiMainly
 
Bool_t DrellYan
 
Bool_t PhotonCollision
 
Bool_t G4only
 
Bool_t setByHand
 
Bool_t Debug
 
Bool_t withEvtGen
 
Bool_t OnlyMuons
 
Bool_t fcharmtarget
 
FairLogger * fLogger
 
Pythia8::Pythia * fPythiaN
 don't make it persistent, magic ROOT command
 
Pythia8::Pythia * fPythiaP
 
EvtGenDecays * evtgenN
 
EvtGenDecays * evtgenP
 
GenieGeneratorfMaterialInvestigator
 
Bool_t withNtuple
 
TNtuple * fNtuple
 special option for Dark Photon physics studies
 
TString targetName
 
TString Option
 
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
 
TFile * fin
 
TNtuple * nTree
 
Float_t n_id
 
Float_t n_px
 
Float_t n_py
 
Float_t n_pz
 
Float_t n_M
 
Float_t n_E
 
Float_t n_mpx
 
Float_t n_mpy
 
Float_t n_mpz
 
Float_t n_mE
 
Float_t n_mid
 
Float_t ck
 
Int_t heartbeat
 

Private Attributes

Pythia8::RndmEngine * fRandomEngine
 

Detailed Description

Definition at line 15 of file FixedTargetGenerator.h.

Constructor & Destructor Documentation

◆ FixedTargetGenerator()

FixedTargetGenerator::FixedTargetGenerator ( )

default constructor

Definition at line 23 of file FixedTargetGenerator.cxx.

24{
25 fUseRandom1 = kFALSE;
26 fUseRandom3 = kTRUE;
27 fSeed = 0;
28 fMom = 400; // proton
29 targetName = "";
30 fcharmtarget = false;
31 xOff=0; yOff=0;
32 tauOnly = false;
33 JpsiMainly = false;
34 DrellYan = false;
35 PhotonCollision = false;
36 G4only = false;
37 OnlyMuons = false;
38 maxCrossSection = 0.;
39 EMax = 0;
40 fBoost = 1.;
41 withEvtGen = kFALSE;
42 withNtuple = kFALSE;
43 chicc=1.7e-3; //prob to produce primary ccbar pair/pot
44 chibb=1.6e-7; //prob to produce primary bbbar pair/pot
45 nrpotspill = 5.E13; // nr of protons / spill
46 setByHand = kFALSE;
47 Debug = kFALSE;
48 Option = "Primary";
49 wspill = 1.; // event weight == 1 for primary events
50 heartbeat = 1000;
51}

◆ ~FixedTargetGenerator()

FixedTargetGenerator::~FixedTargetGenerator ( )
virtual

destructor

Definition at line 296 of file FixedTargetGenerator.cxx.

297{
298}

Member Function Documentation

◆ ClassDef()

FixedTargetGenerator::ClassDef ( FixedTargetGenerator  ,
 
)
protected

◆ GetPotForCharm()

Double_t FixedTargetGenerator::GetPotForCharm ( )
inline

Definition at line 53 of file FixedTargetGenerator.h.

53{return nrpotspill/wspill;}

◆ GetPythia()

Pythia8::Pythia * FixedTargetGenerator::GetPythia ( )
inline

Definition at line 54 of file FixedTargetGenerator.h.

54{return fPythiaP;}
Pythia8::Pythia * fPythiaP

◆ GetPythiaN()

Pythia8::Pythia * FixedTargetGenerator::GetPythiaN ( )
inline

Definition at line 55 of file FixedTargetGenerator.h.

55{return fPythiaN;}
Pythia8::Pythia * fPythiaN
don't make it persistent, magic ROOT command

◆ Init()

Bool_t FixedTargetGenerator::Init ( )
virtual

Definition at line 99 of file FixedTargetGenerator.cxx.

100{
101 fPythiaP = new Pythia8::Pythia(); // pythia needed also for G4only, to copy 2mu BRs
102 if (Option == "Primary" && !G4only){
103 fPythiaN = new Pythia8::Pythia();
104 }else if (Option != "charm" && Option != "beauty" && !G4only) {
105 LOG(INFO) << "FixedTargetGenerator::Init Option not known "<< Option.Data() <<", abort";
106 }
109 std::vector<int> r = { 221, 221, 223, 223, 113, 331, 333};
110 std::vector<int> c = { 6 , 7, 5 , 7, 5, 6, 9}; // decay channel mumu mumuX
111
112 if (Option == "Primary" && !G4only){
113 fPythiaP->settings.mode("Beams:idB", 2212);
114 fPythiaN->settings.mode("Beams:idB", 2112);
115 }else{
116 fPythiaP->readString("ProcessLevel:all = off");
117 }
118 std::array<Pythia8::Pythia*,2> plist = {fPythiaP,fPythiaN};
119 Int_t pcount = 0;
120 for(const auto& fPythia : plist) {
121 if (pcount > 0 && (Option != "Primary" || G4only)){continue;}
122 pcount+=1;
123 fPythia->setRndmEnginePtr(fRandomEngine);
124 fPythia->settings.mode("Random:seed",fSeed);
125 fPythia->settings.mode("Next:numberCount",heartbeat);
126 if (Option == "Primary"){
127 fPythia->settings.mode("Beams:idA", 2212);
128 fPythia->settings.mode("Beams:frameType", 2);
129 fPythia->settings.parm("Beams:eA",fMom);
130 fPythia->settings.parm("Beams:eB",0.);
131 }
132 if (JpsiMainly){
133// use this for all onia productions
134 fPythia->readString("Onia:all(3S1) = on");
135 fPythia->readString("443:new J/psi J/psi 3 0 0 3.09692 0.00009 3.09602 3.09782 0. 1 1 0 1 0");
136 fPythia->readString("443:addChannel = 1 1. 0 -13 13");
137 fPythia->readString("553:new Upsilon Upsilon 3 0 0 9.46030 0.00005 9.45980 9.46080 0.00000e+00 0 1 0 1 0");
138 fPythia->readString("553:addChannel = 1 1. 0 -13 13");
139 }
140 if (DrellYan){
141 fPythia->readString("WeakSingleBoson:ffbar2gmZ = on");
142// Z0/gamma -> mu+ mu-
143 fPythia->readString("23:onMode = off");
144 fPythia->readString("23:onIfAll = 13 13");
145 fPythia->readString("23:mMin = 0.5");
146 fPythia->readString("PhaseSpace:mHatMin = 0.5");
147 fPythia->readString("PartonLevel:FSR = on");
148 fPythia->readString("PDF:pSet = 4");
149 }
150 if (PhotonCollision){
151 fPythia->readString("PhotonCollision:gmgm2mumu = on");
152 }
153 if (!DrellYan && !PhotonCollision){
154 fPythia->readString("SoftQCD:inelastic = on");
155 fPythia->readString("PhotonCollision:gmgm2mumu = on");
156 fPythia->readString("PromptPhoton:all = on");
157 fPythia->readString("WeakBosonExchange:all = on");
158 fPythia->readString("WeakSingleBoson:all = on");
159 }
160 if (tauOnly){
161 fPythia->readString("431:new D_s+ D_s- 1 3 0 1.96849 0.00000 0.00000 0.00000 1.49900e-01 0 1 0 1 0");
162 fPythia->readString("431:addChannel = 1 0.0640000 0 -15 16");
163 }
164
165// find all long lived particles in pythia
166 Int_t n = 1;
167 while(n!=0){
168 n = fPythia->particleData.nextId(n);
169#ifndef PYTHIA8_V
170 ParticleDataEntry* p = fPythia->particleData.particleDataEntryPtr(n);
171#else
172#if PYTHIA8_V<8309
173 ParticleDataEntry* p = fPythia->particleData.particleDataEntryPtr(n);
174#else
175 Pythia8::ParticleDataEntryPtr p = fPythia->particleData.particleDataEntryPtr(n);
176#endif
177#endif
178 if (p->tau0()>1){
179 std::string particle = std::to_string(n)+":mayDecay = false";
180 fPythia->readString(particle);
181 LOG(INFO) << "FixedTargetGenerator::Init Made "<< p->name().c_str() <<" stable for Pythia, should decay in Geant4";
182 }
183 }
184// boost branching fraction of rare di-muon decays
185// eta omega rho0 eta' phi
186 if (fBoost != 1.){
187 LOG(INFO) << "FixedTargetGenerator::Init Rescale BRs of dimuon decays in Pythia: "<<fBoost;
188 for (unsigned int i=0; i<r.size(); ++i) {
189#ifndef PYTHIA8_V
190 Pythia8::ParticleDataEntry* V = fPythia->particleData.particleDataEntryPtr(r[i]);
191#else
192#if PYTHIA8_V<8309
193 Pythia8::ParticleDataEntry* V = fPythia->particleData.particleDataEntryPtr(r[i]);
194#else
195 Pythia8::ParticleDataEntryPtr V = fPythia->particleData.particleDataEntryPtr(r[i]);
196#endif
197#endif
198 Pythia8::DecayChannel ch = V->channel(c[i]);
199 if (TMath::Abs(ch.product(0))!=13 || TMath::Abs(ch.product(1))!=13){
200 LOG(INFO) << "FixedTargetGenerator::Init this is not the right decay channel:"<<r[i] << " "<<c[i];
201 }else{
202 TString tmp="";
203 tmp+=r[i];tmp+=":";tmp+= c[i];
204 tmp+=":bRatio =";
205 tmp+=fBoost*ch.bRatio();
206 fPythia->readString(tmp.Data());
207 }
208 }
209 }
210 fPythia->init();
211 }
212 // Initialize EvtGen.
213 if (withEvtGen){
214 TString SIMPATH=getenv("SIMPATH");
215 TString DecayFile =SIMPATH+"/share/EvtGen/DECAY.DEC";
216 TString ParticleFile =SIMPATH + "/share/EvtGen/evt.pdl";
217 if(SIMPATH == "")
218 {
219 std::cout << "Using $EVTGENDATA "<< getenv("EVTGENDATA") << std::endl;
220 DecayFile =TString(getenv("EVTGENDATA"))+"/DECAY.DEC";
221 ParticleFile =TString(getenv("EVTGENDATA"))+ "/evt.pdl";
222 }
223 EvtAbsRadCorr *fsrPtrIn = 0;
224 EvtExternalGenList *extPtr = new EvtExternalGenList();
225 std::list<EvtDecayBase*> models = extPtr->getListOfModels();
226 // Define the random number generator
227 EvtRandomEngine* eng = new EvtSimpleRandomEngine();
228 EvtRandom::setRandomEngine(eng);
229 EvtGen *myEvtGenPtr = new EvtGen(DecayFile.Data(), ParticleFile.Data(),eng, fsrPtrIn, &models, 1, false);
230 TString UdecayFile = getenv("FAIRSHIP");UdecayFile +="/gconfig/USERDECAY.DEC";
231 evtgenP = new EvtGenDecays(fPythiaP, DecayFile.Data(), ParticleFile.Data(),myEvtGenPtr);
232 evtgenP->readDecayFile(UdecayFile.Data()); // will make update of EvtGen with user decay file
233 // use one instance of EvtGen, requires patch to Pythia8Plugins/EvtGen.h
234 if (Option == "Primary"){
235 evtgenN = new EvtGenDecays(fPythiaN, DecayFile.Data(), ParticleFile.Data(),myEvtGenPtr);
236 }
237 }
238 if (targetName!=""){
240 if (!fcharmtarget){
241 TGeoNavigator* nav = gGeoManager->GetCurrentNavigator();
242 nav->cd(targetName);
243 TGeoNode* target = nav->GetCurrentNode();
244 TObjArray* nodes = target->GetVolume()->GetNodes();
245 TGeoNode* first = (TGeoNode*)nodes->At(0);
246 Int_t ilast = nodes->GetSize()-5; // assumes that the last 5 nodes are for the shielding around the target.
247 TGeoNode* last = (TGeoNode*)nodes->At(ilast);
248 nav->cd(targetName+"/"+first->GetName());
249 TGeoBBox* sha = (TGeoBBox*)first->GetVolume()->GetShape();
250 Double_t dz = sha->GetDZ();
251 Double_t origin[3] = {0,0,-dz};
252 Double_t master[3] = {0,0,0};
253 nav->LocalToMaster(origin,master);
254 startZ = master[2];
255 nav->cd(targetName+"/"+last->GetName());
256 sha = (TGeoBBox*)first->GetVolume()->GetShape();
257 dz = sha->GetDZ();
258 origin[2] = +dz;
259 nav->LocalToMaster(origin,master);
260 endZ = master[2];
261 start[0]=xOff;
262 start[1]=yOff;
263 start[2]=startZ;
264 end[0]=xOff;
265 end[1]=yOff;
266 end[2]=endZ;
267 }
268 else{ //charm geometry uses a different target
269 TGeoVolume* top = gGeoManager->GetTopVolume();
270 TGeoNode* target = top->FindNode(targetName);
271 if (!target){
272 LOG(ERROR) << "FixedTargetGenerator::Init target not found, "<< targetName.Data() <<", program will crash";
273 }
274 Double_t z_middle = target->GetMatrix()->GetTranslation()[2];
275 TGeoBBox* sha = (TGeoBBox*)target->GetVolume()->GetShape();
276 startZ = z_middle - sha->GetDZ();
277 endZ = z_middle + sha->GetDZ();
278 start[0]=xOff;
279 start[1]=yOff;
280 start[2]=startZ;
281 end[0]=xOff;
282 end[1]=yOff;
283 end[2]=endZ;
284 }
285//find maximum interaction length
288 }
289
290 return kTRUE;
291}
GenieGenerator * fMaterialInvestigator
Pythia8::RndmEngine * fRandomEngine
Double_t MeanMaterialBudget(const Double_t *start, const Double_t *end, Double_t *mparam)
int i
Definition ShipAna.py:86
origin(sTree, it)
c
Definition hnl.py:100

◆ InitForCharmOrBeauty()

Bool_t FixedTargetGenerator::InitForCharmOrBeauty ( TString  fInName,
Int_t  nev,
Double_t  npots = 5E13,
Int_t  nStart = 0 
)

Definition at line 52 of file FixedTargetGenerator.cxx.

53{
54 Option = "charm";
55 nEntry = nStart;
56 // open input file with charm or beauty
57 fin = TFile::Open(fInName);
58 nTree = (TNtuple*)fin->FindObjectAny("pythia6"); // old format, simple ntuple
59 nEvents = nTree->GetEntries();
60 nTree->SetBranchAddress("id",&n_id);
61 nTree->SetBranchAddress("px",&n_px);
62 nTree->SetBranchAddress("py",&n_py);
63 nTree->SetBranchAddress("pz",&n_pz);
64 nTree->SetBranchAddress("M",&n_M);
65 nTree->SetBranchAddress("E",&n_E);
66 nTree->SetBranchAddress("mid",&n_mid);
67 nTree->SetBranchAddress("mpx",&n_mpx);
68 nTree->SetBranchAddress("mpy",&n_mpy);
69 nTree->SetBranchAddress("mpz",&n_mpz);
70 nTree->SetBranchAddress("mE",&n_mE);
71 if (nTree->GetBranch("k")){
72 LOG(INFO) << "FixedTargetGenerator::InitForCharmOrBeauty: +++has branch+++ ";
73 nTree->SetBranchAddress("k",&ck);}
74// check if we deal with charm or beauty:
75 nTree->GetEvent(0);
76 if (!setByHand and n_M>5){
77 chicc = chibb;
78 LOG(INFO) << "FixedTargetGenerator::InitForCharmOrBeauty: automatic detection of beauty, configured for beauty";
79 LOG(INFO) << "FixedTargetGenerator::InitForCharmOrBeauty: bb cross section / mbias"<<chicc;
80 }else{
81 LOG(INFO) << "FixedTargetGenerator::InitForCharmOrBeauty: cc cross section / mbias"<<chicc;
82 }
83// convert pot to weight corresponding to one spill of 5e13 pot
84 // get histogram with number of pot to normalise
85 // pot are counted double, i.e. for each signal, i.e. pot/2.
86 Int_t nrcpot=((TH1F*)fin->Get("2"))->GetBinContent(1)/2.; // number of primary interactions
88 LOG(INFO) << "FixedTargetGenerator::InitForCharmOrBeauty: Input file:"<<fInName.Data()<<" with "<<nEvents<< " entries, corresponding to nr-pot="<<nrcpot/chicc;
89 LOG(INFO) << "FixedTargetGenerator::InitForCharmOrBeauty: weight "<< wspill <<", corresponding to "<<nrpotspill<<" p.o.t. per spill for "<< nev <<" events to process";
90
91 pot=0.;
92 //Determine fDs on this file for primaries
93 nDsprim=0;
94 ntotprim=0;
95
96 return kTRUE;
97}
dict nrcpot
Definition makeDecay.py:56

◆ Print()

void FixedTargetGenerator::Print ( )

◆ ReadEvent()

Bool_t FixedTargetGenerator::ReadEvent ( FairPrimaryGenerator *  cpg)

public method ReadEvent

Definition at line 302 of file FixedTargetGenerator.cxx.

303{
304 Double_t zinter=0;
305 Double_t ZoverA = 1.;
306 if (targetName.Data() !=""){
307// calculate primary proton interaction point:
308// loop over trajectory between start and end to pick an interaction point, copied from GenieGenerator and adapted to hadrons
309 Double_t prob2int = -1.;
310 Double_t rndm = 0.;
311 Double_t sigma;
312 Int_t count=0;
313 Double_t zinterStart = start[2];
314 if (Option == "charm" || Option == "beauty"){
315// simulate more downstream interaction points for interactions down in the cascade
316 if (!(nTree->GetBranch("k"))){ck=1;}
317 }
318 else {ck=1;}
319 while (ck>0.5){
320 while (prob2int<rndm) {
321 //place x,y,z uniform along path
322 zinter = gRandom->Uniform(zinterStart,end[2]);
323 Double_t point[3]={xOff,yOff,zinter};
325 Double_t interLength = mparam[8];
326 TGeoNode *node = gGeoManager->FindNode(point[0],point[1],point[2]);
327 TGeoMaterial *mat = 0;
328 if (node && !gGeoManager->IsOutside()) {
329 mat = node->GetVolume()->GetMaterial();
330 Double_t n = mat->GetDensity()/mat->GetA();
331 ZoverA = mat->GetZ() / mat->GetA();
332 sigma = 1./(n*mat->GetIntLen())/mbarn;
333 prob2int = TMath::Exp(-interLength)*sigma/maxCrossSection;
334 }else{
335 prob2int=0.;
336 }
337 rndm = gRandom->Uniform(0.,1.);
338 count+=1;
339 }
340 zinterStart = zinter;
341 ck-=1;
342 }
343 zinter = zinter*cm;
344 }
345 Pythia8::Pythia* fPythia;
346 if (G4only){
347 cpg->AddTrack(2212,0.,0.,fMom,xOff/cm,yOff/cm,start[2],-1,kTRUE,-1.,0.,1.);
348 return kTRUE;
349 }else if (Option == "Primary"){
350 if (gRandom->Uniform(0.,1.) < ZoverA ){
351 fPythiaP->next();
352 if (withEvtGen){evtgenP->decay();}
353 fPythia = fPythiaP;
354 }else{
355 fPythiaN->next();
356 if (withEvtGen){evtgenN->decay();}
357 fPythia = fPythiaN;
358 }
359 }else{
360 if (nEntry==nEvents){
361 LOG(INFO) << "FixedTargetGenerator::readEvent Rewind input file: "<<nEntry;
362 nEntry=0;}
363 nTree->GetEvent(nEntry);
364 nEntry+=1;
365 // sanity check, count number of p.o.t. on input file.
366 Double_t pt=TMath::Sqrt( (n_mpx*n_mpx)+(n_mpy*n_mpy));
367 // every event appears twice, i.e.
368 if (pt<1.e-5 && n_mid==2212){
369 pot+=+0.5;
370 ntotprim+=1;}
371 Int_t idabs=int(TMath::Abs(n_id));
372 if (idabs==431){ nDsprim+=1;}
373 fPythiaP->event.reset();
374 Int_t nID1 = n_id;
375 fPythiaP->event.append(int(n_id),1,0,0,n_px,n_py,n_pz,n_E,n_M,0.,9.);
376 TMCProcess procID = kPTransportation;
377 if (n_mid==2212 && (n_mpx*n_mpx+n_mpy*n_mpy)<1E-5) {procID = kPPrimary;} // probably primary and not from cascade
378 cpg->AddTrack(int(n_mid),n_mpx,n_mpy,n_mpz, xOff/cm,yOff/cm,zinter/cm,-1,kFALSE,n_mE,0.,wspill,procID);
379// second charm hadron in the event
380 nTree->GetEvent(nEntry);
381 if (nID1 * n_id > 0){LOG(INFO) << "FixedTargetGenerator::readEvent same sign charm: "<<nEntry<<" "<< nID1<<" "<< n_id;}
382 nEntry+=1;
383 fPythiaP->event.append(int(n_id),1,0,0,n_px,n_py,n_pz,n_E,n_M,0.,9.);
384 fPythiaP->next();
385 fPythia = fPythiaP;
386 }
387 if (withEvtGen){
388 fPythia->moreDecays();} // let the very short lived resonances decay via Pythia8
389 if(Debug && !G4only){fPythia->event.list();}
390 TMCProcess procID;
391 for(Int_t ii=1; ii<fPythia->event.size(); ii++){
392 Double_t e = fPythia->event[ii].e();
393 Double_t m = fPythia->event[ii].m();
394 Double_t pz = fPythia->event[ii].pz();
395 Int_t id = fPythia->event[ii].id();
396 Bool_t wanttracking=kTRUE;
397 if (e-m<EMax || !fPythia->event[ii].isFinal() || pz<0) {wanttracking=kFALSE;}
399// don't track underlying event
400 if (fabs(id)!=13){wanttracking=kFALSE;}
401 }
402 Double_t z = fPythia->event[ii].zProd()+zinter;
403 Double_t x = fPythia->event[ii].xProd()+xOff;
404 Double_t y = fPythia->event[ii].yProd()+yOff;
405 Double_t tof = fPythia->event[ii].tProd() / (10*c_light) ; // to go from mm to s
406 Double_t px = fPythia->event[ii].px();
407 Double_t py = fPythia->event[ii].py();
408 Int_t im = fPythia->event[ii].mother1()-1;
409 procID = kPPrimary;
410 if (Option != "Primary"){
411 procID = kPDecay;
412 im+=1;
413 if (ii==1) {procID = kPHadronic;}
414 }else{
415 if (ii<3){im=-1;}
416 }
417 cpg->AddTrack(id,px,py,pz,x/cm,y/cm,z/cm,im,wanttracking,e,tof,wspill,procID);
418 if(withNtuple){
419 int idabs = TMath::Abs(id);
420 if (idabs<10) {continue;}
421 if (idabs<18 || idabs==22 || idabs==111 || idabs==221 || idabs==223 || idabs==331
422 || idabs==211 || idabs==113 || idabs==333 || idabs==321 || idabs==2212 ){
423 Int_t moID = fPythia->event[im+1].id();
424 fNtuple->Fill(0,id,px,py,pz,e,x/cm,y/cm,z/cm,moID);
425 }
426 }
427 }
428 return kTRUE;
429}
const Double_t c_light
const Double_t mbarn
Double_t cm
Double_t m
TNtuple * fNtuple
special option for Dark Photon physics studies

◆ SetBoost()

void FixedTargetGenerator::SetBoost ( Double_t  f)
inline

Definition at line 38 of file FixedTargetGenerator.h.

38{ fBoost = f; } // boost factor for rare di-muon decays

◆ SetCharmTarget()

void FixedTargetGenerator::SetCharmTarget ( Bool_t  charmtarget = true)
inline

Definition at line 36 of file FixedTargetGenerator.h.

36{fcharmtarget = charmtarget;}; //charm geometry uses a different target, default one is usual

◆ SetChibb()

void FixedTargetGenerator::SetChibb ( Double_t  x)
inline

Definition at line 46 of file FixedTargetGenerator.h.

46{ chibb = x; } // chibb = bbbar over mbias cross section

◆ SetChicc()

void FixedTargetGenerator::SetChicc ( Double_t  x)
inline

Definition at line 47 of file FixedTargetGenerator.h.

47{ chicc = x; } // chicc = ccbar over mbias cross section

◆ SetDebug()

void FixedTargetGenerator::SetDebug ( Bool_t  x)
inline

Definition at line 51 of file FixedTargetGenerator.h.

51{Debug=x;}

◆ SetDrellYan()

void FixedTargetGenerator::SetDrellYan ( )
inline

Definition at line 43 of file FixedTargetGenerator.h.

43{ DrellYan = true; } // only generate prompt Z0* processes

◆ SetEnergyCut()

void FixedTargetGenerator::SetEnergyCut ( Float_t  emax)
inline

Definition at line 50 of file FixedTargetGenerator.h.

50{EMax=emax;}// min energy to be copied to Geant4

◆ SetG4only()

void FixedTargetGenerator::SetG4only ( )
inline

Definition at line 39 of file FixedTargetGenerator.h.

39{ G4only = true; } // only run Geant4, no pythia primary interaction

◆ SetHeartBeat()

void FixedTargetGenerator::SetHeartBeat ( Int_t  x)
inline

Definition at line 49 of file FixedTargetGenerator.h.

49{heartbeat=x;}

◆ SetJpsiMainly()

void FixedTargetGenerator::SetJpsiMainly ( )
inline

Definition at line 41 of file FixedTargetGenerator.h.

41{ JpsiMainly = true; } // let all Jpsi decay to mumu

◆ SetMom()

void FixedTargetGenerator::SetMom ( Double_t  mom)
inline

Definition at line 33 of file FixedTargetGenerator.h.

33{ fMom = mom; };

◆ SetOnlyMuons()

void FixedTargetGenerator::SetOnlyMuons ( )
inline

Definition at line 42 of file FixedTargetGenerator.h.

42{ OnlyMuons = true; } // only transport muons

◆ SetOpt4DP()

void FixedTargetGenerator::SetOpt4DP ( TNtuple *  t)
inline

Definition at line 52 of file FixedTargetGenerator.h.

52{withNtuple=kTRUE; fNtuple = t ; }

◆ SetParameters()

void FixedTargetGenerator::SetParameters ( char *  )

◆ SetPhotonCollision()

void FixedTargetGenerator::SetPhotonCollision ( )
inline

Definition at line 44 of file FixedTargetGenerator.h.

44{ PhotonCollision = true; } // only generate prompt photon processes

◆ SetSeed()

void FixedTargetGenerator::SetSeed ( Double_t  seed)
inline

Definition at line 48 of file FixedTargetGenerator.h.

◆ SetTarget()

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

Definition at line 37 of file FixedTargetGenerator.h.

◆ SetTauOnly()

void FixedTargetGenerator::SetTauOnly ( )
inline

Definition at line 40 of file FixedTargetGenerator.h.

40{ tauOnly = true; } // only have Ds decay to tau

◆ UseRandom1()

void FixedTargetGenerator::UseRandom1 ( )
inline

Definition at line 34 of file FixedTargetGenerator.h.

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

◆ UseRandom3()

void FixedTargetGenerator::UseRandom3 ( )
inline

Definition at line 35 of file FixedTargetGenerator.h.

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

◆ WithEvtGen()

void FixedTargetGenerator::WithEvtGen ( )
inline

Definition at line 45 of file FixedTargetGenerator.h.

45{ withEvtGen = true;} // use EvtGen as external decayer to Pythia, experimental phase, only works for one Pythia instance

Member Data Documentation

◆ bparam

Double_t FixedTargetGenerator::bparam
protected

Definition at line 82 of file FixedTargetGenerator.h.

◆ chibb

Double_t FixedTargetGenerator::chibb
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ chicc

Double_t FixedTargetGenerator::chicc
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ ck

Float_t FixedTargetGenerator::ck
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ Debug

Bool_t FixedTargetGenerator::Debug
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ DrellYan

Bool_t FixedTargetGenerator::DrellYan
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ EMax

Double_t FixedTargetGenerator::EMax
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ end

Double_t FixedTargetGenerator::end[3]
protected

Definition at line 81 of file FixedTargetGenerator.h.

◆ endZ

Double_t FixedTargetGenerator::endZ
protected

Definition at line 85 of file FixedTargetGenerator.h.

◆ evtgenN

EvtGenDecays* FixedTargetGenerator::evtgenN
protected

Definition at line 72 of file FixedTargetGenerator.h.

◆ evtgenP

EvtGenDecays* FixedTargetGenerator::evtgenP
protected

Definition at line 73 of file FixedTargetGenerator.h.

◆ fBoost

Double_t FixedTargetGenerator::fBoost
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ fcharmtarget

Bool_t FixedTargetGenerator::fcharmtarget
protected

Definition at line 68 of file FixedTargetGenerator.h.

◆ fin

TFile* FixedTargetGenerator::fin
protected

Definition at line 87 of file FixedTargetGenerator.h.

◆ fLogger

FairLogger* FixedTargetGenerator::fLogger
protected

Definition at line 69 of file FixedTargetGenerator.h.

◆ fMaterialInvestigator

GenieGenerator* FixedTargetGenerator::fMaterialInvestigator
protected

Definition at line 74 of file FixedTargetGenerator.h.

◆ fMom

Double_t FixedTargetGenerator::fMom
protected

Definition at line 62 of file FixedTargetGenerator.h.

◆ fNtuple

TNtuple* FixedTargetGenerator::fNtuple
protected

special option for Dark Photon physics studies

Definition at line 76 of file FixedTargetGenerator.h.

◆ fPythiaN

Pythia8::Pythia* FixedTargetGenerator::fPythiaN
protected

don't make it persistent, magic ROOT command

Definition at line 70 of file FixedTargetGenerator.h.

◆ fPythiaP

Pythia8::Pythia* FixedTargetGenerator::fPythiaP
protected

Definition at line 71 of file FixedTargetGenerator.h.

◆ fRandomEngine

Pythia8::RndmEngine* FixedTargetGenerator::fRandomEngine
private

Definition at line 58 of file FixedTargetGenerator.h.

◆ fSeed

Double_t FixedTargetGenerator::fSeed
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ fUseRandom1

Bool_t FixedTargetGenerator::fUseRandom1
protected

Definition at line 63 of file FixedTargetGenerator.h.

◆ fUseRandom3

Bool_t FixedTargetGenerator::fUseRandom3
protected

Definition at line 64 of file FixedTargetGenerator.h.

◆ G4only

Bool_t FixedTargetGenerator::G4only
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ heartbeat

Int_t FixedTargetGenerator::heartbeat
protected

Definition at line 90 of file FixedTargetGenerator.h.

◆ JpsiMainly

Bool_t FixedTargetGenerator::JpsiMainly
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ maxCrossSection

Double_t FixedTargetGenerator::maxCrossSection
protected

Definition at line 86 of file FixedTargetGenerator.h.

◆ mparam

Double_t FixedTargetGenerator::mparam[10]
protected

Definition at line 83 of file FixedTargetGenerator.h.

◆ n_E

Float_t FixedTargetGenerator::n_E
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_id

Float_t FixedTargetGenerator::n_id
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_M

Float_t FixedTargetGenerator::n_M
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_mE

Float_t FixedTargetGenerator::n_mE
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_mid

Float_t FixedTargetGenerator::n_mid
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_mpx

Float_t FixedTargetGenerator::n_mpx
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_mpy

Float_t FixedTargetGenerator::n_mpy
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_mpz

Float_t FixedTargetGenerator::n_mpz
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_px

Float_t FixedTargetGenerator::n_px
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_py

Float_t FixedTargetGenerator::n_py
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ n_pz

Float_t FixedTargetGenerator::n_pz
protected

Definition at line 89 of file FixedTargetGenerator.h.

◆ nDsprim

Int_t FixedTargetGenerator::nDsprim
protected

Definition at line 66 of file FixedTargetGenerator.h.

◆ nEntry

Int_t FixedTargetGenerator::nEntry
protected

Definition at line 66 of file FixedTargetGenerator.h.

◆ nEvents

Int_t FixedTargetGenerator::nEvents
protected

Definition at line 66 of file FixedTargetGenerator.h.

◆ nrpotspill

Double_t FixedTargetGenerator::nrpotspill
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ ntotprim

Int_t FixedTargetGenerator::ntotprim
protected

Definition at line 66 of file FixedTargetGenerator.h.

◆ nTree

TNtuple* FixedTargetGenerator::nTree
protected

Definition at line 88 of file FixedTargetGenerator.h.

◆ OnlyMuons

Bool_t FixedTargetGenerator::OnlyMuons
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ Option

TString FixedTargetGenerator::Option
protected

Definition at line 77 of file FixedTargetGenerator.h.

◆ PhotonCollision

Bool_t FixedTargetGenerator::PhotonCollision
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ pot

Int_t FixedTargetGenerator::pot
protected

Definition at line 66 of file FixedTargetGenerator.h.

◆ setByHand

Bool_t FixedTargetGenerator::setByHand
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ start

Double_t FixedTargetGenerator::start[3]
protected

Definition at line 80 of file FixedTargetGenerator.h.

◆ startZ

Double_t FixedTargetGenerator::startZ
protected

Definition at line 84 of file FixedTargetGenerator.h.

◆ targetName

TString FixedTargetGenerator::targetName
protected

Definition at line 77 of file FixedTargetGenerator.h.

◆ tauOnly

Bool_t FixedTargetGenerator::tauOnly
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ withEvtGen

Bool_t FixedTargetGenerator::withEvtGen
protected

Definition at line 67 of file FixedTargetGenerator.h.

◆ withNtuple

Bool_t FixedTargetGenerator::withNtuple
protected

Definition at line 75 of file FixedTargetGenerator.h.

◆ wspill

Double_t FixedTargetGenerator::wspill
protected

Definition at line 65 of file FixedTargetGenerator.h.

◆ xOff

Double_t FixedTargetGenerator::xOff
protected

Definition at line 78 of file FixedTargetGenerator.h.

◆ yOff

Double_t FixedTargetGenerator::yOff
protected

Definition at line 79 of file FixedTargetGenerator.h.


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