BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSParallelWorldCurvilinear.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 "BDSAuxiliaryNavigator.hh"
21#include "BDSDebug.hh"
22#include "BDSDetectorConstruction.hh"
23#include "BDSGlobalConstants.hh"
24#include "BDSParallelWorldCurvilinear.hh"
25
26#include "G4String.hh"
27#include "G4LogicalVolume.hh"
28#include "G4VisAttributes.hh"
29#include "G4VPhysicalVolume.hh"
30
31class BDSBeamline;
32
34 G4VUserParallelWorld("CurvilinearWorld_" + name),
35 suffix(name),
36 clWorldVis(nullptr)
37{;}
38
39BDSParallelWorldCurvilinear::~BDSParallelWorldCurvilinear()
40{
41 delete clWorldVis;
42}
43
45{// general abbreviation: 'cl' == curvilinear
46#ifdef BDSDEBUG
47 G4cout << __METHOD_NAME__ << G4endl;
48#endif
49
50 G4VPhysicalVolume* clWorld = GetWorld();
51
52 // TODO - only register main one for now
53 // register read out world PV with our auxiliary navigator.
54 if (suffix == "main")
56
57 G4LogicalVolume* clWorldLV = clWorld->GetLogicalVolume();
58
59 // visualisation
61 clWorldVis = new G4VisAttributes(*(globals->VisibleDebugVisAttr()));
62 clWorldVis->SetForceWireframe(true);//just wireframe so we can see inside it
63 clWorldLV->SetVisAttributes(clWorldVis);
64
65 BDSBeamlineSet blSet = BDSAcceleratorModel::Instance()->BeamlineSet(suffix);
66
67 BDSDetectorConstruction::PlaceBeamlineInWorld(blSet.curvilinearWorld, clWorld,
68 globals->CheckOverlaps(), false, true, true);
69}
const BDSBeamlineSet & BeamlineSet(const G4String &name) const
Accessor.
static void AttachWorldVolumeToNavigatorCL(G4VPhysicalVolume *curvilinearWorldPVIn)
Simple struct to return a beamline plus associated beam lines.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
static void PlaceBeamlineInWorld(BDSBeamline *beamline, G4VPhysicalVolume *containerPV, G4bool checkOverlaps=false, G4bool setRegions=false, G4bool registerInfo=false, G4bool useCLPlacementTransform=false, G4bool useIncrementalCopyNumbers=false)
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
BDSParallelWorldCurvilinear()=delete
No default constructor.
G4String suffix
Just the input part of the name.
G4VisAttributes * clWorldVis
Visualisation attributes for the world volume.