BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSDetectorConstruction.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 "BDSAcceleratorComponent.hh"
20#include "BDSAcceleratorComponentRegistry.hh"
21#include "BDSAcceleratorModel.hh"
22#include "BDSApertureInfo.hh"
23#include "BDSAuxiliaryNavigator.hh"
24#include "BDSBeamline.hh"
25#include "BDSBeamlineBLMBuilder.hh"
26#include "BDSBeamlineEndPieceBuilder.hh"
27#include "BDSBeamlineElement.hh"
28#include "BDSBeamlinePlacementBuilder.hh"
29#include "BDSBeamlineSet.hh"
30#include "BDSBeamPipeInfo.hh"
31#include "BDSBLM.hh"
32#include "BDSBLMRegistry.hh"
33#include "BDSBOptrMultiParticleChangeCrossSection.hh"
34#include "BDSComponentFactory.hh"
35#include "BDSComponentFactoryUser.hh"
36#include "BDSCurvilinearBuilder.hh"
37#include "BDSDebug.hh"
38#include "BDSDetectorConstruction.hh"
39#include "BDSException.hh"
40#include "BDSExtent.hh"
41#include "BDSFieldBuilder.hh"
42#include "BDSFieldObjects.hh"
43#include "BDSFieldQuery.hh"
44#include "BDSFieldQueryInfo.hh"
45#include "BDSFieldLoaderQueryPoints.hh"
46#include "BDSGap.hh"
47#include "BDSGeometryComponent.hh"
48#include "BDSGeometryExternal.hh"
49#include "BDSGeometryFactory.hh"
50#include "BDSGlobalConstants.hh"
51#include "BDSHistBinMapper.hh"
52#include "BDSIntegratorSet.hh"
53#include "BDSLine.hh"
54#include "BDSMaterials.hh"
55#include "BDSParser.hh"
56#include "BDSPhysicalVolumeInfo.hh"
57#include "BDSPhysicalVolumeInfoRegistry.hh"
58#include "BDSRegion.hh"
59#include "BDSSamplerInfo.hh"
60#include "BDSSamplerType.hh"
61#include "BDSScorerFactory.hh"
62#include "BDSScorerInfo.hh"
63#include "BDSScorerMeshInfo.hh"
64#include "BDSScoringMeshBox.hh"
65#include "BDSScoringMeshCylinder.hh"
66#include "BDSSDEnergyDeposition.hh"
67#include "BDSSDManager.hh"
68#include "BDSSDType.hh"
69#include "BDSSurvey.hh"
70#include "BDSTeleporter.hh"
71#include "BDSTrajectoryPoint.hh"
72#include "BDSTunnelBuilder.hh"
73#include "BDSUtilities.hh"
74#include "BDSWarning.hh"
75
76#include "parser/blmplacement.h"
77#include "parser/element.h"
78#include "parser/fastlist.h"
79#include "parser/options.h"
80#include "parser/physicsbiasing.h"
81#include "parser/placement.h"
82#include "parser/samplerplacement.h"
83#include "parser/scorermesh.h"
84
85#include "globals.hh"
86#include "G4AffineTransform.hh"
87#include "G4Box.hh"
88#include "G4LogicalVolume.hh"
89#include "G4Material.hh"
90#include "G4ProductionCuts.hh"
91#include "G4PVPlacement.hh"
92#include "G4VPrimitiveScorer.hh"
93#include "G4Region.hh"
94#include "G4ScoringManager.hh"
95#include "G4String.hh"
96#include "G4Transform3D.hh"
97#include "G4Version.hh"
98#include "G4VisAttributes.hh"
99#include "G4VPhysicalVolume.hh"
100#if G4VERSION_NUMBER > 1039
101#include "G4ChannelingOptrMultiParticleChangeCrossSection.hh"
102#endif
103
104#ifdef BDSCHECKUSERLIMITS
105#include "G4UserLimits.hh"
106#endif
107
108#include "CLHEP/Units/SystemOfUnits.h"
109#include "CLHEP/Vector/EulerAngles.h"
110
111#include <iterator>
112#include <limits>
113#include <list>
114#include <map>
115#include <set>
116#include <sstream>
117#include <string>
118#include <utility>
119#include <vector>
120
121
122BDSDetectorConstruction::BDSDetectorConstruction(BDSComponentFactoryUser* userComponentFactoryIn):
123 placementBL(nullptr),
124 designParticle(nullptr),
125 userComponentFactory(userComponentFactoryIn),
126 nSamplers(0),
127 buildPlacementFieldsWorld(false),
128 worldLogicalVolume(nullptr)
129{
131 verbose = globals->Verbose();
132 checkOverlaps = globals->CheckOverlaps();
133 circular = globals->Circular();
134
135 BDSTrajectoryPoint::dEThresholdForScattering = globals->DEThresholdForScattering();
136
137 // instantiate the accelerator model holding class
138 acceleratorModel = BDSAcceleratorModel::Instance();
139 canSampleAngledFaces = true;
140 BDSIntegratorSetType integratorSetType = globals->IntegratorSet();
141 if ( (integratorSetType == BDSIntegratorSetType::bdsimtwo)
142 || (integratorSetType == BDSIntegratorSetType::geant4)
143#if G4VERSION_NUMBER > 1039
144 || (integratorSetType == BDSIntegratorSetType::geant4dp)
145#endif
146 )
147 { // set to be value of option, default is false.
148 canSampleAngledFaces = globals->SampleElementsWithPoleface();
149 }
150
151 UpdateSamplerDiameterAndCountSamplers();
152 PrepareExtraSamplerSDs();
153 CountPlacementFields();
154}
155
157{
158 nSamplers = 0;
159 auto beamline = BDSParser::Instance()->GetBeamline(); // main beam line
160 G4double maxBendingRatio = 1e-9;
161 for (const auto& blElement : beamline)
162 {
163 // count number of samplers
164 auto st = BDS::DetermineSamplerType(blElement.samplerType);
165 if (st != BDSSamplerType::none)
166 {nSamplers++;}
167
168 G4double length = blElement.l;
169 G4double angle = blElement.angle;
170 if (!BDS::IsFinite(length))
171 {continue;} // avoid divide by 0
172 G4double ratio = angle / length;
173 maxBendingRatio = std::max(maxBendingRatio, ratio);
174 }
175
177 G4double curvilinearRadius = 0.5*globals->CurvilinearDiameter();
178 G4double tolerance = 0.9; // 10% tolerance -> factor of 0.9
179 if (maxBendingRatio > 0.4*tolerance) // max ratio for a 2.5m sampler diameter
180 {
181 G4double curvilinearRadiusBends = (tolerance / maxBendingRatio)*CLHEP::m;
182 if (curvilinearRadiusBends < curvilinearRadius)
183 {
184 G4cout << __METHOD_NAME__ << "Reducing curvilinear diameter from " << 2*curvilinearRadius / CLHEP::m
185 << "m to " << 2*curvilinearRadiusBends / CLHEP::m << "m" << G4endl;
186 globals->SetCurvilinearDiameter(2*curvilinearRadiusBends);
188 }
189 G4double sd = globals->SamplerDiameter();
190 if (curvilinearRadius*2 < sd)
191 {
192 G4cout << __METHOD_NAME__ << "Reducing sampler diameter from " << sd / CLHEP::m << "m to the same" << G4endl;
193 globals->SetSamplerDiameter(2*curvilinearRadius);
194 }
195 }
196
197 // add number of sampler placements to count of samplers
199}
200
202{
203 const auto& samplerFilterIDtoPDGSet = BDSParser::Instance()->GetSamplerFilterIDToSet();
204 BDSSDManager::Instance()->ConstructSamplerSDsForParticleSets(samplerFilterIDtoPDGSet);
205}
206
208{
209 G4int nFields = 0;
210 const auto& placements = BDSParser::Instance()->GetPlacements();
211 for (const auto& placement : placements)
212 {// here we assume if a bdsim element is used at all that it's active even though it may not be
213 if (!placement.fieldAll.empty() || !placement.bdsimElement.empty())
214 {nFields++;}
215 }
216 buildPlacementFieldsWorld = nFields > 0;
217}
218
220{
221 if (verbose || debug)
222 {G4cout << __METHOD_NAME__ << "starting accelerator geometry construction\n" << G4endl;}
223
224 // construct all parser defined regions
226
227 // construct all parser defined aperture objects
229
230 // construct the main beam line and any other secondary beam lines
232
233 // construct placement geometry from parser
234 BDSBeamline* mainBeamLine = BDSAcceleratorModel::Instance()->BeamlineSetMain().massWorld;
235 auto componentFactory = new BDSComponentFactory(designParticle, userComponentFactory);
237 mainBeamLine,
238 componentFactory);
239 BDSAcceleratorModel::Instance()->RegisterPlacementBeamline(placementBL); // Acc model owns it
240 delete componentFactory;
241
243 mainBeamLine);
244 if (blms)
245 {BDSAcceleratorModel::Instance()->RegisterBLMs(blms);} // Acc model owns it
246
247 // build the tunnel and supports
249 {BuildTunnel();}
250
251 // build world and calculate coordinates
252 auto worldPV = BuildWorld();
253
254 // placement procedure - put everything in the world
255 ComponentPlacement(worldPV);
256
257 if (verbose || debug)
258 {G4cout << __METHOD_NAME__ << "detector Construction done" << G4endl;}
259
260 fieldQueries = BDSDetectorConstruction::PrepareFieldQueries(mainBeamLine);
261
262#ifdef BDSDEBUG
263 G4cout << G4endl << __METHOD_NAME__ << "printing material table" << G4endl;
264 G4cout << *(G4Material::GetMaterialTable()) << G4endl << G4endl;
265 if (verbose || debug) {G4cout << "Finished listing materials, returning physiWorld" << G4endl;}
266#endif
267#ifdef BDSCHECKUSERLIMITS
268 PrintUserLimitsSummary(worldPV);
269#endif
270 return worldPV;
271}
272
273BDSDetectorConstruction::~BDSDetectorConstruction()
274{
275#if G4VERSION_NUMBER > 1009
276 // delete bias objects
277 for (auto i : biasObjects)
278 {delete i;}
279#endif
280 for (auto q : fieldQueries)
281 {delete q;}
282}
283
285{
286 BDSRegion* defaultRegion = new BDSRegion("default");
287 for (const GMAD::Region& r : BDSParser::Instance()->GetRegions())
288 {
289 BDSRegion* reg = new BDSRegion(r, defaultRegion);
290 G4cout << "New region defined: " << G4endl << *reg << G4endl;
292 }
293 delete defaultRegion;
294}
295
297{
298 std::map<G4String, BDSApertureInfo*> apertures;
299 for (const GMAD::Aperture& a : BDSParser::Instance()->GetApertures())
300 {
301 BDSApertureInfo* ap = new BDSApertureInfo(a.apertureType,
302 a.aper1 * CLHEP::m,
303 a.aper2 * CLHEP::m,
304 a.aper3 * CLHEP::m,
305 a.aper4 * CLHEP::m,
306 a.name);
307 apertures[a.name] = ap;
308 }
310}
311
313{
314 // build main beam line
315 if (verbose || debug)
316 {G4cout << "parsing the beamline element list..."<< G4endl;}
317 G4Transform3D initialTransform = BDSGlobalConstants::Instance()->BeamlineTransform();
318 G4double initialS = BDSGlobalConstants::Instance()->BeamlineS();
319
320 BDSBeamlineSet mainBeamline = BuildBeamline(BDSParser::Instance()->GetBeamline(),
321 "main beam line",
322 initialTransform,
323 initialS,
324 circular);
325
326#ifdef BDSDEBUG
327 G4cout << "Registry size "
329 G4cout << "Parser beam line size "
330 << BDSParser::Instance()->GetBeamline().size() << G4endl;
332#endif
333
334 // print warning if beam line is approximately circular but flag isn't specified
335 if (!circular && mainBeamline.massWorld)
336 {
337 if (mainBeamline.massWorld->ElementAnglesSumToCircle())
338 {BDS::Warning("Total sum of all element angles is approximately 2*pi but the circular option was not specified, this simulation may run indefinitely");}
339 }
340
341 // register the beam line in the holder class for the full model
343
344 // build secondary beam lines
345 // loop over placements and check if any are beam lines (have sequences specified)
346 auto placements = BDSParser::Instance()->GetPlacements();
347 for (const auto& placement : placements)
348 {
349 if (placement.sequence.empty())
350 {continue;} // no sequence specified -> just a placement
351 auto parserLine = BDSParser::Instance()->GetSequence(placement.sequence);
352
353 // determine offset in world for extra beam line
354 const BDSBeamline* mbl = mainBeamline.massWorld;
355 // TODO - so by default if placement.s is finite, it'll be made w.r.t. the main beam line
356 // but this could be any beam line in future if we find the right beam line to pass in.
357 G4Transform3D startTransform = CreatePlacementTransform(placement, mbl);
358 G4double startS = mbl ? mbl->back()->GetSPositionEnd() : 0;
359
360 // aux beam line must be non-circular by definition to branch off of beam line (for now)
361 // TODO - the naming convention here is repeated in BDSParallelWorldInfo which is registered
362 // beforehand separately - fix by making the information originate in one place despite
363 // the parallel world instantiated first before Construct() in this class is called.
364 G4String beamlineName = placement.name + "_" + placement.sequence;
365 BDSBeamlineSet extraBeamline = BuildBeamline(parserLine,
366 beamlineName,
367 startTransform,
368 startS,
369 false, // circular
370 true); // is placement
371
372 acceleratorModel->RegisterBeamlineSetExtra(beamlineName, extraBeamline);
373 }
374}
375
377{
378 BDSSamplerType sType = BDS::DetermineSamplerType(element->samplerType);
379 if (sType == BDSSamplerType::none)
380 {return nullptr;}
381 BDSSamplerInfo* result = new BDSSamplerInfo(element->samplerName,
382 sType,
383 (G4int)element->samplerParticleSetID);
384 return result;
385}
386
388 const G4String& name,
389 const G4Transform3D& initialTransform,
390 G4double initialS,
391 G4bool beamlineIsCircular,
392 G4bool isPlacementBeamline)
393{
394 if (beamLine.empty()) // note a line always has a 'line' element first so an empty line will not be 'empty'
395 {return BDSBeamlineSet();}
396
397 if (userComponentFactory)
398 {userComponentFactory->SetDesignParticle(designParticle);}
399 BDSComponentFactory* theComponentFactory = new BDSComponentFactory(designParticle, userComponentFactory);
400 BDSBeamline* massWorld = new BDSBeamline(initialTransform, initialS);
401
402 if (beamlineIsCircular)
403 {
404 G4bool unsuitable = UnsuitableFirstElement(beamLine.begin());
405 if (unsuitable)
406 {
407 G4cerr << "The first element in the beam line is unsuitable for a circular "
408 << "model as the first element will " << G4endl << "overlap with the "
409 << "teleporter and terminator - the necessary mechanics for a circular "
410 << "model in Geant4" << G4endl;
411 throw BDSException(__METHOD_NAME__, "check construction for circular machine");
412 }
413 }
414
415 if (beamLine.size() <= 1) // if an empty LINE it still has 1 item in it
416 {throw BDSException(__METHOD_NAME__, "BDSIM requires the sequence defined with the use command to have at least one element.");}
417
418 for (auto elementIt = beamLine.begin(); elementIt != beamLine.end(); ++elementIt)
419 {
420 // find next and previous element, but ignore special elements or thin elements.
421 const GMAD::Element* prevElement = nullptr;
422 auto prevIt = elementIt;
423 while (prevIt != beamLine.begin())
424 {
425 --prevIt;
426 if (prevIt->isSpecial() == false && prevIt->l > BDSGlobalConstants::Instance()->ThinElementLength())
427 {
428 prevElement = &(*prevIt);
429 break;
430 }
431 }
432
433 const GMAD::Element* nextElement = nullptr;
434 auto nextIt = elementIt;
435 ++nextIt;
436 G4double nextElementInputFace = 0; // get poleface angle for next element whilst testing if next element exists
437 while (nextIt != beamLine.end())
438 {
439 if (nextIt->isSpecial() == false && nextIt->l > BDSGlobalConstants::Instance()->ThinElementLength())
440 {
441 nextElement = &(*nextIt);
442 //rotated entrance face of the next element may modify the exit face of the current element.
443 nextElementInputFace = nextElement->e1;
444 break;
445 }
446 ++nextIt;
447 }
448 G4double currentArcLength = massWorld->GetTotalArcLength();
449 BDSAcceleratorComponent* temp = theComponentFactory->CreateComponent(&(*elementIt),
450 prevElement,
451 nextElement,
452 currentArcLength);
453 if(temp)
454 {
455 G4bool forceNoSamplerOnThisElement = false;
456 if ((!canSampleAngledFaces) && (BDS::IsFinite((*elementIt).e2)))
457 {forceNoSamplerOnThisElement = true;}
458 if ((!canSampleAngledFaces) && (BDS::IsFinite(nextElementInputFace)))
459 {forceNoSamplerOnThisElement = true;}
460 if (temp->GetType() == "dump") // don't sample after a dump as there'll be nothing
461 {forceNoSamplerOnThisElement = true;}
462 BDSSamplerInfo* samplerInfo = forceNoSamplerOnThisElement ? nullptr : BuildSamplerInfo(&(*elementIt));
463 BDSTiltOffset* tiltOffset = BDSComponentFactory::CreateTiltOffset(&(*elementIt));
464 massWorld->AddComponent(temp, tiltOffset, samplerInfo);
465 }
466 }
467
468 // Special circular machine bits
469 // Add terminator to do ring turn counting logic and kill particles
470 // Add teleporter to account for slight ring offset
471 if (beamlineIsCircular && !massWorld->empty())
472 {
473#ifdef BDSDEBUG
474 G4cout << __METHOD_NAME__ << "Circular machine - creating terminator & teleporter" << G4endl;
475#endif
476 G4double teleporterLength = 0;
477 G4Transform3D teleporterTransform = BDS::CalculateTeleporterDelta(massWorld, teleporterLength);
478
479 auto hasBeamPipeInfo = [](BDSBeamlineElement* ble) {return ble->GetBeamPipeInfo() != nullptr;};
480 auto firstElementWithBPInfo = std::find_if(massWorld->begin(), massWorld->end(), hasBeamPipeInfo);
481 auto lastElementWithBPInfo = std::find_if(massWorld->rbegin(), massWorld->rend(), hasBeamPipeInfo);
482
483 G4double firstbeamPipeMaxExtent = (*firstElementWithBPInfo)->GetBeamPipeInfo()->Extent().MaximumAbsTransverse();
484 G4double lastbeamPipeMaxExtent = (*lastElementWithBPInfo)->GetBeamPipeInfo()->Extent().MaximumAbsTransverse();
485
486 // the extent is a half width, so we double it - also the terminator width.
487 G4double teleporterHorizontalWidth = 2 * std::max(firstbeamPipeMaxExtent, lastbeamPipeMaxExtent);
488
489 BDSAcceleratorComponent* terminator = theComponentFactory->CreateTerminator(teleporterHorizontalWidth);
490 if (terminator)
491 {
492 terminator->Initialise();
493 massWorld->AddComponent(terminator);
494 }
495
496 BDSAcceleratorComponent* teleporter = theComponentFactory->CreateTeleporter(teleporterLength,
497 teleporterHorizontalWidth,
498 teleporterTransform);
499 if (teleporter)
500 {
501 teleporter->Initialise();
502 massWorld->AddComponent(teleporter);
503 }
504 }
505
506 if (BDSGlobalConstants::Instance()->Survey())
507 {
508 G4String surveyFileName = BDSGlobalConstants::Instance()->SurveyFileName() + ".dat";
509 if (isPlacementBeamline)
510 {surveyFileName = BDSGlobalConstants::Instance()->SurveyFileName() + "_" + name + ".dat";}
511 BDSSurvey* survey = new BDSSurvey(surveyFileName);
512 survey->Write(massWorld);
513 delete survey;
514 }
515 delete theComponentFactory;
516
517 // print summary
518 G4cout << __METHOD_NAME__ << "\"" << name << "\" " << G4endl << *massWorld;
519
520 // Build curvilinear geometry w.r.t. beam line.
522 BDSBeamline* clBeamline = clBuilder->BuildCurvilinearBeamLine1To1(massWorld, beamlineIsCircular);
523 BDSBeamline* clBridgeBeamline = clBuilder->BuildCurvilinearBridgeBeamLine(clBeamline);
524 delete clBuilder;
525
526 // construct beamline of end pieces
527 BDSBeamline* endPieces = BDS::BuildEndPieceBeamline(massWorld, circular);
528 BDSBeamlineSet beamlineSet;
529 beamlineSet.massWorld = massWorld;
530 beamlineSet.curvilinearWorld = clBeamline;
531 beamlineSet.curvilinearBridgeWorld = clBridgeBeamline;
532 beamlineSet.endPieces = endPieces;
533 return beamlineSet;
534}
535
537{
538 const BDSBeamline* mainBeamLine = acceleratorModel->BeamlineMain();
539 BDSBeamline* tunnelBeamline;
540 BDSTunnelBuilder* tb = new BDSTunnelBuilder(BDSGlobalConstants::Instance()->TunnelMaxSegmentLength());
541 tunnelBeamline = tb->BuildTunnelSections(mainBeamLine);
542 delete tb;
543
545}
546
548{
549 // calculate extents of everything we need to place in the world first
550 std::vector<BDSExtentGlobal> extents;
551
552 // these beam lines should always exist so are safe to access.
553 const auto& blMain = acceleratorModel->BeamlineSetMain();
554 blMain.GetExtentGlobals(extents);
555
556 // check optional placement beam line (like vector of placements)
558 if (plBeamline) // optional - may be nullptr
559 {extents.push_back(plBeamline->GetExtentGlobal());}
560
561 // check tunnel beam line
562 BDSBeamline* tunnelBeamline = acceleratorModel->TunnelBeamline();
563 if (tunnelBeamline)
564 {extents.push_back(tunnelBeamline->GetExtentGlobal());}
565
566 // check extra beam lines
567 const auto& extras = BDSAcceleratorModel::Instance()->ExtraBeamlines();
568 // extras is a map, so iterator has first and second for key and value
569 for (const auto& bl : extras)
570 {bl.second.GetExtentGlobals(extents);}
571
572 // inspect any sampler placements and calculate their extent without constructing them
573 extents.push_back(CalculateExtentOfSamplerPlacements(blMain.massWorld));
574
575 // inspect any scoring meshes and calculate their extent without constructing them
576 extents.push_back(CalculateExtentOfScorerMeshes(blMain.massWorld));
577
578 // Expand to maximum extents of each beam line.
579 G4ThreeVector worldR;
580 // loop over all extents from all beam lines
581 for (const auto& ext : extents)
582 {
583 for (G4int i = 0; i < 3; i++)
584 {worldR[i] = std::max(worldR[i], ext.GetMaximumExtentAbsolute()[i]);} // expand with the maximum
585 }
586
587 G4String worldName = "World";
588 G4VSolid* worldSolid = nullptr;
589 G4LogicalVolume* worldLV = nullptr;
590
591 G4String worldGeometryFile = BDSGlobalConstants::Instance()->WorldGeometryFile();
592 if (!worldGeometryFile.empty())
593 {
594 if (BDSGlobalConstants::Instance()->WorldMaterialSet())
595 {BDS::Warning(__METHOD_NAME__, "conflicting options - world material option specified but material will be taken from world GDML file");}
596 G4bool ac = BDSGlobalConstants::Instance()->AutoColourWorldGeometryFile();
597
598 std::vector<G4String> namedWorldVacuumVolumes = BDS::SplitOnWhiteSpace(BDSGlobalConstants::Instance()->WorldVacuumVolumeNames());
599
601 worldGeometryFile,
602 nullptr,
603 ac,
604 0, 0,
605 &namedWorldVacuumVolumes,
606 true,
607 BDSSDType::energydepworldcontents);
608
609 // get list of 'material' and 'vacuum' volumes for possible biasing of this geometry
610 worldVacuumLogicalVolumes = geom->VacuumVolumes();
611 auto allWorldVolumes = geom->GetAllLogicalVolumes();
612 allWorldVolumes.erase(geom->GetContainerLogicalVolume());
613 for (auto volume : worldVacuumLogicalVolumes)
614 {allWorldVolumes.erase(volume);}
615 worldContentsLogicalVolumes = allWorldVolumes; // cache volumes for biasing
616
617 worldExtent = geom->GetExtent();
618 BDSExtentGlobal worldExtentGlobal = BDSExtentGlobal(worldExtent, G4Transform3D());
619 G4bool worldContainsAllBeamlines = worldExtentGlobal.Encompasses(extents);
620
621 G4cout << "External world geometry: \"" << worldGeometryFile << "\"" << G4endl;
622 G4cout << "Loaded world extent: \n" << worldExtent << G4endl;
623
624 // cannot construct world if any beamline extent is greater than the world extents
625 if (!worldContainsAllBeamlines)
626 {
627 G4String message = "Beamlines cannot be constructed, beamline extents are larger than \n";
628 message += "the extents of the external world";
629 throw BDSException(__METHOD_NAME__, message);
630 }
631
632 worldSolid = geom->GetContainerSolid();
633 worldLV = geom->GetContainerLogicalVolume();
634
635 // make the world sensitive to energy deposition with its own unique hits collection
636 // this will be a nullptr depending on the options.
637 // make world sensitive if importance sampling is needed
638 if (BDSGlobalConstants::Instance()->StoreELossWorld()
639 || BDSGlobalConstants::Instance()->UseImportanceSampling()
640 || BDSGlobalConstants::Instance()->StoreELossWorldContents())
641 {
643 // override the logical volume itself with a specific SD
644 worldLV->SetSensitiveDetector(BDSSDManager::Instance()->WorldComplete());
645 }
646 }
647 else
648 {
649 // add on margin for constructed world volume
650#ifdef BDSDEBUG
651 G4cout << __METHOD_NAME__ << "world extent absolute: " << worldR << G4endl;
652#endif
653 G4double margin = BDSGlobalConstants::Instance()->WorldVolumeMargin();
654 margin = std::max(margin, 2*CLHEP::m); // minimum margin of 2m.
655 worldR += G4ThreeVector(margin,margin,margin); //add 5m extra in every dimension
656#ifdef BDSDEBUG
657 G4cout << __METHOD_NAME__ << "with " << margin << "m margin, it becomes in all dimensions: " << worldR << G4endl;
658#else
659 G4cout << __METHOD_NAME__ << "World dimensions: " << worldR / CLHEP::m << " m" << G4endl;
660#endif
661 G4String worldMaterialName = BDSGlobalConstants::Instance()->WorldMaterial();
662 G4Material* worldMaterial = BDSMaterials::Instance()->GetMaterial(worldMaterialName);
663 worldExtent = BDSExtent(worldR);
664 worldSolid = new G4Box(worldName + "_solid", worldR.x(), worldR.y(), worldR.z());
665
666
667 worldLV = new G4LogicalVolume(worldSolid, // solid
668 worldMaterial, // material
669 worldName + "_lv"); // name
670
671 // make the world sensitive to energy deposition with its own unique hits collection
672 if (BDSGlobalConstants::Instance()->StoreELossWorld())
673 {worldLV->SetSensitiveDetector(BDSSDManager::Instance()->WorldComplete());}
674 }
675
676 // visual attributes
677 // copy the debug vis attributes but change to force wireframe
678 G4VisAttributes* debugWorldVis = new G4VisAttributes(*(BDSGlobalConstants::Instance()->ContainerVisAttr()));
679 debugWorldVis->SetForceWireframe(true);//just wireframe so we can see inside it
680 worldLV->SetVisAttributes(debugWorldVis);
681
682 // set limits
683 worldLV->SetUserLimits(BDSGlobalConstants::Instance()->DefaultUserLimits());
684
685 // place the world
686 G4VPhysicalVolume* worldPV = new G4PVPlacement(nullptr, // no rotation
687 G4ThreeVector(), // at (0,0,0)
688 worldLV, // its logical volume
689 worldName, // its name
690 nullptr, // its mother volume
691 false, // no boolean operation
692 0, // copy number
693 checkOverlaps); // overlap checking
694
695 // Register the lv & pvs to the our holder class for the model
699
700 // Register world PV with our auxiliary navigator so steppers and magnetic
701 // fields know which geometry to navigate to get local / global transforms.
702 // This is the regular world used as a backup to the curvilinear world.
705
708
709 worldLogicalVolume = worldLV; // set for possible biasing
710 return worldPV;
711}
712
713void BDSDetectorConstruction::ComponentPlacement(G4VPhysicalVolume* worldPV)
714{
715 // We musn't place parallel world geometry here - their world is produced by
716 // Geant4 at the right time, so we have a separate placement call for them
717 BDSBeamlineSet mainBL = BDSAcceleratorModel::Instance()->BeamlineSetMain();
718 PlaceBeamlineInWorld(mainBL.massWorld,
719 worldPV, checkOverlaps, true);
720 PlaceBeamlineInWorld(mainBL.endPieces,
721 worldPV, checkOverlaps);
723 {
725 worldPV, checkOverlaps);
726 }
727 // No energy counter SD added here as individual placements have that attached
728 // during construction time
730
731 // Place BLMs. Similarly no sensitivity set here - done at construction time.
732 PlaceBeamlineInWorld(BDSAcceleratorModel::Instance()->BLMsBeamline(),
733 worldPV,
735 false,
736 false,
737 false,
738 true); // use incremental copy numbers
739
740 const auto& extras = BDSAcceleratorModel::Instance()->ExtraBeamlines();
741 for (auto const& bl : extras)
742 {// extras is a map so iterator has first and second for key and value
743 // note these are currently not sensitive as there's no CL frame for them
744 PlaceBeamlineInWorld(bl.second.massWorld, worldPV, checkOverlaps);
745 PlaceBeamlineInWorld(bl.second.endPieces, worldPV, checkOverlaps);
746 }
747}
748
750 G4VPhysicalVolume* containerPV,
751 G4bool checkOverlaps,
752 G4bool setRegions,
753 G4bool registerInfo,
754 G4bool useCLPlacementTransform,
755 G4bool useIncrementalCopyNumbers)
756{
757 if (!beamline)
758 {return;}
759
760 G4int i = 0;
761 for (auto element : *beamline)
762 {
763 // Do nothing for gap element
764 if (dynamic_cast<BDSGap*>(element->GetAcceleratorComponent()))
765 {continue;}
766
767 if (setRegions)
768 {
769 auto accComp = element->GetAcceleratorComponent();
770 const G4String regionName = accComp->GetRegion();
771 if(!regionName.empty()) // ie string is defined so we should attach region
772 {
773 G4Region* region = BDSAcceleratorModel::Instance()->Region(regionName);
774 auto contLV = accComp->GetContainerLogicalVolume();
775 contLV->SetRegion(region);
776 region->AddRootLogicalVolume(contLV);
777 }
778 }
779
780 // setup the sensitivity
781 element->GetAcceleratorComponent()->AttachSensitiveDetectors();
782
783 // make the placement
784 G4int copyNumber = useIncrementalCopyNumbers ? i : element->GetCopyNo();
785 G4String placementName = element->GetPlacementName() + "_pv";
786 std::set<G4VPhysicalVolume*> pvs = element->PlaceElement(placementName, containerPV, useCLPlacementTransform,
787 copyNumber, checkOverlaps);
788
789 if (registerInfo)
790 {
791 BDSPhysicalVolumeInfo* theinfo = new BDSPhysicalVolumeInfo(element->GetName(),
792 placementName,
793 element->GetSPositionMiddle(),
794 element->GetIndex(),
795 beamline);
796
798 }
799 i++; // for incremental copy numbers
800 }
801}
802
804 const BDSBeamline* beamLine,
805 G4double* S,
806 BDSExtent* placementExtent,
807 const G4String& objectTypeForErrorMsg)
808{
809 G4Transform3D result;
810
811 // 3 scenarios
812 // 1) global placement X,Y,Z + rotation
813 // 2) w.r.t. beam line placement x,y,S + rotation
814 // 3) w.r.t. element in beam line placement elementName + x,y,s + rotation
815
816 // in all cases, need the rotation
817 G4RotationMatrix rm = G4RotationMatrix();
818 if (placement.axisAngle)
819 {
820 G4ThreeVector axis = G4ThreeVector(placement.axisX,
821 placement.axisY,
822 placement.axisZ);
823 rm = G4RotationMatrix(axis, placement.angle*CLHEP::rad);
824 }
825 else
826 {
827 if (BDS::IsFinite(placement.phi) ||
828 BDS::IsFinite(placement.theta) ||
829 BDS::IsFinite(placement.psi))
830 {// only build if finite
831 CLHEP::HepEulerAngles ang = CLHEP::HepEulerAngles(placement.phi*CLHEP::rad,
832 placement.theta*CLHEP::rad,
833 placement.psi*CLHEP::rad);
834 rm = G4RotationMatrix(ang);
835 }
836 }
837
838 // create a transform w.r.t. the beam line if s is finite and it's not w.r.t a
839 // particular element. If it's w.r.t. a particular element, treat s as local curvilinear
840 // s and use as local 'z' in the transform.
841 if (!placement.referenceElement.empty())
842 {// scenario 3
843 if (!beamLine)
844 {throw BDSException(__METHOD_NAME__, "no valid beam line yet for " + objectTypeForErrorMsg + " w.r.t. a beam line.");}
845 const BDSBeamlineElement* element = beamLine->GetElement(placement.referenceElement,
846 placement.referenceElementNumber);
847 if (!element)
848 {
849 G4cerr << __METHOD_NAME__ << "No element named \"" << placement.referenceElement << "\" (instance #"
850 << placement.referenceElementNumber << ") in " << objectTypeForErrorMsg << " \""
851 << placement.name << "\"" << G4endl;
852 G4cout << "Note, this may be because the element is a bend and split into " << G4endl;
853 G4cout << "multiple sections with unique names. Run the visualiser to get " << G4endl;
854 G4cout << "the name of the segment, or place w.r.t. the element before / after." << G4endl;
855 throw BDSException("Error in "+objectTypeForErrorMsg+".");
856 }
857 // in this case we should use s for longitudinal offset - warn user if mistakenly using z
858 if (BDS::IsFinite(placement.z))
859 {
860 G4String message = objectTypeForErrorMsg + " \"" + placement.name + "\" is placed using";
861 message += " a referenceElement but the z offset is\n non zero. Note, s should be used to";
862 message += " offset the placement in this case and z will\n have no effect.";
863 BDS::Warning(message);
864 }
865 G4double sCoordinate = element->GetSPositionMiddle(); // start from middle of element
866 sCoordinate += placement.s * CLHEP::m; // add on (what's considered) 'local' s from the placement
867 if (S)
868 {*S = sCoordinate;}
869
870 G4ThreeVector offset = G4ThreeVector();
871 if (placementExtent)
872 {offset = SideToLocalOffset(placement, beamLine, *placementExtent);}
873
874 G4Transform3D beamlinePart = beamLine->GetGlobalEuclideanTransform(sCoordinate,
875 placement.x * CLHEP::m + offset.x(),
876 placement.y * CLHEP::m + offset.y());
877 G4Transform3D localRotation(rm, G4ThreeVector());
878 result = beamlinePart * localRotation;
879 }
880 else if (BDS::IsFinite(placement.s))
881 {// scenario 2
882 if (!beamLine)
883 {throw BDSException(__METHOD_NAME__, "no valid beam line yet placement w.r.t. a beam line.");}
884 G4ThreeVector offset = G4ThreeVector();
885 if (placementExtent)
886 {offset = SideToLocalOffset(placement, beamLine, *placementExtent);}
887
888 G4Transform3D beamlinePart = beamLine->GetGlobalEuclideanTransform(placement.s * CLHEP::m,
889 placement.x * CLHEP::m + offset.x(),
890 placement.y * CLHEP::m + offset.y());
891
892 G4Transform3D localRotation(rm, G4ThreeVector());
893 result = beamlinePart * localRotation;
894 if (S)
895 {*S = placement.s*CLHEP::m;}
896 }
897 else
898 {// scenario 1
899 G4ThreeVector translation = G4ThreeVector(placement.x*CLHEP::m,
900 placement.y*CLHEP::m,
901 placement.z*CLHEP::m);
902
903
904 result = G4Transform3D(rm, translation);
905 if (S)
906 {*S = -1000;} // default
907 }
908
909 return result;
910}
911
913 const BDSBeamline* beamLine,
914 G4double* S)
915{
916 // convert a scorermesh to a general placement for generation of the transform only.
917 GMAD::Placement convertedPlacement(scorerMesh);
918 return CreatePlacementTransform(convertedPlacement, beamLine, S, nullptr, "scorermesh");
919}
920
922 const BDSBeamline* beamLine,
923 G4double* S)
924{
925 // convert a sampler placement to a general placement for generation of the transform only.
926 GMAD::Placement convertedPlacement(samplerPlacement);
927 return CreatePlacementTransform(convertedPlacement, beamLine, S, nullptr, "samplerplacement");
928}
929
931 const BDSBeamline* beamLine,
932 G4double* S,
933 BDSExtent* blmExtent)
934{
935 // convert a sampler placement to a general placement for generation of the transform.
936 GMAD::Placement convertedPlacement(blmPlacement);
937 return CreatePlacementTransform(convertedPlacement, beamLine, S, blmExtent, "blm");
938}
939
941 const BDSBeamline* beamLine,
942 G4double* S)
943{
944 GMAD::Placement convertedPlacement(queryPlacement);
945 return CreatePlacementTransform(convertedPlacement, beamLine, S, nullptr, "query");
946}
947
949 const BDSBeamline* beamLine,
950 const BDSExtent& placementExtent)
951{
952 G4ThreeVector result = G4ThreeVector();
953 G4String side = G4String(placement.side);
954 side = BDS::LowerCase(side);
955
956 // Get the iterators pointing to the first and last elements
957 // that the placement lines up with.
958 G4double pathLength = placement.s*CLHEP::m;
959 std::pair<G4double, G4double> extentZ = placementExtent.ExtentZ();
960 G4double sLow = pathLength + extentZ.first;
961 G4double sHigh = pathLength + extentZ.second;
962 // iterator pointing to lower bound
963 auto start = beamLine->FindFromS(sLow);
964 auto end = beamLine->FindFromS(sHigh);
965 if (end != beamLine->end())
966 {end++;}
967
968 // Fold across the extents returning the greatest extent. The transverse extents
969 // will give be the transverse extents of the beamline section.
970 BDSExtent sectionMaxExtent = BDSExtent();
971 for (auto iter = start; iter != end; ++iter)
972 {sectionMaxExtent = BDS::MaximumCombinedExtent((*iter)->GetExtent(), sectionMaxExtent);}
973
974 // Multiplied by 5 because it works...
975 G4double ls = 5 * BDSGlobalConstants::Instance()->LengthSafetyLarge();
976
977 if (BDS::IsFinite(placement.sideOffset))
978 {ls = placement.sideOffset * CLHEP::m;}
979
980 if (side == "top")
981 {
982 result.setY(sectionMaxExtent.YPos() + placementExtent.YPos() + ls);
983 G4double xOffset = sectionMaxExtent.XPos() - 0.5*sectionMaxExtent.DX();
984 result.setX(xOffset);
985 }
986 else if (side == "bottom")
987 {
988 result.setY(sectionMaxExtent.YNeg() + placementExtent.YNeg() - ls);
989 G4double xOffset = sectionMaxExtent.XPos() - 0.5*sectionMaxExtent.DX();
990 result.setX(xOffset);
991 }
992 else if (side == "left")
993 {result.setX(sectionMaxExtent.XPos() + placementExtent.XPos() + ls);}
994 else if (side == "right")
995 {result.setX(sectionMaxExtent.XNeg() + placementExtent.XNeg() - ls);}
996 else if (side != "")
997 {throw BDSException(std::string("Unknown side in placement: " + side));}
998 return result;
999}
1000
1002{
1003 BDSExtent apertureExtent;
1004 if (sp.apertureModel.empty())
1005 {
1006 if (sp.samplerType == "plane")
1007 {
1008 BDSApertureInfo aperture = BDSApertureInfo(sp.shape,
1009 sp.aper1 * CLHEP::m,
1010 sp.aper2 * CLHEP::m,
1011 sp.aper3 * CLHEP::m,
1012 sp.aper4 * CLHEP::m,
1013 sp.name);
1014 apertureExtent = aperture.Extent();
1015 }
1016 else if (sp.samplerType == "cylinder" || sp.samplerType == "cylinderforward") // we ignore the possibility of only a part-cylinder
1017 {apertureExtent = BDSExtent(sp.aper1*CLHEP::m, sp.aper1*CLHEP::m, sp.aper2*CLHEP::m);}
1018 else if (sp.samplerType == "sphere" || sp.samplerType == "sphereforward") // we ignore the possibility of only a part-sphere
1019 {apertureExtent = BDSExtent(sp.aper1*CLHEP::m, sp.aper1*CLHEP::m, sp.aper1*CLHEP::m);}
1020 else
1021 {throw BDSException(__METHOD_NAME__, "unknown samplerType \"" + sp.samplerType + "\" in samplerplacement \"" + sp.name + "\"");}
1022 }
1023 else
1024 {
1025 BDSApertureInfo* aperture = BDSAcceleratorModel::Instance()->Aperture(sp.apertureModel);
1026 apertureExtent = aperture->Extent();
1027 }
1028
1029 // aperture is only transverse - fiddle z
1030 BDSExtent result = BDSExtent(apertureExtent.XNeg(), apertureExtent.XPos(),
1031 apertureExtent.YNeg(), apertureExtent.YPos(),
1032 1*CLHEP::um, 1*CLHEP::um);
1033 return result;
1034}
1035
1037{
1038 BDSExtentGlobal result;
1039 std::vector<GMAD::SamplerPlacement> samplerPlacements = BDSParser::Instance()->GetSamplerPlacements();
1040 for (const auto& samplerPlacement : samplerPlacements)
1041 {
1042 BDSExtent samplerExtent = CalculateExtentOfSamplerPlacement(samplerPlacement);
1043 G4Transform3D placementTransform = CreatePlacementTransform(samplerPlacement, beamLine);
1044 BDSExtentGlobal samplerExtentG = BDSExtentGlobal(samplerExtent, placementTransform);
1045 result = result.ExpandToEncompass(samplerExtentG);
1046 }
1047 return result;
1048}
1049
1051{
1052 BDSScorerMeshInfo recipe(sm);
1053 return recipe.Extent();
1054}
1055
1057{
1058 BDSExtentGlobal result;
1059 std::vector<GMAD::ScorerMesh> scorerMeshes = BDSParser::Instance()->GetScorerMesh();
1060 for (const auto& mesh : scorerMeshes)
1061 {
1062 BDSExtent meshExtent = CalculateExtentOfScorerMesh(mesh);
1063 G4Transform3D placementTransform = CreatePlacementTransform(mesh, beamLine);
1064 BDSExtentGlobal meshExtentG = BDSExtentGlobal(meshExtent, placementTransform);
1065 result = result.ExpandToEncompass(meshExtentG);
1066 }
1067 return result;
1068}
1069
1070#if G4VERSION_NUMBER > 1009
1072BDSDetectorConstruction::BuildCrossSectionBias(const std::list<std::string>& biasList,
1073 const std::list<std::string>& defaultBias,
1074 const G4String& elementName)
1075{
1076 // no accelerator components to bias
1077 if (biasList.empty() && defaultBias.empty())
1078 {return nullptr;}
1079
1080 std::list<std::string> biasesAll = biasList.empty() ? defaultBias : biasList;
1081
1082 // build a unique 'key' as the sorted set of bias names
1083 std::set<std::string> biasNamesSorted = {biasesAll.begin(), biasesAll.end()};
1084 G4String biasSetKey;
1085 G4String biasSetPrintOut;
1086 for (const auto& n : biasNamesSorted)
1087 {
1088 biasSetKey += n + "_";
1089 biasSetPrintOut += n + " ";
1090 }
1091 biasSetKey = BDS::StrStrip(biasSetKey, '_', BDS::StringStripType::trailing);
1092
1093 auto exists = biasSetObjects.find(biasSetKey);
1094 if (exists != biasSetObjects.end())
1095 {return exists->second;}
1096
1097 // loop over all physics biasing
1098 G4cout << "Bias> Creating unique set of bias objects ( " << biasSetPrintOut << ")" << G4endl;
1100
1101 const auto& biasObjectList = BDSParser::Instance()->GetBiasing();
1102 for (std::string const & bs : biasesAll)
1103 {
1104 auto result = biasObjectList.find(bs);
1105 if (result == biasObjectList.end())
1106 {throw BDSException("Error: bias named \"" + bs + "\" not found for element named \"" + elementName + "\"");}
1107 const GMAD::PhysicsBiasing& pb = *result;
1108
1109 if(debug)
1110 {G4cout << __METHOD_NAME__ << "bias loop : " << bs << " " << pb.particle << " " << pb.process << G4endl;}
1111
1112 eg->AddParticle(pb.particle);
1113
1114 // loop through all processes
1115 if (pb.flag.size() != pb.processList.size())
1116 {throw BDSException(__METHOD_NAME__, "number of flag entries in \"" + pb.name + "\" doesn't match number of processes");}
1117 if (pb.factor.size() != pb.processList.size())
1118 {throw BDSException(__METHOD_NAME__, "number of factor entries in \"" + pb.name + "\" doesn't match number of processes");}
1119 for (unsigned int p = 0; p < pb.processList.size(); ++p)
1120 {eg->SetBias(bs, pb.particle,pb.processList[p],pb.factor[p],(G4int)pb.flag[p]);}
1121 }
1122
1123 biasObjects.push_back(eg);
1124 biasSetObjects[biasSetKey] = eg; // cache it
1125 return eg;
1126}
1127#endif
1128
1130{
1131#if G4VERSION_NUMBER > 1009
1133 if (debug)
1134 {G4cout << __METHOD_NAME__ << "registry=" << registry << G4endl;}
1135
1136#if G4VERSION_NUMBER > 1039
1137 // only available in 4.10.4 onwards
1138 // crystal biasing necessary for implementation of variable density
1139 std::set<G4LogicalVolume*>* crystals = BDSAcceleratorModel::Instance()->VolumeSet("crystals");
1140 if (!crystals->empty())
1141 {
1142 G4cout << __METHOD_NAME__ << "Using crystal biasing: true" << G4endl; // to match print out style further down
1143 auto crystalBiasing = new G4ChannelingOptrMultiParticleChangeCrossSection();
1144 for (auto crystal : *crystals)
1145 {crystalBiasing->AttachTo(crystal);}
1146 }
1147#endif
1148
1149 auto blmSet = BDSBLMRegistry::Instance()->BLMs();
1150 for (auto blm : blmSet)
1151 {
1152 G4String biasNamesS = blm->Bias();
1153 if (biasNamesS.empty())
1154 {continue;}
1155 auto biasNamesV = BDS::SplitOnWhiteSpace(biasNamesS);
1156 std::list<std::string> biasNames = {biasNamesV.begin(), biasNamesV.end()};
1157 std::list<std::string> emptyDefaultBias;
1158 auto biasForBLM = BuildCrossSectionBias(biasNames, emptyDefaultBias, blm->GetName());
1159 for (auto lv : blm->GetAllLogicalVolumes())
1160 {biasForBLM->AttachTo(lv);}
1161 biasForBLM->AttachTo(blm->GetContainerLogicalVolume()); // in some cases it's just a single volume
1162 }
1163
1165 G4String defaultBiasVacuum = BDSParser::Instance()->GetOptions().defaultBiasVacuum;
1166 auto defaultBiasVacuumVector = BDS::SplitOnWhiteSpace(defaultBiasVacuum);
1167 auto defaultBiasVacuumList = std::list<std::string>(defaultBiasVacuumVector.begin(), defaultBiasVacuumVector.end());
1168 G4String defaultBiasMaterial = BDSParser::Instance()->GetOptions().defaultBiasMaterial;
1169 auto defaultBiasMaterialVector = BDS::SplitOnWhiteSpace(defaultBiasMaterial);
1170 auto defaultBiasMaterialList = std::list<std::string>(defaultBiasMaterialVector.begin(), defaultBiasMaterialVector.end());
1171 G4String biasForWorldVolume = g->BiasForWorldVolume();
1172 auto biasForWorldVolumeVector = BDS::SplitOnWhiteSpace(biasForWorldVolume);
1173 auto biasForWorldVolumeList = std::list<std::string>(biasForWorldVolumeVector.begin(), biasForWorldVolumeVector.end());
1174 G4String biasForWorldContents = g->BiasForWorldContents();
1175 auto biasForWorldContentsVector = BDS::SplitOnWhiteSpace(biasForWorldContents);
1176 auto biasForWorldContentsList = std::list<std::string>(biasForWorldContentsVector.begin(), biasForWorldContentsVector.end());
1177 G4String biasForWorldVacuum = g->BiasForWorldVacuum();
1178 auto biasForWorldVacuumVector = BDS::SplitOnWhiteSpace(biasForWorldVacuum);
1179 auto biasForWorldVacuumList = std::list<std::string>(biasForWorldVacuumVector.begin(), biasForWorldVacuumVector.end());
1180
1181 G4bool useDefaultBiasVacuum = !defaultBiasVacuum.empty();
1182 G4bool useDefaultBiasMaterial = !defaultBiasMaterial.empty();
1183 const auto& biasObjectList = BDSParser::Instance()->GetBiasing();
1184 G4bool useBiasForWorldVolume = !biasForWorldVolume.empty();
1185 G4bool useBiasForWorldContents = !biasForWorldContents.empty();
1186 G4bool useBiasForWorldVacuum = !biasForWorldVacuum.empty();
1187 G4bool biasesDefined = !biasObjectList.empty();
1188
1189 G4bool overallUseBiasing = useDefaultBiasVacuum || useDefaultBiasMaterial || biasesDefined || useBiasForWorldVolume || useBiasForWorldContents;
1190 G4cout << __METHOD_NAME__ << "Using generic biasing: " << BDS::BoolToString(overallUseBiasing) << G4endl;
1191 if (!overallUseBiasing)
1192 {return;} // no biasing used -> dont attach as just overhead for no reason
1193
1194 // apply per element biases
1195 std::map<G4String, BDSAcceleratorComponent*> allAcceleratorComponents = registry->AllComponentsIncludingUnique();
1196 for (auto const & item : allAcceleratorComponents)
1197 {
1198 if (debug)
1199 {G4cout << __METHOD_NAME__ << "checking component named: " << item.first << G4endl;}
1200 BDSAcceleratorComponent* accCom = item.second;
1201 BDSLine* l = dynamic_cast<BDSLine*>(accCom);
1202 if (l)
1203 {continue;} // do nothing for lines because each sub-component already has all definitions
1204 G4String accName = accCom->GetName();
1205 std::set<G4LogicalVolume*> vacuumLVs = accCom->GetAcceleratorVacuumLogicalVolumes();
1206
1207 // Build vacuum bias object based on vacuum bias list in the component
1208 auto vacuumBiasList = accCom->GetBiasVacuumList();
1209 if (!vacuumBiasList.empty() && vacuumLVs.empty())
1210 {BDS::Warning("biasVacuum set for component \"" + accName + "\" but there's no 'vacuum' volume for it and it can't be biased.\nRemove biasVacuum or name it with the namedVacuumVolumes parameter.");}
1211 if ((!vacuumBiasList.empty() || useDefaultBiasVacuum) && !vacuumLVs.empty())
1212 {
1213 auto egVacuum = BuildCrossSectionBias(vacuumBiasList, defaultBiasVacuumList, accName);
1214 for (auto lv : vacuumLVs)
1215 {
1216 if (debug)
1217 {G4cout << __METHOD_NAME__ << "vacuum volume name: " << lv << " " << lv->GetName() << G4endl;}
1218 egVacuum->AttachTo(lv);
1219 }
1220 }
1221
1222 // Build material bias object based on material bias list in the component
1223 auto materialBiasList = accCom->GetBiasMaterialList();
1224 if (!materialBiasList.empty() || useDefaultBiasMaterial)
1225 {
1226 auto egMaterial = BuildCrossSectionBias(materialBiasList, defaultBiasMaterialList, accName);
1227 auto allLVs = accCom->GetAcceleratorMaterialLogicalVolumes();
1228 if (debug)
1229 {G4cout << __METHOD_NAME__ << "# of logical volumes for biasing under 'material': " << allLVs.size() << G4endl;}
1230 for (auto lv : allLVs)
1231 {// BDSAcceleratorComponent automatically removes 'vacuum' volumes from all so we don't need to check
1232 if (debug)
1233 {G4cout << __METHOD_NAME__ << "Biasing 'material' logical volume: " << lv << " " << lv->GetName() << G4endl;}
1234 egMaterial->AttachTo(lv);
1235 }
1236 }
1237 }
1238
1239 if (useBiasForWorldContents)
1240 {
1241 std::list<std::string> emptyList;
1242 auto egWC = BuildCrossSectionBias(emptyList, biasForWorldContentsList, "world_contents_bias");
1243 for (auto lv : worldContentsLogicalVolumes)
1244 {egWC->AttachTo(lv);}
1245 }
1246 if (useBiasForWorldVolume)
1247 {
1248 std::list<std::string> emptyList;
1249 auto egWV = BuildCrossSectionBias(emptyList, biasForWorldVolumeList, "world_volume_bias");
1250 egWV->AttachTo(worldLogicalVolume);
1251 }
1252 if (useBiasForWorldVacuum)
1253 {
1254 std::list<std::string> emptyList;
1255 auto egWVac = BuildCrossSectionBias(emptyList, biasForWorldVacuumList, "world_vacuum_bias");
1256 for (auto lv : worldVacuumLogicalVolumes)
1257 {egWVac->AttachTo(lv);}
1258 }
1259#endif
1260}
1261
1263{
1264 auto flds = BDSFieldBuilder::Instance()->CreateAndAttachAll(); // avoid shadowing 'fields'
1266
1268}
1269
1270G4bool BDSDetectorConstruction::UnsuitableFirstElement(GMAD::FastList<GMAD::Element>::FastListConstIterator element)
1271{
1272 // skip past any line elements in parser to find first non-line element
1273 while ((*element).type == GMAD::ElementType::_LINE)
1274 {element++;}
1275
1276 if ((*element).type == GMAD::ElementType::_RBEND)
1277 {return true;} // unsuitable
1278 else if (BDS::IsFinite((*element).e1))
1279 {return true;} // unsuitable
1280 else
1281 {return false;} // suitable
1282}
1283
1285{
1286 // needed for filtering
1287 G4LogicalVolume* worldLV = acceleratorModel->WorldLV();
1288
1289 std::vector<GMAD::ScorerMesh> scoringMeshes = BDSParser::Instance()->GetScorerMesh();
1290 std::vector<GMAD::Scorer> scorers = BDSParser::Instance()->GetScorers();
1291
1292 if (scoringMeshes.empty())
1293 {return;}
1294
1295 G4ScoringManager* scManager = G4ScoringManager::GetScoringManager();
1296 scManager->SetVerboseLevel(1);
1297
1298 // convert all the parser scorer definitions into recipes (including parameter checking)
1299 std::map<G4String, BDSScorerInfo> scorerRecipes;
1300 for (const auto& scorer : scorers)
1301 {
1302 BDSScorerInfo si = BDSScorerInfo(scorer, true); // true = upgrade to 3d as required for 3d mesh here
1303 scorerRecipes.insert(std::make_pair(si.name, si));
1304 }
1305
1306 // construct meshes
1307 BDSScorerFactory scorerFactory;
1308 for (const auto& mesh : scoringMeshes)
1309 {
1310 // convert to recipe class as this checks parameters
1311 BDSScorerMeshInfo meshRecipe = BDSScorerMeshInfo(mesh);
1312
1313 // name we'll use for the mesh
1314 G4String meshName = meshRecipe.name;
1315
1316 // TBC - could be any beam line in future - just w.r.t. main beam line just now
1317 const BDSBeamline* mbl = BDSAcceleratorModel::Instance()->BeamlineMain();
1318 G4Transform3D placement = CreatePlacementTransform(mesh, mbl);
1319
1320 BDSScoringMeshBox* scorerBox = nullptr;
1321 BDSScoringMeshCylinder* scorerCylindrical = nullptr;
1322 const BDSHistBinMapper* mapper = nullptr;
1323
1324 G4String geometryType = BDS::LowerCase(G4String(mesh.geometryType));
1325
1326 if (geometryType == "box")
1327 {// create a scoring box
1328 scorerBox = new BDSScoringMeshBox(meshName, meshRecipe, placement);
1329 mapper = scorerBox->Mapper();
1330 }
1331 else if (geometryType == "cylindrical")
1332 {// create a scoring cylinder
1333 scorerCylindrical = new BDSScoringMeshCylinder(meshName, meshRecipe, placement);
1334 mapper = scorerCylindrical->Mapper();
1335 }
1336 else
1337 {
1338 G4String msg = "mesh geometry type \"" + geometryType + "\" is not correct. The possible options are \"box\" and \"cylindrical\"";
1339 throw BDSException(__METHOD_NAME__, msg);
1340 }
1341
1342 // add the scorer(s) to the scoring mesh
1343 std::vector<G4String> meshPrimitiveScorerNames; // final vector of unique mesh + ps names
1344 std::vector<G4double> meshPrimitiveScorerUnits;
1345 std::vector<G4String> scorerNames;
1346
1347 std::vector<G4String> words = BDS::SplitOnWhiteSpace(mesh.scoreQuantity);
1348 for (const auto& word : words)
1349 {
1350 auto search = scorerRecipes.find(word);
1351 if (search == scorerRecipes.end())
1352 {throw BDSException(__METHOD_NAME__, "scorerQuantity \"" + word + "\" for mesh \"" + meshName + "\" not found.");}
1353
1354 G4double psUnit = 1.0;
1355 G4VPrimitiveScorer* ps = scorerFactory.CreateScorer(&(search->second), mapper, &psUnit, worldLV);
1356 // The mesh internally creates a multifunctional detector which is an SD and has
1357 // the name of the mesh. Any primitive scorer attached is added to the mfd. To get
1358 // the hits map we need the full name of the unique primitive scorer so we build that
1359 // name here and store it.
1360 G4String uniqueName = meshName + "/" + ps->GetName();
1361 meshPrimitiveScorerNames.push_back(uniqueName);
1362 meshPrimitiveScorerUnits.push_back(psUnit);
1363
1364 // sets the current ps but appends to list of multiple
1365 if (geometryType == "box")
1366 {scorerBox->SetPrimitiveScorer(ps);}
1367 else if (geometryType == "cylindrical")
1368 {scorerCylindrical->SetPrimitiveScorer(ps);}
1369
1370 BDSScorerHistogramDef outputHistogram(meshRecipe, uniqueName, ps->GetName(), psUnit, *mapper);
1371 BDSAcceleratorModel::Instance()->RegisterScorerHistogramDefinition(outputHistogram);
1372 BDSAcceleratorModel::Instance()->RegisterScorerPlacement(meshName, placement);
1373 }
1374
1375 if (geometryType == "box")
1376 {scManager->RegisterScoringMesh(scorerBox);} // sets the current ps but appends to list of multiple
1377 else if (geometryType == "cylindrical")
1378 {scManager->RegisterScoringMesh(scorerCylindrical);}// sets the current ps but appends to list of multiple
1379
1380 // register it with the sd manager as this is where we get all collection IDs from
1381 // in the end of event action. This must come from the mesh as it creates the
1382 // multifunctionaldetector and therefore has the complete name of the scorer collection
1383 BDSSDManager::Instance()->RegisterPrimitiveScorerNames(meshPrimitiveScorerNames, &meshPrimitiveScorerUnits);
1384 }
1385}
1386
1387std::vector<BDSFieldQueryInfo*> BDSDetectorConstruction::PrepareFieldQueries(const BDSBeamline* mainBeamline)
1388{
1389 std::vector<BDSFieldQueryInfo*> result;
1390 const std::vector<GMAD::Query>& parserQueries = BDSParser::Instance()->GetQuery();
1391 for (const auto& def : parserQueries)
1392 {
1393 if (!def.queryMagneticField && !def.queryElectricField)
1394 {throw BDSException(__METHOD_NAME__, "neither \"queryMagneticField\" nor \"queryElectricField\" are true (=1) - one must be turned on.");}
1395
1396 if (!def.pointsFile.empty())
1397 {
1398 std::vector<G4String> columnNames;
1399 auto points = BDS::LoadFieldQueryPoints(G4String(def.pointsFile), &columnNames);
1400 result.emplace_back(new BDSFieldQueryInfo(G4String(def.name),
1401 G4String(def.outfileMagnetic),
1402 G4String(def.outfileElectric),
1403 G4bool(def.queryMagneticField),
1404 G4bool(def.queryElectricField),
1405 points,
1406 columnNames,
1407 G4bool(def.overwriteExistingFiles),
1408 G4String(def.fieldObject),
1409 def.checkParameters,
1410 def.drawArrows,
1411 def.drawZeroValuePoints,
1412 def.drawBoxes,
1413 def.boxAlpha));
1414 }
1415 else
1416 {
1417 G4Transform3D globalTransform3D = CreatePlacementTransform(def, mainBeamline);
1418 auto rot = globalTransform3D.getRotation();
1419 rot = rot.inverse();
1420 G4AffineTransform globalTransform(rot, globalTransform3D.getTranslation());
1421
1422 result.emplace_back(new BDSFieldQueryInfo(G4String(def.name),
1423 G4String(def.outfileMagnetic),
1424 G4String(def.outfileElectric),
1425 G4bool(def.queryMagneticField),
1426 G4bool(def.queryElectricField),
1427 {def.nx, def.xmin*CLHEP::m, def.xmax*CLHEP::m},
1428 {def.ny, def.ymin*CLHEP::m, def.ymax*CLHEP::m},
1429 {def.nz, def.zmin*CLHEP::m, def.zmax*CLHEP::m},
1430 {def.nt, def.tmin*CLHEP::ns, def.tmax*CLHEP::ns},
1431 globalTransform,
1432 G4bool(def.overwriteExistingFiles),
1433 G4String(def.fieldObject),
1434 G4bool(def.printTransform),
1435 def.checkParameters,
1436 def.drawArrows,
1437 def.drawZeroValuePoints,
1438 def.drawBoxes,
1439 def.boxAlpha));
1440 }
1441 }
1442 return result;
1443}
1444
1445#ifdef BDSCHECKUSERLIMITS
1446void BDSDetectorConstruction::PrintUserLimitsSummary(const G4VPhysicalVolume* world) const
1447{
1448 G4cout << "USERLIMITS START" << G4endl;
1449 G4double globalMinEK = BDSGlobalConstants::Instance()->MinimumKineticEnergy();
1450 PrintUserLimitsPV(world, globalMinEK);
1451 G4cout << "USERLIMITS END" << G4endl;
1452}
1453
1454void BDSDetectorConstruction::PrintUserLimitsPV(const G4VPhysicalVolume* aPV, G4double globalMinEK) const
1455{
1456 const auto lv = aPV->GetLogicalVolume();
1457 const auto ul = lv->GetUserLimits();
1458 if (ul)
1459 {
1460 G4Track dummyTrack;
1461 G4double ekUL = ul->GetUserMinEkine(dummyTrack);
1462 //G4cout << lv->GetName() << " Ek Min: " << ekUL << G4endl;
1463 if (ekUL < globalMinEK)
1464 {G4cout << lv->GetName() << " Ek Min: " << ekUL << " MeV < global: " << globalMinEK << " MeV" << G4endl;}
1465 }
1466 else
1467 {G4cout << lv->GetName() << " no G4UserLimits" << G4endl;}
1468 for (G4int i = 0; i < (G4int)lv->GetNoDaughters(); i++)
1469 {PrintUserLimitsPV(lv->GetDaughter(i), globalMinEK);}
1470}
1471#endif
A registry of constructed BDSAcceleratorComponent instances that can be searched.
static BDSAcceleratorComponentRegistry * Instance()
Singleton accessor.
std::map< G4String, BDSAcceleratorComponent * > AllComponentsIncludingUnique() const
void PrintNumberOfEachType() const
Print out the number of each type of component registered.
Abstract class that represents a component of an accelerator.
virtual G4String GetName() const
The name of the component without modification.
virtual std::set< G4LogicalVolume * > GetAcceleratorVacuumLogicalVolumes() const
Access the 'vacuum' volume(s) in this component if any. Default is empty set.
virtual std::set< G4LogicalVolume * > GetAcceleratorMaterialLogicalVolumes() const
Return a set of logical volumes excluding the ones in the 'vacuum' set.
std::list< std::string > GetBiasVacuumList() const
Access the bias list copied from parser.
G4String GetType() const
Get a string describing the type of the component.
std::list< std::string > GetBiasMaterialList() const
Access the bias list copied from parser.
void RegisterFields(std::vector< BDSFieldObjects * > &fieldsIn)
Register all field objects.
void RegisterBeamlineSetExtra(const G4String &name, const BDSBeamlineSet &setIn)
Register a set of beam lines to be managed and cleared up at the end of the simulation.
void RegisterWorldLV(G4LogicalVolume *worldIn)
Register constituent of world.
std::set< G4LogicalVolume * > * VolumeSet(const G4String &name)
Returns pointer to a set of logical volumes. If no set by that name exits, create it.
BDSApertureInfo * Aperture(G4String name) const
const std::map< G4String, BDSBeamlineSet > & ExtraBeamlines() const
Accessor.
BDSBeamline * TunnelBeamline() const
Access the beam line containing all the tunnel segments.
void RegisterApertures(const std::map< G4String, BDSApertureInfo * > &aperturesIn)
Register a map of apertures with associated names.
void RegisterScorerPlacement(const G4String &meshName, const G4Transform3D &placement)
Register a copy of the scorer placement so it can be used in the output.
const BDSBeamlineSet & BeamlineSetMain() const
Accessor.
const BDSBeamline * BeamlineMain() const
Accessor.
void RegisterTunnelBeamline(BDSBeamline *beamlineIn)
Register the beam line containing all the tunnel segments.
BDSBeamline * PlacementBeamline() const
Access the beam line of arbitrary placements.
void RegisterWorldPV(G4VPhysicalVolume *worldIn)
Register constituent of world.
void RegisterPlacementBeamline(BDSBeamline *placementBLIn)
Register the beam line of arbitrary placements.
void RegisterBeamlineSetMain(const BDSBeamlineSet &setIn)
Register the main beam line set.
void RegisterScorerHistogramDefinition(const BDSScorerHistogramDef &def)
G4LogicalVolume * WorldLV() const
Access the logical volume of the world.
void RegisterBLMs(BDSBeamline *blmsIn)
Register a beam line of blm placements.
G4Region * Region(const G4String &name) const
Access region information. Will exit if not found.
void RegisterWorldSolid(G4VSolid *worldIn)
Register constituent of world.
void RegisterRegion(BDSRegion *region)
Register a region.
Holder class for all information required to describe an aperture.
BDSExtent Extent() const
static void AttachWorldVolumeToNavigator(G4VPhysicalVolume *worldPVIn)
Setup the navigator w.r.t. to a world volume - typically real world.
static BDSBLMRegistry * Instance()
Accessor for registry.
A class that holds a fully constructed BDSAcceleratorComponent as well as any information relevant to...
G4double GetSPositionEnd() const
Accessor.
G4double GetSPositionMiddle() const
Accessor.
Simple struct to return a beamline plus associated beam lines.
void GetExtentGlobals(std::vector< BDSExtentGlobal > &extents) const
Append global extents of all beamlines to supplied vector.
A vector of BDSBeamlineElement instances - a beamline.
Definition: BDSBeamline.hh:61
const BDSBeamlineElement * GetElement(G4String acceleratorComponentName, G4int i=0) const
Get the ith placement of an element in the beam line. Returns null pointer if not found.
Definition: BDSBeamline.cc:736
BDSBeamlineElement * back() const
Return a reference to the last element.
Definition: BDSBeamline.hh:211
G4bool ElementAnglesSumToCircle() const
Check if the sum of the angle of all elements is two pi.
Definition: BDSBeamline.cc:919
reverse_iterator rend()
Iterator mechanics.
Definition: BDSBeamline.hh:172
const_iterator FindFromS(G4double S) const
Returns an iterator to the beamline element at s.
Definition: BDSBeamline.cc:678
G4bool empty() const
Iterator mechanics.
Definition: BDSBeamline.hh:175
G4Transform3D GetGlobalEuclideanTransform(G4double s, G4double x=0, G4double y=0, G4int *indexOfFoundElement=nullptr) const
Definition: BDSBeamline.cc:598
BDSExtentGlobal GetExtentGlobal() const
Get the global extents for this beamline.
Definition: BDSBeamline.cc:924
reverse_iterator rbegin()
Iterator mechanics.
Definition: BDSBeamline.hh:171
iterator begin()
Iterator mechanics.
Definition: BDSBeamline.hh:167
void AddComponent(BDSAcceleratorComponent *component, BDSTiltOffset *tiltOffset=nullptr, BDSSamplerInfo *samplerInfo=nullptr)
Definition: BDSBeamline.cc:122
iterator end()
Iterator mechanics.
Definition: BDSBeamline.hh:168
G4double GetTotalArcLength() const
Get the total ARC length for the beamline - ie the maximum s position.
Definition: BDSBeamline.hh:139
Factory for user specified accelerator components.
void SetDesignParticle(const BDSParticleDefinition *designParticleIn)
Update values for nominal rigidity and relativisitic beta of the beam particle.
Factory to produce all types of BDSAcceleratorComponents.
BDSAcceleratorComponent * CreateTeleporter(const G4double teleporterLength, const G4double teleporterWidth, const G4Transform3D &transformIn)
BDSAcceleratorComponent * CreateTerminator(G4double witdth)
BDSAcceleratorComponent * CreateComponent(GMAD::Element const *elementIn, GMAD::Element const *prevElementIn, GMAD::Element const *nextElementIn, G4double currentArcLength=0)
static BDSTiltOffset * CreateTiltOffset(GMAD::Element const *el)
Create the tilt and offset information object by inspecting the parser element.
Factory for simple parallel geometry for curvilinear coordinates.
BDSBeamline * BuildCurvilinearBeamLine1To1(BDSBeamline const *const beamline, const G4bool circular)
Build a beam line of curvilinear geometry based on another beam line.
BDSBeamline * BuildCurvilinearBridgeBeamLine(BDSBeamline const *const beamline)
Build bridging volumes to join the curvilinear ones.
G4bool circular
Whether or not we're building a circular machine.
BDSBeamlineSet BuildBeamline(const GMAD::FastList< GMAD::Element > &beamLine, const G4String &name, const G4Transform3D &initialTransform=G4Transform3D(), G4double initialS=0.0, G4bool beamlineIsCircular=false, G4bool isPlacementBeamline=false)
BDSBOptrMultiParticleChangeCrossSection * BuildCrossSectionBias(const std::list< std::string > &biasList, const std::list< std::string > &defaultBias, const G4String &elementName)
Function that creates physics biasing cross section.
void InitialiseRegions()
Create and set parameters for various G4Regions.
std::vector< BDSBOptrMultiParticleChangeCrossSection * > biasObjects
List of bias objects - for memory management.
BDSExtent CalculateExtentOfScorerMesh(const GMAD::ScorerMesh &sm) const
Calculate local extent of scorer mesh in 3D.
static void PlaceBeamlineInWorld(BDSBeamline *beamline, G4VPhysicalVolume *containerPV, G4bool checkOverlaps=false, G4bool setRegions=false, G4bool registerInfo=false, G4bool useCLPlacementTransform=false, G4bool useIncrementalCopyNumbers=false)
BDSAcceleratorModel * acceleratorModel
Accelerator model pointer.
static G4ThreeVector SideToLocalOffset(const GMAD::Placement &placement, const BDSBeamline *beamLine, const BDSExtent &placementExtent)
Attach component with extent2 to component with extent1 with placement.
G4bool verbose
Variable copied from global constants.
static std::vector< BDSFieldQueryInfo * > PrepareFieldQueries(const BDSBeamline *mainBeamline)
Prepare field queries from parser information.
G4VPhysicalVolume * BuildWorld()
void InitialiseApertures()
Create all aperture definitions from parser and store in BDSAcceleratorModel.
void CountPlacementFields()
Count number of fields required for placements.
G4bool canSampleAngledFaces
Whether the integrator set permits sampling elements with angled faces.
BDSExtent worldExtent
Record of the world extent.
std::set< G4LogicalVolume * > worldContentsLogicalVolumes
Cache of possibly loaded logical volumes from a world geometry file - used for biasing.
static BDSSamplerInfo * BuildSamplerInfo(const GMAD::Element *element)
BDSExtent CalculateExtentOfSamplerPlacement(const GMAD::SamplerPlacement &sp) const
Calculate local extent of custom user sampler.
void ConstructScoringMeshes()
Construct scoring meshes.
BDSExtentGlobal CalculateExtentOfScorerMeshes(const BDSBeamline *bl) const
G4bool checkOverlaps
Variable copied from global constants.
G4int nSamplers
Count of number of samplers to be built.
void ComponentPlacement(G4VPhysicalVolume *worldPV)
Place beam line, tunnel beam line, end pieces and placements in world.
void BuildBeamlines()
Build the main beam line and then any other required beam lines.
const BDSParticleDefinition * designParticle
Particle definition all components are built w.r.t. Includes rigidity etc.
void BuildTunnel()
Build the tunnel around the already constructed flat beam line.
virtual void ConstructSDandField()
Construct sensitive detectors and fields.
virtual G4VPhysicalVolume * Construct()
static G4Transform3D CreatePlacementTransform(const GMAD::Placement &placement, const BDSBeamline *beamLine, G4double *S=nullptr, BDSExtent *placementExtent=nullptr, const G4String &objectTypeForErrorMsg="placement")
G4bool UnsuitableFirstElement(std::list< GMAD::Element >::const_iterator element)
BDSExtentGlobal CalculateExtentOfSamplerPlacements(const BDSBeamline *beamLine) const
General exception with possible name of object and message.
Definition: BDSException.hh:35
Holder for +- extents in 3 dimensions with a rotation and translation.
BDSExtentGlobal ExpandToEncompass(const BDSExtentGlobal &other) const
Return a copy of this extent but expanded to encompass another global extent.
G4bool Encompasses(const BDSExtentGlobal &otherExtent)
Return whether the extent encompasses another extent.
Holder for +- extents in 3 dimensions.
Definition: BDSExtent.hh:39
G4double XPos() const
Accessor.
Definition: BDSExtent.hh:66
G4double XNeg() const
Accessor.
Definition: BDSExtent.hh:65
G4double YNeg() const
Accessor.
Definition: BDSExtent.hh:67
std::pair< G4double, G4double > ExtentZ() const
Accessor.
Definition: BDSExtent.hh:63
G4double DX() const
The difference in a dimension.
Definition: BDSExtent.hh:83
G4double YPos() const
Accessor.
Definition: BDSExtent.hh:68
static BDSFieldBuilder * Instance()
Singleton pattern accessor.
Holder class for all information required for a field query.
static void AttachWorldVolumeToNavigator(G4VPhysicalVolume *worldPVIn)
Setup the navigator w.r.t. to a world volume.
Gap in accelerator beamline.
Definition: BDSGap.hh:32
G4LogicalVolume * GetContainerLogicalVolume() const
Accessor - see member for more info.
BDSExtent GetExtent() const
Accessor - see member for more info.
virtual void AttachSensitiveDetectors()
Attach a sensitive detector class to all registered sensitive volumes in this component.
virtual std::set< G4LogicalVolume * > GetAllLogicalVolumes() const
Access all logical volumes belonging to this component.
G4VSolid * GetContainerSolid() const
Accessor - see member for more info.
A loaded piece of externally provided geometry.
const std::set< G4LogicalVolume * > & VacuumVolumes() const
Access the vacuum volumes.
static BDSGeometryFactory * Instance()
Singleton accessor.
BDSGeometryExternal * BuildGeometry(G4String componentName, const G4String &formatAndFilePath, std::map< G4String, G4Colour * > *colourMapping=nullptr, G4bool autoColour=true, G4double suggestedLength=0, G4double suggestedHorizontalWidth=0, std::vector< G4String > *namedVacuumVolumes=nullptr, G4bool makeSensitive=true, BDSSDType sensitivityType=BDSSDType::energydep, G4bool stripOuterVolumeAndMakeAssembly=false, G4UserLimits *userLimitsToAttachToAllLVs=nullptr, G4bool dontReloadGeometry=false)
A class that holds global options and constants.
void SetCurvilinearDiameterShrunkForBends()
Setter.
static BDSGlobalConstants * Instance()
Access method.
void SetCurvilinearDiameter(G4double curvilinearDiameterIn)
Setter.
void SetSamplerDiameter(G4double samplerDiameterIn)
Setter.
Mapping from axis indices to 1D index.
A class that hold multiple accelerator components.
Definition: BDSLine.hh:38
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
G4Material * GetMaterial(G4String material) const
Get material by name.
const GMAD::FastList< GMAD::Element > & GetSequence(const std::string &name)
Return sequence.
Definition: BDSParser.hh:81
std::vector< GMAD::Scorer > GetScorers() const
Return scorer list.
Definition: BDSParser.hh:126
static BDSParser * Instance()
Access method.
Definition: BDSParser.cc:28
const GMAD::Options & GetOptions() const
Return options.
Definition: BDSParser.hh:50
std::vector< GMAD::Query > GetQuery() const
Query list.
Definition: BDSParser.hh:117
std::vector< GMAD::SamplerPlacement > GetSamplerPlacements() const
Return sampler placement list.
Definition: BDSParser.hh:123
const GMAD::FastList< GMAD::PhysicsBiasing > & GetBiasing() const
Return biasing list.
Definition: BDSParser.hh:95
std::vector< GMAD::ScorerMesh > GetScorerMesh() const
Return scorermesh list.
Definition: BDSParser.hh:129
std::vector< GMAD::Placement > GetPlacements() const
Return the vector of placement objects.
Definition: BDSParser.hh:114
void RegisterExcludedPV(G4VPhysicalVolume *physicalVolume)
void RegisterInfo(G4VPhysicalVolume *physicalVolume, BDSPhysicalVolumeInfo *info, G4bool isReadOutVolume=false, G4bool isTunnel=false)
static BDSPhysicalVolumeInfoRegistry * Instance()
Singleton accessor.
A class holding any information pertaining to a particular physical volume in a BDSIM geant4 model.
Range cuts for a region. Help with defaults.
Definition: BDSRegion.hh:41
void ConstructSamplerSDsForParticleSets(const std::map< int, std::set< int > > &samplerFilterIDtoPDGSetIn)
void RegisterPrimitiveScorerNames(const std::vector< G4String > &namesIn, const std::vector< G4double > *units=nullptr)
Loop over a vector and apply above single name function.
All info required to build a sampler but not place it.
Create primitive scorers on demand.
G4VPrimitiveScorer * CreateScorer(const BDSScorerInfo *info, const BDSHistBinMapper *mapper, G4double *unit=nullptr, G4LogicalVolume *worldLV=nullptr)
Main function to create a scorer.
Definition for a scorer histogram.
Recipe class for scorer. Checks values.
G4String name
Scorer name.
Recipe class for a scoring mesh.
Wrapper for G4ScoringBox to allow full access to placement.
Wrapper for G4ScoringCylinder to allow full access to placement.
A class of functions to output Geant4/Mokka/BDSIM parameters for the beamline.
Definition: BDSSurvey.hh:39
void Write(BDSBeamlineElement *beamlineElement)
write line
Definition: BDSSurvey.cc:109
A holder for any placement offsets and rotations for a BDSAcceleratorComponent.
static G4double dEThresholdForScattering
A class that constructs tunnel segments around a beamline.
BDSBeamline * BuildTunnelSections(const BDSBeamline *flatBeamLine) const
Aperture class.
Definition: aperture.h:38
blm for parser.
Definition: blmplacement.h:39
FastListIterator begin()
Definition: fastlist.h:215
FastListIterator end()
Definition: fastlist.h:220
bool empty() const
Whether the list is empty.
Definition: fastlist.h:68
FastListConstIterator find(std::string name, unsigned int count=1) const
Definition: fastlist.h:245
int size() const
size of list
Definition: fastlist.h:164
const FastList< Element > & GetBeamline() const
Definition: parser.cc:896
Physics biasing class for parser.
std::vector< double > factor
factors corresponding to process
std::vector< PhysicsBiasingType > flag
flag which particles are biased
std::string particle
particle name
std::string process
geant4 process: single string, but can have multiple processes separated with a space
Placement class for parser.
Definition: placement.h:41
bool axisAngle
Flag to use the axis angle construction of rotation.
Definition: placement.h:65
double axisZ
Definition: placement.h:61
std::string name
Name of this placement.
Definition: placement.h:43
std::string side
which side to attach to: top, bottom, left, right.
Definition: placement.h:66
double theta
Euler angle for rotation.
Definition: placement.h:55
double sideOffset
Gap between side and component.
Definition: placement.h:67
double psi
Euler angle for rotation.
Definition: placement.h:56
double y
Offset in y.
Definition: placement.h:51
double s
Curvilinear s position to place w.r.t..
Definition: placement.h:49
double axisY
Definition: placement.h:60
double angle
Definition: placement.h:62
double x
Offset in x.
Definition: placement.h:50
int referenceElementNumber
Index of repetition of element if there are multiple uses.
Definition: placement.h:48
double phi
Euler angle for rotation.
Definition: placement.h:54
std::string referenceElement
Name of reference element w.r.t. to place to.
Definition: placement.h:47
double axisX
Definition: placement.h:59
Query structure class for parser.
Definition: query.h:37
Region class for parser.
Definition: region.h:36
Sampler placement class for parser.
std::string name
Name of this samplerplacement.
std::string samplerType
Plane, Cylinder, Sphere.
ScorerMesh class for parser.
Definition: scorermesh.h:38
BDSBeamline * BuildEndPieceBeamline(const BDSBeamline *beamline, const G4bool circularMachine)
G4Transform3D CalculateTeleporterDelta(const BDSBeamline *beamline, G4double &teleporterLength)
std::vector< BDSFourVector< G4double > > LoadFieldQueryPoints(const G4String &fileName, std::vector< G4String > *columnNamesIn)
BDSExtent MaximumCombinedExtent(const BDSExtent &first, const BDSExtent &second)
Returns the extent which is the greatest extent in all six directions.
Definition: BDSExtent.cc:248
BDSBeamline * BuildPlacementGeometry(const std::vector< GMAD::Placement > &placements, const BDSBeamline *parentBeamLine, BDSComponentFactory *componentFactory)
G4String LowerCase(const G4String &str)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
std::vector< G4String > SplitOnWhiteSpace(const G4String &input)
Split a string on whitespace and return a vector of these 'words'.
BDSBeamline * BuildBLMs(const std::vector< GMAD::BLMPlacement > &blmPlacements, const BDSBeamline *parentBeamLine)
G4bool IsFinite(G4double value, G4double tolerance=std::numeric_limits< double >::epsilon())
G4String StrStrip(const G4String &str, char ch, StringStripType stripType=StringStripType::both)
Utility function to simplify lots of syntax changes for pedantic g4 changes.
Definition: BDSUtilities.cc:98
Element class.
Definition: element.h:43
std::string samplerType
element has a sampler of this type (default "none")
Definition: element.h:200
int samplerParticleSetID
Definition: element.h:205
double e1
input pole face rotation for bends
Definition: element.h:60
std::string samplerName
name of sampler (default empty)
Definition: element.h:199