19#include "BDSAcceleratorComponent.hh"
21#include "BDSDegrader.hh"
22#include "BDSException.hh"
23#include "BDSSDType.hh"
26#include "G4ExtrudedSolid.hh"
27#include "G4LogicalVolume.hh"
28#include "G4PVPlacement.hh"
29#include "G4RotationMatrix.hh"
30#include "G4TwoVector.hh"
31#include "G4VisAttributes.hh"
38#include "CLHEP/Units/SystemOfUnits.h"
42 G4double horizontalWidthIn,
44 G4double wedgeLengthIn,
45 G4double degraderHeightIn,
46 G4double degraderOffsetIn,
48 G4Material* degraderMaterialIn,
49 G4Material* vacuumMaterialIn,
52 horizontalWidth(horizontalWidthIn),
53 numberWedges(numberWedgesIn),
54 wedgeLength(wedgeLengthIn),
55 degraderHeight(degraderHeightIn),
56 degraderOffset(degraderOffsetIn),
57 baseWidth(baseWidthIn),
58 degraderMaterial(degraderMaterialIn),
59 vacuumMaterial(vacuumMaterialIn),
63BDSDegrader::~BDSDegrader()
69 if (horizontalWidth <= 0)
70 {
throw BDSException(__METHOD_NAME__,
"option \"horizontalWidth\" is not defined or must be greater than 0 for element \"" +
name +
"\"");}
73 {
throw BDSException(__METHOD_NAME__,
"option \"numberWedges\" is not defined or must be greater than 0 for element \"" +
name +
"\"");}
76 {
throw BDSException(__METHOD_NAME__,
"option \"wedgeLength\" is not defined or must be greater than 0 for element \"" +
name +
"\"");}
78 if (degraderHeight <= 0)
79 {
throw BDSException(__METHOD_NAME__,
"option \"degraderHeight\" is not defined or must be greater than 0 for element \"" +
name +
"\"");}
81 if (degraderHeight > (0.5*horizontalWidth))
82 {
throw BDSException(__METHOD_NAME__,
"option \"degraderHeight\" must be less than 0.5 times \"horizontalWidth\" for element \"" +
name +
"\"");}
84 if ((degraderOffset + baseWidth) > 0.5*horizontalWidth)
85 {
throw BDSException(__METHOD_NAME__,
"option \"degraderOffset\" must be less than (0.5 times \"horizontalWidth\", minus aper1) for element \"" +
name +
"\"");}
88 degraderOffset += baseWidth;
91 G4double minimumWidth = wedgeLength;
92 G4double containerWidth = degraderOffset;
93 if (degraderOffset < minimumWidth)
94 {containerWidth = minimumWidth;}
96 containerSolid =
new G4Box(
name +
"_container_solid",
101 containerLogicalVolume =
new G4LogicalVolume(containerSolid,
103 name +
"_container_lv");
115 G4double addedBasewidth = baseWidth* wedgeWidth/(2.0*wedgeLength) - 10*
lengthSafety;
117 std::vector<G4TwoVector> fullWedgeVertices;
118 std::vector<G4TwoVector> halfWedgeVertices;
121 fullWedgeVertices.push_back(G4TwoVector(0.5*wedgeWidth + addedBasewidth, 0));
122 fullWedgeVertices.push_back(G4TwoVector(0, wedgeLength+baseWidth));
123 fullWedgeVertices.push_back(G4TwoVector(-0.5*wedgeWidth - addedBasewidth, 0));
126 halfWedgeVertices.push_back(G4TwoVector(0, 0));
127 halfWedgeVertices.push_back(G4TwoVector(0.5*wedgeWidth + addedBasewidth, 0));
128 halfWedgeVertices.push_back(G4TwoVector(0, wedgeLength+baseWidth));
131 G4ExtrudedSolid* fullWedge =
new G4ExtrudedSolid(
name +
"_fullwedge_solid",
134 G4TwoVector(),1, G4TwoVector(), 1);
136 G4LogicalVolume* fullWedgeLV =
new G4LogicalVolume(fullWedge,
138 name +
"_fullwedge_lv");
142 G4ExtrudedSolid* halfWedge =
new G4ExtrudedSolid(
name +
"_halfwedge_solid",
145 G4TwoVector(),1, G4TwoVector(), 1);
147 G4LogicalVolume* halfWedgeLV =
new G4LogicalVolume(halfWedge,
149 name +
"_halfwedge_lv");
153 G4VisAttributes* degraderVisAttr =
new G4VisAttributes(colour);
154 fullWedgeLV->SetVisAttributes(degraderVisAttr);
155 halfWedgeLV->SetVisAttributes(degraderVisAttr);
159 G4RotationMatrix* rightRot =
new G4RotationMatrix();
160 rightRot->rotateY(-CLHEP::halfpi);
161 rightRot->rotateX(-CLHEP::halfpi);
164 G4RotationMatrix* leftRot =
new G4RotationMatrix();
165 leftRot->rotateY(CLHEP::halfpi);
166 leftRot->rotateX(-CLHEP::halfpi);
170 G4RotationMatrix* altRightRot =
new G4RotationMatrix();
171 altRightRot->rotateY(CLHEP::halfpi);
172 altRightRot->rotateX(CLHEP::halfpi);
176 for (G4int i = 0; i < (numberWedges + 1); i++)
178 G4String wedgeName =
name +
"_" + std::to_string(i)+
"_pv";
179 if (IsEven(numberWedges))
183 else if (i==numberWedges)
186 {
PlaceWedge(
true, i*wedgeWidth, wedgeName, fullWedgeLV, rightRot);}
188 {
PlaceWedge(
false, i*wedgeWidth, wedgeName, fullWedgeLV, leftRot);}
190 else if (IsOdd(numberWedges))
194 else if (i == numberWedges)
197 {
PlaceWedge(
false, i*wedgeWidth, wedgeName, fullWedgeLV, leftRot);}
199 {
PlaceWedge(
true, i*wedgeWidth, wedgeName, fullWedgeLV, rightRot);}
212 const G4String& placementName,
214 G4RotationMatrix* rotation)
216 G4double xOffset = -degraderOffset;
218 {xOffset = degraderOffset;}
220 G4VPhysicalVolume* wedgePV =
new G4PVPlacement(rotation,
221 G4ThreeVector(xOffset, 0, -0.5*
arcLength + zOffset),
224 containerLogicalVolume,
Abstract class that represents a component of an accelerator.
const G4double arcLength
Const protected member variable that may not be changed by derived classes.
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4bool sensitiveOuter
Useful variable often used in construction.
static G4double lengthSafety
Useful variable often used in construction.
static G4bool checkOverlaps
Useful variable often used in construction.
G4double chordLength
Protected member variable that can be modified by derived classes.
void PlaceWedge(G4bool placeRight, G4double zOffset, const G4String &placementName, G4LogicalVolume *lv, G4RotationMatrix *rotation)
BDSDegrader()=delete
Private default constructor to force the use of the supplied one.
virtual void BuildContainerLogicalVolume()
General exception with possible name of object and message.
void RegisterRotationMatrix(G4RotationMatrix *rotationMatrix)
void RegisterLogicalVolume(G4LogicalVolume *logicalVolume)
void RegisterPhysicalVolume(G4VPhysicalVolume *physicalVolume)
void RegisterVisAttributes(G4VisAttributes *visAttribute)
void RegisterSolid(G4VSolid *solid)
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)