BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSBeamPipe.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2023.
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 "BDSBeamPipe.hh"
20#include "BDSUtilities.hh"
21
22#include "globals.hh" // geant4 globals / types
23#include "G4VSolid.hh"
24
25BDSBeamPipe::BDSBeamPipe(G4VSolid* containerSolidIn,
26 G4LogicalVolume* containerLVIn,
27 BDSExtent extentIn,
28 G4VSolid* containerSubtractionSolidIn,
29 G4LogicalVolume* vacuumLVIn,
30 G4bool containerIsCircularIn,
31 G4double containerRadiusIn,
32 G4ThreeVector inputFaceNormalIn,
33 G4ThreeVector outputFaceNormalIn):
34 BDSGeometryComponent(containerSolidIn, containerLVIn, extentIn),
35 containerSubtractionSolid(containerSubtractionSolidIn),
36 vacuumLogicalVolume(vacuumLVIn),
37 containerIsCircular(containerIsCircularIn),
38 containerRadius(containerRadiusIn),
39 inputFaceNormal(inputFaceNormalIn.unit()),
40 outputFaceNormal(outputFaceNormalIn.unit())
41{;}
42
43BDSBeamPipe::~BDSBeamPipe()
44{
45 // we leak containerSubtractionSolid as it's slow to delete lots
46 // of them because of g4 deregistering them - g4 will clean up
47 //delete containerSubtractionSolid;
48}
49
50std::set<G4LogicalVolume*> BDSBeamPipe::GetVolumesForField() const
51{
52 std::set<G4LogicalVolume*> result = GetAllLogicalVolumes(); // from base class
53 // We avoid setting the field on a tight-fitting container volume to avoid strong
54 // field in small gaps. However, in the case of only a container volume, we should
55 // return that - for when there's no beam pipe geometry.
56 if (result.size() > 1)
57 {result.erase(GetContainerLogicalVolume());}
58 else if (!BDS::IsFinite((G4double)result.size()))
59 {result = {GetContainerLogicalVolume()};}
60
61 return result;
62}
BDSBeamPipe(G4VSolid *containerSolidIn, G4LogicalVolume *containerLVIn, BDSExtent extentIn, G4VSolid *containerSubtractionSolidIn, G4LogicalVolume *vacuumLVIn, G4bool containerIsCircularIn=false, G4double containerRadiusIn=0.0, G4ThreeVector inputFaceNormalIn=G4ThreeVector(0, 0,-1), G4ThreeVector outputFaceNormalIn=G4ThreeVector(0, 0, 1))
Definition: BDSBeamPipe.cc:25
std::set< G4LogicalVolume * > GetVolumesForField() const
Definition: BDSBeamPipe.cc:50
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
A generic geometry component for a bdsim model.
G4LogicalVolume * GetContainerLogicalVolume() const
Accessor - see member for more info.
virtual std::set< G4LogicalVolume * > GetAllLogicalVolumes() const
Access all logical volumes belonging to this component.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())