BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSCollimatorCrystal.hh
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
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#ifndef BDSCOLLIMATORCRYSTAL_H
20#define BDSCOLLIMATORCRYSTAL_H
21
22#include "globals.hh"
23
24#include "BDSAcceleratorComponent.hh"
25
26class BDSBeamPipeInfo;
27class BDSCrystal;
28class BDSCrystalInfo;
29
37{
38public:
39 BDSCollimatorCrystal(const G4String& name,
40 G4double length,
41 BDSBeamPipeInfo* beamPipeInfoIn,
42 BDSCrystalInfo* crystalInfoLeftIn,
43 BDSCrystalInfo* crystalInfoRightIn = nullptr,
44 G4double halfGapLeftIn = 0,
45 G4double halfGapRightIn = 0,
46 G4double angleYAxisLeftIn = 0,
47 G4double angleYAxisRightIn = 0);
48 virtual ~BDSCollimatorCrystal();
49
51 virtual G4String Material() const;
52
53protected:
55 virtual void Build();
56
57private:
60
63
66 void LongitudinalOverlap(const BDSExtent& extCrystal,
67 G4double crystalAngle,
68 const G4String& side) const;
69
71 void RegisterCrystalLVs(const BDSCrystal* crystal) const;
72
75 G4double TransverseOffsetToEdge(const BDSCrystal* crystal,
76 G4double placementAngle,
77 G4bool left) const;
78
81 const G4double halfGapLeft;
82 const G4double halfGapRight;
83 const G4double angleYAxisLeft;
84 const G4double angleYAxisRight;
85 BDSCrystal* crystalLeft;
86 BDSCrystal* crystalRight;
87};
88
89#endif
Abstract class that represents a component of an accelerator.
const G4String name
Const protected member variable that may not be changed by derived classes.
Holder class for all information required to describe a beam pipe model.
A piece of vacuum beam pipe with a crystal for channelling.
void RegisterCrystalLVs(const BDSCrystal *crystal) const
Register logical volumes in crystals and collimator sets in accelerator model.
virtual void BuildContainerLogicalVolume()
Void function to fulfill BDSAcceleratorComponent requirements.
BDSCrystalInfo * crystalInfoRight
Model associated with right crystal.
BDSCollimatorCrystal()=delete
No default constructor.
virtual G4String Material() const
Override base class version and return crystal material.
BDSCrystalInfo * crystalInfoLeft
Model associated with left crystal.
virtual void Build()
Construct geometry.
G4double TransverseOffsetToEdge(const BDSCrystal *crystal, G4double placementAngle, G4bool left) const
void LongitudinalOverlap(const BDSExtent &extCrystal, G4double crystalAngle, const G4String &side) const
Holder for all information required to create a crystal.
An object for both a crystal from a factory.
Definition: BDSCrystal.hh:40
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39