Skip to content

Commit 3204b13

Browse files
committed
feat(r3bgen):Added new function to R3BINCLRootGenerator and new CI test
feat(input):Added input/p_U238_500.root
1 parent 76109ec commit 3204b13

6 files changed

Lines changed: 187 additions & 29 deletions

File tree

input/p_U238_500.root

25.7 KB
Binary file not shown.

r3bgen/R3BINCLRootGenerator.cxx

Lines changed: 14 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -159,6 +159,20 @@ bool R3BINCLRootGenerator::ReadEvent(FairPrimaryGenerator* primGen)
159159
px = pt * TMath::Cos(fPhi[j] * TMath::DegToRad());
160160
py = pt * TMath::Sin(fPhi[j] * TMath::DegToRad());
161161
}
162+
163+
// Rotation according to the beam direction
164+
if (fRotXBeam != 0 || fRotYBeam != 0)
165+
{
166+
TVector3 ptotal(px, py, pz);
167+
auto rotx = fRotXBeam * TMath::DegToRad();
168+
auto roty = fRotYBeam * TMath::DegToRad();
169+
ptotal.RotateY(roty);
170+
ptotal.RotateX(rotx);
171+
px = ptotal.X();
172+
py = ptotal.Y();
173+
pz = ptotal.Z();
174+
}
175+
162176
R3BLOG(debug, "PDG:Px:Py:Pz " << pdg << " " << px << " " << py << " " << pz);
163177
if (fOnlyFragments)
164178
{

r3bgen/R3BINCLRootGenerator.h

Lines changed: 8 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -81,6 +81,12 @@ class R3BINCLRootGenerator : public FairGenerator
8181
void SetOnlyfragments(bool Opt = true) { fOnlyFragments = Opt; }
8282
void SetMinPdgCode(int Opt) { fPdgCodeMin = Opt; }
8383

84+
/**
85+
** Methods to rotate the particles according to the beam direction (deg)
86+
**/
87+
void SetRotationX(double rot) { fRotXBeam = rot; }
88+
void SetRotationY(double rot) { fRotYBeam = rot; }
89+
8490
private:
8591
TString fFileName; // Input file name
8692
TFile* fInput;
@@ -89,6 +95,8 @@ class R3BINCLRootGenerator : public FairGenerator
8995
bool fOnlySpallation = false; // True if we want to simulate only spallation events
9096
bool fOnlyFragments = false; // True if we want to simulate only fragments
9197
int fPdgCodeMin = 1000050070; // Limit in Boro-7
98+
double fRotXBeam = 0.;
99+
double fRotYBeam = 0.;
92100

93101
/** Private method RegisterIons. Goes through the input file and registers
94102
** any ion needed. TODO: Should not be needed by FairRoot. **/

r3bgen/test/CMakeLists.txt

Lines changed: 18 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -1,4 +1,22 @@
1+
##############################################################################
2+
# Copyright (C) 2020 GSI Helmholtzzentrum für Schwerionenforschung GmbH #
3+
# Copyright (C) 2020-2026 Members of R3B Collaboration #
4+
# #
5+
# This software is distributed under the terms of the #
6+
# GNU General Public Licence (GPL) version 3, #
7+
# copied verbatim in the file "LICENSE". #
8+
# #
9+
# In applying this license GSI does not waive the privileges and immunities #
10+
# granted to it by virtue of its status as an Intergovernmental Organization #
11+
# or submit itself to any jurisdiction. #
12+
##############################################################################
13+
114
generate_root_test_script(${R3BROOT_SOURCE_DIR}/r3bgen/test/testR3BPhaseSpaceGeneratorIntegration.C)
215
add_test(NAME testR3BPhaseSpaceGeneratorIntegration COMMAND ${R3BROOT_BINARY_DIR}/r3bgen/test/testR3BPhaseSpaceGeneratorIntegration.sh)
316
set_tests_properties(testR3BPhaseSpaceGeneratorIntegration PROPERTIES TIMEOUT "100")
417
set_tests_properties(testR3BPhaseSpaceGeneratorIntegration PROPERTIES PASS_REGULAR_EXPRESSION "Macro finished successfully.")
18+
19+
generate_root_test_script(${R3BROOT_SOURCE_DIR}/r3bgen/test/testINCLRootGenerator.C)
20+
add_test(NAME testINCLRootGenerator COMMAND ${R3BROOT_BINARY_DIR}/r3bgen/test/testINCLRootGenerator.sh)
21+
set_tests_properties(testINCLRootGenerator PROPERTIES TIMEOUT "100")
22+
set_tests_properties(testINCLRootGenerator PROPERTIES PASS_REGULAR_EXPRESSION "Macro finished successfully.")
Lines changed: 110 additions & 0 deletions
Original file line numberDiff line numberDiff line change
@@ -0,0 +1,110 @@
1+
/******************************************************************************
2+
* Copyright (C) 2025 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3+
* Copyright (C) 2025-2026 Members of R3B Collaboration *
4+
* *
5+
* This software is distributed under the terms of the *
6+
* GNU General Public Licence (GPL) version 3, *
7+
* copied verbatim in the file "LICENSE". *
8+
* *
9+
* In applying this license GSI does not waive the privileges and immunities *
10+
* granted to it by virtue of its status as an Intergovernmental Organization *
11+
* or submit itself to any jurisdiction. *
12+
******************************************************************************/
13+
14+
#include <TStopwatch.h>
15+
#include <TString.h>
16+
#include <TSystem.h>
17+
#include <memory>
18+
19+
void testINCLRootGenerator(const int nbevents = 4)
20+
{
21+
// Timer
22+
TStopwatch timer;
23+
timer.Start();
24+
25+
// Logging
26+
auto logger = FairLogger::GetLogger();
27+
logger->SetLogVerbosityLevel("low");
28+
logger->SetLogScreenLevel("warn");
29+
logger->SetColoredLog(true);
30+
31+
// System paths
32+
const TString workDirectory = getenv("VMCWORKDIR");
33+
gSystem->Setenv("GEOMPATH", workDirectory + "/geometry");
34+
gSystem->Setenv("CONFIG_DIR", workDirectory + "/gconfig");
35+
36+
// Input file
37+
const TString INCLinput = workDirectory + "/input/p_U238_500.root";
38+
39+
// Output files
40+
const TString simufile = "incl.simu.root";
41+
const TString parafile = "incl.para.root";
42+
43+
// Input GLAD geometry
44+
const TString fGladGeo = "glad_v2025.1.geo.root";
45+
46+
// Basic simulation setup
47+
auto run = std::make_unique<FairRunSim>();
48+
run->SetName("TGeant4");
49+
run->SetStoreTraj(false);
50+
run->SetMaterials("media_r3b.geo");
51+
52+
auto config = std::make_unique<FairGenericVMCConfig>();
53+
run->SetSimulationConfig(std::move(config));
54+
run->SetSink(std::make_unique<FairRootFileSink>(simufile.Data()));
55+
56+
// ----- Runtime data base --------------------------------------------
57+
auto* rtdb = run->GetRuntimeDb();
58+
UInt_t runId = 1;
59+
rtdb->initContainers(runId);
60+
61+
// Primary particle generator
62+
auto inclGen = std::make_unique<R3BINCLRootGenerator>(INCLinput.Data());
63+
inclGen->SetOnlyFission();
64+
inclGen->SetOnlyfragments();
65+
inclGen->SetRotationY(1.); // deg
66+
inclGen->SetXYZ(0., 0., -138.); // cm
67+
inclGen->SetDxDyDz(0.5, 0.5, 0.5); // cm
68+
69+
auto primGen = std::make_unique<FairPrimaryGenerator>();
70+
primGen->AddGenerator(inclGen.release());
71+
run->SetGenerator(primGen.release());
72+
73+
// Geometry: Cave
74+
auto cave = std::make_unique<R3BCave>("CAVE");
75+
cave->SetGeometryFileName("r3b_cave_vacuum.geo");
76+
run->AddModule(cave.release());
77+
78+
// Geometry: GLAD
79+
run->AddModule(new R3BGladMagnet(fGladGeo.Data()));
80+
81+
// GLAD Filed
82+
auto* GladField = new R3BGladFieldMap("R3BGladMap");
83+
GladField->SetFieldfromCurrent(2600.); // Current in Amperes
84+
run->SetField(GladField);
85+
86+
// Init
87+
run->Init();
88+
89+
// Save field parameters
90+
auto* fieldPar = dynamic_cast<R3BFieldPar*>(rtdb->getContainer("R3BFieldPar"));
91+
fieldPar->SetParameters(GladField);
92+
fieldPar->setChanged();
93+
94+
// Output file with parameters
95+
auto parOut = std::make_unique<FairParRootFileIo>(true);
96+
parOut->open(parafile.Data());
97+
rtdb->setOutput(parOut.release());
98+
rtdb->saveOutput();
99+
rtdb->print();
100+
101+
// Simulate
102+
if (nbevents > 0)
103+
run->Run(nbevents);
104+
105+
// Report
106+
timer.Stop();
107+
std::cout << "Real time: " << timer.RealTime() << "s, CPU time: " << timer.CpuTime() << "s" << std::endl;
108+
std::cout << "Macro finished successfully." << std::endl;
109+
gApplication->Terminate();
110+
}
Lines changed: 37 additions & 29 deletions
Original file line numberDiff line numberDiff line change
@@ -1,6 +1,6 @@
11
/******************************************************************************
2-
* Copyright (C) 2019 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3-
* Copyright (C) 2019-2025 Members of R3B Collaboration *
2+
* Copyright (C) 2020 GSI Helmholtzzentrum für Schwerionenforschung GmbH *
3+
* Copyright (C) 2020-2026 Members of R3B Collaboration *
44
* *
55
* This software is distributed under the terms of the *
66
* GNU General Public Licence (GPL) version 3, *
@@ -17,7 +17,7 @@
1717
#include <iostream>
1818
#include <memory>
1919

20-
void testR3BPhaseSpaceGeneratorIntegration()
20+
void testR3BPhaseSpaceGeneratorIntegration(const int nbevents = 10)
2121
{
2222
// Timer
2323
TStopwatch timer;
@@ -34,58 +34,65 @@ void testR3BPhaseSpaceGeneratorIntegration()
3434
gSystem->Setenv("GEOMPATH", workDirectory + "/geometry");
3535
gSystem->Setenv("CONFIG_DIR", workDirectory + "/gconfig");
3636

37-
// Simulation Base
37+
// Output files
38+
const TString simufile = "PhaseSpace.simu.root";
39+
const TString parafile = "PhaseSpace.para.root";
40+
41+
// Basic simulation setup
3842
FairRunSim run;
3943
run.SetName("TGeant4");
40-
run.SetOutputFile("testR3BPhaseSpaceGeneratorIntegration.root");
44+
run.SetSink(std::make_unique<FairRootFileSink>(simufile.Data()));
4145
run.SetMaterials("media_r3b.geo");
4246

43-
// World
44-
auto cave = new R3BCave("Cave");
47+
// ----- Runtime data base --------------------------------------------
48+
auto* rtdb = run.GetRuntimeDb();
49+
UInt_t runId = 1;
50+
rtdb->initContainers(runId);
51+
52+
// Geometry: Cave
53+
auto cave = std::make_unique<R3BCave>("CAVE");
4554
cave->SetGeometryFileName("r3b_cave_vacuum.geo");
46-
run.AddModule(cave);
55+
run.AddModule(cave.release());
4756

4857
// Magnet
49-
// run.AddModule(new R3BGladMagnet("glad_v2023.1.geo.root"));
50-
auto magField = new R3BGladFieldMap("R3BGladMap");
51-
magField->SetScale(-0.6);
52-
run.SetField(magField);
58+
auto GladField = new R3BGladFieldMap("R3BGladMap");
59+
GladField->SetScale(-0.6);
60+
run.SetField(GladField);
5361

5462
// Primaries
55-
auto primGen = new FairPrimaryGenerator();
63+
auto primGen = std::make_unique<FairPrimaryGenerator>();
5664
auto gen = new R3BPhaseSpaceGenerator();
5765
gen->GetBeam().SetEnergyDistribution(R3BDistribution1D::Delta(600));
5866
gen->SetErelDistribution(R3BDistribution1D::Delta(100));
5967
gen->AddParticle(5, 17);
6068
gen->AddNeutron();
6169
gen->AddProton();
6270
primGen->AddGenerator(gen);
63-
run.SetGenerator(primGen);
71+
run.SetGenerator(primGen.release());
6472

6573
// Logging
6674
run.SetStoreTraj(false);
67-
FairLogger::GetLogger()->SetLogVerbosityLevel("low");
68-
FairLogger::GetLogger()->SetLogScreenLevel("info");
6975

7076
// Init & Special MC Settings
7177
run.Init();
7278

73-
// Database
74-
auto rtdb = run.GetRuntimeDb();
75-
auto fieldPar = (R3BFieldPar*)rtdb->getContainer("R3BFieldPar");
76-
fieldPar->SetParameters(magField);
79+
// Save field parameters
80+
auto* fieldPar = dynamic_cast<R3BFieldPar*>(rtdb->getContainer("R3BFieldPar"));
81+
fieldPar->SetParameters(GladField);
7782
fieldPar->setChanged();
78-
auto parOut = new FairParRootFileIo(true);
79-
parOut->open("testR3BPhaseSpaceGeneratorIntegration.para.root");
80-
rtdb->setOutput(parOut);
83+
84+
// Output file with parameters
85+
auto parOut = std::make_unique<FairParRootFileIo>(true);
86+
parOut->open(parafile.Data());
87+
rtdb->setOutput(parOut.release());
8188
rtdb->saveOutput();
8289
rtdb->print();
8390

84-
// Go
85-
run.Run(10);
91+
// Simulate
92+
run.Run(nbevents);
8693

87-
// Test Output
88-
auto file = TFile::Open("testR3BPhaseSpaceGeneratorIntegration.root");
94+
// Report
95+
auto file = TFile::Open("PhaseSpace.simu.root");
8996
auto tree = (TTree*)file->Get("evt");
9097
auto mctc = new TClonesArray("R3BMCTrack");
9198
tree->SetBranchAddress("MCTrack", &mctc);
@@ -97,14 +104,14 @@ void testR3BPhaseSpaceGeneratorIntegration()
97104
return;
98105
}
99106

100-
auto track = (R3BMCTrack*)mctc->At(1);
107+
auto track = dynamic_cast<R3BMCTrack*>(mctc->At(1));
101108
if (track->GetPdgCode() != 2112 || track->GetMotherId() != -1)
102109
{
103110
std::cout << "Not the correct primary particle" << std::endl;
104111
return;
105112
}
106113

107-
const Double_t ekin = track->GetEnergy() - track->GetMass();
114+
const auto ekin = track->GetEnergy() - track->GetMass();
108115
std::cout << "Ekin of primary neutron:" << ekin << std::endl;
109116
if (abs(ekin - 0.6) > 0.02)
110117
{
@@ -116,4 +123,5 @@ void testR3BPhaseSpaceGeneratorIntegration()
116123
timer.Stop();
117124
std::cout << "Real time: " << timer.RealTime() << "s, CPU time: " << timer.CpuTime() << "s" << std::endl;
118125
std::cout << "Macro finished successfully." << std::endl;
126+
gApplication->Terminate();
119127
}

0 commit comments

Comments
 (0)