BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSParallelWorldUtilities.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSAcceleratorModel.hh"
20#include "BDSDebug.hh"
21#include "BDSException.hh"
22#include "BDSGlobalConstants.hh"
23#include "BDSLinkDetectorConstruction.hh"
24#include "BDSParallelWorldCurvilinear.hh"
25#include "BDSParallelWorldCurvilinearBridge.hh"
26#include "BDSParallelWorldImportance.hh"
27#include "BDSParallelWorldInfo.hh"
28#include "BDSParallelWorldPlacementFields.hh"
29#include "BDSParallelWorldSampler.hh"
30#include "BDSParallelWorldUtilities.hh"
31#include "BDSParser.hh"
32#include "BDSSamplerType.hh"
33#include "BDSUtilities.hh"
34
35#include "parser/element.h"
36#include "parser/fastlist.h"
37#include "parser/placement.h"
38
39#include "G4GeometrySampler.hh"
40#include "G4ImportanceBiasing.hh"
41#include "G4IStore.hh"
42#include "G4ParallelWorldPhysics.hh"
43#include "G4VModularPhysicsList.hh"
44#include "G4VUserDetectorConstruction.hh"
45#include "G4VUserParallelWorld.hh"
46
47#include "globals.hh"
48
49#include <map>
50#include <vector>
51
52std::vector<BDSParallelWorldInfo> BDS::NumberOfExtraWorldsRequired()
53{
54 std::vector<BDSParallelWorldInfo> worlds;
55
56 const std::vector<GMAD::Placement> placements = BDSParser::Instance()->GetPlacements();
57 for (const auto& pl : placements)
58 {// loop over all placement definitions and see if any are beam line placements
59 if (!pl.sequence.empty()) // it's a beam line placement
60 {// curvilinear world required at least
61 G4bool samplerWorldRequired = false;
62 for (const auto& el : BDSParser::Instance()->GetSequence(pl.sequence))
63 {// do any elements require a sampler -> sampler world required
64 BDSSamplerType sType = BDS::DetermineSamplerType(el.samplerType);
65 if (sType != BDSSamplerType::none)
66 {
67 samplerWorldRequired = true;
68 break; // no need to loop over rest of sequence
69 }
70 }
71 BDSParallelWorldInfo info(pl.sequence, true, samplerWorldRequired, pl.name);
72 worlds.push_back(info);
73 }
74 }
75 return worlds;
76}
77
78std::vector<G4VUserParallelWorld*> BDS::ConstructAndRegisterParallelWorlds(G4VUserDetectorConstruction* massWorld,
79 G4bool buildSamplerWorld,
80 G4bool buildPlacementFieldsWorld)
81{
82 BDSAcceleratorModel* acceleratorModel = BDSAcceleratorModel::Instance();
83
84 // registry of all created worlds that require the physics process so
85 // that their boundaries affect tracking
86 std::vector<G4VUserParallelWorld*> worldsRequiringPhysics;
87
88 // standard worlds
89 if (buildSamplerWorld) // optional
90 {
91 auto samplerWorld = new BDSParallelWorldSampler("main");
92 massWorld->RegisterParallelWorld(samplerWorld);
93 auto massWorldBDS = dynamic_cast<BDSLinkDetectorConstruction*>(massWorld);
94 if (massWorldBDS)
95 {
96 G4int samplerWorldID = massWorld->GetNumberOfParallelWorld() - 1;
97 massWorldBDS->SetSamplerWorldID(samplerWorldID);
98 }
99 acceleratorModel->RegisterParallelWorld(samplerWorld);
100 worldsRequiringPhysics.push_back(dynamic_cast<G4VUserParallelWorld*>(samplerWorld));
101 }
102
103 auto curvilinearWorld = new BDSParallelWorldCurvilinear("main");
104 auto curvilinearBridgeWorld = new BDSParallelWorldCurvilinearBridge("main");
105 massWorld->RegisterParallelWorld(curvilinearWorld);
106 massWorld->RegisterParallelWorld(curvilinearBridgeWorld);
107
108 // G4VUserDetectorConstruction doesn't delete parallel worlds so we should
109 acceleratorModel->RegisterParallelWorld(curvilinearWorld);
110 acceleratorModel->RegisterParallelWorld(curvilinearBridgeWorld);
111
112 // extra worlds for additional beam line placements
113 std::vector<BDSParallelWorldInfo> worldInfos = BDS::NumberOfExtraWorldsRequired();
114
115 for (const auto& info : worldInfos)
116 {
117 if (info.curvilinearWorld)
118 {
119 auto cLWorld = new BDSParallelWorldCurvilinear(info.curvilinearWorldName);
120 auto cLBridgeWorld = new BDSParallelWorldCurvilinearBridge(info.curvilinearWorldName);
121 massWorld->RegisterParallelWorld(cLWorld);
122 massWorld->RegisterParallelWorld(cLBridgeWorld);
123 acceleratorModel->RegisterParallelWorld(cLWorld);
124 acceleratorModel->RegisterParallelWorld(cLBridgeWorld);
125 }
126 if (info.samplerWorld)
127 {
128 auto sWorld = new BDSParallelWorldSampler(info.sequenceName);
129 acceleratorModel->RegisterParallelWorld(sWorld);
130 worldsRequiringPhysics.push_back(dynamic_cast<G4VUserParallelWorld*>(sWorld));
131 massWorld->RegisterParallelWorld(sWorld);
132 }
133 }
134
135 // only create the importance parallel world if the file is specified
136 if (BDSGlobalConstants::Instance()->UseImportanceSampling())
137 {
138 G4String importanceWorldGeometryFile = BDSGlobalConstants::Instance()->ImportanceWorldGeometryFile();
139 G4String importanceVolumeMapFile = BDSGlobalConstants::Instance()->ImportanceVolumeMapFile();
140 auto importanceWorld = new BDSParallelWorldImportance("main",
141 importanceWorldGeometryFile,
142 importanceVolumeMapFile);
143 acceleratorModel->RegisterParallelWorld(importanceWorld);
144 massWorld->RegisterParallelWorld(importanceWorld);
145 worldsRequiringPhysics.push_back(dynamic_cast<G4VUserParallelWorld*>(importanceWorld));
146 }
147
148 // optional parallel world for coordinate transforms for fields attached to placements of geometry
149 if (buildPlacementFieldsWorld)
150 {
151 auto placementFieldsPW = new BDSParallelWorldPlacementFields("placement_fields");
152 acceleratorModel->RegisterParallelWorld(placementFieldsPW);
153 massWorld->RegisterParallelWorld(placementFieldsPW);
154 }
155
156 return worldsRequiringPhysics;
157}
158
159std::vector<G4ParallelWorldPhysics*> BDS::ConstructParallelWorldPhysics(const std::vector<G4VUserParallelWorld*>& worlds)
160{
161 std::vector<G4ParallelWorldPhysics*> result;
162 for (auto world : worlds)
163 {result.push_back(new G4ParallelWorldPhysics(world->GetName()));}
164 return result;
165}
166
167void BDS::RegisterSamplerPhysics(const std::vector<G4ParallelWorldPhysics*>& processes,
168 G4VModularPhysicsList* physicsList)
169{
170 for (auto process : processes)
171 {physicsList->RegisterPhysics(process);}
172}
173
174void BDS::AddIStore(const std::vector<G4VUserParallelWorld*>& worlds)
175 {
177 //only add importance store if the world exists
178 if (importanceWorld)
179 {importanceWorld->AddIStore();}
180 else
181 {throw BDSException(__METHOD_NAME__, "Importance sampling world not found.");}
182 }
183
184void BDS::RegisterImportanceBiasing(const std::vector<G4VUserParallelWorld*>& worlds,
185 G4VModularPhysicsList* physList)
186{
188
189 // create world geometry sampler
190 G4GeometrySampler* pgs = new G4GeometrySampler(importanceWorld->GetWorldVolume(), "neutron");
191 pgs->SetParallel(true);
192 physList->RegisterPhysics(new G4ImportanceBiasing(pgs,importanceWorld->GetName()));
193}
194
195BDSParallelWorldImportance* BDS::GetImportanceSamplingWorld(const std::vector<G4VUserParallelWorld*>& worlds)
196{
197 // get importance world
198 G4String importanceWorldName = "importanceWorld_main";
199 G4VUserParallelWorld* importanceWorld = nullptr;
200 for (auto world : worlds)
201 {
202 if (std::strcmp(world->GetName(),importanceWorldName) == 0)
203 {importanceWorld = world; break;}
204 }
205 BDSParallelWorldImportance* iworld = dynamic_cast<BDSParallelWorldImportance*>(importanceWorld);
206 return iworld;
207}
A holder class for all representations of the accelerator model created in BDSIM.
General exception with possible name of object and message.
Definition: BDSException.hh:35
static BDSGlobalConstants * Instance()
Access method.
Construction of the geometry in the case of a link model.
A parallel world for bridiging curvilinear volumes.
A parallel world for curvilinear coordinates.
Class that constructs a parallel importance world.
void AddIStore()
Create IStore for all importance sampling geometry cells.
G4VPhysicalVolume * GetWorldVolume()
World volume getter required in parallel world utilities.
A simple class to hold what parallel geometry is required for each sequence.
A parallel world for placement field transform shapes.
A parallel world for sampler planes.
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
std::vector< GMAD::Placement > GetPlacements() const
Return the vector of placement objects.
Definition: BDSParser.hh:114
std::vector< G4ParallelWorldPhysics * > ConstructParallelWorldPhysics(const std::vector< G4VUserParallelWorld * > &worlds)
Construct the parallel physics process for each sampler world.
void RegisterSamplerPhysics(const std::vector< G4ParallelWorldPhysics * > &processes, G4VModularPhysicsList *physicsList)
Register each parallel physics process to the main physics list.
std::vector< BDSParallelWorldInfo > NumberOfExtraWorldsRequired()
void RegisterImportanceBiasing(const std::vector< G4VUserParallelWorld * > &worlds, G4VModularPhysicsList *physicsList)
Create importance geometry sampler and register importance biasing with physics list.
std::vector< G4VUserParallelWorld * > ConstructAndRegisterParallelWorlds(G4VUserDetectorConstruction *massWorld, G4bool buildSamplerWorld, G4bool buildPlacementFieldsWorld)
BDSParallelWorldImportance * GetImportanceSamplingWorld(const std::vector< G4VUserParallelWorld * > &worlds)
Get importance sampling world from list of all parallel worlds.
void AddIStore(const std::vector< G4VUserParallelWorld * > &worlds)
Get store, and prepare importance sampling for importance geometry sampler.