BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Public Member Functions | Protected Member Functions | Protected Attributes | Private Member Functions | Private Attributes
BDSBeamPipeFactoryPoints Class Referenceabstract

Factory for beam pipes defined by a series of x,y points that are extruded. More...

#include <BDSBeamPipeFactoryPoints.hh>

Inheritance diagram for BDSBeamPipeFactoryPoints:
Inheritance graph
Collaboration diagram for BDSBeamPipeFactoryPoints:
Collaboration graph

Public Member Functions

virtual BDSBeamPipeCreateBeamPipe (const G4String &nameIn, G4double lengthIn, G4double aper1In, G4double aper2In, G4double aper3In, G4double aper4In, G4Material *vacuumMaterialIn, G4double beamPipeThicknessIn, G4Material *beamPipeMaterialIn, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")
 Required overloaded method from BDSBeamPipeFactoryBase to build straight piece of beam pipe. More...
 
virtual BDSBeamPipeCreateBeamPipe (const G4String &nameIn, G4double lengthIn, const G4ThreeVector &inputFaceNormalIn, const G4ThreeVector &outputFaceNormalIn, G4double aper1In, G4double aper2In, G4double aper3In, G4double aper4In, G4Material *vacuumMaterialIn, G4double beamPipeThicknessIn, G4Material *beamPipeMaterialIn, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")
 Required overloaded method from BDSBeamPipeFactoryBase to build angled piece of beam pipe. More...
 
- Public Member Functions inherited from BDSBeamPipeFactoryBase
virtual BDSBeamPipeCreateBeamPipe (const G4String &nameIn, G4double lengthIn, G4double aper1=0, G4double aper2=0, G4double aper3=0, G4double aper4=0, G4Material *vacuumMaterialIn=nullptr, G4double beamPipeThicknessIn=0, G4Material *beamPipeMaterialIn=nullptr, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")=0
 create a flat ended beampipe More...
 
virtual BDSBeamPipeCreateBeamPipe (const G4String &nameIn, G4double lengthIn, const G4ThreeVector &inputFaceNormalIn, const G4ThreeVector &outputFaceNormalIn, G4double aper1=0, G4double aper2=0, G4double aper3=0, G4double aper4=0, G4Material *vacuumMaterialIn=nullptr, G4double beamPipeThicknessIn=0, G4Material *beamPipeMaterialIn=nullptr, const G4String &pointsFileIn="", const G4String &pointsUnitIn="")=0
 
virtual ~BDSBeamPipeFactoryBase ()
 Virtual base destructor. More...
 
- Public Member Functions inherited from BDSFactoryBase
virtual void FactoryBaseCleanUp ()
 Empty containers for next use - factories are never deleted so can't rely on scope. More...
 

Protected Member Functions

virtual void GeneratePoints (G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness, G4int pointsPerTwoPi=40)=0
 
virtual G4double CalculateIntersectionRadius (G4double aper1, G4double aper2, G4double aper3, G4double aper4, G4double beamPipeThickness)=0
 
virtual void CleanUp ()
 Clear member vectors and run base class clean up to clear pointers between runs. More...
 
void AppendPoint (std::vector< G4TwoVector > &vec, G4double x, G4double y)
 
void AppendAngle (std::vector< G4TwoVector > &vec, G4double startAngle, G4double finishAngle, G4double radius, G4int nPoints=10, G4double xOffset=0, G4double yOffset=0)
 
void AppendAngleEllipse (std::vector< G4TwoVector > &vec, G4double startAngle, G4double finishAngle, G4double radiusA, G4double radiusB, G4int nPoints=10, G4double xOffset=0, G4double yOffset=0)
 Generate 2-vector points (and append them) about an ellipse. More...
 
- Protected Member Functions inherited from BDSBeamPipeFactoryBase
 BDSBeamPipeFactoryBase ()
 base constructor More...
 
void CleanUpBase ()
 
virtual void CleanUp ()
 
void CommonConstruction (const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double length)
 finalise beampipe construction More...
 
BDSBeamPipeBuildBeamPipeAndRegisterVolumes (BDSExtent extent, G4double containerRadius, G4bool containerIsCircular=false)
 build beampipe and register logical volumes More...
 
virtual void BuildLogicalVolumes (const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn)
 build logical volumes More...
 
virtual void SetVisAttributes (G4Material *beamPipeMaterialIn)
 Set visual attributes. More...
 
virtual void SetUserLimits (G4double length)
 Set user limits. More...
 
virtual void PlaceComponents (const G4String &nameIn)
 Place volumes. More...
 

Protected Attributes

std::vector< G4TwoVector > vacuumEdge
 Vector of x,y coordinates for vacuum extruded solid edge. More...
 
std::vector< G4TwoVector > beamPipeInnerEdge
 Vector of x,y coordinates for beam pipe inner edge. More...
 
std::vector< G4TwoVector > beamPipeOuterEdge
 Vector of x,y coordinates for beam pipe outer edge. More...
 
std::vector< G4TwoVector > containerEdge
 Vector of x,y coordinates for the container volume. More...
 
std::vector< G4TwoVector > containerSubtractionEdge
 Vector of x,y coordinates for the container subtraction volume. More...
 
G4double extentX
 Extents prepared by GeneratePoints function as only it knows the extents. More...
 
G4double extentY
 Extents prepared by GeneratePoints function as only it knows the extents. More...
 
G4String pointsFile
 Temporary copy of input variables. More...
 
G4String pointsUnit
 Temporary copy of input variables. More...
 
- Protected Attributes inherited from BDSBeamPipeFactoryBase
G4bool sensitiveBeamPipe
 Whether the beam pipe will record energy deposition. More...
 
G4bool sensitiveVacuum
 Whether the vacuum will record any energy deposition. More...
 
G4bool storeApertureImpacts
 Whether to store aperture impacts. More...
 
G4VSolid * vacuumSolid
 
G4VSolid * beamPipeSolid
 
G4VSolid * containerSolid
 
G4VSolid * containerSubtractionSolid
 Longer (in length) version of container solid for unambiguous subtraction. More...
 
G4LogicalVolume * vacuumLV
 
G4LogicalVolume * beamPipeLV
 
G4LogicalVolume * containerLV
 
G4PVPlacement * vacuumPV
 
G4PVPlacement * beamPipePV
 
G4ThreeVector inputFaceNormal
 For recording the face normals in the finished pipe component. More...
 
G4ThreeVector outputFaceNormal
 For recording the face normals in the finished pipe component. More...
 
- Protected Attributes inherited from BDSFactoryBase
std::set< G4LogicalVolume * > allLogicalVolumes
 
std::set< G4VPhysicalVolume * > allPhysicalVolumes
 
std::set< G4RotationMatrix * > allRotationMatrices
 
std::set< G4UserLimits * > allUserLimits
 
std::set< G4VSolid * > allSolids
 
std::set< G4VisAttributes * > allVisAttributes
 
G4double lengthSafety
 Cache of global constants variable. More...
 
G4double lengthSafetyLarge
 Cache of global constants variable. More...
 
G4bool checkOverlaps
 Cache of global constants variable. More...
 
G4bool visDebug
 Cache of global constants variable. More...
 
G4double nSegmentsPerCircle
 Cache of global constants variable. More...
 
G4VisAttributes * containerVisAttr
 Cache of global constants variable. More...
 
G4UserLimits * defaultUserLimits
 Cache of global constants variable. More...
 

Private Member Functions

void CleanUpPoints ()
 Actual clean up here for this class only. More...
 
void CreateSolids (const G4String &name, G4double length, G4bool buildLongForIntersection=false)
 
void CreateSolidsAngled (const G4String &name, G4double length, const G4ThreeVector &inputFace, const G4ThreeVector &outputFace)
 Create angled solids. Uses CreateSolids() method. More...
 
virtual BDSBeamPipeCommonFinalConstruction (const G4String &nameIn, G4Material *vacuumMaterialIn, G4Material *beamPipeMaterialIn, G4double lengthIn)
 

Private Attributes

G4double intersectionRadius
 
G4VSolid * beamPipeInnerSolid
 Solid for inner surface of beam pipe. More...
 
G4VSolid * beamPipeOuterSolid
 Solid for outer surface of beam pipe. More...
 

Additional Inherited Members

- Static Protected Member Functions inherited from BDSBeamPipeFactoryBase
static void CheckAngledVolumeCanBeBuilt (G4double length, const G4ThreeVector &inputfaceAngle, const G4ThreeVector &outputfaceAngle, G4double horizontalWidth, const G4String &name)
 check if a beam pipe volume with angled faces can be constructed More...
 

Detailed Description

Factory for beam pipes defined by a series of x,y points that are extruded.

This is a pure virtual class that provides common functionality for aperture models derived from a series of points transversely. The derived class need therefore only generate the vector of points that will make the edges of the vacuum and inner & outer beam pipe surfaces. The solids are of type G4ExtrudedSolid, or in the case of angled faces, a G4ExtrudedSolid intersected with a G4CutTubs (cylinder with angled faces).

Author
Laurie Nevay

Definition at line 47 of file BDSBeamPipeFactoryPoints.hh.

Constructor & Destructor Documentation

◆ BDSBeamPipeFactoryPoints()

BDSBeamPipeFactoryPoints::BDSBeamPipeFactoryPoints ( )

Definition at line 36 of file BDSBeamPipeFactoryPoints.cc.

◆ ~BDSBeamPipeFactoryPoints()

BDSBeamPipeFactoryPoints::~BDSBeamPipeFactoryPoints ( )
virtual

Definition at line 41 of file BDSBeamPipeFactoryPoints.cc.

Member Function Documentation

◆ AppendAngle()

void BDSBeamPipeFactoryPoints::AppendAngle ( std::vector< G4TwoVector > &  vec,
G4double  startAngle,
G4double  finishAngle,
G4double  radius,
G4int  nPoints = 10,
G4double  xOffset = 0,
G4double  yOffset = 0 
)
protected

Generate 2-vector points (and append them) about a circle. Uses ellipse code with equal radii.

Definition at line 76 of file BDSBeamPipeFactoryPoints.cc.

References AppendAngleEllipse().

Here is the call graph for this function:

◆ AppendAngleEllipse()

void BDSBeamPipeFactoryPoints::AppendAngleEllipse ( std::vector< G4TwoVector > &  vec,
G4double  startAngle,
G4double  finishAngle,
G4double  radiusA,
G4double  radiusB,
G4int  nPoints = 10,
G4double  xOffset = 0,
G4double  yOffset = 0 
)
protected

Generate 2-vector points (and append them) about an ellipse.

Definition at line 87 of file BDSBeamPipeFactoryPoints.cc.

Referenced by AppendAngle().

Here is the caller graph for this function:

◆ AppendPoint()

void BDSBeamPipeFactoryPoints::AppendPoint ( std::vector< G4TwoVector > &  vec,
G4double  x,
G4double  y 
)
protected

Definition at line 69 of file BDSBeamPipeFactoryPoints.cc.

◆ CalculateIntersectionRadius()

virtual G4double BDSBeamPipeFactoryPoints::CalculateIntersectionRadius ( G4double  aper1,
G4double  aper2,
G4double  aper3,
G4double  aper4,
G4double  beamPipeThickness 
)
protectedpure virtual

Function derived class must provide to calculated the intersection radius, which is the radius of the intersection cylinder to create the aperture faces. This should be big enough to contain the cross-section of the beam pipe but not overly large, so as to avoid intersecting plane errors (and therefore exits) from large angles on short elements.

Implemented in BDSBeamPipeFactoryClicPCL, BDSBeamPipeFactoryOctagonal, BDSBeamPipeFactoryPointsFile, and BDSBeamPipeFactoryRaceTrack.

Referenced by CreateBeamPipe().

Here is the caller graph for this function:

◆ CleanUp()

void BDSBeamPipeFactoryPoints::CleanUp ( )
protectedvirtual

Clear member vectors and run base class clean up to clear pointers between runs.

Reimplemented from BDSBeamPipeFactoryBase.

Reimplemented in BDSBeamPipeFactoryClicPCL.

Definition at line 44 of file BDSBeamPipeFactoryPoints.cc.

References CleanUpPoints().

Referenced by BDSBeamPipeFactoryClicPCL::CleanUp(), and CreateBeamPipe().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CleanUpPoints()

void BDSBeamPipeFactoryPoints::CleanUpPoints ( )
private

Actual clean up here for this class only.

Definition at line 50 of file BDSBeamPipeFactoryPoints.cc.

References beamPipeInnerEdge, beamPipeInnerSolid, beamPipeOuterEdge, beamPipeOuterSolid, containerEdge, containerSubtractionEdge, extentX, extentY, intersectionRadius, pointsFile, pointsUnit, and vacuumEdge.

Referenced by CleanUp().

Here is the caller graph for this function:

◆ CommonFinalConstruction()

BDSBeamPipe * BDSBeamPipeFactoryPoints::CommonFinalConstruction ( const G4String &  nameIn,
G4Material *  vacuumMaterialIn,
G4Material *  beamPipeMaterialIn,
G4double  lengthIn 
)
privatevirtual

Overloads BDSBeamPipeFactoryBase method to do common construction tasks and constructs BDSBeamPipe instance to return.

Reimplemented in BDSBeamPipeFactoryClicPCL.

Definition at line 284 of file BDSBeamPipeFactoryPoints.cc.

References BDSBeamPipeFactoryBase::BuildBeamPipeAndRegisterVolumes(), BDSBeamPipeFactoryBase::CommonConstruction(), extentX, extentY, and intersectionRadius.

Referenced by CreateBeamPipe().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ CreateBeamPipe() [1/2]

BDSBeamPipe * BDSBeamPipeFactoryPoints::CreateBeamPipe ( const G4String &  nameIn,
G4double  lengthIn,
const G4ThreeVector &  inputFaceNormalIn,
const G4ThreeVector &  outputFaceNormalIn,
G4double  aper1In,
G4double  aper2In,
G4double  aper3In,
G4double  aper4In,
G4Material *  vacuumMaterialIn,
G4double  beamPipeThicknessIn,
G4Material *  beamPipeMaterialIn,
const G4String &  pointsFileIn = "",
const G4String &  pointsUnitIn = "" 
)
virtual

◆ CreateBeamPipe() [2/2]

BDSBeamPipe * BDSBeamPipeFactoryPoints::CreateBeamPipe ( const G4String &  nameIn,
G4double  lengthIn,
G4double  aper1In,
G4double  aper2In,
G4double  aper3In,
G4double  aper4In,
G4Material *  vacuumMaterialIn,
G4double  beamPipeThicknessIn,
G4Material *  beamPipeMaterialIn,
const G4String &  pointsFileIn = "",
const G4String &  pointsUnitIn = "" 
)
virtual

Required overloaded method from BDSBeamPipeFactoryBase to build straight piece of beam pipe.

Implements BDSBeamPipeFactoryBase.

Definition at line 219 of file BDSBeamPipeFactoryPoints.cc.

References CalculateIntersectionRadius(), CleanUp(), CommonFinalConstruction(), CreateSolids(), GeneratePoints(), intersectionRadius, pointsFile, and pointsUnit.

Here is the call graph for this function:

◆ CreateSolids()

void BDSBeamPipeFactoryPoints::CreateSolids ( const G4String &  name,
G4double  length,
G4bool  buildLongForIntersection = false 
)
private

Create the vacuum, inner & outer, container and container subtraction solids. Optional flag buildLongForIntersection builds the solids longer by 1.5x length in preparation for them to be intersected for angled faces, saving code duplication.

Make the solids all longer for intersection

Definition at line 108 of file BDSBeamPipeFactoryPoints.cc.

References beamPipeInnerEdge, beamPipeInnerSolid, beamPipeOuterEdge, beamPipeOuterSolid, containerEdge, containerSubtractionEdge, BDSBeamPipeFactoryBase::containerSubtractionSolid, BDSFactoryBase::lengthSafety, and vacuumEdge.

Referenced by CreateBeamPipe(), and CreateSolidsAngled().

Here is the caller graph for this function:

◆ CreateSolidsAngled()

void BDSBeamPipeFactoryPoints::CreateSolidsAngled ( const G4String &  name,
G4double  length,
const G4ThreeVector &  inputFace,
const G4ThreeVector &  outputFace 
)
private

Create angled solids. Uses CreateSolids() method.

Definition at line 157 of file BDSBeamPipeFactoryPoints.cc.

References BDS::CalculateSafeAngledVolumeLength(), BDSBeamPipeFactoryBase::CheckAngledVolumeCanBeBuilt(), CreateSolids(), intersectionRadius, and BDSFactoryBase::lengthSafety.

Referenced by CreateBeamPipe().

Here is the call graph for this function:
Here is the caller graph for this function:

◆ GeneratePoints()

virtual void BDSBeamPipeFactoryPoints::GeneratePoints ( G4double  aper1,
G4double  aper2,
G4double  aper3,
G4double  aper4,
G4double  beamPipeThickness,
G4int  pointsPerTwoPi = 40 
)
protectedpure virtual

Function to generate transverse points. Should set member variables vacuumEdge, beamPipeInnerEdge, beamPipeOuterEdge, containerEdge and containerSubtractionEdge. This must also set extentX and extentY member variables for +ve symmetric extents.

Implemented in BDSBeamPipeFactoryClicPCL, BDSBeamPipeFactoryOctagonal, BDSBeamPipeFactoryRaceTrack, and BDSBeamPipeFactoryPointsFile.

Referenced by CreateBeamPipe().

Here is the caller graph for this function:

Field Documentation

◆ beamPipeInnerEdge

std::vector<G4TwoVector> BDSBeamPipeFactoryPoints::beamPipeInnerEdge
protected

◆ beamPipeInnerSolid

G4VSolid* BDSBeamPipeFactoryPoints::beamPipeInnerSolid
private

Solid for inner surface of beam pipe.

Definition at line 164 of file BDSBeamPipeFactoryPoints.hh.

Referenced by CleanUpPoints(), and CreateSolids().

◆ beamPipeOuterEdge

std::vector<G4TwoVector> BDSBeamPipeFactoryPoints::beamPipeOuterEdge
protected

◆ beamPipeOuterSolid

G4VSolid* BDSBeamPipeFactoryPoints::beamPipeOuterSolid
private

Solid for outer surface of beam pipe.

Definition at line 167 of file BDSBeamPipeFactoryPoints.hh.

Referenced by CleanUpPoints(), and CreateSolids().

◆ containerEdge

std::vector<G4TwoVector> BDSBeamPipeFactoryPoints::containerEdge
protected

◆ containerSubtractionEdge

std::vector<G4TwoVector> BDSBeamPipeFactoryPoints::containerSubtractionEdge
protected

◆ extentX

G4double BDSBeamPipeFactoryPoints::extentX
protected

◆ extentY

G4double BDSBeamPipeFactoryPoints::extentY
protected

Extents prepared by GeneratePoints function as only it knows the extents.

Definition at line 147 of file BDSBeamPipeFactoryPoints.hh.

Referenced by CleanUpPoints(), CommonFinalConstruction(), BDSBeamPipeFactoryOctagonal::GeneratePoints(), and BDSBeamPipeFactoryRaceTrack::GeneratePoints().

◆ intersectionRadius

G4double BDSBeamPipeFactoryPoints::intersectionRadius
private

Radius of intersection cylinder for angled faced beam pipe pieces. Should be determined and set in GeneratePoints() by derived class. Should be big enough to enclose the beam pipe section transversely but not too large that overlapping z planes errors might occur with short pieces of beam pipe with large angled faces.

Definition at line 161 of file BDSBeamPipeFactoryPoints.hh.

Referenced by CleanUpPoints(), CommonFinalConstruction(), CreateBeamPipe(), and CreateSolidsAngled().

◆ pointsFile

G4String BDSBeamPipeFactoryPoints::pointsFile
protected

Temporary copy of input variables.

Definition at line 151 of file BDSBeamPipeFactoryPoints.hh.

Referenced by CleanUpPoints(), CreateBeamPipe(), and BDSBeamPipeFactoryPointsFile::GeneratePoints().

◆ pointsUnit

G4String BDSBeamPipeFactoryPoints::pointsUnit
protected

Temporary copy of input variables.

Definition at line 152 of file BDSBeamPipeFactoryPoints.hh.

Referenced by CleanUpPoints(), CreateBeamPipe(), and BDSBeamPipeFactoryPointsFile::GeneratePoints().

◆ vacuumEdge

std::vector<G4TwoVector> BDSBeamPipeFactoryPoints::vacuumEdge
protected

The documentation for this class was generated from the following files: