FairShip
Loading...
Searching...
No Matches
FixedTargetGenerator Class Reference

#include <FixedTargetGenerator.h>

Inheritance diagram for FixedTargetGenerator:
Collaboration diagram for FixedTargetGenerator:

Public Member Functions

 FixedTargetGenerator ()
 
 ~FixedTargetGenerator () override
 
Bool_t ReadEvent (FairPrimaryGenerator *) override
 
void SetParameters (char *)
 
void Print ()
 
Bool_t Init (const char *inFile) override
 
Bool_t Init (const char *inFile, int startEvent) override
 
Bool_t Init () override
 
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 SetTarget (TString s, Double_t x, Double_t y)
 
void SetZoffset (Double_t z)
 
void SetXoffset (Double_t x)
 
void SetYoffset (Double_t y)
 
void SetSmearBeam (Double_t sb)
 
void SetPaintRadius (Double_t r)
 
void SetTargetCoordinates (Double_t start_z, Double_t end_z)
 
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 ()
 
- Public Member Functions inherited from SHiP::Generator
 Generator ()
 
virtual ~Generator ()
 
virtual Bool_t Init (const char *, int)=0
 
virtual Bool_t Init (const char *)=0
 
virtual Bool_t Init (const std::vector< std::string > &inFiles, int startNumber)
 
virtual Bool_t Init (const std::vector< std::string > &inFiles)
 
virtual void UseExternalFile (std::string x, Int_t i)
 
virtual void UseExternalFile (std::vector< std::string > &inFiles, Int_t i)
 

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
 
FairLogger * fLogger
 
Pythia8::Pythia * fPythiaN
 don't make it persistent, magic ROOT command
 
Pythia8::Pythia * fPythiaP
 
Pythia8::EvtGenDecays * evtgenN
 
Pythia8::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 zOff
 
Double_t fsmearBeam
 
Double_t fPaintBeam
 
Double_t start [3]
 
Double_t end [3]
 
Double_t bparam
 
Double_t mparam [10]
 
Double_t startZ
 
Double_t endZ
 
Bool_t targetFromGeometry
 
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
 
- Protected Attributes inherited from SHiP::Generator
std::optional< std::string > fextFile
 
Int_t firstEvent = 0
 

Private Attributes

std::shared_ptr< Pythia8::RndmEngine > fRandomEngine
 

Detailed Description

Definition at line 22 of file FixedTargetGenerator.h.

Constructor & Destructor Documentation

◆ FixedTargetGenerator()

FixedTargetGenerator::FixedTargetGenerator ( )

default constructor

Definition at line 33 of file FixedTargetGenerator.cxx.

33 {
34 fUseRandom1 = kFALSE;
35 fUseRandom3 = kTRUE;
36 fSeed = 0;
37 fMom = 400; // proton
38 fLogger = FairLogger::GetLogger();
39 targetName = "cave_1/target_vacuum_box_1/TargetArea_1/HeVolume_1";
40 xOff = 0;
41 yOff = 0;
42 zOff = 0;
43 fsmearBeam = 8 * mm; // default value for smearing beam (8 mm)
44 fPaintBeam = 5 * cm; // default value for painting beam (5 cm)
45 tauOnly = false;
46 JpsiMainly = false;
47 DrellYan = false;
48 PhotonCollision = false;
49 G4only = false;
50 OnlyMuons = false;
51 maxCrossSection = 0.;
52 EMax = 0;
53 fBoost = 1.;
54 withEvtGen = kFALSE;
55 withNtuple = kFALSE;
56 chicc = 1.7e-3; // prob to produce primary ccbar pair/pot
57 chibb = 1.6e-7; // prob to produce primary bbbar pair/pot
58 nrpotspill = 5.E13; // nr of protons / spill
59 setByHand = kFALSE;
60 Debug = kFALSE;
61 targetFromGeometry = kFALSE;
62 Option = "Primary";
63 wspill = 1.; // event weight == 1 for primary events
64 heartbeat = 1000;
65}
Double_t cm
Double_t mm

◆ ~FixedTargetGenerator()

FixedTargetGenerator::~FixedTargetGenerator ( )
override

destructor

Definition at line 325 of file FixedTargetGenerator.cxx.

325{}

Member Function Documentation

◆ GetPotForCharm()

Double_t FixedTargetGenerator::GetPotForCharm ( )
inline

Definition at line 104 of file FixedTargetGenerator.h.

104{ return nrpotspill / wspill; }

◆ GetPythia()

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

Definition at line 105 of file FixedTargetGenerator.h.

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

◆ GetPythiaN()

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

Definition at line 106 of file FixedTargetGenerator.h.

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

◆ Init() [1/3]

Bool_t FixedTargetGenerator::Init ( )
override

Definition at line 127 of file FixedTargetGenerator.cxx.

127 {
128 fPythiaP =
129 new Pythia8::Pythia(); // pythia needed also for G4only, to copy 2mu BRs
130 if (Option == "Primary" && !G4only) {
131 fPythiaN = new Pythia8::Pythia();
132 } else if (Option != "charm" && Option != "beauty" && !G4only) {
133 LOG(error) << "Option not known " << Option.Data() << ", abort";
134 }
135 if (fUseRandom1) fRandomEngine = std::make_shared<PyTr1Rng>();
136 if (fUseRandom3) fRandomEngine = std::make_shared<PyTr3Rng>();
137 std::vector<int> r = {221, 221, 223, 223, 113, 331, 333};
138 std::vector<int> c = {6, 7, 5, 7, 5, 6, 9}; // decay channel mumu mumuX
139
140 if (Option == "Primary" && !G4only) {
141 fPythiaP->settings.mode("Beams:idB", 2212);
142 fPythiaN->settings.mode("Beams:idB", 2112);
143 } else {
144 fPythiaP->readString("ProcessLevel:all = off");
145 }
146 std::array<Pythia8::Pythia*, 2> plist = {fPythiaP, fPythiaN};
147 Int_t pcount = 0;
148 for (const auto& fPythia : plist) {
149 if (pcount > 0 && (Option != "Primary" || G4only)) {
150 continue;
151 }
152 pcount += 1;
153 fPythia->setRndmEnginePtr(fRandomEngine);
154 fPythia->settings.mode("Random:seed", fSeed);
155 fPythia->settings.mode("Next:numberCount", heartbeat);
156 if (Option == "Primary") {
157 fPythia->settings.mode("Beams:idA", 2212);
158 fPythia->settings.mode("Beams:frameType", 2);
159 fPythia->settings.parm("Beams:eA", fMom); // codespell:ignore parm
160 fPythia->settings.parm("Beams:eB", 0.); // codespell:ignore parm
161 }
162 if (JpsiMainly) {
163 // use this for all onia productions
164 fPythia->readString("Onia:all(3S1) = on");
165 fPythia->readString(
166 "443:new J/psi J/psi 3 0 0 3.09692 0.00009 3.09602 "
167 " 3.09782 0. 1 1 0 1 0");
168 fPythia->readString("443:addChannel = 1 1. 0 -13 13");
169 fPythia->readString(
170 "553:new Upsilon Upsilon 3 0 0 9.46030 0.00005 "
171 "9.45980 9.46080 0.00000e+00 0 1 0 1 0");
172 fPythia->readString("553:addChannel = 1 1. 0 -13 13");
173 }
174 if (DrellYan) {
175 fPythia->readString("WeakSingleBoson:ffbar2gmZ = on");
176 // Z0/gamma -> mu+ mu-
177 fPythia->readString("23:onMode = off");
178 fPythia->readString("23:onIfAll = 13 13");
179 fPythia->readString("23:mMin = 0.5");
180 fPythia->readString("PhaseSpace:mHatMin = 0.5");
181 fPythia->readString("PartonLevel:FSR = on");
182 fPythia->readString("PDF:pSet = 4");
183 }
184 if (PhotonCollision) {
185 fPythia->readString("PhotonCollision:gmgm2mumu = on");
186 }
187 if (!DrellYan && !PhotonCollision) {
188 fPythia->readString("SoftQCD:inelastic = on");
189 fPythia->readString("PhotonCollision:gmgm2mumu = on");
190 fPythia->readString("PromptPhoton:all = on");
191 fPythia->readString("WeakBosonExchange:all = on");
192 fPythia->readString("WeakSingleBoson:all = on");
193 }
194 if (tauOnly) {
195 fPythia->readString(
196 "431:new D_s+ D_s- 1 3 0 1.96849 0.00000 0.00000 "
197 "0.00000 1.49900e-01 0 1 0 1 0");
198 fPythia->readString(
199 "431:addChannel = 1 0.0640000 0 -15 16");
200 }
201
202 // find all long lived particles in pythia
203 Int_t n = 1;
204 while (n != 0) {
205 n = fPythia->particleData.nextId(n);
206 std::shared_ptr<Pythia8::ParticleDataEntry> p =
207 fPythia->particleData.particleDataEntryPtr(n);
208 if (p->tau0() > 1) {
209 std::string particle = std::to_string(n) + ":mayDecay = false";
210 fPythia->readString(particle);
211 LOG(info) << "Made " << p->name().c_str()
212 << " stable for Pythia, should decay in Geant4";
213 }
214 }
215 // boost branching fraction of rare di-muon decays
216 // eta omega rho0 eta' phi
217 if (fBoost != 1.) {
218 LOG(info) << "Rescale BRs of dimuon decays in Pythia: " << fBoost;
219 for (unsigned int i = 0; i < r.size(); ++i) {
220 std::shared_ptr<Pythia8::ParticleDataEntry> V =
221 fPythia->particleData.particleDataEntryPtr(r[i]);
222 Pythia8::DecayChannel ch = V->channel(c[i]);
223 if (TMath::Abs(ch.product(0)) != 13 ||
224 TMath::Abs(ch.product(1)) != 13) {
225 LOG(info) << "this is not the right decay channel: " << r[i] << " "
226 << c[i];
227 } else {
228 TString tmp = "";
229 tmp += r[i];
230 tmp += ":";
231 tmp += c[i];
232 tmp += ":bRatio =";
233 tmp += fBoost * ch.bRatio();
234 fPythia->readString(tmp.Data());
235 }
236 }
237 }
238 fPythia->init();
239 }
240 // Initialize EvtGen.
241 if (withEvtGen) {
242 // Use EVTGENDATA environment variable to find EvtGen data files
243 const char* evtgendata = getenv("EVTGENDATA");
244 if (!evtgendata) {
245 LOG(fatal) << "EVTGENDATA environment variable not set";
246 }
247 TString DecayFile = TString(evtgendata) + "/DECAY.DEC";
248
249 TString ParticleFile = TString(evtgendata) + "/evt.pdl";
250 std::cout << "Using $EVTGENDATA " << evtgendata << std::endl;
251 EvtExternalGenList* extPtr = new EvtExternalGenList();
252 std::list<EvtDecayBase*> models = extPtr->getListOfModels();
253 TString UdecayFile = getenv("FAIRSHIP");
254 UdecayFile += "/gconfig/USERDECAY.DEC";
255 evtgenP = new Pythia8::EvtGenDecays(fPythiaP, DecayFile.Data(),
256 ParticleFile.Data(), extPtr);
257 evtgenP->readDecayFile(
258 UdecayFile.Data()); // will make update of EvtGen with user decay file
259 // use one instance of EvtGen, requires patch to Pythia8Plugins/EvtGen.h
260 if (Option == "Primary") {
261 evtgenN =
262 new Pythia8::EvtGenDecays(fPythiaN, DecayFile.Data(), "", extPtr);
263 }
264 }
265 if (targetFromGeometry) {
266 // Use geometry-based coordinates passed from run_simScript.py
268 start[0] = xOff;
269 start[1] = yOff;
270 start[2] = startZ + zOff;
271 end[0] = xOff;
272 end[1] = yOff;
273 end[2] = endZ;
274 LOG(info) << "FixedTargetGenerator: Using geometry-based target "
275 "coordinates startZ="
276 << startZ << " endZ=" << endZ;
277 // find maximum interaction length
280 } else if (targetName != "") {
281 // Fallback to fragile TGeo navigation for backward compatibility
283 TGeoNavigator* nav = gGeoManager->GetCurrentNavigator();
284 if (nav->CheckPath(targetName)) {
285 nav->cd(targetName);
286 } else {
287 LOG(fatal) << "Invalid target volume specified";
288 }
289 TGeoNode* target = nav->GetCurrentNode();
290 TObjArray* nodes = target->GetVolume()->GetNodes();
291 // Get the first and last node of the target to calculate the material seen
292 TGeoNode* first = static_cast<TGeoNode*>(nodes->At(0));
293 TGeoNode* last = static_cast<TGeoNode*>(nodes->At(nodes->GetSize() - 1));
294 nav->cd(targetName + "/" + first->GetName());
295 TGeoBBox* sha = static_cast<TGeoBBox*>(first->GetVolume()->GetShape());
296 Double_t dz = sha->GetDZ();
297 Double_t origin[3] = {0, 0, -dz};
298 Double_t master[3] = {0, 0, 0};
299 nav->LocalToMaster(origin, master);
300 startZ = master[2];
301 nav->cd(targetName + "/" + last->GetName());
302 sha = static_cast<TGeoBBox*>(first->GetVolume()->GetShape());
303 dz = sha->GetDZ();
304 origin[2] = +dz;
305 nav->LocalToMaster(origin, master);
306 endZ = master[2];
307 start[0] = xOff;
308 start[1] = yOff;
309 start[2] = startZ + zOff;
310 end[0] = xOff;
311 end[1] = yOff;
312 end[2] = endZ;
313 // find maximum interaction length
316 } else {
317 LOG(fatal) << "No target set.";
318 }
319
320 return kTRUE;
321}
GenieGenerator * fMaterialInvestigator
Pythia8::EvtGenDecays * evtgenN
Pythia8::EvtGenDecays * evtgenP
std::shared_ptr< Pythia8::RndmEngine > fRandomEngine
int i
Definition: ShipAna.py:97
None origin(sTree, it)
constants c
Definition: hnl.py:115
ROOT p
Definition: makeDecay.py:76
double MeanMaterialBudget(std::span< const double, 3 > start, std::span< const double, 3 > end, std::span< double, 10 > mparam)

◆ Init() [2/3]

Bool_t FixedTargetGenerator::Init ( const char *  inFile)
inlineoverridevirtual

Implements SHiP::Generator.

Definition at line 34 of file FixedTargetGenerator.h.

34{ return Init(inFile, 0); };

◆ Init() [3/3]

Bool_t FixedTargetGenerator::Init ( const char *  inFile,
int  startEvent 
)
inlineoverridevirtual

Implements SHiP::Generator.

Definition at line 36 of file FixedTargetGenerator.h.

36 {
37 LOG(warning) << "Init with files not implemented for FixedTargetGenerator. "
38 "Using default Init() instead";
39 return Init();
40 };

◆ InitForCharmOrBeauty()

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

Definition at line 66 of file FixedTargetGenerator.cxx.

67 {
68 Option = "charm";
69 nEntry = nStart;
70 // open input file with charm or beauty
71 fin = TFile::Open(fInName);
72 if (!fin) {
73 LOG(fatal) << "Could not open input file: " << fInName.Data();
74 return kFALSE;
75 }
76 nTree = dynamic_cast<TNtuple*>(
77 fin->FindObjectAny("pythia6")); // old format, simple ntuple
78 nEvents = nTree->GetEntries();
79 nTree->SetBranchAddress("id", &n_id);
80 nTree->SetBranchAddress("px", &n_px);
81 nTree->SetBranchAddress("py", &n_py);
82 nTree->SetBranchAddress("pz", &n_pz);
83 nTree->SetBranchAddress("M", &n_M);
84 nTree->SetBranchAddress("E", &n_E);
85 nTree->SetBranchAddress("mid", &n_mid);
86 nTree->SetBranchAddress("mpx", &n_mpx);
87 nTree->SetBranchAddress("mpy", &n_mpy);
88 nTree->SetBranchAddress("mpz", &n_mpz);
89 nTree->SetBranchAddress("mE", &n_mE);
90 if (nTree->GetBranch("k")) {
91 LOG(info) << "+++has branch+++";
92 nTree->SetBranchAddress("k", &ck);
93 }
94 // check if we deal with charm or beauty:
95 nTree->GetEvent(0);
96 if (!setByHand && n_M > 5) {
97 chicc = chibb;
98 LOG(info) << "automatic detection of beauty, configured for beauty";
99 LOG(info) << "bb cross section / mbias " << chicc;
100 } else {
101 LOG(info) << "cc cross section / mbias " << chicc;
102 }
103 // convert pot to weight corresponding to one spill of 5e13 pot
104 // get histogram with number of pot to normalise
105 // pot are counted double, i.e. for each signal, i.e. pot/2.
106 auto* potHist = dynamic_cast<TH1F*>(fin->Get("2"));
107 if (!potHist) {
108 LOG(fatal) << "Histogram '2' not found in input file";
109 return kFALSE;
110 }
111 Int_t nrcpot =
112 potHist->GetBinContent(1) / 2.; // number of primary interactions
114 LOG(info) << "Input file: " << fInName.Data() << " with " << nEvents
115 << " entries, corresponding to nr-pot=" << (nrcpot / chicc);
116 LOG(info) << "weight " << wspill << " corresponding to " << nrpotspill
117 << " p.o.t. per spill for " << nev << " events to process";
118
119 pot = 0.;
120 // Determine fDs on this file for primaries
121 nDsprim = 0;
122 ntotprim = 0;
123
124 return kTRUE;
125}
dict nrcpot
Definition: makeDecay.py:63

◆ Print()

void FixedTargetGenerator::Print ( )

◆ ReadEvent()

Bool_t FixedTargetGenerator::ReadEvent ( FairPrimaryGenerator *  cpg)
override

public method ReadEvent

Definition at line 329 of file FixedTargetGenerator.cxx.

329 {
330 // Calculate beam smearing and painting
332
333 Double_t zinter = 0;
334 Double_t ZoverA = 1.;
335 if (!targetName.IsNull()) {
336 // calculate primary proton interaction point:
337 // loop over trajectory between start and end to pick an interaction point,
338 // copied from GenieGenerator and adapted to hadrons
339 Double_t prob2int = -1.;
340 Double_t rndm = 0.;
341 Double_t sigma;
342 Double_t zinterStart = start[2];
343 if (Option == "charm" || Option == "beauty") {
344 // simulate more downstream interaction points for interactions down in
345 // the cascade
346 if (!(nTree->GetBranch("k"))) {
347 ck = 1;
348 }
349 } else {
350 ck = 1;
351 }
352 while (ck > 0.5) {
353 while (prob2int < rndm) {
354 // place x,y,z uniform along path
355 zinter = gRandom->Uniform(zinterStart, end[2]);
356 Double_t point[3] = {xOff, yOff, zinter};
358 Double_t interLength = mparam[8];
359 TGeoNode* node = gGeoManager->FindNode(point[0], point[1], point[2]);
360 TGeoMaterial* mat = nullptr;
361 if (node && !gGeoManager->IsOutside()) {
362 mat = node->GetVolume()->GetMaterial();
363 Double_t n = mat->GetDensity() / mat->GetA();
364 ZoverA = mat->GetZ() / mat->GetA();
365 sigma = 1. / (n * mat->GetIntLen()) / mbarn;
366 prob2int = TMath::Exp(-interLength) * sigma / maxCrossSection;
367 } else {
368 prob2int = 0.;
369 }
370 rndm = gRandom->Uniform(0., 1.);
371 }
372 zinterStart = zinter;
373 ck -= 1;
374 }
375 zinter = zinter * cm;
376 }
377 Pythia8::Pythia* fPythia;
378 if (G4only) {
379 cpg->AddTrack(2212, 0., 0., fMom, (xOff + dx) / cm, (yOff + dy) / cm,
380 start[2], -1, kTRUE, -1., 0., 1.);
381 return kTRUE;
382 } else if (Option == "Primary") {
383 if (gRandom->Uniform(0., 1.) < ZoverA) {
384 fPythiaP->next();
385 if (withEvtGen) {
386 evtgenP->decay();
387 }
388 fPythia = fPythiaP;
389 } else {
390 fPythiaN->next();
391 if (withEvtGen) {
392 evtgenN->decay();
393 }
394 fPythia = fPythiaN;
395 }
396 } else {
397 if (nEntry == nEvents) {
398 LOG(info) << "Rewind input file: " << nEntry;
399 nEntry = 0;
400 }
401 nTree->GetEvent(nEntry);
402 nEntry += 1;
403 // sanity check, count number of p.o.t. on input file.
404 Double_t pt = TMath::Sqrt((n_mpx * n_mpx) + (n_mpy * n_mpy));
405 // every event appears twice, i.e.
406 if (pt < 1.e-5 && n_mid == 2212) {
407 pot += +0.5;
408 ntotprim += 1;
409 }
410 Int_t idabs = static_cast<int>(TMath::Abs(n_id));
411 if (idabs == 431) {
412 nDsprim += 1;
413 }
414 fPythiaP->event.reset();
415 Int_t nID1 = n_id;
416 fPythiaP->event.append(static_cast<int>(n_id), 1, 0, 0, n_px, n_py, n_pz,
417 n_E, n_M, 0., 9.);
418 TMCProcess procID = kPTransportation;
419 if (n_mid == 2212 && (n_mpx * n_mpx + n_mpy * n_mpy) < 1E-5) {
420 procID = kPPrimary;
421 } // probably primary and not from cascade
422 cpg->AddTrack(static_cast<int>(n_mid), n_mpx, n_mpy, n_mpz,
423 (xOff + dx) * cm, (yOff + dy) * cm, zinter * cm, -1, kFALSE,
424 n_mE, 0., wspill, procID);
425 // second charm hadron in the event
426 nTree->GetEvent(nEntry);
427 if (nID1 * n_id > 0) {
428 LOG(info) << "same sign charm: " << nEntry << ", " << nID1 << ", "
429 << n_id;
430 }
431 nEntry += 1;
432 fPythiaP->event.append(static_cast<int>(n_id), 1, 0, 0, n_px, n_py, n_pz,
433 n_E, n_M, 0., 9.);
434 fPythiaP->next();
435 fPythia = fPythiaP;
436 }
437 if (withEvtGen) {
438 fPythia->moreDecays();
439 } // let the very short lived resonances decay via Pythia8
440 if (Debug && !G4only) {
441 fPythia->event.list();
442 }
443 TMCProcess procID;
444 for (Int_t ii = 1; ii < fPythia->event.size(); ii++) {
445 Double_t e = fPythia->event[ii].e();
446 Double_t m = fPythia->event[ii].m();
447 Double_t pz = fPythia->event[ii].pz();
448 Int_t id = fPythia->event[ii].id();
449 Bool_t wanttracking = kTRUE;
450 if (e - m < EMax || !fPythia->event[ii].isFinal() || pz < 0) {
451 wanttracking = kFALSE;
452 }
454 // don't track underlying event
455 if (fabs(id) != 13) {
456 wanttracking = kFALSE;
457 }
458 }
459 Double_t z = fPythia->event[ii].zProd() * mm + zinter * cm;
460 Double_t x = fPythia->event[ii].xProd() * mm + xOff * cm + dx * cm;
461 Double_t y = fPythia->event[ii].yProd() * mm + yOff * cm + dy * cm;
462 Double_t tof =
463 fPythia->event[ii].tProd() / (10 * c_light); // to go from mm to s
464 Double_t px = fPythia->event[ii].px();
465 Double_t py = fPythia->event[ii].py();
466 Int_t im = fPythia->event[ii].mother1() - 1;
467 procID = kPPrimary;
468 if (Option != "Primary") {
469 procID = kPDecay;
470 im += 1;
471 if (ii == 1) {
472 procID = kPHadronic;
473 }
474 } else {
475 if (ii < 3) {
476 im = -1;
477 }
478 }
479 cpg->AddTrack(id, px, py, pz, x * cm, y * cm, z * cm, im, wanttracking, e,
480 tof, wspill, procID);
481 if (withNtuple) {
482 int idabs = TMath::Abs(id);
483 if (idabs < 10) {
484 continue;
485 }
486 if (idabs < 18 || idabs == 22 || idabs == 111 || idabs == 221 ||
487 idabs == 223 || idabs == 331 || idabs == 211 || idabs == 113 ||
488 idabs == 333 || idabs == 321 || idabs == 2212) {
489 Int_t moID = fPythia->event[im + 1].id();
490 fNtuple->Fill(0, id, px, py, pz, e, x * cm, y * cm, z * cm, moID);
491 }
492 }
493 }
494 return kTRUE;
495}
std::pair< Double_t, Double_t > CalculateBeamOffset(Double_t smearBeam, Double_t paintBeam)
const Double_t c_light
const Double_t mbarn
Double_t m
Definition: diagrams_e.h:4
TNtuple * fNtuple
special option for Dark Photon physics studies
ROOT pt
Definition: makeDecay.py:125
dict sigma
Definition: runPythia8.py:187

◆ SetBoost()

void FixedTargetGenerator::SetBoost ( Double_t  f)
inline

Definition at line 71 of file FixedTargetGenerator.h.

71 {
72 fBoost = f;
73 } // boost factor for rare di-muon decays

◆ SetChibb()

void FixedTargetGenerator::SetChibb ( Double_t  x)
inline

Definition at line 88 of file FixedTargetGenerator.h.

88 {
89 chibb = x;
90 } // chibb = bbbar over mbias cross section

◆ SetChicc()

void FixedTargetGenerator::SetChicc ( Double_t  x)
inline

Definition at line 91 of file FixedTargetGenerator.h.

91 {
92 chicc = x;
93 } // chicc = ccbar over mbias cross section

◆ SetDebug()

void FixedTargetGenerator::SetDebug ( Bool_t  x)
inline

Definition at line 99 of file FixedTargetGenerator.h.

99{ Debug = x; }

◆ SetDrellYan()

void FixedTargetGenerator::SetDrellYan ( )
inline

Definition at line 80 of file FixedTargetGenerator.h.

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

◆ SetEnergyCut()

void FixedTargetGenerator::SetEnergyCut ( Float_t  emax)
inline

Definition at line 96 of file FixedTargetGenerator.h.

96 {
97 EMax = emax;
98 } // min energy to be copied to Geant4

◆ SetG4only()

void FixedTargetGenerator::SetG4only ( )
inline

Definition at line 74 of file FixedTargetGenerator.h.

74 {
75 G4only = true;
76 } // only run Geant4, no pythia primary interaction

◆ SetHeartBeat()

void FixedTargetGenerator::SetHeartBeat ( Int_t  x)
inline

Definition at line 95 of file FixedTargetGenerator.h.

95{ heartbeat = x; }

◆ SetJpsiMainly()

void FixedTargetGenerator::SetJpsiMainly ( )
inline

Definition at line 78 of file FixedTargetGenerator.h.

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

◆ SetMom()

void FixedTargetGenerator::SetMom ( Double_t  mom)
inline

Definition at line 47 of file FixedTargetGenerator.h.

47{ fMom = mom; };

◆ SetOnlyMuons()

void FixedTargetGenerator::SetOnlyMuons ( )
inline

Definition at line 79 of file FixedTargetGenerator.h.

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

◆ SetOpt4DP()

void FixedTargetGenerator::SetOpt4DP ( TNtuple *  t)
inline

Definition at line 100 of file FixedTargetGenerator.h.

100 {
101 withNtuple = kTRUE;
102 fNtuple = t;
103 }

◆ SetPaintRadius()

void FixedTargetGenerator::SetPaintRadius ( Double_t  r)
inline

Definition at line 65 of file FixedTargetGenerator.h.

65{ fPaintBeam = r; };

◆ SetParameters()

void FixedTargetGenerator::SetParameters ( char *  )

◆ SetPhotonCollision()

void FixedTargetGenerator::SetPhotonCollision ( )
inline

Definition at line 81 of file FixedTargetGenerator.h.

81 {
82 PhotonCollision = true;
83 } // only generate prompt photon processes

◆ SetSeed()

void FixedTargetGenerator::SetSeed ( Double_t  seed)
inline

Definition at line 94 of file FixedTargetGenerator.h.

94{ fSeed = seed; }

◆ SetSmearBeam()

void FixedTargetGenerator::SetSmearBeam ( Double_t  sb)
inline

Definition at line 64 of file FixedTargetGenerator.h.

64{ fsmearBeam = sb; };

◆ SetTarget()

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

Definition at line 56 of file FixedTargetGenerator.h.

56 {
57 targetName = s;
58 xOff = x;
59 yOff = y;
60 };

◆ SetTargetCoordinates()

void FixedTargetGenerator::SetTargetCoordinates ( Double_t  start_z,
Double_t  end_z 
)
inline

Definition at line 66 of file FixedTargetGenerator.h.

66 {
67 startZ = start_z;
68 endZ = end_z;
69 targetFromGeometry = true;
70 };

◆ SetTauOnly()

void FixedTargetGenerator::SetTauOnly ( )
inline

Definition at line 77 of file FixedTargetGenerator.h.

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

◆ SetXoffset()

void FixedTargetGenerator::SetXoffset ( Double_t  x)
inline

Definition at line 62 of file FixedTargetGenerator.h.

62{ xOff = x; };

◆ SetYoffset()

void FixedTargetGenerator::SetYoffset ( Double_t  y)
inline

Definition at line 63 of file FixedTargetGenerator.h.

63{ yOff = y; };

◆ SetZoffset()

void FixedTargetGenerator::SetZoffset ( Double_t  z)
inline

Definition at line 61 of file FixedTargetGenerator.h.

61{ zOff = z; };

◆ UseRandom1()

void FixedTargetGenerator::UseRandom1 ( )
inline

Definition at line 48 of file FixedTargetGenerator.h.

48 {
49 fUseRandom1 = kTRUE;
50 fUseRandom3 = kFALSE;
51 };

◆ UseRandom3()

void FixedTargetGenerator::UseRandom3 ( )
inline

Definition at line 52 of file FixedTargetGenerator.h.

52 {
53 fUseRandom1 = kFALSE;
54 fUseRandom3 = kTRUE;
55 };

◆ WithEvtGen()

void FixedTargetGenerator::WithEvtGen ( )
inline

Definition at line 84 of file FixedTargetGenerator.h.

84 {
85 withEvtGen = true;
86 } // use EvtGen as external decayer to Pythia, experimental phase, only works

Member Data Documentation

◆ bparam

Double_t FixedTargetGenerator::bparam
protected

Definition at line 135 of file FixedTargetGenerator.h.

◆ chibb

Double_t FixedTargetGenerator::chibb
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ chicc

Double_t FixedTargetGenerator::chicc
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ ck

Float_t FixedTargetGenerator::ck
protected

Definition at line 144 of file FixedTargetGenerator.h.

◆ Debug

Bool_t FixedTargetGenerator::Debug
protected

Definition at line 118 of file FixedTargetGenerator.h.

◆ DrellYan

Bool_t FixedTargetGenerator::DrellYan
protected

Definition at line 117 of file FixedTargetGenerator.h.

◆ EMax

Double_t FixedTargetGenerator::EMax
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ end

Double_t FixedTargetGenerator::end[3]
protected

Definition at line 134 of file FixedTargetGenerator.h.

◆ endZ

Double_t FixedTargetGenerator::endZ
protected

Definition at line 138 of file FixedTargetGenerator.h.

◆ evtgenN

Pythia8::EvtGenDecays* FixedTargetGenerator::evtgenN
protected

Definition at line 122 of file FixedTargetGenerator.h.

◆ evtgenP

Pythia8::EvtGenDecays* FixedTargetGenerator::evtgenP
protected

Definition at line 123 of file FixedTargetGenerator.h.

◆ fBoost

Double_t FixedTargetGenerator::fBoost
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ fin

TFile* FixedTargetGenerator::fin
protected

Definition at line 141 of file FixedTargetGenerator.h.

◆ fLogger

FairLogger* FixedTargetGenerator::fLogger
protected

Definition at line 119 of file FixedTargetGenerator.h.

◆ fMaterialInvestigator

GenieGenerator* FixedTargetGenerator::fMaterialInvestigator
protected

Definition at line 124 of file FixedTargetGenerator.h.

◆ fMom

Double_t FixedTargetGenerator::fMom
protected

Definition at line 112 of file FixedTargetGenerator.h.

◆ fNtuple

TNtuple* FixedTargetGenerator::fNtuple
protected

special option for Dark Photon physics studies

Definition at line 126 of file FixedTargetGenerator.h.

◆ fPaintBeam

Double_t FixedTargetGenerator::fPaintBeam
protected

Definition at line 132 of file FixedTargetGenerator.h.

◆ fPythiaN

Pythia8::Pythia* FixedTargetGenerator::fPythiaN
protected

don't make it persistent, magic ROOT command

Definition at line 120 of file FixedTargetGenerator.h.

◆ fPythiaP

Pythia8::Pythia* FixedTargetGenerator::fPythiaP
protected

Definition at line 121 of file FixedTargetGenerator.h.

◆ fRandomEngine

std::shared_ptr<Pythia8::RndmEngine> FixedTargetGenerator::fRandomEngine
private

Definition at line 109 of file FixedTargetGenerator.h.

◆ fSeed

Double_t FixedTargetGenerator::fSeed
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ fsmearBeam

Double_t FixedTargetGenerator::fsmearBeam
protected

Definition at line 131 of file FixedTargetGenerator.h.

◆ fUseRandom1

Bool_t FixedTargetGenerator::fUseRandom1
protected

Definition at line 113 of file FixedTargetGenerator.h.

◆ fUseRandom3

Bool_t FixedTargetGenerator::fUseRandom3
protected

Definition at line 114 of file FixedTargetGenerator.h.

◆ G4only

Bool_t FixedTargetGenerator::G4only
protected

Definition at line 117 of file FixedTargetGenerator.h.

◆ heartbeat

Int_t FixedTargetGenerator::heartbeat
protected

Definition at line 145 of file FixedTargetGenerator.h.

◆ JpsiMainly

Bool_t FixedTargetGenerator::JpsiMainly
protected

Definition at line 117 of file FixedTargetGenerator.h.

◆ maxCrossSection

Double_t FixedTargetGenerator::maxCrossSection
protected

Definition at line 140 of file FixedTargetGenerator.h.

◆ mparam

Double_t FixedTargetGenerator::mparam[10]
protected

Definition at line 136 of file FixedTargetGenerator.h.

◆ n_E

Float_t FixedTargetGenerator::n_E
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_id

Float_t FixedTargetGenerator::n_id
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_M

Float_t FixedTargetGenerator::n_M
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_mE

Float_t FixedTargetGenerator::n_mE
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_mid

Float_t FixedTargetGenerator::n_mid
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_mpx

Float_t FixedTargetGenerator::n_mpx
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_mpy

Float_t FixedTargetGenerator::n_mpy
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_mpz

Float_t FixedTargetGenerator::n_mpz
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_px

Float_t FixedTargetGenerator::n_px
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_py

Float_t FixedTargetGenerator::n_py
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ n_pz

Float_t FixedTargetGenerator::n_pz
protected

Definition at line 143 of file FixedTargetGenerator.h.

◆ nDsprim

Int_t FixedTargetGenerator::nDsprim
protected

Definition at line 116 of file FixedTargetGenerator.h.

◆ nEntry

Int_t FixedTargetGenerator::nEntry
protected

Definition at line 116 of file FixedTargetGenerator.h.

◆ nEvents

Int_t FixedTargetGenerator::nEvents
protected

Definition at line 116 of file FixedTargetGenerator.h.

◆ nrpotspill

Double_t FixedTargetGenerator::nrpotspill
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ ntotprim

Int_t FixedTargetGenerator::ntotprim
protected

Definition at line 116 of file FixedTargetGenerator.h.

◆ nTree

TNtuple* FixedTargetGenerator::nTree
protected

Definition at line 142 of file FixedTargetGenerator.h.

◆ OnlyMuons

Bool_t FixedTargetGenerator::OnlyMuons
protected

Definition at line 118 of file FixedTargetGenerator.h.

◆ Option

TString FixedTargetGenerator::Option
protected

Definition at line 127 of file FixedTargetGenerator.h.

◆ PhotonCollision

Bool_t FixedTargetGenerator::PhotonCollision
protected

Definition at line 117 of file FixedTargetGenerator.h.

◆ pot

Int_t FixedTargetGenerator::pot
protected

Definition at line 116 of file FixedTargetGenerator.h.

◆ setByHand

Bool_t FixedTargetGenerator::setByHand
protected

Definition at line 117 of file FixedTargetGenerator.h.

◆ start

Double_t FixedTargetGenerator::start[3]
protected

Definition at line 133 of file FixedTargetGenerator.h.

◆ startZ

Double_t FixedTargetGenerator::startZ
protected

Definition at line 137 of file FixedTargetGenerator.h.

◆ targetFromGeometry

Bool_t FixedTargetGenerator::targetFromGeometry
protected

Definition at line 139 of file FixedTargetGenerator.h.

◆ targetName

TString FixedTargetGenerator::targetName
protected

Definition at line 127 of file FixedTargetGenerator.h.

◆ tauOnly

Bool_t FixedTargetGenerator::tauOnly
protected

Definition at line 117 of file FixedTargetGenerator.h.

◆ withEvtGen

Bool_t FixedTargetGenerator::withEvtGen
protected

Definition at line 118 of file FixedTargetGenerator.h.

◆ withNtuple

Bool_t FixedTargetGenerator::withNtuple
protected

Definition at line 125 of file FixedTargetGenerator.h.

◆ wspill

Double_t FixedTargetGenerator::wspill
protected

Definition at line 115 of file FixedTargetGenerator.h.

◆ xOff

Double_t FixedTargetGenerator::xOff
protected

Definition at line 128 of file FixedTargetGenerator.h.

◆ yOff

Double_t FixedTargetGenerator::yOff
protected

Definition at line 129 of file FixedTargetGenerator.h.

◆ zOff

Double_t FixedTargetGenerator::zOff
protected

Definition at line 130 of file FixedTargetGenerator.h.


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