BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSParallelWorldSampler.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 "BDSAcceleratorModel.hh"
20#include "BDSApertureInfo.hh"
21#include "BDSBeamline.hh"
22#include "BDSBeamlineElement.hh"
23#include "BDSBeamPipeInfo.hh"
24#include "BDSDebug.hh"
25#include "BDSDetectorConstruction.hh"
26#include "BDSExtent.hh"
27#include "BDSException.hh"
28#include "BDSGlobalConstants.hh"
29#include "BDSOutput.hh"
30#include "BDSParallelWorldSampler.hh"
31#include "BDSParser.hh"
32#include "BDSSampler.hh"
33#include "BDSSamplerCustom.hh"
34#include "BDSSamplerCylinder.hh"
35#include "BDSSamplerInfo.hh"
36#include "BDSSamplerPlane.hh"
37#include "BDSSamplerRegistry.hh"
38#include "BDSSamplerSphere.hh"
39#include "BDSSamplerType.hh"
40#include "BDSSDSampler.hh"
41#include "BDSUtilities.hh"
42#include "BDSWarning.hh"
43
44#include "parser/samplerplacement.h"
45
46#include "globals.hh" // geant4 types / globals
47#include "G4LogicalVolume.hh"
48#include "G4PVPlacement.hh"
49#include "G4Transform3D.hh"
50#include "G4VisAttributes.hh"
51#include "G4VPhysicalVolume.hh"
52#include "G4VUserParallelWorld.hh"
53
54#include <map>
55#include <set>
56#include <vector>
57
58
60 G4VUserParallelWorld("SamplerWorld_" + name),
61 suffix(name),
62 samplerWorldVis(nullptr),
63 generalPlane(nullptr),
64 samplerWorldLV(nullptr)
65{;}
66
67BDSParallelWorldSampler::~BDSParallelWorldSampler()
68{
69 // we leak all the placements as g4 is incredibly slow to delete
70 // them as it deregisters them - g4 will clean up
71 //for (auto placement : placements)
72 // {delete placement;}
73 delete samplerWorldVis;
74 // samplers are deleted by the sampler registry
75}
76
78{
80 G4VPhysicalVolume* samplerWorld = GetWorld();
81 samplerWorldLV = samplerWorld->GetLogicalVolume();
82
83 samplerWorldVis = new G4VisAttributes(*(globals->VisibleDebugVisAttr()));
84 samplerWorldVis->SetForceWireframe(true);//just wireframe so we can see inside it
85 samplerWorldLV->SetVisAttributes(samplerWorldVis);
86
87 const G4double samplerRadius = 0.5*globals->SamplerDiameter();
88 const BDSBeamline* beamline = BDSAcceleratorModel::Instance()->BeamlineMain();
89 const G4bool checkOverlaps = globals->CheckOverlaps();
90
91 // Construct the one sampler typically used for a general sampler.
92 generalPlane = new BDSSamplerPlane("Plane_sampler", samplerRadius);
93 samplerInstances[-1] = generalPlane;
94 // Construct any others we need.
95 std::map<int, std::set<int>> filterIDToSet = BDSParser::Instance()->GetSamplerFilterIDToSet();
96 for (const auto& kv : filterIDToSet)
97 {
98 BDSSamplerPlane* samp = new BDSSamplerPlane("Plane_sampler_w_filter_set_" + std::to_string(kv.first),
99 samplerRadius,
100 (G4int)kv.first);
101 samplerInstances[(G4int)kv.first] = samp;
102 }
103
104 // For each element in the beamline (if a beamline has been defined) construct and place the
105 // appropriate type of sampler if required. Info encoded in BDSBeamlineElement instance
106 if (beamline != nullptr)
107 {
108 for (auto element : *beamline)
109 {Place(element, samplerRadius);}
110 }
111
112 // Now user customised samplers
113 std::vector<GMAD::SamplerPlacement> samplerPlacements = BDSParser::Instance()->GetSamplerPlacements();
114
115 for (const auto& samplerPlacement : samplerPlacements)
116 {
117 G4String sn = samplerPlacement.name;
118 G4cout << "User placed sampler: \"" << sn << "\"" << G4endl;
119
120 if (BDSOutput::InvalidSamplerName(samplerPlacement.name))
121 {
122 G4cerr << __METHOD_NAME__ << "invalid sampler name \"" << sn << "\"" << G4endl;
124 throw BDSException(__METHOD_NAME__, "");
125 }
126
127 // use main beamline - in future, multiple beam lines
128 G4Transform3D transform = BDSDetectorConstruction::CreatePlacementTransform(samplerPlacement, beamline);
129
130 BDSSamplerType st = BDS::DetermineSamplerType(samplerPlacement.samplerType);
131 AdjustTransform(transform, st); // for 'forward' samplers we add an extra rotation
132 G4double radius = 0;
133 BDSSampler* sampler = BuildSampler(samplerPlacement, st, radius);
134 G4int samplerID = BDSSamplerRegistry::Instance()->RegisterSampler(sn,
135 sampler,
136 transform,
137 -1000,
138 nullptr,
139 st,
140 radius);
141
142 G4String uniqueName = BDSSamplerRegistry::Instance()->GetNameUnique(samplerID);
143 if (uniqueName != sn)
144 {BDS::Warning("Sampler placement with name \"" + sn + "\" will be named \"" + uniqueName + "\" in the output");}
145
146 G4PVPlacement* pl = new G4PVPlacement(transform,
147 sampler->GetContainerLogicalVolume(),
148 sn + "_pv",
149 samplerWorldLV,
150 false,
151 samplerID,
152 checkOverlaps);
153 placements.push_back(pl);
154 }
155}
156
158 const G4String& variableName,
159 const G4String& objectName) const
160{
161 if (value <= 0)
162 {throw BDSException("\"" + variableName + "\" is <= 0 and must be > 0 in samplerplacement \"" + objectName + "\"");}
163}
164
167 G4double& radius) const
168{
169 BDSSampler* result = nullptr;
170 G4String samplerName = G4String(samplerPlacement.name);
171
172 switch (st.underlying())
173 {
174 case BDSSamplerType::plane:
175 {
176 BDSApertureInfo* shape;
177 if (samplerPlacement.apertureModel.empty())
178 {
179 shape = new BDSApertureInfo(samplerPlacement.shape,
180 samplerPlacement.aper1 * CLHEP::m,
181 samplerPlacement.aper2 * CLHEP::m,
182 samplerPlacement.aper3 * CLHEP::m,
183 samplerPlacement.aper4 * CLHEP::m,
184 samplerName);
185 }
186 else
187 {shape = BDSAcceleratorModel::Instance()->Aperture(samplerPlacement.apertureModel);}
188
189 result = new BDSSamplerCustom(samplerName, *shape, samplerPlacement.partIDSetID);
190 break;
191 }
192 case BDSSamplerType::cylinder:
193 {
194 ErrorNonPositive(samplerPlacement.aper1, "aper1", samplerName);
195 ErrorNonPositive(samplerPlacement.aper2, "aper2", samplerName);
196 G4double startAnglePhi = samplerPlacement.startAnglePhi * CLHEP::rad;
197 G4double sweepAnglePhi = samplerPlacement.sweepAnglePhi * CLHEP::rad;
198 if (sweepAnglePhi <= 0) // default in parser is -1 to flag we should use 2pi now we have units
199 {sweepAnglePhi = CLHEP::twopi;}
200 else if (sweepAnglePhi > CLHEP::twopi + 1e-6)
201 {throw BDSException(__METHOD_NAME__, "\"sweepAnglePhi\" must be in range (0 to pi] in samplerplacement \"" + samplerName + "\"");}
202
203 result = new BDSSamplerCylinder(samplerName,
204 samplerPlacement.aper1 * CLHEP::m,
205 2 * samplerPlacement.aper2 * CLHEP::m,
206 startAnglePhi,
207 sweepAnglePhi,
208 samplerPlacement.partIDSetID);
209 radius = samplerPlacement.aper1 * CLHEP::m;
210 break;
211 }
212 case BDSSamplerType::cylinderforward:
213 {
214 ErrorNonPositive(samplerPlacement.aper1, "aper1", samplerName);
215 ErrorNonPositive(samplerPlacement.aper2, "aper2", samplerName);
216 if (BDS::IsFinite(samplerPlacement.startAnglePhi))
217 {BDS::Warning("\"startAnglePhi\" in samplerplacement \""+samplerName+"\" != 0 -> this has no effect for a cylinderforward sampler");}
218
219 G4double sweepAnglePhi = samplerPlacement.sweepAnglePhi * CLHEP::rad;
220 if (sweepAnglePhi <= 0) // default in parser is -1 to flag we should use 2pi now we have units
221 {sweepAnglePhi = CLHEP::twopi;}
222 else if (sweepAnglePhi > CLHEP::twopi + 1e-6)
223 {throw BDSException(__METHOD_NAME__, "\"sweepAnglePhi\" must be in range (0 to pi] in samplerplacement \"" + samplerName + "\"");}
224 G4double startAnglePhi = -0.5*sweepAnglePhi;
225
226 result = new BDSSamplerCylinder(samplerName,
227 samplerPlacement.aper1 * CLHEP::m,
228 2 * samplerPlacement.aper2 * CLHEP::m,
229 startAnglePhi,
230 sweepAnglePhi,
231 samplerPlacement.partIDSetID);
232 radius = samplerPlacement.aper1 * CLHEP::m;
233 break;
234 }
235 case BDSSamplerType::sphere:
236 {
237 ErrorNonPositive(samplerPlacement.aper1, "aper1", samplerName);
238 G4double startAnglePhi = samplerPlacement.startAnglePhi * CLHEP::rad;
239 G4double sweepAnglePhi = samplerPlacement.sweepAnglePhi * CLHEP::rad;
240 if (sweepAnglePhi <= 0) // default in parser is -1 to flag we should use 2pi now we have units
241 {sweepAnglePhi = CLHEP::twopi;}
242 else if (sweepAnglePhi > CLHEP::twopi + 1e-6)
243 {throw BDSException(__METHOD_NAME__, "\"sweepAnglePhi\" must be in range (0 to 2 pi] in samplerplacement \"" + samplerName + "\"");}
244 G4double startAngleTheta = samplerPlacement.startAngleTheta * CLHEP::rad;
245 if (startAngleTheta < 0)
246 {throw BDSException(__METHOD_NAME__, "\"startAngleTheta\" must be in the range [0 to pi] in samplerplacement \"" + samplerName + "\"");}
247 G4double sweepAngleTheta = samplerPlacement.sweepAngleTheta * CLHEP::rad;
248 if (sweepAngleTheta <= 0) // default in parser is -1 to flag we should use 2pi now we have units
249 {sweepAngleTheta = CLHEP::pi;}
250 else if (sweepAngleTheta > CLHEP::pi + 1e-6)
251 {throw BDSException(__METHOD_NAME__, "\"sweepAngleTheta\" must be in range (0 to pi] in samplerplacement \"" + samplerName + "\"");}
252 else if (sweepAngleTheta - startAngleTheta > CLHEP::pi + 1e-6)
253 {throw BDSException(__METHOD_NAME__, "\"startAngleTheta\" + \"sweepAngleTheta\" must be in range (0 to pi] in samplerplacement \"" + samplerName + "\"");}
254
255 result = new BDSSamplerSphere(samplerName,
256 samplerPlacement.aper1 * CLHEP::m,
257 startAnglePhi,
258 sweepAnglePhi,
259 startAngleTheta,
260 sweepAngleTheta,
261 samplerPlacement.partIDSetID);
262 radius = samplerPlacement.aper1 * CLHEP::m;
263 break;
264 }
265 case BDSSamplerType::sphereforward:
266 {
267 ErrorNonPositive(samplerPlacement.aper1, "aper1", samplerName);
268 if (BDS::IsFinite(samplerPlacement.startAnglePhi))
269 {BDS::Warning("\"startAnglePhi\" in samplerplacement \""+samplerName+"\" != 0 -> this has no effect for a sphereforward sampler");}
270 if (BDS::IsFinite(samplerPlacement.startAngleTheta))
271 {BDS::Warning("\"startAngleTheta\" in samplerplacement \""+samplerName+"\" != 0 -> this has no effect for a sphereforward sampler");}
272
273 G4double sweepAnglePhi = samplerPlacement.sweepAnglePhi * CLHEP::rad;
274 if (sweepAnglePhi <= 0) // default in parser is -1 to flag we should use 2pi now we have units
275 {sweepAnglePhi = CLHEP::twopi;}
276 else if (sweepAnglePhi > CLHEP::twopi + 1e-6)
277 {throw BDSException(__METHOD_NAME__, "\"sweepAnglePhi\" must be in range (0 to 2 pi] in samplerplacement \"" + samplerName + "\"");}
278 G4double sweepAngleTheta = samplerPlacement.sweepAngleTheta * CLHEP::rad;
279 if (sweepAngleTheta <= 0) // default in parser is -1 to flag we should use 2pi now we have units
280 {sweepAngleTheta = CLHEP::pi;}
281 else if (sweepAngleTheta > CLHEP::pi + 1e-6)
282 {throw BDSException(__METHOD_NAME__, "\"sweepAngleTheta\" must be in range (0 to pi] in samplerplacement \"" + samplerName + "\"");}
283
284 G4double startAngleTheta = CLHEP::halfpi -0.5*sweepAngleTheta;
285 G4double startAnglePhi = -0.5*sweepAnglePhi;
286
287 result = new BDSSamplerSphere(samplerName,
288 samplerPlacement.aper1 * CLHEP::m,
289 startAnglePhi,
290 sweepAnglePhi,
291 startAngleTheta,
292 sweepAngleTheta,
293 samplerPlacement.partIDSetID);
294 radius = samplerPlacement.aper1 * CLHEP::m;
295 break;
296 }
297 }
298 return result;
299}
300
302 BDSSamplerType st) const
303{
304 switch (st.underlying())
305 {
306 case BDSSamplerType::cylinderforward:
307 {
308 auto rmExisting = trans.getRotation();
309 G4RotationMatrix rm;
310 rm.rotate(-CLHEP::halfpi, {0,0,1});
311 rm.rotate(-CLHEP::halfpi, {1,0,0});
312 rmExisting *= rm;
313 trans = G4Transform3D(rmExisting, trans.getTranslation());
314 break;
315 }
316 case BDSSamplerType::sphereforward:
317 {
318 auto rmExisting = trans.getRotation();
319 G4RotationMatrix rm;
320 rm.rotate(-CLHEP::halfpi, {0,1,0});
321 rmExisting *= rm;
322 trans = G4Transform3D(rmExisting, trans.getTranslation());
323 break;
324 }
325 default:
326 {break;}
327 }
328}
329
331 G4double samplerRadius)
332{
333 BDSSamplerInfo* samplerInfo = element->GetSamplerInfo();
334 if (!samplerInfo)
335 {return;}
336 BDSSamplerType samplerType = samplerInfo->samplerType;
337 if (samplerType == BDSSamplerType::none)
338 {return;}
339 // else must be a valid sampler
340
341 G4String name = samplerInfo->name;
342 G4double sEnd = element->GetSPositionEnd();
343
344 BDSSampler* sampler = nullptr;
345 switch (samplerType.underlying())
346 {
347 case BDSSamplerType::plane:
348 {sampler = samplerInstances[samplerInfo->pdgSetID]; break;}
349 case BDSSamplerType::cylinder:
350 {// these are built uniquely for each instance so the length matches exactly
351 // for a cylindrical sampler we make it 'fit' the component
352 G4double boundingRadius = element->GetExtent().TransverseBoundingRadius();
353 if (boundingRadius > 0) // protect against bad extent, i.e. 0
354 {samplerRadius = boundingRadius;}
355
356 G4double length = element->GetChordLength();
357 G4double angle = element->GetAngle();
358 if (samplerInfo->startElement && samplerInfo->finishElement)
359 {// to cover a range of split components, e.g. an sbend
360 const auto startElement = samplerInfo->startElement;
361 const auto finishElement = samplerInfo->finishElement;
362 G4ThreeVector connectingVector = finishElement->GetReferencePositionEnd() - startElement->GetReferencePositionStart();
363 G4double chordLength = connectingVector.mag();
364
365 G4ThreeVector directionFrameStart = G4ThreeVector(0,0,1).transform(*(startElement->GetReferenceRotationStart()));
366 G4ThreeVector directionFrameFinish = G4ThreeVector(0,0,1).transform(*(finishElement->GetReferenceRotationEnd()));
367
368 G4ThreeVector dChordFrameStart = connectingVector.unit() - directionFrameStart;
369 G4ThreeVector dChordFrameFinish = connectingVector.unit() - directionFrameFinish;
370
371 // we purposively ignore possible pole face rotations or a different input normal
372 // this is convention for cylindrical samplers where we do the same for drifts
373 G4ThreeVector ipfnCSampler = G4ThreeVector(0,0,-1) + dChordFrameStart;
374 G4ThreeVector opfnCSampler = G4ThreeVector(0,0,1) - dChordFrameFinish;
375
376 sampler = new BDSSamplerCylinder(name,
377 samplerRadius,
378 chordLength,
379 ipfnCSampler,
380 opfnCSampler,
381 0, CLHEP::twopi,
382 samplerInfo->pdgSetID);
383
384 }
385 else if (BDS::IsFinite(angle))
386 {
387 const auto ipfn = element->GetAcceleratorComponent()->InputFaceNormal();
388 const auto opfn = element->GetAcceleratorComponent()->OutputFaceNormal();
389 sampler = new BDSSamplerCylinder(name,
390 samplerRadius,
391 length,
392 ipfn,
393 opfn,
394 0, CLHEP::twopi,
395 samplerInfo->pdgSetID);
396 }
397 else
398 {
399 sampler = new BDSSamplerCylinder(name,
400 samplerRadius,
401 length,
402 0, CLHEP::twopi,
403 samplerInfo->pdgSetID);
404 }
405 break;
406 }
407 default: // no spherical samplers for beam line samplers
408 {break;} // leave as nullptr - shouldn't occur due to if at top
409 }
410
411 if (sampler)
412 {
413 G4Transform3D* pt = new G4Transform3D(*element->GetSamplerPlacementTransform());
414
415 G4int samplerID = BDSSamplerRegistry::Instance()->RegisterSampler(name,
416 sampler,
417 *pt,
418 sEnd,
419 element,
420 samplerType,
421 samplerRadius);
422
423 G4VPhysicalVolume* samplerWorld = GetWorld();
424 samplerWorldLV = samplerWorld->GetLogicalVolume();
425 const G4bool checkOverlaps = BDSGlobalConstants::Instance()->CheckOverlaps();
426 // record placements for cleaning up at destruction.
427 G4PVPlacement* pl = new G4PVPlacement(*pt, // placement transform
428 sampler->GetContainerLogicalVolume(), // logical volume
429 name + "_pv", // name of placement
430 samplerWorldLV, // mother volume
431 false, // no boolean operation
432 samplerID, // copy number
433 checkOverlaps);
434
435 placements.push_back(pl);
436 }
437}
G4ThreeVector OutputFaceNormal() const
G4ThreeVector InputFaceNormal() const
BDSApertureInfo * Aperture(G4String name) const
const BDSBeamline * BeamlineMain() const
Accessor.
Holder class for all information required to describe an aperture.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4double GetChordLength() const
Accessor.
G4double GetAngle() const
Accessor.
G4ThreeVector GetReferencePositionEnd() const
Accessor.
BDSExtent GetExtent() const
Accessor.
BDSSamplerInfo * GetSamplerInfo() const
Accessor.
G4double GetSPositionEnd() const
Accessor.
BDSAcceleratorComponent * GetAcceleratorComponent() const
Accessor.
G4Transform3D * GetSamplerPlacementTransform() const
Accessor.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
static G4Transform3D CreatePlacementTransform(const GMAD::Placement &placement, const BDSBeamline *beamLine, G4double *S=nullptr, BDSExtent *placementExtent=nullptr, const G4String &objectTypeForErrorMsg="placement")
General exception with possible name of object and message.
Definition: BDSException.hh:35
G4double TransverseBoundingRadius() const
Return a radius that would encompass the maximum x,y extent.
Definition: BDSExtent.cc:183
G4LogicalVolume * GetContainerLogicalVolume() const
Accessor - see member for more info.
A class that holds global options and constants.
static BDSGlobalConstants * Instance()
Access method.
static void PrintProtectedNames(std::ostream &out)
Feedback for protected names.
Definition: BDSOutput.cc:383
static G4bool InvalidSamplerName(const G4String &samplerName)
Test whether a sampler name is invalid or not.
Definition: BDSOutput.cc:378
BDSParallelWorldSampler()=delete
No default constructor.
BDSSamplerPlane * generalPlane
General single sampler we use for plane samplers.
void AdjustTransform(G4Transform3D &trans, BDSSamplerType st) const
void ErrorNonPositive(G4double value, const G4String &variableName, const G4String &objectName) const
Utility function to reduce code.
void Place(const BDSBeamlineElement *element, G4double samplerRadius)
Place a sampler from a single element.
std::vector< G4VPhysicalVolume * > placements
Cache of the placements to clean up at the end.
BDSSampler * BuildSampler(const GMAD::SamplerPlacement &samplerPlacement, BDSSamplerType st, G4double &radius) const
Construct the geometry for a sampler. Update 'radius' by reference if applicable.
G4VisAttributes * samplerWorldVis
Visualisation attributes for the sampler world.
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
std::vector< GMAD::SamplerPlacement > GetSamplerPlacements() const
Return sampler placement list.
Definition: BDSParser.hh:123
Custom shaped sampler with fixed thickness.
Cylindrical sampler around an object.
All info required to build a sampler but not place it.
Rectangular sampler with fixed thickness but variable x,y.
static BDSSamplerRegistry * Instance()
Accessor for registry.
G4String GetNameUnique(G4int index) const
Accessor.
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)
Spherical sampler around an object.
Base class and registry of sampler instances.
Definition: BDSSampler.hh:36
type underlying() const
return underlying value (can be used in switch statement)
Sampler placement class for parser.
std::string name
Name of this samplerplacement.
int partIDSetID
The unique ID of the particle set given by the parser.
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())