BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSCollimatorElliptical.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 "BDSCollimatorElliptical.hh"
20#include "BDSDebug.hh"
21#include "BDSUtilities.hh"
22#include "BDSWarning.hh"
23
24#include "globals.hh" // geant4 globals / types
25#include "G4EllipticalTube.hh"
26#include "G4EllipticalCone.hh"
27#include "G4VSolid.hh"
28
30 G4double lengthIn,
31 G4double horizontalWidthIn,
32 G4Material* collimatorMaterialIn,
33 G4Material* vacuumMaterialIn,
34 G4double xApertureIn,
35 G4double yApertureIn,
36 G4double xApertureOutIn,
37 G4double yApertureOutIn,
38 G4Colour* colourIn,
39 G4bool circularOuterIn):
40 BDSCollimator(nameIn, lengthIn, horizontalWidthIn, "ecol",
41 collimatorMaterialIn, vacuumMaterialIn,
42 xApertureIn, yApertureIn, xApertureOutIn, yApertureOutIn, colourIn, circularOuterIn)
43{;}
44
46{
50 {
52 {BDS::Warning(__METHOD_NAME__, "element: \"" + name + "\": X/Y half axes ratio at entrance and exit apertures are not equal");}
53 }
54}
55
57{
58 if (tapered)
59 {
60 // we use the notes from G4EllipticalCone.hh
61 // 1. half length in Z = zTopCut
62 // 2. Dx and Dy = half length of ellipse axis at z = -zTopCut
63 // 3. dx and dy = half length of ellipse axis at z = zTopCut
64 // ! Attention : dx/dy=Dx/Dy
65 // xSemiAxis = (Dx-dx)/(2*zTopCut)
66 // ySemiAxis = (Dy-dy)/(2*zTopCut)
67 // z height = (Dx+dx)/(2*xSemiAxis)
68 G4double xSemiAxis = (xAperture - xApertureOut) / chordLength;
69 G4double ySemiAxis = (yAperture - yApertureOut) / chordLength;
70 G4double zHeight = (xAperture + xApertureOut) / (2*xSemiAxis);
71
72 innerSolid = new G4EllipticalCone(name + "_inner_solid", // name
73 xSemiAxis, // major axis of largest ellipse
74 ySemiAxis, // minor axis of largest ellipse
75 zHeight, // height of cone
76 1.02*chordLength*0.5); // cut
77
78 vacuumSolid = new G4EllipticalCone(name + "_vacuum_solid",// name
79 xSemiAxis, // major axis of largest ellipse
80 ySemiAxis, // minor axis of largest ellipse
81 zHeight, // height of cone
82 chordLength*0.5 - lengthSafety); // cut
83 }
84 else
85 {
86 innerSolid = new G4EllipticalTube(name + "_inner_solid", // name
87 xAperture, // x half width
88 yAperture, // y half width
89 chordLength); // z half length
90 // z half length long for unambiguous subtraction
91
92 vacuumSolid = new G4EllipticalTube(name + "_inner_solid", // name
93 xAperture - lengthSafety, // x half width
94 yAperture - lengthSafety, // y half width
95 chordLength * 0.5); // z half length
96 }
97
100}
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4double lengthSafety
Useful variable often used in construction.
G4double chordLength
Protected member variable that can be modified by derived classes.
BDSCollimatorElliptical()
Private default constructor to force the use of the supplied one.
Base class for collimators with common construction.
G4double xAperture
Aperture at entrance in x dimension.
G4bool tapered
Flag for tapered collimator.
G4double yApertureOut
Aperture at exit in y dimension.
G4VSolid * innerSolid
Geometrical objects:
virtual void CheckParameters()
G4double xApertureOut
Aperture at exit in x dimension.
G4VSolid * vacuumSolid
Geometrical objects:
G4double yAperture
Aperture at entrance in y dimension.
void RegisterSolid(G4VSolid *solid)
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())