BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSParallelWorldPlacementFields.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 "BDSDetectorConstruction.hh"
21#include "BDSException.hh"
22#include "BDSGlobalConstants.hh"
23#include "BDSNavigatorPlacements.hh"
24#include "BDSParallelWorldPlacementFields.hh"
25
26#include "G4LogicalVolume.hh"
27#include "G4PVPlacement.hh"
28#include "G4String.hh"
29#include "G4Types.hh"
30#include "G4VisAttributes.hh"
31#include "G4VPhysicalVolume.hh"
32
33#include <vector>
34
35
36BDSParallelWorldPlacementFields::BDSParallelWorldPlacementFields(const G4String& name):
37 G4VUserParallelWorld("PlacementFieldsWorld_" + name),
38 pfWorldVis(nullptr)
39{;}
40
41BDSParallelWorldPlacementFields::~BDSParallelWorldPlacementFields()
42{
43 // we leak all the placements as g4 is incredibly slow to delete
44 // them as it deregisters them - g4 will clean up
45 //for (auto placement : placements)
46 // {delete placement;}
47 delete pfWorldVis;
48}
49
51{
54 G4VPhysicalVolume* pfWorld = GetWorld();
55 G4LogicalVolume* pfWorldLV = pfWorld->GetLogicalVolume();
56
59
60 pfWorldVis = new G4VisAttributes(*(globals->VisibleDebugVisAttr()));
61 pfWorldVis->SetForceWireframe(true);//just wireframe so we can see inside it
62 pfWorldLV->SetVisAttributes(pfWorldVis);
63
64 G4bool checkOverlaps = globals->CheckOverlaps();
65 const auto& placementsToMake = BDSAcceleratorModel::Instance()->PlacementFieldPWPlacements();
66 for (const auto& placement : placementsToMake)
67 {
68 // record placements for cleaning up at destruction.
69 G4PVPlacement* pl = new G4PVPlacement(placement.transform,
70 placement.name,
71 placement.lvToPlace,
72 pfWorld,
73 false,
74 placement.copyNumber,
75 checkOverlaps);
76 placements.push_back(pl);
77 }
78}
const std::vector< BDSPlacementToMake > & PlacementFieldPWPlacements() const
Access field volumes for parallel world.
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
static void AttachWorldVolumeToNavigator(G4VPhysicalVolume *worldPVIn)
Setup the navigator w.r.t. to a world volume - typically real world.
std::vector< G4VPhysicalVolume * > placements
Cache of the placements to clean up at the end.