CbmRoot
Loading...
Searching...
No Matches
CbmTransport.cxx
Go to the documentation of this file.
1/* Copyright (C) 2019-2020 GSI Helmholtzzentrum fuer Schwerionenforschung, Darmstadt
2 SPDX-License-Identifier: GPL-3.0-only
3 Authors: Volker Friese, Florian Uhlig [committer] */
4
10#include "CbmTransport.h"
11
12#include "CbmBeamProfile.h"
13#include "CbmEventGenerator.h"
14#include "CbmFieldMap.h"
15#include "CbmFieldPar.h"
16#include "CbmFileUtils.h"
17#include "CbmGeant3Settings.h"
18#include "CbmGeant4Settings.h"
19#include "CbmPlutoGenerator.h"
20#include "CbmSetup.h"
21#include "CbmStack.h"
22#include "CbmTarget.h"
23#include "CbmUnigenGenerator.h"
24
25#include <FairMonitor.h>
26#include <FairParRootFileIo.h>
27#include <FairRootFileSink.h>
28#include <FairRunSim.h>
29#include <FairRuntimeDb.h>
30#include <FairSystemInfo.h>
31#include <FairUrqmdGenerator.h>
32#include <Logger.h>
33
34#include <TDatabasePDG.h>
35#include <TG4RunConfiguration.h>
36#include <TGeant3.h>
37#include <TGeant3TGeo.h>
38#include <TGeant4.h>
39#include <TGeoManager.h>
40#include <TPythia6Decayer.h>
41#include <TROOT.h>
42#include <TRandom.h>
43#include <TStopwatch.h>
44#include <TString.h>
45#include <TSystem.h>
46#include <TVector3.h>
47#include <TVirtualMC.h>
48
49#include <boost/filesystem.hpp>
50
51#include <array>
52#include <cassert>
53#include <climits>
54#include <iostream>
55#include <sstream>
56#include <string>
57
58using std::stringstream;
59
60
61// ----- Constructor ----------------------------------------------------
63 : TNamed("CbmTransport", "Transport Run")
64 , fSetup(CbmSetup::Instance())
65 , fField(nullptr)
66 , fTarget()
67 , fEventGen(new CbmEventGenerator())
68 , fEventFilter(new CbmMCEventFilter())
69 , fRun(new FairRunSim())
70 , fOutFileName()
71 , fParFileName()
72 , fGeoFileName()
73 , fGenerators()
74 , fRealTimeInit(0.)
75 , fRealTimeRun(0.)
76 , fCpuTime(0.)
77 , fEngine(kGeant3)
78 , fStackFilter(new CbmStackFilter())
79 , fGenerateRunInfo(kFALSE)
80 , fStoreTrajectories(kFALSE)
81{
82 // TODO: I do not like instantiating FairRunSim from this constructor;
83 // It should be done in Run(). However, the presence of a FairRunSim
84 // is required by CbmUnigenGenerator. Not a good construction; should
85 // be done better.
86
87 // Initialisation of the TDatabasePDG. This is done here because in the
88 // course of the run, new particles may be added. The method
89 // ReadPDGTable, however, is not executed from the constructor, but only
90 // from GetParticle(), if the particle list is not empty.
91 // So, if one adds particles before the first call to GetParticle(),
92 // the particle table is never loaded.
93 // TDatabasePDG is a singleton, but there is no way to check whether
94 // it has already read the particle table file, nor to see if there are
95 // any contents, nor to clean the particle list. A truly remarkable
96 // implementation.
97 auto pdgdb = TDatabasePDG::Instance();
98 pdgdb->ReadPDGTable();
99
100 // By default, vertex smearing along the beam is activated
101 fEventGen->SmearVertexZ(kTRUE);
102 fRun->SetGenerator(fEventGen); // has to be available for some generators
103}
104// --------------------------------------------------------------------------
105
106
107// ----- Destructor -----------------------------------------------------
109// --------------------------------------------------------------------------
110
111
112// ----- Add a file-based input -----------------------------------------
113void CbmTransport::AddInput(const char* fileName, ECbmGenerator genType)
114{
115
116 FairGenerator* generator = NULL;
117
118 if (genType == kUrqmd) {
119 if (gSystem->AccessPathName(fileName)) {
120 LOG(fatal) << GetName() << ": Input file " << fileName << " not found!";
121 return;
122 }
123 }
124 else {
125 if (!Cbm::File::IsRootFile(fileName)) {
126 LOG(fatal) << GetName() << ": Input file " << fileName << " not found!";
127 return;
128 }
129 }
130
131 switch (genType) {
132 case kUnigen: generator = new CbmUnigenGenerator(TString(fileName)); break;
133 case kUrqmd: generator = new FairUrqmdGenerator(fileName); break;
134 case kPluto: generator = new CbmPlutoGenerator(fileName); break;
135 }
136
137 assert(generator);
138 fEventGen->AddGenerator(generator);
139}
140// --------------------------------------------------------------------------
141
142
143// ----- Add a generator-based input ------------------------------------
144void CbmTransport::AddInput(FairGenerator* generator)
145{
146 assert(generator);
147 fEventGen->AddGenerator(generator);
148}
149// --------------------------------------------------------------------------
150
151
152// ----- Configure the TVirtualMC ---------------------------------------
154{
155
156 std::cout << std::endl;
157 LOG(info) << GetName() << ": Configuring VMC...";
158 TVirtualMC* vmc = nullptr;
159
160 if (fEngine == kGeant3) {
161 TString* gModel = fRun->GetGeoModel();
162 if (strncmp(gModel->Data(), "TGeo", 4) == 0) {
163 LOG(info) << GetName() << ": Create TGeant3TGeo";
164 vmc = new TGeant3TGeo("C++ Interface to Geant3 with TGeo");
165 } //? Geant3 with TGeo
166 else {
167 LOG(info) << GetName() << ": Create TGeant3";
168 vmc = new TGeant3("C++ Interface to Geant3");
169 } //? Native Geant3
171 fGeant3Settings->Init(vmc);
172 } //? Geant3
173
174 else if (fEngine == kGeant4) {
175 LOG(info) << GetName() << ": Create TGeant4";
177 std::array<std::string, 3> runConfigSettings = fGeant4Settings->GetG4RunConfig();
178 TG4RunConfiguration* runConfig =
179 new TG4RunConfiguration(runConfigSettings[0], runConfigSettings[1], runConfigSettings[2]);
180 vmc = new TGeant4("TGeant4", "C++ Interface to Geant4", runConfig);
181 fGeant4Settings->Init(vmc);
182 } //? Geant4
183
184 else
185 LOG(fatal) << GetName() << ": unknown transport engine!";
186
187 // Create stack
188 std::unique_ptr<CbmStack> stack(new CbmStack());
189 stack->SetFilter(fStackFilter);
190 if (vmc) vmc->SetStack(stack.release());
191}
192// --------------------------------------------------------------------------
193
194
195// ----- Force creation of event vertex at a given z position -----------
196void CbmTransport::ForceVertexAtZ(Double_t zVertex)
197{
198 assert(fEventGen);
199 fEventGen->ForceVertexAtZ(zVertex);
200}
201// --------------------------------------------------------------------------
202
203
204// ----- Force creation of event vertex in the target -------------------
206{
207 assert(fEventGen);
209}
210// --------------------------------------------------------------------------
211
212
213// ----- Force user-defined single-mode decays --------------------------
215{
216
217 assert(gMC);
218 auto pdgdb = TDatabasePDG::Instance();
219
220 // --- Setting user decays does not work with TGeant4
221 if ((!fDecayModes.empty()) && fEngine == kGeant4)
222 LOG(fatal) << GetName() << ": Forcing decay modes is not possible with TGeant4!";
223
224 for (auto& decay : fDecayModes) {
225
226 Int_t pdg = decay.first;
227 UInt_t nDaughters = decay.second.size();
228 stringstream log;
229 log << GetName() << ": Force decay " << pdgdb->GetParticle(pdg)->GetName() << " -> ";
230
231 // First check whether VMC knows the particle. Not all particles
232 // in TDatabasePDG are necessarily defined in VMC.
233 // This check is there because the call to TVirtualMC::SetUserDecay
234 // has no return value signifying success. If the particle is not
235 // found, just an error message is printed, which most likely
236 // goes unnoticed by the user.
237 // The access to ParticleMCTpye seems to me the only way to check.
238 // No method like Bool_t CheckParticle(Int_t) is there, which any
239 // sensible programmer would have put.
240 if (gMC->ParticleMCType(pdg) == kPTUndefined) { // At least TGeant3 delivers that
241 LOG(info) << log.str();
242 LOG(fatal) << GetName() << ": PDG " << pdg << " not in VMC!";
243 continue;
244 }
245
246 // For up to three daughters, the native decayer is used
247 if (nDaughters <= 3) {
248 Float_t branch[6] = {100., 0., 0., 0., 0., 0.}; // branching ratios
249 Int_t mode[6][3] = {{0, 0, 0}, {0, 0, 0}, {0, 0, 0}, {0, 0, 0}, {0, 0, 0}, {0, 0, 0}}; // decay modes
250 for (UInt_t iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
251 mode[0][iDaughter] = decay.second[iDaughter];
252 log << pdgdb->GetParticle(decay.second[iDaughter])->GetName() << " ";
253 }
254 Bool_t success = gMC->SetDecayMode(pdg, branch, mode);
255 if (!success) {
256 LOG(info) << log.str();
257 LOG(fatal) << GetName() << ": Setting decay mode failed!";
258 }
259 log << ", using native decayer.";
260 } //? not more than three daughters
261
262 // For more than three daughters, we must use TPythia6 as external decayer
263 else {
264 auto p6decayer = TPythia6Decayer::Instance();
265 Int_t daughterPdg[nDaughters];
266 Int_t multiplicity[nDaughters];
267 for (UInt_t iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
268 daughterPdg[iDaughter] = decay.second[iDaughter];
269 multiplicity[iDaughter] = 1;
270 log << pdgdb->GetParticle(decay.second[iDaughter])->GetName() << " ";
271 } //# daughters
272 p6decayer->ForceParticleDecay(pdg, daughterPdg, multiplicity, nDaughters);
273 // We have to tell the VMC to use the Pythia decayer for this particle
274 gMC->SetUserDecay(pdg);
275 } //? more than three daughters
276
277 LOG(info) << log.str();
278 } //# user-defined decay modes
279}
280// --------------------------------------------------------------------------
281
282// ----- Initialisation of event generator ------------------------------
284{
285
286 // --- Set the target properties to the event generator
288
289 // --- Log output
290 std::cout << std::endl;
291 LOG(info) << "----- Settings for event generator";
292 fEventGen->Print();
293 LOG(info) << "----- End settings for event generator";
294 std::cout << std::endl;
295
296 // --- If a target is specified, check its consistency with the beam
297 // --- profile. The average beam must hit the target surface.
298 if (fTarget) {
300 LOG(fatal) << GetName() << ": Beam profile is not consistent with target!";
301 } //? Target not consistent with beam
302 } //? Target specified
303}
304// --------------------------------------------------------------------------
305
306
307// ----- Set the source the setup will be loaded from -------------------
309// --------------------------------------------------------------------------
310
311
312// ----- Load a standard setup ------------------------------------------
313void CbmTransport::LoadSetup(const char* setupName) { fSetup->LoadSetup(setupName); }
314// --------------------------------------------------------------------------
315
316
317// ----- Register ions to the TDatabsePDG -------------------------------
319{
320
321 // TODO: Better would be loading the additional particles from a text file.
322 // TDatabasePDG reads the particle definitions from pdg_table.txt
323 // (in $SIMPATH/share/root/etc). The method TDatabase::ReadPDGTable()
324 // is triggered on first call to TDatabasePDG::GetParticle(Int_t), if the
325 // database is still empty.
326 // We could call TDatabasePDG::ReadPDGTable(name_of_cbm_file) after the
327 // first initialisation of the database; the there defined particles would
328 // be added on top of the existing ones.
329
330 // Particle database and variables
331 TDatabasePDG* pdgdb = TDatabasePDG::Instance();
332 const char* name = "";
333 Int_t code = 0;
334 Double_t mass = 0.;
335 Bool_t stable = kTRUE;
336 Double_t charge = 0.;
337
338 // --- deuteron and anti-deuteron
339 name = "d+";
340 code = 1000010020;
341 mass = 1.876124;
342 stable = kTRUE;
343 charge = 1.;
344 pdgdb->AddParticle(name, name, mass, stable, 0., charge, "Ion", code);
345 pdgdb->AddAntiParticle("d-", -1 * code);
346
347 // --- tritium and anti-tritium
348 name = "t+";
349 code = 1000010030;
350 mass = 2.809432;
351 stable = kTRUE;
352 charge = 1.;
353 pdgdb->AddParticle(name, name, mass, stable, 0., charge, "Ion", code);
354 pdgdb->AddAntiParticle("t-", -1 * code);
355
356 // --- Helium_3 and its anti-nucleus
357 name = "He3+";
358 code = 1000020030;
359 mass = 2.809413;
360 stable = kTRUE;
361 charge = 2.;
362 pdgdb->AddParticle(name, name, mass, stable, 0., charge, "Ion", code);
363 pdgdb->AddAntiParticle("He3-", -1 * code);
364
365 // --- Helium_4 and its anti-nucleus
366 name = "He4+";
367 code = 1000020040;
368 mass = 3.7284;
369 stable = kTRUE;
370 charge = 2.;
371 pdgdb->AddParticle(name, name, mass, stable, 0., charge, "Ion", code);
372 pdgdb->AddAntiParticle("He3-", -1 * code);
373}
374// --------------------------------------------------------------------------
375
376
377// ----- Register radiation length --------------------------------------
379{
380 assert(fRun);
381 fRun->SetRadLenRegister(choice);
382 LOG(info) << GetName() << ": Radiation length register is enabled";
383}
384// --------------------------------------------------------------------------
385
386
387// ----- Create and register the setup modules --------------------------
389// --------------------------------------------------------------------------
390
391
392// ----- Set correct decay modes for pi0 and eta ------------------------
393void CbmTransport::PiAndEtaDecay(TVirtualMC* vmc)
394{
395
396 assert(vmc);
397 LOG(info) << GetName() << ": Set decay modes for pi0 and eta";
398
399 if (fEngine == kGeant4) {
400 std::cout << std::endl << std::endl;
401 LOG(warn) << "***********************************";
402 LOG(warn) << "***********************************";
403 LOG(warn) << GetName() << ": User decay modes cannot be set with TGeant4!";
404 LOG(warn) << GetName() << ": Built-in decay modes for pi0 and eta will be used!";
405 LOG(warn) << "***********************************";
406 LOG(warn) << "***********************************";
407 std::cout << std::endl << std::endl;
408 return;
409 }
410
411 TGeant3* gMC3 = static_cast<TGeant3*>(vmc);
412 assert(vmc);
413
414 // Decay modes for eta mesons (PDG 2016)
415 Int_t modeEta[6]; // decay modes
416 Float_t brEta[6]; // branching ratios in %
417
418 // --- eta -> gamma gamma
419 modeEta[0] = 101;
420 brEta[0] = 39.41;
421
422 // --- eta -> pi0 pi0 pi0
423 modeEta[1] = 70707;
424 brEta[1] = 32.68;
425
426 // --- eta -> pi+ pi- pi0
427 modeEta[2] = 80907;
428 brEta[2] = 22.92;
429
430 // --- eta -> pi+ pi- gamma
431 modeEta[3] = 80901;
432 brEta[3] = 4.22;
433
434 // --- eta -> e+ e- gamma
435 modeEta[4] = 30201;
436 brEta[4] = 0.69;
437
438 // --- eta -> pi0 gamma gamma
439 modeEta[5] = 10107;
440 brEta[5] = 2.56e-2;
441
442 // --- Set the eta decays
443 gMC3->Gsdk(17, brEta, modeEta);
444
445 // --- Decay modes for pi0
446 Int_t modePi[6]; // decay modes
447 Float_t brPi[6]; // branching ratios in %
448
449 // --- pi0 -> gamma gamma
450 modePi[0] = 101;
451 brPi[0] = 98.823;
452
453 // --- pi0 -> e+ e- gamma
454 modePi[1] = 30201;
455 brPi[1] = 1.174;
456
457 // --- No other channels for pi0
458 for (Int_t iMode = 2; iMode < 6; iMode++) {
459 modePi[iMode] = 0;
460 brPi[iMode] = 0.;
461 }
462
463 // --- Set the pi0 decays
464 gMC3->Gsdk(7, brPi, modePi);
465}
466// --------------------------------------------------------------------------
467
468
469// ----- Execute transport run ------------------------------------------
470void CbmTransport::Run(Int_t nEvents)
471{
472
473 // Get the minimum number of events from all file based generators
474 // Set the number of events to process to this minimum number of events
475 Int_t numAvailEvents {0};
476 Int_t numMinAvailEvents {INT_MAX};
477 TObjArray* genList = fEventGen->GetListOfGenerators();
478 for (Int_t i = 0; i < genList->GetEntries(); i++) {
479 CbmUnigenGenerator* gen = dynamic_cast<CbmUnigenGenerator*>(genList->At(i));
480 if (gen) {
481 numAvailEvents = gen->GetNumAvailableEvents();
482 if (nEvents > numAvailEvents) {
483 if (numAvailEvents < numMinAvailEvents) { numMinAvailEvents = numAvailEvents; }
484 }
485 }
486 CbmPlutoGenerator* pgen = dynamic_cast<CbmPlutoGenerator*>(genList->At(i));
487 if (pgen) {
488 numAvailEvents = pgen->GetNumAvailableEvents();
489 if (nEvents > numAvailEvents) {
490 if (numAvailEvents < numMinAvailEvents) { numMinAvailEvents = numAvailEvents; }
491 }
492 }
493 }
494 if (nEvents > numMinAvailEvents) {
495 LOG(warning) << "";
496 LOG(warning) << "The number of requested events (" << nEvents << ") is larger than the number of available events ("
497 << numMinAvailEvents << ")";
498 LOG(warning) << "Set the number of events to process to " << numMinAvailEvents;
499 LOG(warning) << "";
500 nEvents = numMinAvailEvents;
501 }
502
503 // --- Timer
504 TStopwatch timer;
505
506 // --- Set the global random seed
507 gRandom->SetSeed(fRandomSeed);
508
509 // --- Check presence of required requisites
510 if (fOutFileName.IsNull()) LOG(fatal) << GetName() << ": No output file specified!";
511 if (fParFileName.IsNull()) LOG(fatal) << GetName() << ": No parameter file specified!";
512 std::cout << std::endl << std::endl;
513
514
515 // --- Add some required particles to the TDatabasePDG
516 RegisterIons();
517
518
519 // --- Set transport engine
520 const char* engineName = "";
521 switch (fEngine) {
522 case kGeant3: engineName = "TGeant3"; break;
523 case kGeant4: engineName = "TGeant4"; break;
524 default: {
525 LOG(warn) << GetName() << ": Unknown transport engine ";
526 engineName = "TGeant3";
527 break;
528 }
529 } //? engine
530 LOG(info) << GetName() << ": Using engine " << engineName;
531 fRun->SetName(engineName);
532
533
534 // --- Create file sink using output file name
535 // TODO: remove release after switching to FairRoot v18.8
536 // fRun->SetSink(std::make_unique<FairRootFileSink>(fOutFileName));
537 fRun->SetSink(std::make_unique<FairRootFileSink>(fOutFileName).release());
538
539 // --- Create and register the setup modules, field and media with FairRoot
541
542 // --- Create and register the target
543 if (fTarget) {
544 LOG(info) << fTarget->ToString();
545 fRun->AddModule(fTarget.get());
546 }
547 else
548 LOG(warn) << GetName() << ": No target defined!";
549
550
551 // --- Create the magnetic field
552 LOG(info) << GetName() << ": Register magnetic field";
553 LOG(info) << fField;
555 fField->Print("");
556 fRun->SetField(fField);
557
558
559 // --- Initialise the event generator
561
562 // --- Trigger generation of run info
563 fRun->SetGenerateRunInfo(fGenerateRunInfo);
564
565
566 // --- Trigger storage of trajectories, if chosen
567 fRun->SetStoreTraj(fStoreTrajectories);
568
569
570 // --- Set VMC configuration
571 std::function<void()> f = std::bind(&CbmTransport::ConfigureVMC, this);
572 fRun->SetSimSetup(f);
573
574
575 // --- Set event filter task
576 fRun->AddTask(fEventFilter.get());
577
578
579 // --- Initialise run
580 fRun->Init();
581
582
583 // --- Force user-defined decays. This has to happen after FairRunSim::Init()
584 // because otherwise there seem to be no particles in GEANT.
586
587
588 // --- Set correct decay modes for pi0 and eta
589 // --- This is needed only when using Geant3
590 if (fEngine == kGeant3) PiAndEtaDecay(gMC);
591
592
593 // --- Runtime database
594 FairRuntimeDb* rtdb = fRun->GetRuntimeDb();
595 CbmFieldPar* fieldPar = static_cast<CbmFieldPar*>(rtdb->getContainer("CbmFieldPar"));
596 fieldPar->SetParameters(fField);
597 fieldPar->setChanged();
598 FairParRootFileIo* parOut = new FairParRootFileIo(kTRUE);
599 parOut->open(fParFileName.Data());
600 rtdb->setOutput(parOut);
601 rtdb->saveOutput();
602 rtdb->print();
603
604
605 // --- Measure time for initialisation
606 timer.Stop();
607 fRealTimeInit = timer.RealTime();
608
609 TFile* old = gFile;
610 TDirectory* oldDir = gDirectory;
611
612 // Write Transport Settings to the output file
613 auto sink = fRun->GetSink();
614 assert(sink->GetSinkType() == kFILESINK);
615 auto rootFileSink = static_cast<FairRootFileSink*>(sink);
616 TFile* outfile = rootFileSink->GetRootFile();
617 ;
618 outfile->cd();
619
620 LOG(info) << "Here I am";
621 if (fEngine == kGeant3) {
622 LOG(info) << "Write Geant3Settings";
623 fGeant3Settings->Write();
624 }
625 else if (fEngine == kGeant4) {
626 LOG(info) << "Write Geant4Settings";
627 fGeant4Settings->Write();
628 }
629 gFile = old;
630 gDirectory = oldDir;
631
632 // --- Start run
633 timer.Start(kFALSE); // without reset
634 fRun->Run(nEvents);
635 timer.Stop();
636 fRealTimeRun = timer.RealTime() - fRealTimeInit;
637 fCpuTime = timer.CpuTime();
638
639 // --- Create a geometry file if required
640 if (!fGeoFileName.IsNull()) fRun->CreateGeometryFile(fGeoFileName);
641
642
643 // --- Screen log
644 std::cout << std::endl;
645 LOG(info) << GetName() << ": Run finished successfully.";
646 LOG(info) << GetName() << ": Wall time for Init : " << fRealTimeInit << " s ";
647 LOG(info) << GetName() << ": Wall time for Run : " << fRealTimeRun << " s ("
648 << fRealTimeRun / fEventFilter->GetNofInputEvents() << " s / event)";
649 LOG(info) << GetName() << ": Output file : " << fOutFileName;
650 LOG(info) << GetName() << ": Parameter file : " << fParFileName;
651 if (!fGeoFileName.IsNull()) LOG(info) << GetName() << ": Geometry file : " << fGeoFileName;
652 std::cout << std::endl << std::endl;
653
654
655 // --- Remove TGeoManager
656 // To avoid crashes when exiting. Reason for this behaviour is unknown.
657 if (gGeoManager) {
658 if (gROOT->GetVersionInt() >= 60602) {
659 gGeoManager->GetListOfVolumes()->Delete();
660 gGeoManager->GetListOfShapes()->Delete();
661 }
662 delete gGeoManager;
663 }
664}
665// --------------------------------------------------------------------------
666
667
668// ----- Set the beam angle distribution --------------------------------
669void CbmTransport::SetBeamAngle(Double_t x0, Double_t y0, Double_t sigmaX, Double_t sigmaY)
670{
671 assert(fEventGen);
672 fEventGen->SetBeamAngle(x0, y0, sigmaX, sigmaY);
673}
674// --------------------------------------------------------------------------
675
676
677// ----- Set the beam position ------------------------------------------
678void CbmTransport::SetBeamPosition(Double_t x0, Double_t y0, Double_t sigmaX, Double_t sigmaY, Double_t z)
679{
680 assert(fEventGen);
681 fEventGen->SetBeamPosition(x0, y0, sigmaX, sigmaY, z);
682}
683// --------------------------------------------------------------------------
684
685
686// ----- Set a decay mode -----------------------------------------------
687void CbmTransport::SetDecayMode(Int_t pdg, UInt_t nDaughters, Int_t* daughterPdg)
688{
689
690 if (fDecayModes.count(pdg) != 0) {
691 LOG(fatal) << GetName() << ": User decay mode for PDG " << pdg << " is already defined!";
692 return;
693 }
694
695 for (UInt_t iDaughter = 0; iDaughter < nDaughters; iDaughter++) {
696 fDecayModes[pdg].push_back(daughterPdg[iDaughter]);
697 }
698}
699// --------------------------------------------------------------------------
700
701
702// ----- Set geometry file name -----------------------------------------
703void CbmTransport::SetGeoFileName(TString fileName)
704{
705
706 // Check for the directory
707 std::string name = fileName.Data();
708 Int_t found = name.find_last_of("/");
709 if (found >= 0) {
710 TString geoDir = name.substr(0, found);
711 if (gSystem->AccessPathName(geoDir.Data())) {
712 LOG(error) << GetName() << ": Directory for geometry file " << geoDir
713 << " does not exist; the file will not be created.";
714 return;
715 } //? Directory of geometry file does not exist
716 } //? File name contains directory path
717
718 fGeoFileName = fileName;
719}
720// --------------------------------------------------------------------------
721
722
723// ----- Set parameter file name ----------------------------------------
724void CbmTransport::SetParFileName(TString fileName)
725{
726
727 // --- If file does not exist, check the directory
728 if (gSystem->AccessPathName(fileName)) {
729 std::string name = fileName.Data();
730 Int_t found = name.find_last_of("/");
731 if (found >= 0) {
732 TString parDir = name.substr(0, found);
733 if (gSystem->AccessPathName(parDir.Data())) {
734 LOG(fatal) << GetName() << ": Parameter directory " << parDir << " does not exist!";
735 return;
736 } //? Directory of parameter file does not exist
737 } //? File name contains directory path
738 } //? Parameter file does not exist
739
740 fParFileName = fileName;
741}
742// --------------------------------------------------------------------------
743
744
745// ----- Set random event plane generation ------------------------------
746void CbmTransport::SetRandomEventPlane(Double_t phiMin, Double_t phiMax) { fEventGen->SetEventPlane(phiMin, phiMax); }
747// --------------------------------------------------------------------------
748
749
750// ----- Set output file name -------------------------------------------
751void CbmTransport::SetOutFileName(TString fileName, Bool_t overwrite)
752{
753
754 // --- Protect against overwriting an existing file
755 if ((!gSystem->AccessPathName(fileName.Data())) && (!overwrite)) {
756 LOG(fatal) << fName << ": output file " << fileName << " already exists!";
757 return;
758 }
759
760 // --- If the directory does not yet exist, create it
761 const char* directory = gSystem->DirName(fileName.Data());
762 if (gSystem->AccessPathName(directory)) {
763 Int_t success = gSystem->mkdir(directory, kTRUE);
764 if (success == -1)
765 LOG(fatal) << fName << ": output directory " << directory << " does not exist and cannot be created!";
766 else
767 LOG(info) << fName << ": created directory " << directory;
768 }
769
770 fOutFileName = fileName;
771}
772// --------------------------------------------------------------------------
773
774
775// ----- Define the target ----------------------------------------------
776void CbmTransport::SetTarget(const char* medium, Double_t thickness, Double_t diameter, Double_t x, Double_t y,
777 Double_t z, Double_t angle, Double_t density)
778{
779 if (density >= 0.) {
780 fTarget.reset(new CbmTarget(medium, thickness, diameter, density));
781 }
782 else {
783 fTarget.reset(new CbmTarget(medium, thickness, diameter));
784 }
785 fTarget->SetPosition(x, y, z);
786 fTarget->SetRotation(angle);
787}
788// --------------------------------------------------------------------------
789
790
791// ----- Enable vertex distribution in x and y --------------------------
793{
794 assert(fEventGen);
795 fEventGen->SmearGausVertexXY(choice);
796}
797// --------------------------------------------------------------------------
798
799
800// ----- Enable vertex distribution z -----------------------------------
802{
803 assert(fEventGen);
804 fEventGen->SmearVertexZ(choice);
805}
806// --------------------------------------------------------------------------
807
ECbmSetupSource
Definition CbmSetup.h:34
ClassImp(CbmTransport)
ECbmGenerator
@ kPluto
@ kUnigen
@ kUrqmd
@ kGeant3
@ kGeant4
friend fvec log(const fvec &a)
Bool_t CheckWithTarget(const CbmTarget &target) const
Check consistency with a target.
CbmEventGenerator.
virtual void Print(Option_t *opt="") const
Log output.
void SetBeamAngle(Double_t meanThetaX, Double_t meanThetaY, Double_t sigmaThetaX=-1., Double_t sigmaThetaY=-1.)
Set the beam angle in the focal plane.
void SetTarget(std::shared_ptr< const CbmTarget > target)
Set target properties.
void ForceVertexAtZ(Double_t zVertex)
Force event vertex to be at a given z.
const CbmBeamProfile & GetBeamProfile()
Beam profile.
void SetBeamPosition(Double_t meanX, Double_t meanY, Double_t sigmaX=-1., Double_t sigmaY=-1., Double_t zF=0.)
Set the beam position in the focal plane.
void ForceVertexInTarget(Bool_t choice=kTRUE)
Enable or disable forcing the vertex to be in the target.
void SetParameters(FairField *field)
User interface class to define the Geant3 simulation settings.
void Init(TVirtualMC *)
Set all parameters defined in this class.
User interface class to define the Geant4 simulation settings.
std::array< std::string, 3 > GetG4RunConfig()
Get the Geant4 run configuration.
void Init(TVirtualMC *)
Set all parameters defined in this class.
Class deciding whether to store an MC event.
Int_t GetNumAvailableEvents()
Get the maximum number of events available in the input file.
void SetSetupSource(ECbmSetupSource setupSource)
Set the source the setup will be loaded from.
Definition CbmSetup.cxx:263
CbmFieldMap * CreateFieldMap()
Definition CbmSetup.cxx:66
void LoadSetup(const char *setupName)
Definition CbmSetup.h:64
void RegisterSetup()
Definition CbmSetup.h:69
Class for filtering the stack before writing.
Class for constructing the geometry of the CBM target.
Definition CbmTarget.h:37
User interface class for transport simulation.
void SetDecayMode(Int_t pdg, UInt_t nDaughters, Int_t *daughterPdg)
Set a decay mode for a particle.
void SetTarget(const char *medium="Gold", Double_t thickness=0.025, Double_t diameter=2.5, Double_t x=0., Double_t y=0., Double_t z=-44.0, Double_t rot=0., Double_t density=-1)
Define the target.
std::shared_ptr< CbmTarget > fTarget
void RegisterRadLength(Bool_t choice=kTRUE)
Enable registration of radiation length.
virtual ~CbmTransport()
Destructor
void SetOutFileName(TString name, Bool_t overwrite=kFALSE)
Set the output file name.
CbmTransport()
Constructor.
CbmGeant3Settings * fGeant3Settings
std::unique_ptr< CbmMCEventFilter > fEventFilter
FairRunSim * fRun
Double_t fRealTimeInit
void SetGeoFileName(TString name)
Define geometry file name (output)
void LoadSetup(const char *setupName)
Use a standard setup.
void PiAndEtaDecay(TVirtualMC *vmc)
Correct decay modes for pi0 and eta.
void SetBeamPosition(Double_t x0, Double_t y0, Double_t sigmaX=-1., Double_t sigmaY=-1., Double_t zF=0.)
Set the beam position.
void SetRandomEventPlane(Double_t phiMin=0., Double_t phiMax=2. *TMath::Pi())
Activate random event plane.
ULong_t fRandomSeed
Double_t fCpuTime
void AddInput(const char *fileName, ECbmGenerator generator=kUnigen)
Add an input by file name and generator type.
void SetVertexSmearZ(Bool_t choice=kTRUE)
Enable smearing of event vertex in z.
void ConfigureVMC()
Set the parameters for the TVirtualMC.
void ForceVertexAtZ(Double_t zVertex)
Force the event vertex to be at a given z position.
void RegisterIons()
Register ions.
void SetParFileName(TString name)
Define parameter file name.
std::unique_ptr< CbmStackFilter > fStackFilter
TString fParFileName
TString fOutFileName
CbmSetup * fSetup
void SetVertexSmearXY(Bool_t choice=kTRUE)
Enable smearing of event vertex in x and y.
void ForceUserDecays()
Force user-defined single-mode decays.
Bool_t fStoreTrajectories
ECbmEngine fEngine
FairField * fField
void InitEventGenerator()
Event generator initialisation.
Double_t fRealTimeRun
void SetSetupSource(ECbmSetupSource setupSource)
Set the source the setup will be loaded from.
void RegisterSetup()
Create and register the setup modules.
CbmGeant4Settings * fGeant4Settings
void Run(Int_t nEvents)
Execute transport run.
TString fGeoFileName
Bool_t fGenerateRunInfo
std::map< Int_t, std::vector< Int_t > > fDecayModes
CbmEventGenerator * fEventGen
void SetBeamAngle(Double_t x0, Double_t y0, Double_t sigmaX=-1., Double_t sigmaY=-1.)
Set the beam angle (emittency at the beam position)
void ForceVertexInTarget(Bool_t choice=kTRUE)
Enable or disable forcing the vertex to be in the target.
Generates input to transport simulation from files in Unigen format.
Int_t GetNumAvailableEvents()
Get the maximum number of events available in the input file.
bool IsRootFile(std::string filename)
Check if a the file exist and is a root file.