BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSScreenLayer.cc
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#include "BDSDebug.hh"
20#include "BDSException.hh"
21#include "BDSGlobalConstants.hh"
22#include "BDSMaterials.hh"
23#include "BDSSamplerRegistry.hh"
24#include "BDSScreenLayer.hh"
25#include "BDSSDManager.hh"
26#include "BDSSDSampler.hh"
27#include "BDSUtilities.hh"
28
29#include "globals.hh"
30#include "G4Box.hh"
31#include "G4LogicalBorderSurface.hh"
32#include "G4LogicalVolume.hh"
33#include "G4MaterialPropertiesTable.hh"
34#include "G4OpticalSurface.hh"
35#include "G4PVPlacement.hh"
36#include "G4ThreeVector.hh"
37#include "G4VisAttributes.hh"
38
40{}
41
42BDSScreenLayer::BDSScreenLayer(G4ThreeVector sizeIn,
43 G4String nameIn,
44 G4String materialIn,
45 G4double grooveWidthIn,
46 G4double grooveSpatialFrequencyIn):
47 size(sizeIn),
48 name(nameIn),
49 material(materialIn),
50 grooveWidth(grooveWidthIn),
51 grooveSpatialFrequency(grooveSpatialFrequencyIn)
52{
53 if (!BDS::IsFinite(size.z()))
54 {throw BDSException(__METHOD_NAME__, "Insufficent length for screen layer \"" + name + "\"");}
55 colour=G4Colour(0.1,0.8,0.1,0.3);
56 Build();
57}
58
59void BDSScreenLayer::Build()
60{
61 BuildGroove();
62 BuildScreen();
63 SetVisAttributes();
64}
65
66void BDSScreenLayer::BuildGroove()
67{
68 //There may or may not be grooves in the screen layer.
69 if (!BDS::IsFinite(grooveWidth))
70 {return;}
71
72 grooveSolid = new G4Box((name+"_grooveSolid").c_str(),
73 grooveWidth/2.0,
74 size.y()/2.0,
75 size.z()/2.0);
76 grooveLog = new G4LogicalVolume(grooveSolid,
77 BDSMaterials::Instance()->GetMaterial("air"),
78 (name+"_grooveLog").c_str());
79}
80
81void BDSScreenLayer::BuildScreen()
82{
83 //A small safety length to separate the screen layers.
84 G4double tinyLenSaf = 1e-3*CLHEP::nm;
85 solid = new G4Box((name+"_solid").c_str(),
86 size.x()/2.0,
87 size.y()/2.0,
88 size.z()/2.0-tinyLenSaf);
89 log = new G4LogicalVolume(solid,
90 BDSMaterials::Instance()->GetMaterial(material),
91 (name+"_log").c_str());
92 CutGrooves();
93
94 log->SetUserLimits(BDSGlobalConstants::Instance()->DefaultUserLimits());
95}
96
97void BDSScreenLayer::CutGrooves()
98{
99 if(grooveWidth>0)
100 {
101 for(G4double xPosition=-size.x()/2.0+grooveWidth/2.0;
102 xPosition<((size.x()/2.0)-grooveWidth/2.0);
103 xPosition+=grooveSpatialFrequency)
104 { CutGroove(xPosition);}
105 }
106}
107
108void BDSScreenLayer::CutGroove(G4double xPosition)
109{
110 if (!grooveLog)
111 {return;}
112 G4ThreeVector pos(xPosition, 0, 0);
113
114 //Create a new physical volume placement for each groove in the screen.
115 new G4PVPlacement(nullptr,
116 pos,
117 grooveLog,
118 (G4String)(name+"_groove"),
119 log,
120 true,
121 nGrooves,
122 false);
123
124 nGrooves++; //Increment the counter for the number of grooves in the screen.
125}
126
127void BDSScreenLayer::SetVisAttributes()
128{
129 G4VisAttributes* visAtt=new G4VisAttributes(colour);
130 visAtt->SetForceSolid(true);
131 G4VisAttributes* visAttGroove=new G4VisAttributes(G4Colour(0.0,0.0,0.0));
132 visAttGroove->SetForceSolid(true);
133 log->SetVisAttributes(visAtt);
134 if(grooveLog)
135 {grooveLog->SetVisAttributes(visAttGroove);}
136}
137
138void BDSScreenLayer::SetPhys(G4PVPlacement* physIn)
139{
140 if (physIn->GetLogicalVolume() != GetLog())
141 {throw BDSException(__METHOD_NAME__, "physical volume placement does not match logical volume.");}
142 phys = physIn;
143}
144
145void BDSScreenLayer::SetColour(G4Colour col)
146{
147 colour=col;
148 SetVisAttributes();
149}
150
151
153{
154 G4String samplerName = name;
155 log->SetSensitiveDetector(BDSSDManager::Instance()->SamplerPlane());
156 samplerID=BDSSamplerRegistry::Instance()->RegisterSampler(samplerName,nullptr);
157 log->SetUserLimits(BDSGlobalConstants::Instance()->DefaultUserLimits());
158}
159
160/*
161void BDSScreenLayer::backInternalMirror()
162{
163 internalMirror = new InternalMirror(InternalMirror::BACK, size,material,log,phys);
164}
165
166void BDSScreenLayer::frontInternalMirror()
167{
168 internalMirror = new InternalMirror(InternalMirror::FRONT,size,material,log,phys);
169}
170
171BDSScreenLayer::InternalMirror::~InternalMirror()
172{;}
173
174BDSScreenLayer::InternalMirror::InternalMirror(G4int varsideIn,
175 G4ThreeVector sizeIn,
176 G4String materialIn,
177 G4LogicalVolume* motherLogIn,
178 G4PVPlacement* motherPhysIn):
179 side(varsideIn),
180 motherSize(sizeIn),
181 motherMaterial(materialIn),
182 motherLog(motherLogIn),
183 motherPhys(motherPhysIn)
184{
185 thickness=1e-9*CLHEP::m;
186 compute();
187 geom();
188 place();
189 optical();
190}
191
192void BDSScreenLayer::InternalMirror::geom()
193{
194 G4double tinyLenSaf = 1e-6*CLHEP::nm;
195 solid = new G4Box("internalMirrorSolid",
196 motherSize.x()/2.0,
197 motherSize.y()/2.0,
198 thickness - tinyLenSaf);
199 log = new G4LogicalVolume(solid,
200 BDSMaterials::Instance()->GetMaterial("titanium"),
201 "internalMirrorLog");
202}
203
204void BDSScreenLayer::InternalMirror::place()
205{
206 phys=new G4PVPlacement(nullptr, //Create a new physical volume placement for each groove in the screen.
207 G4ThreeVector(0,0,pos),
208 log,
209 "internalMirrorPhys",
210 motherLog,
211 false,
212 0,
213 true);
214}
215
216void BDSScreenLayer::InternalMirror::optical()
217{
218 G4OpticalSurface* OpSurface=new G4OpticalSurface("OpSurface");
219 //G4LogicalBorderSurface* LogSurface =
220 new G4LogicalBorderSurface("LogSurface", motherPhys, phys, OpSurface);
221 // G4LogicalSkinSurface* LogSurface = new G4LogicalSkinSurface("LogSurface",screenLayer(1)->log(),OpSurface);
222 OpSurface->SetType(dielectric_metal);
223 OpSurface->SetModel(unified);
224 OpSurface->SetFinish(polished);
225 G4MaterialPropertiesTable* SMPT = new G4MaterialPropertiesTable();
226 SMPT->AddConstProperty("REFLECTIVITY",0.8);
227 // SMPT->AddConstProperty("RINDEX",1.5750);
228 OpSurface->SetMaterialPropertiesTable(SMPT);
229}
230
231void BDSScreenLayer::InternalMirror::compute()
232{
233 G4double sign=0;
234 try{
235 if(side==BACK){
236 sign = 1;
237 }else if(side==FRONT){
238 sign = -1;
239 }else{
240 throw 1;
241 }
242 }catch(int e){
243 G4cerr<< "BDSScreenLayer::computInternalMirror - exception number " << e << " occurred. Exiting." << G4endl;
244 exit(e);
245 }
246
247 pos = sign*(motherSize.z()/2.0-thickness/2.0);
248}
249*/
General exception with possible name of object and message.
Definition: BDSException.hh:35
static BDSGlobalConstants * Instance()
Access method.
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
static BDSSamplerRegistry * Instance()
Accessor for registry.
G4int RegisterSampler(const G4String &name, BDSSampler *sampler, const G4Transform3D &transform=G4Transform3D(), G4double S=-1000, const BDSBeamlineElement *element=nullptr, BDSSamplerType type=BDSSamplerType::plane, G4double radius=0)
G4LogicalVolume * GetLog() const
Accessor.
void AssignSampler()
Make this plane a sampling plane.
G4int nGrooves
Counter for the number of grooves etched into the screen.
BDSScreenLayer()
Used in modules/AWAKE/include/BDSMultiFacetLayer.hh.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())