BDSIM
BDSIM is a Geant4 extension toolkit for simulation of particle transport in accelerator beamlines.
BDSGeometryFactoryGDML.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#ifdef USE_GDML
20#include "BDSAcceleratorModel.hh"
21#include "BDSColourFromMaterial.hh"
22#include "BDSColours.hh"
23#include "BDSDebug.hh"
24#include "BDSException.hh"
25#include "BDSGeometryExternal.hh"
26#include "BDSGeometryFactoryGDML.hh"
27#include "BDSGeometryInspector.hh"
28#include "BDSGDMLPreprocessor.hh" // also only available with USE_GDML
29#include "BDSGlobalConstants.hh"
30#include "BDSMaterials.hh"
31#include "BDSWarning.hh"
32
33#include "globals.hh"
34#include "G4Colour.hh"
35#include "G4GDMLParser.hh"
36#include "G4LogicalVolume.hh"
37#include "G4UserLimits.hh"
38#include "G4VisAttributes.hh"
39#include "G4VPhysicalVolume.hh"
40
41#include <fstream>
42#include <iostream>
43#include <iterator>
44#include <set>
45#include <sstream>
46#include <string>
47#include <utility>
48#include <vector>
49
50class G4VisAttributes;
51class G4VSolid;
52
53BDSGeometryFactoryGDML::BDSGeometryFactoryGDML()
54{;}
55
57 G4String fileName,
58 std::map<G4String, G4Colour*>* mapping,
59 G4bool autoColour,
60 G4double /*suggestedLength*/,
61 G4double /*suggestedHorizontalWidth*/,
62 std::vector<G4String>* namedVacuumVolumes,
63 G4UserLimits* userLimitsToAttachToAllLVs)
64{
65 CleanUp();
66
67 // Compensate for G4GDMLParser deficiency in loading more than one file with similar names
68 // in objects. Prepend all names with component name.
69 G4String processedFile;
70 G4bool preprocessGDML = BDSGlobalConstants::Instance()->PreprocessGDML();
71 G4bool preprocessGDMLSchema = BDSGlobalConstants::Instance()->PreprocessGDMLSchema();
72 if (preprocessGDML)
73 {processedFile = BDS::PreprocessGDML(fileName, componentName, preprocessGDMLSchema);} // use all in one method
74 else if (preprocessGDMLSchema) // generally don't process the file but process the schema to local copy only
75 {processedFile = BDS::PreprocessGDMLSchemaOnly(fileName);} // use schema only method
76 else // no processing
77 {processedFile = fileName;}
78 G4String preprocessNameToStrip = preprocessGDML ? componentName+"_" : "";
79
80 G4GDMLParser* parser = new G4GDMLParser();
81 parser->SetOverlapCheck(BDSGlobalConstants::Instance()->CheckOverlaps());
82 parser->Read(processedFile, /*validate=*/true);
83
84 G4VPhysicalVolume* containerPV = parser->GetWorldVolume();
85 G4LogicalVolume* containerLV = containerPV->GetLogicalVolume();
86 G4VSolid* containerSolid = containerLV->GetSolid();
87 G4ThreeVector gdmlWorldOrigin = G4ThreeVector();
88 if (containerPV->GetName() == "world_volume_lv_PV")
89 {
90 gdmlWorldOrigin = parser->GetPosition("PygdmlOrigin"); // TODO check if Pygdml geometry
91 gdmlWorldOrigin[2] = 0.0;
92 }
93
94 // record all pvs and lvs used in this loaded geometry
95 std::set<G4VPhysicalVolume*> pvsGDML;
96 std::set<G4LogicalVolume*> lvsGDML;
97 std::map<G4String, G4Material*> materialsGDML;
98 GetAllLogicalPhysicalAndMaterials(containerPV, pvsGDML, lvsGDML, materialsGDML);
99 BDSMaterials::Instance()->CacheMaterialsFromGDML(materialsGDML, componentName, preprocessGDML);
100
101 // load possible colours in auxiliary tags
102 std::map<G4String, G4Colour*> gdmlColours;
103 G4int iColour = 0;
104 for (auto lv : lvsGDML)
105 {
106 auto auxInfo = parser->GetVolumeAuxiliaryInformation(lv);
107 for (const auto& af : auxInfo)
108 {
109 if (af.type == "colour")
110 {
111 std::stringstream ss(af.value);
112 std::vector<G4String> colVals((std::istream_iterator<G4String>(ss)), std::istream_iterator<G4String>());
113 if (colVals.size() != 4)
114 {BDS::Warning(__METHOD_NAME__, "invalid number of colour values for logical volume " + lv->GetName());}
115 G4String colourName = componentName + "_colour_"+std::to_string(iColour);
116 iColour++;
117 G4String colourString = colourName + ":";
118 for (const auto& c : colVals)
119 {colourString += " " + c;}
120 // false = don't normalise to 255 as already done so
121 G4Colour* colour = BDSColours::Instance()->GetColour(colourString, false);
122 gdmlColours[lv->GetName()] = colour;
123 }
124 }
125 }
126
127 G4cout << "Loaded GDML file \"" << fileName << "\" containing:" << G4endl;
128 G4cout << pvsGDML.size() << " physical volumes, and " << lvsGDML.size() << " logical volumes" << G4endl;
129
130 // resolve loaded map with possible external map with minimal copying
131 std::map<G4String, G4Colour*>* mappingToUse = nullptr;
132 G4bool deleteMap = false;
133 if (!gdmlColours.empty())
134 {
135 if (mapping)
136 {// copy and extend the map
137 mappingToUse = new std::map<G4String, G4Colour*>(*mapping);
138 mappingToUse->insert(gdmlColours.begin(), gdmlColours.end());
139 deleteMap = true;
140 }
141 else
142 {mappingToUse = &gdmlColours;}
143 }
144 else
145 {mappingToUse = mapping;}
146
147 auto visesGDML = ApplyColourMapping(lvsGDML, mappingToUse, autoColour, preprocessNameToStrip);
148 if (deleteMap)
149 {delete mappingToUse;}
150
151 G4UserLimits* ul = userLimitsToAttachToAllLVs ? userLimitsToAttachToAllLVs : BDSGlobalConstants::Instance()->DefaultUserLimits();
152 ApplyUserLimits(lvsGDML, ul);
153
154 // make sure container is visible - Geant4 always makes the container invisible.
155 G4Colour* c = BDSColourFromMaterial::Instance()->GetColour(containerLV->GetMaterial());
156 G4VisAttributes* vis = new G4VisAttributes(*c);
157 vis->SetVisibility(true);
158 visesGDML.insert(vis);
159 containerLV->SetVisAttributes(vis);
160
161 std::pair<BDSExtent, BDSExtent> outerInner = BDS::DetermineExtents(containerSolid);
162
163 BDSGeometryExternal* result = new BDSGeometryExternal(containerSolid, containerLV,
164 outerInner.first, /*outer*/
165 outerInner.second, /*inner*/
166 gdmlWorldOrigin);
167 result->RegisterLogicalVolume(lvsGDML);
168 result->RegisterPhysicalVolume(pvsGDML);
169 result->RegisterVisAttributes(visesGDML);
170 result->RegisterVacuumVolumes(GetVolumes(lvsGDML, namedVacuumVolumes, preprocessGDML, componentName));
171
172 delete parser;
173 return result;
174}
175
176G4String BDSGeometryFactoryGDML::PreprocessedName(const G4String& objectName,
177 const G4String& acceleratorComponentName) const
178{return BDSGDMLPreprocessor::ProcessedNodeName(objectName, acceleratorComponentName);}
179
181 std::set<G4VPhysicalVolume*>& pvsIn,
182 std::set<G4LogicalVolume*>& lvsIn,
183 std::map<G4String, G4Material*>& materialsGDML)
184{
185 const auto& lv = volume->GetLogicalVolume();
186 lvsIn.insert(lv);
187 G4Material* mat = lv->GetMaterial();
188 materialsGDML[mat->GetName()] = mat;
189 for (G4int i = 0; i < (G4int)lv->GetNoDaughters(); i++)
190 {
191 const auto& pv = lv->GetDaughter(i);
192 pvsIn.insert(pv);
193 GetAllLogicalPhysicalAndMaterials(pv, pvsIn, lvsIn, materialsGDML); // recurse into daughter
194 }
195}
196
198 const G4String& outputFileName,
199 const G4String& key,
200 const G4String& replacement)
201{
202 // open input file in read mode
203 std::ifstream ifs(fileName);
204
205 // verify file open.
206 if (!ifs.is_open())
207 {throw BDSException(__METHOD_NAME__, "Cannot open file \"" + fileName + "\"");}
208
209 std::ofstream fout(outputFileName);
210#ifdef BDSDEBUG
211 G4cout << __METHOD_NAME__ << "Original file: " << fileName << G4endl;
212 G4cout << __METHOD_NAME__ << "Temporary file: " << outputFileName << G4endl;
213#endif
214
215 int lenOfKey = key.size();
216
217 // loop over and replace
218 std::string buffer;
219 while (std::getline(ifs, buffer))
220 {// if we find key, replace it
221 int f = buffer.find(key);
222 if (f != -1)
223 {
224 std::string outputString = std::string(buffer);
225 outputString.replace(f, lenOfKey, replacement);
226 fout << outputString << "\n"; // getline strips \n
227 }
228 else // copy line to temp file as is
229 {fout << buffer << "\n";}
230 }
231
232 // clean up
233 ifs.close();
234 fout.close();
235}
236
237#else
238// insert empty function to avoid no symbols warning
239void _SymbolToPreventWarningGeomFacGDML(){;}
240#endif
static BDSColourFromMaterial * Instance()
Singleton pattern.
G4Colour * GetColour(const G4Material *material, const G4String &prefixToStripFromName="")
Get colour from name.
static BDSColours * Instance()
singleton pattern
Definition: BDSColours.cc:33
G4Colour * GetColour(const G4String &type, G4bool normaliseTo255=true)
Get colour from name.
Definition: BDSColours.cc:202
General exception with possible name of object and message.
Definition: BDSException.hh:35
static G4String ProcessedNodeName(const G4String &nodeName, const G4String &prefix)
void RegisterLogicalVolume(G4LogicalVolume *logicalVolume)
void RegisterPhysicalVolume(G4VPhysicalVolume *physicalVolume)
void RegisterVisAttributes(G4VisAttributes *visAttribute)
A loaded piece of externally provided geometry.
void RegisterVacuumVolumes(const std::set< G4LogicalVolume * > &vacuumVolumesIn)
Register a set of volumes to be identified as vacuum volumes for the BDSAcceleratorComponent.
std::set< G4LogicalVolume * > GetVolumes(const std::set< G4LogicalVolume * > &allLVs, std::vector< G4String > *volumeNames, G4bool preprocessGDML, const G4String &componentName) const
Get the volumes that match the name. Volume names are matched exactly and are case sensitive.
virtual void CleanUp()
Virtual clean up that derived classes can override that calls CleanUpBase().
virtual void ApplyUserLimits(const std::set< G4LogicalVolume * > &lvsIn, G4UserLimits *userLimits)
Attach a set of user limits to every logical volume supplied.
virtual std::set< G4VisAttributes * > ApplyColourMapping(std::set< G4LogicalVolume * > &lvs, std::map< G4String, G4Colour * > *mapping, G4bool autoColour, const G4String &preprocessPrefix="")
void GetAllLogicalPhysicalAndMaterials(const G4VPhysicalVolume *volume, std::set< G4VPhysicalVolume * > &pvs, std::set< G4LogicalVolume * > &lvs, std::map< G4String, G4Material * > &materialsGDML)
virtual G4String PreprocessedName(const G4String &objectName, const G4String &acceleratorComponentName) const
Use the GDML preprocessing scheme to prepare the preprocesseed volume names.
void ReplaceStringInFile(const G4String &filename, const G4String &outputFileName, const G4String &key, const G4String &replacement)
virtual BDSGeometryExternal * Build(G4String componentName, G4String fileName, std::map< G4String, G4Colour * > *colourMapping=nullptr, G4bool autoColour=true, G4double suggestedLength=0, G4double suggestedHorizontalWidth=0, std::vector< G4String > *namedVacuumVolumes=nullptr, G4UserLimits *userLimitsToAttachToAllLVs=nullptr)
static BDSGlobalConstants * Instance()
Access method.
static BDSMaterials * Instance()
Singleton pattern access.
Definition: BDSMaterials.cc:38
void CacheMaterialsFromGDML(const std::map< G4String, G4Material * > &materialsGDML, const G4String &prepend, G4bool prependWasUsed)
std::pair< BDSExtent, BDSExtent > DetermineExtents(const G4VSolid *solid)