BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
Loading...
Searching...
No Matches
BDSScreen.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 "BDSMultilayerScreen.hh"
20#include "BDSScreen.hh"
21#include "BDSScreenLayer.hh"
22#include "BDSSDType.hh"
23
24#include "globals.hh"
25#include "G4Colour.hh"
26#include "G4RotationMatrix.hh"
27#include "G4VisAttributes.hh"
28
29#include <list>
30#include <sstream>
31
32BDSScreen::BDSScreen(G4String nameIn,
33 G4double chordLengthIn,
34 BDSBeamPipeInfo* beamPipeInfoIn,
35 G4TwoVector sizeIn, //X Y dimensions of screen
36 G4double screenAngleIn):
37 BDSDrift(nameIn, chordLengthIn, beamPipeInfoIn),
38 size(sizeIn),
39 screenAngle(screenAngleIn),
40 screenPos(G4ThreeVector()),
41 nLayers(0)
42{
43 mlScreen = new BDSMultilayerScreen(size, nameIn+"_mlscreen");
44
45 screenRot = new G4RotationMatrix();
46 screenRot->rotateY(screenAngle);
47}
48
49BDSScreen::~BDSScreen()
50{
51 delete screenRot;
52 delete mlScreen;
53}
54
56{
57 BDSDrift::Build(); //Build the beam pipe geometry.
58
59 // Make the beam pipe wireframe
60 G4VisAttributes* VisAtt1 = new G4VisAttributes(G4Colour(0.4, 0.4, 0.4, 0.3));
61 VisAtt1->SetForceWireframe(true);
62 VisAtt1->SetVisibility(true);
63 containerLogicalVolume->SetVisAttributes(VisAtt1);
64
65 PlaceScreen(); //Place the screen in the beam pipe
66}
67
68void BDSScreen::AddScreenLayer(G4double thickness, G4String material, G4int isSampler)
69{
70 std::stringstream ss;
71 ss << nLayers;
72 G4String lNum = ss.str();
73 G4String lName = name+"_"+lNum;
74 mlScreen->AddScreenLayer(thickness,material,lName, isSampler);
75 if(!isSampler && sensitiveOuter)
76 {RegisterSensitiveVolume(mlScreen->LastLayer()->GetLog(), BDSSDType::energydep);}
77 nLayers++;
78}
79
80void BDSScreen::PlaceScreen()
81{
82 mlScreen->Build();//Build the screen.
83 G4LogicalVolume* vacuumLV = *(GetAcceleratorVacuumLogicalVolumes().begin());
84 mlScreen->Place(screenRot, screenPos, vacuumLV); //Place the screen in the beampipe centre.
85}
const G4String name
Const protected member variable that may not be changed by derived classes.
static G4bool sensitiveOuter
Useful variable often used in construction.
virtual std::set< G4LogicalVolume * > GetAcceleratorVacuumLogicalVolumes() const
Access the 'vacuum' volume(s) in this component if any. Default is empty set.
Holder class for all information required to describe a beam pipe model.
A piece of vacuum beam pipe.
Definition: BDSDrift.hh:39
virtual void Build()
Construct geometry.
Definition: BDSDrift.cc:37
void RegisterSensitiveVolume(G4LogicalVolume *sensitiveVolume, BDSSDType sensitivityType)
An accelerator component for diagnostics screens e.g. OTR. Screen inside beam pipe.
void AddScreenLayer(G4double thickness, G4String material, G4String name, G4int isSampler=0, G4double grooveWidth=0, G4double grooveSpatialFrequency=0)
Construct and add a screen layer.
virtual void Place(G4RotationMatrix *rot, const G4ThreeVector &pos, G4LogicalVolume *motherVol)
Place the member logical volume 'log' with a given transform in a given mother volume.
void Build()
Construct container, compute dimensions then place layers.
BDSScreenLayer * LastLayer() const
Get the last layer of the screen.
G4LogicalVolume * GetLog() const
Accessor.
virtual void Build()
Construct geometry.
Definition: BDSScreen.cc:55
void AddScreenLayer(G4double thickness, G4String material, G4int isSampler=0)
Add a screen layer.
Definition: BDSScreen.cc:68