BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSFieldBuilder.cc
1/*
2Beam Delivery Simulation (BDSIM) Copyright (C) Royal Holloway,
3University of London 2001 - 2022.
4
5This file is part of BDSIM.
6
7BDSIM is free software: you can redistribute it and/or modify
8it under the terms of the GNU General Public License as published
9by the Free Software Foundation version 3 of the License.
10
11BDSIM is distributed in the hope that it will be useful, but
12WITHOUT ANY WARRANTY; without even the implied warranty of
13MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14GNU General Public License for more details.
15
16You should have received a copy of the GNU General Public License
17along with BDSIM. If not, see <http://www.gnu.org/licenses/>.
18*/
19#include "BDSDebug.hh"
20#include "BDSFieldBuilder.hh"
21#include "BDSFieldFactory.hh"
22#include "BDSFieldInfo.hh"
23#include "BDSFieldObjects.hh"
24
25#include "G4LogicalVolume.hh"
26
27#include <set>
28#include <vector>
29
31
33{
34 if (!instance)
35 {instance = new BDSFieldBuilder();}
36 return instance;
37}
38
39BDSFieldBuilder::~BDSFieldBuilder()
40{
41 instance = nullptr;
42}
43
45{
46 size_t defaultSize = 30;
47 infos.reserve(defaultSize);
48 lvs.reserve(defaultSize);
49 propagators.reserve(defaultSize);
50}
51
53 const std::vector<G4LogicalVolume*>& logicalVolumes,
54 const G4bool propagateToDaughters,
55 const BDSMagnetStrength* magnetStrengthForScaling,
56 const G4String& scalingKey)
57{
58 if (info)
59 {
60#ifdef BDSDEBUG
61 G4cout << __METHOD_NAME__ << "Registering info: " << info
62 << " to volume(s): ";
63 for (auto vol : logicalVolumes)
64 {G4cout << vol->GetName() << " ";}
65 G4cout << G4endl;
66#endif
67 infos.push_back(info);
68 lvs.push_back(logicalVolumes);
69 propagators.push_back(propagateToDaughters);
70 if (info->AutoScale())
71 {// only store if we're going to use for autoscaling
72 G4int index = (G4int)infos.size() - 1;
73 scalingStrengths[index] = magnetStrengthForScaling;
74 scalingKeys[index] = scalingKey;
75 }
76 }
77}
78
80 G4LogicalVolume* logicalVolume,
81 const G4bool propagateToDaughters,
82 const BDSMagnetStrength* magnetStrengthForScaling,
83 const G4String& scalingKey)
84{
85 std::vector<G4LogicalVolume*> lvsForThisInfo = {logicalVolume};
87 lvsForThisInfo,
88 propagateToDaughters,
89 magnetStrengthForScaling,
90 scalingKey);
91}
92
94 const std::set<G4LogicalVolume*>& logicalVolumes,
95 const G4bool propagateToDaughters,
96 const BDSMagnetStrength* magnetStrengthForScaling,
97 const G4String& scalingKey)
98{
99 // copy into vector for this interface
100 std::vector<G4LogicalVolume*> lvsForThisInfo;
101 for (auto lv : logicalVolumes)
102 {lvsForThisInfo.push_back(lv);}
104 lvsForThisInfo,
105 propagateToDaughters,
106 magnetStrengthForScaling,
107 scalingKey);
108}
109
110std::vector<BDSFieldObjects*> BDSFieldBuilder::CreateAndAttachAll()
111{
112 std::vector<BDSFieldObjects*> fields;
113 fields.reserve(infos.size());
114 for (G4int i = 0; i < (G4int)infos.size(); i++)
115 {
116 BDSFieldObjects* field = nullptr;
117 const BDSFieldInfo* currentInf = infos[i];
118 if (currentInf->AutoScale())
119 {
122 scalingKeys[i]);
123 }
124 else
125 {field = BDSFieldFactory::Instance()->CreateField(*(infos[i]));}
126 if (field)
127 {
128 fields.push_back(field);
129 field->AttachToVolume(lvs[i], propagators[i]); // works with vector of LVs*
130 }
131 }
132 return fields;
133}
Register for all fields to be built and volumes to be attached to.
std::map< G4int, G4String > scalingKeys
Optional register of scaling strengths and keys.
std::vector< const BDSFieldInfo * > infos
Register of components to build.
std::map< G4int, const BDSMagnetStrength * > scalingStrengths
Optional register of scaling strengths and keys.
std::vector< G4bool > propagators
Register of components to build.
void RegisterFieldForConstruction(const BDSFieldInfo *info, G4LogicalVolume *logicalVolume, const G4bool propagateToDaughters=false, const BDSMagnetStrength *magnetStrengthForScaling=nullptr, const G4String &scalingKey="none")
std::vector< std::vector< G4LogicalVolume * > > lvs
Register of components to build.
static BDSFieldBuilder * instance
Singleton instance.
BDSFieldBuilder()
Private default constructor to enforce singleton pattern.
static BDSFieldBuilder * Instance()
Singleton pattern accessor.
static BDSFieldFactory * Instance()
Public accessor method for singleton pattern.
BDSFieldObjects * CreateField(const BDSFieldInfo &info, const BDSMagnetStrength *scalingStrength=nullptr, const G4String &scalingKey="none")
Main interface to field factory.
All info required to build complete field of any type.
Definition: BDSFieldInfo.hh:65
G4bool AutoScale() const
Accessor.
A holder for all the Geant4 field related objects.
void AttachToVolume(G4LogicalVolume *volume, G4bool penetrateToDaughterVolumes=true) const
Interface to easily attach to logical volume.
Efficient storage of magnet strengths.