BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSSamplerCustom.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 "BDSApertureInfo.hh"
20#include "BDSDebug.hh"
21#include "BDSException.hh"
22#include "BDSGlobalConstants.hh"
23#include "BDSSamplerCustom.hh"
24
25#include "globals.hh" // geant4 types / globals
26#include "G4Box.hh"
27#include "G4LogicalVolume.hh"
28#include "G4Tubs.hh"
29
30G4double BDSSamplerCustom::chordLength = -1;
31
32BDSSamplerCustom::BDSSamplerCustom(const G4String& nameIn,
33 const BDSApertureInfo& shape,
34 G4int filterSetIDIn):
35 BDSSampler(nameIn, filterSetIDIn)
36{
37 /*
38 BDSApertureFactory fac;
39 containerSolid = fac.CreateAperture(name + "_aperture",
40 BDSSamplerPlane::chordLength,
41 shape);
42 */
43 // We make the sampler 10x bigger than normal as it's still really small
44 // but less likely to cause overlap problems. The original sampler width
45 // is designed to be functional but as small as possible to avoid introducing
46 // extra length for optical tracking.
47 switch (shape.apertureType.underlying())
48 {
49 case BDSApertureType::circular:
50 case BDSApertureType::circularvacuum:
51 {
52 containerSolid = new G4Tubs(name + "_solid",
53 0,
54 shape.aper1,
55 0.5*chordLength,
56 0,
57 CLHEP::twopi);
58 break;
59 }
60 case BDSApertureType::rectangular:
61 {
62 containerSolid = new G4Box(name + "_solid",
63 shape.aper1,
64 shape.aper2,
65 0.5*chordLength);
66 break;
67 }
68 default:
69 {
70 std::string msg = "Shape \"" + shape.apertureType.ToString() + "\" is not currently supported.\n";
71 msg += "Please use circular or rectangular.";
72 throw BDSException(__METHOD_NAME__, msg);
73 break;
74 }
75 }
76
77 BDSExtent ae = shape.Extent();
78 G4double dz = 0.5*chordLength;
79 SetExtent(BDSExtent(ae.XNeg(), ae.XPos(), ae.YNeg(), ae.YPos(), -dz, dz));
80
81 CommonConstruction();
82}
Holder class for all information required to describe an aperture.
G4double aper1
Public member for direct access.
G4double aper2
Public member for direct access.
BDSExtent Extent() const
BDSApertureType apertureType
Public member for direct access.
General exception with possible name of object and message.
Definition: BDSException.hh:35
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double XPos() const
Accessor.
Definition: BDSExtent.hh:66
G4double XNeg() const
Accessor.
Definition: BDSExtent.hh:65
G4double YNeg() const
Accessor.
Definition: BDSExtent.hh:67
G4double YPos() const
Accessor.
Definition: BDSExtent.hh:68
BDSSamplerCustom()=delete
No default constructor.
Base class and registry of sampler instances.
Definition: BDSSampler.hh:36
type underlying() const
return underlying value (can be used in switch statement)