Geant4-11
Public Member Functions | Private Member Functions | Private Attributes
G4GDMLParser Class Reference

#include <G4GDMLParser.hh>

Public Member Functions

void AddAuxiliary (G4GDMLAuxStructType myaux)
 
void AddModule (const G4int depth)
 
void AddModule (const G4VPhysicalVolume *const physvol)
 
void AddVolumeAuxiliary (G4GDMLAuxStructType myaux, const G4LogicalVolume *const lvol)
 
void Clear ()
 
 G4GDMLParser ()
 
 G4GDMLParser (G4GDMLReadStructure *)
 
 G4GDMLParser (G4GDMLReadStructure *, G4GDMLWriteStructure *)
 
const G4GDMLAuxListTypeGetAuxList () const
 
const G4GDMLAuxMapTypeGetAuxMap () const
 
G4double GetConstant (const G4String &name) const
 
G4GDMLMatrix GetMatrix (const G4String &name) const
 
G4int GetMaxExportLevel () const
 
G4VPhysicalVolumeGetPhysVolume (const G4String &name) const
 
G4ThreeVector GetPosition (const G4String &name) const
 
G4double GetQuantity (const G4String &name) const
 
G4ThreeVector GetRotation (const G4String &name) const
 
G4ThreeVector GetScale (const G4String &name) const
 
G4double GetVariable (const G4String &name) const
 
G4LogicalVolumeGetVolume (const G4String &name) const
 
G4GDMLAuxListType GetVolumeAuxiliaryInformation (G4LogicalVolume *lvol) const
 
G4VPhysicalVolumeGetWorldVolume (const G4String &setupName="Default") const
 
G4bool IsValid (const G4String &name) const
 
G4LogicalVolumeParseST (const G4String &name, G4Material *medium, G4Material *solid)
 
void Read (const G4String &filename, G4bool Validate=true)
 
void ReadModule (const G4String &filename, G4bool Validate=true)
 
void SetAddPointerToName (G4bool set)
 
void SetEnergyCutsExport (G4bool)
 
void SetMaxExportLevel (G4int)
 
void SetOutputFileOverwrite (G4bool flag)
 
void SetOverlapCheck (G4bool)
 
void SetRegionExport (G4bool)
 
void SetReverseSearch (G4bool)
 
void SetSDExport (G4bool)
 
void SetStripFlag (G4bool)
 
void StripNamePointers () const
 
void Write (const G4String &filename, const G4LogicalVolume *lvol, G4bool storeReferences=true, const G4String &SchemaLocation=G4GDML_DEFAULT_SCHEMALOCATION)
 
void Write (const G4String &filename, const G4VPhysicalVolume *pvol=0, G4bool storeReferences=true, const G4String &SchemaLocation=G4GDML_DEFAULT_SCHEMALOCATION)
 
 ~G4GDMLParser ()
 

Private Member Functions

void ExportRegions (G4bool storeReferences=true)
 
void ImportRegions ()
 

Private Attributes

G4GDMLEvaluator eval
 
G4GDMLMessengermessenger = nullptr
 
G4GDMLReadStructurereader = nullptr
 
G4bool rexp = false
 
G4GDMLAuxListTyperlist = nullptr
 
G4bool strip = false
 
G4GDMLAuxListTypeullist = nullptr
 
G4bool urcode = false
 
G4bool uwcode = false
 
G4GDMLWriteStructurewriter = nullptr
 

Detailed Description

Definition at line 51 of file G4GDMLParser.hh.

Constructor & Destructor Documentation

◆ G4GDMLParser() [1/3]

G4GDMLParser::G4GDMLParser ( )

Definition at line 42 of file G4GDMLParser.cc.

43 : strip(true)
44{
47 messenger = new G4GDMLMessenger(this);
48
49 xercesc::XMLPlatformUtils::Initialize();
50}
G4GDMLWriteStructure * writer
G4GDMLMessenger * messenger
G4GDMLReadStructure * reader

References messenger, reader, and writer.

◆ G4GDMLParser() [2/3]

G4GDMLParser::G4GDMLParser ( G4GDMLReadStructure extr)

Definition at line 53 of file G4GDMLParser.cc.

54 : urcode(true), strip(true)
55{
56 reader = extr;
58 messenger = new G4GDMLMessenger(this);
59
60 xercesc::XMLPlatformUtils::Initialize();
61}

References messenger, reader, and writer.

◆ G4GDMLParser() [3/3]

G4GDMLParser::G4GDMLParser ( G4GDMLReadStructure extr,
G4GDMLWriteStructure extw 
)

Definition at line 64 of file G4GDMLParser.cc.

66 : urcode(true), uwcode(true), strip(true)
67{
68 reader = extr;
69 writer = extw;
70 messenger = new G4GDMLMessenger(this);
71
72 xercesc::XMLPlatformUtils::Initialize();
73}

References messenger, reader, and writer.

◆ ~G4GDMLParser()

G4GDMLParser::~G4GDMLParser ( )

Definition at line 76 of file G4GDMLParser.cc.

77{
78 xercesc::XMLPlatformUtils::Terminate();
79 if(!urcode)
80 {
81 delete reader;
82 }
83 if(!uwcode)
84 {
85 delete writer;
86 }
87 delete ullist;
88 delete rlist;
89
90 delete messenger;
91}
G4GDMLAuxListType * rlist
G4GDMLAuxListType * ullist

References messenger, reader, rlist, ullist, urcode, uwcode, and writer.

Member Function Documentation

◆ AddAuxiliary()

void G4GDMLParser::AddAuxiliary ( G4GDMLAuxStructType  myaux)
inline

Referenced by ExportRegions().

◆ AddModule() [1/2]

void G4GDMLParser::AddModule ( const G4int  depth)
inline

◆ AddModule() [2/2]

void G4GDMLParser::AddModule ( const G4VPhysicalVolume *const  physvol)
inline

◆ AddVolumeAuxiliary()

void G4GDMLParser::AddVolumeAuxiliary ( G4GDMLAuxStructType  myaux,
const G4LogicalVolume *const  lvol 
)
inline

◆ Clear()

void G4GDMLParser::Clear ( )
inline

◆ ExportRegions()

void G4GDMLParser::ExportRegions ( G4bool  storeReferences = true)
private

Definition at line 272 of file G4GDMLParser.cc.

273{
276 for(std::size_t i = 0; i < rstore->size(); ++i)
277 // Skip default regions associated to worlds
278 {
279 const G4String& tname = (*rstore)[i]->GetName();
280 if(G4StrUtil::contains(tname, "DefaultRegionForParallelWorld"))
281 continue;
282 const G4String& rname = writer->GenerateName(tname, (*rstore)[i]);
283 rlist = new G4GDMLAuxListType();
284 G4GDMLAuxStructType raux = { "Region", rname, "", rlist };
285 auto rlvol_iter = (*rstore)[i]->GetRootLogicalVolumeIterator();
286 for(std::size_t j = 0; j < (*rstore)[i]->GetNumberOfRootVolumes(); ++j)
287 {
288 G4LogicalVolume* rlvol = *rlvol_iter;
289 if(reflFactory->IsReflected(rlvol))
290 continue;
291 G4String vname = writer->GenerateName(rlvol->GetName(), rlvol);
292 if(!storeReferences)
293 {
294 reader->StripName(vname);
295 }
296 G4GDMLAuxStructType rsubaux = { "volume", vname, "", 0 };
297 rlist->push_back(rsubaux);
298 ++rlvol_iter;
299 }
300 G4double gam_cut
301 = (*rstore)[i]->GetProductionCuts()->GetProductionCut("gamma");
303 = { "gamcut", eval.ConvertToString(gam_cut), "mm", 0 };
304 rlist->push_back(caux1);
305 G4double e_cut = (*rstore)[i]->GetProductionCuts()->GetProductionCut("e-");
307 = { "ecut", eval.ConvertToString(e_cut), "mm", 0 };
308 rlist->push_back(caux2);
309 G4double pos_cut
310 = (*rstore)[i]->GetProductionCuts()->GetProductionCut("e+");
312 = { "poscut", eval.ConvertToString(pos_cut), "mm", 0 };
313 rlist->push_back(caux3);
314 G4double p_cut
315 = (*rstore)[i]->GetProductionCuts()->GetProductionCut("proton");
317 = { "pcut", eval.ConvertToString(p_cut), "mm", 0 };
318 rlist->push_back(caux4);
319 if((*rstore)[i]->GetUserLimits())
320 {
321 const G4Track fake_trk;
323 const G4String& utype = (*rstore)[i]->GetUserLimits()->GetType();
324 G4GDMLAuxStructType uaux = { "ulimits", utype, "mm", ullist };
325 G4double max_step
326 = (*rstore)[i]->GetUserLimits()->GetMaxAllowedStep(fake_trk);
328 = { "ustepMax", eval.ConvertToString(max_step), "mm", 0 };
329 ullist->push_back(ulaux1);
330 G4double max_trk
331 = (*rstore)[i]->GetUserLimits()->GetUserMaxTrackLength(fake_trk);
333 = { "utrakMax", eval.ConvertToString(max_trk), "mm", 0 };
334 ullist->push_back(ulaux2);
335 G4double max_time
336 = (*rstore)[i]->GetUserLimits()->GetUserMaxTime(fake_trk);
338 = { "utimeMax", eval.ConvertToString(max_time), "mm", 0 };
339 ullist->push_back(ulaux3);
340 G4double min_ekin
341 = (*rstore)[i]->GetUserLimits()->GetUserMinEkine(fake_trk);
343 = { "uekinMin", eval.ConvertToString(min_ekin), "mm", 0 };
344 ullist->push_back(ulaux4);
345 G4double min_rng
346 = (*rstore)[i]->GetUserLimits()->GetUserMinRange(fake_trk);
348 = { "urangMin", eval.ConvertToString(min_rng), "mm", 0 };
349 ullist->push_back(ulaux5);
350 rlist->push_back(uaux);
351 }
352 AddAuxiliary(raux);
353 }
354}
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
double G4double
Definition: G4Types.hh:83
G4String ConvertToString(G4int ival)
G4GDMLEvaluator eval
void AddAuxiliary(G4GDMLAuxStructType myaux)
void StripName(G4String &) const
Definition: G4GDMLRead.cc:112
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:132
const G4String & GetName() const
G4bool IsReflected(G4LogicalVolume *lv) const
static G4ReflectionFactory * Instance()
static G4RegionStore * GetInstance()
G4bool contains(const G4String &str, std::string_view ss)
Check if a string contains a given substring.

References AddAuxiliary(), G4StrUtil::contains(), G4GDMLEvaluator::ConvertToString(), eval, G4GDMLWrite::GenerateName(), G4RegionStore::GetInstance(), G4LogicalVolume::GetName(), G4ReflectionFactory::Instance(), G4ReflectionFactory::IsReflected(), reader, rlist, G4GDMLRead::StripName(), ullist, and writer.

◆ GetAuxList()

const G4GDMLAuxListType * G4GDMLParser::GetAuxList ( ) const
inline

Referenced by ImportRegions().

◆ GetAuxMap()

const G4GDMLAuxMapType * G4GDMLParser::GetAuxMap ( ) const
inline

◆ GetConstant()

G4double G4GDMLParser::GetConstant ( const G4String name) const
inline

◆ GetMatrix()

G4GDMLMatrix G4GDMLParser::GetMatrix ( const G4String name) const
inline

◆ GetMaxExportLevel()

G4int G4GDMLParser::GetMaxExportLevel ( ) const
inline

◆ GetPhysVolume()

G4VPhysicalVolume * G4GDMLParser::GetPhysVolume ( const G4String name) const
inline

◆ GetPosition()

G4ThreeVector G4GDMLParser::GetPosition ( const G4String name) const
inline

◆ GetQuantity()

G4double G4GDMLParser::GetQuantity ( const G4String name) const
inline

◆ GetRotation()

G4ThreeVector G4GDMLParser::GetRotation ( const G4String name) const
inline

◆ GetScale()

G4ThreeVector G4GDMLParser::GetScale ( const G4String name) const
inline

◆ GetVariable()

G4double G4GDMLParser::GetVariable ( const G4String name) const
inline

◆ GetVolume()

G4LogicalVolume * G4GDMLParser::GetVolume ( const G4String name) const
inline

◆ GetVolumeAuxiliaryInformation()

G4GDMLAuxListType G4GDMLParser::GetVolumeAuxiliaryInformation ( G4LogicalVolume lvol) const
inline

◆ GetWorldVolume()

G4VPhysicalVolume * G4GDMLParser::GetWorldVolume ( const G4String setupName = "Default") const
inline

◆ ImportRegions()

void G4GDMLParser::ImportRegions ( )
private

Definition at line 94 of file G4GDMLParser.cc.

95{
97 const G4GDMLAuxListType* auxInfoList = GetAuxList();
98 for(auto iaux = auxInfoList->cbegin(); iaux != auxInfoList->cend(); ++iaux)
99 {
100 if(iaux->type != "Region")
101 continue;
102
103 G4String name = iaux->value;
104 if(strip)
105 {
107 }
108 if(G4StrUtil::contains(name, "DefaultRegionForTheWorld"))
109 continue;
110
111 if(!iaux->auxList)
112 {
113 G4Exception("G4GDMLParser::ImportRegions()", "ReadError", FatalException,
114 "Invalid definition of geometrical region!");
115 }
116 else // Create region and loop over all region attributes
117 {
118 G4Region* aRegion = new G4Region(name);
119 G4ProductionCuts* pcuts = new G4ProductionCuts();
120 aRegion->SetProductionCuts(pcuts);
121 for(auto raux = iaux->auxList->cbegin();
122 raux != iaux->auxList->cend(); ++raux)
123 {
124 const G4String& tag = raux->type;
125 if(tag == "volume")
126 {
127 G4String volname = raux->value;
128 if(strip)
129 {
130 reader->StripName(volname);
131 }
133 auto pos = store->GetMap().find(volname);
134 if(pos != store->GetMap().cend())
135 {
136 // Scan for all possible volumes with same name
137 // and set them as root logical volumes for the region...
138 // Issue a notification in case more than one logical volume
139 // with same name exist and get set.
140 //
141 if (pos->second.size()>1)
142 {
143 std::ostringstream message;
144 message << "There exists more than ONE logical volume "
145 << "in store named: " << volname << "." << G4endl
146 << "NOTE: assigning all such volumes as root logical "
147 << "volumes for region: " << name << "!";
148 G4Exception("G4GDMLParser::ImportRegions()",
149 "Notification", JustWarning, message);
150 }
151 for (auto vpos = pos->second.cbegin();
152 vpos != pos->second.cend(); ++vpos)
153 {
154 aRegion->AddRootLogicalVolume(*vpos);
155 if(reflFactory->IsConstituent(*vpos))
156 aRegion->AddRootLogicalVolume(reflFactory->GetReflectedLV(*vpos));
157 }
158 }
159 else
160 {
161 std::ostringstream message;
162 message << "Volume NOT found in store !" << G4endl
163 << " Volume " << volname << " NOT found in store !"
164 << G4endl
165 << " No region is being set.";
166 G4Exception("G4GDMLParser::ImportRegions()",
167 "InvalidSetup", JustWarning, message);
168 }
169 }
170 else if(tag == "pcut")
171 {
172 const G4String& cvalue = raux->value;
173 const G4String& cunit = raux->unit;
174 if(G4UnitDefinition::GetCategory(cunit) != "Length")
175 {
176 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
177 FatalException, "Invalid unit for length!");
178 }
179 G4double cut =
181 pcuts->SetProductionCut(cut, "proton");
182 }
183 else if(tag == "ecut")
184 {
185 const G4String& cvalue = raux->value;
186 const G4String& cunit = raux->unit;
187 if(G4UnitDefinition::GetCategory(cunit) != "Length")
188 {
189 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
190 FatalException, "Invalid unit for length!");
191 }
192 G4double cut =
194 pcuts->SetProductionCut(cut, "e-");
195 }
196 else if(tag == "poscut")
197 {
198 const G4String& cvalue = raux->value;
199 const G4String& cunit = raux->unit;
200 if(G4UnitDefinition::GetCategory(cunit) != "Length")
201 {
202 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
203 FatalException, "Invalid unit for length!");
204 }
205 G4double cut =
207 pcuts->SetProductionCut(cut, "e+");
208 }
209 else if(tag == "gamcut")
210 {
211 const G4String& cvalue = raux->value;
212 const G4String& cunit = raux->unit;
213 if(G4UnitDefinition::GetCategory(cunit) != "Length")
214 {
215 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
216 FatalException, "Invalid unit for length!");
217 }
218 G4double cut =
220 pcuts->SetProductionCut(cut, "gamma");
221 }
222 else if(tag == "ulimits")
223 {
224 G4double ustepMax = DBL_MAX, utrakMax = DBL_MAX, utimeMax = DBL_MAX;
225 G4double uekinMin = 0., urangMin = 0.;
226 const G4String& ulname = raux->value;
227 for(auto uaux = raux->auxList->cbegin();
228 uaux != raux->auxList->cend(); ++uaux)
229 {
230 const G4String& ultag = uaux->type;
231 const G4String& uvalue = uaux->value;
232 const G4String& uunit = uaux->unit;
233 G4double ulvalue = eval.Evaluate(uvalue) * eval.Evaluate(uunit);
234 if(ultag == "ustepMax")
235 {
236 ustepMax = ulvalue;
237 }
238 else if(ultag == "utrakMax")
239 {
240 utrakMax = ulvalue;
241 }
242 else if(ultag == "utimeMax")
243 {
244 utimeMax = ulvalue;
245 }
246 else if(ultag == "uekinMin")
247 {
248 uekinMin = ulvalue;
249 }
250 else if(ultag == "urangMin")
251 {
252 urangMin = ulvalue;
253 }
254 else
255 {
256 G4Exception("G4GDMLParser::ImportRegions()", "ReadError",
257 FatalException, "Invalid definition of user-limits!");
258 }
259 }
260 G4UserLimits* ulimits = new G4UserLimits(
261 ulname, ustepMax, utrakMax, utimeMax, uekinMin, urangMin);
262 aRegion->SetUserLimits(ulimits);
263 }
264 else
265 continue; // Ignore unknown tags
266 }
267 }
268 }
269}
static const G4double pos
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
#define G4endl
Definition: G4ios.hh:57
G4double Evaluate(const G4String &)
const G4GDMLAuxListType * GetAuxList() const
const std::map< G4String, std::vector< G4LogicalVolume * > > & GetMap() const
static G4LogicalVolumeStore * GetInstance()
void SetProductionCut(G4double cut, G4int index=-1)
G4LogicalVolume * GetReflectedLV(G4LogicalVolume *lv) const
G4bool IsConstituent(G4LogicalVolume *lv) const
void SetProductionCuts(G4ProductionCuts *cut)
void SetUserLimits(G4UserLimits *ul)
void AddRootLogicalVolume(G4LogicalVolume *lv, G4bool search=true)
Definition: G4Region.cc:293
static G4double GetValueOf(const G4String &)
static G4String GetCategory(const G4String &)
const char * name(G4int ptype)
Definition: xmlparse.cc:187
#define DBL_MAX
Definition: templates.hh:62

References G4Region::AddRootLogicalVolume(), G4StrUtil::contains(), DBL_MAX, eval, G4GDMLEvaluator::Evaluate(), FatalException, G4endl, G4Exception(), GetAuxList(), G4UnitDefinition::GetCategory(), G4LogicalVolumeStore::GetInstance(), G4LogicalVolumeStore::GetMap(), G4ReflectionFactory::GetReflectedLV(), G4UnitDefinition::GetValueOf(), G4ReflectionFactory::Instance(), G4ReflectionFactory::IsConstituent(), JustWarning, G4InuclParticleNames::name(), pos, reader, G4ProductionCuts::SetProductionCut(), G4Region::SetProductionCuts(), G4Region::SetUserLimits(), strip, and G4GDMLRead::StripName().

◆ IsValid()

G4bool G4GDMLParser::IsValid ( const G4String name) const
inline

◆ ParseST()

G4LogicalVolume * G4GDMLParser::ParseST ( const G4String name,
G4Material medium,
G4Material solid 
)
inline

Referenced by export_G4GDMLParser().

◆ Read()

void G4GDMLParser::Read ( const G4String filename,
G4bool  Validate = true 
)
inline

◆ ReadModule()

void G4GDMLParser::ReadModule ( const G4String filename,
G4bool  Validate = true 
)
inline

Referenced by export_G4GDMLParser().

◆ SetAddPointerToName()

void G4GDMLParser::SetAddPointerToName ( G4bool  set)
inline

◆ SetEnergyCutsExport()

void G4GDMLParser::SetEnergyCutsExport ( G4bool  )
inline

◆ SetMaxExportLevel()

void G4GDMLParser::SetMaxExportLevel ( G4int  )
inline

◆ SetOutputFileOverwrite()

void G4GDMLParser::SetOutputFileOverwrite ( G4bool  flag)
inline

◆ SetOverlapCheck()

void G4GDMLParser::SetOverlapCheck ( G4bool  )
inline

◆ SetRegionExport()

void G4GDMLParser::SetRegionExport ( G4bool  )
inline

◆ SetReverseSearch()

void G4GDMLParser::SetReverseSearch ( G4bool  )
inline

◆ SetSDExport()

void G4GDMLParser::SetSDExport ( G4bool  )
inline

◆ SetStripFlag()

void G4GDMLParser::SetStripFlag ( G4bool  )
inline

◆ StripNamePointers()

void G4GDMLParser::StripNamePointers ( ) const
inline

◆ Write() [1/2]

void G4GDMLParser::Write ( const G4String filename,
const G4LogicalVolume lvol,
G4bool  storeReferences = true,
const G4String SchemaLocation = G4GDML_DEFAULT_SCHEMALOCATION 
)
inline

◆ Write() [2/2]

void G4GDMLParser::Write ( const G4String filename,
const G4VPhysicalVolume pvol = 0,
G4bool  storeReferences = true,
const G4String SchemaLocation = G4GDML_DEFAULT_SCHEMALOCATION 
)
inline

Field Documentation

◆ eval

G4GDMLEvaluator G4GDMLParser::eval
private

Definition at line 150 of file G4GDMLParser.hh.

Referenced by ExportRegions(), and ImportRegions().

◆ messenger

G4GDMLMessenger* G4GDMLParser::messenger = nullptr
private

Definition at line 154 of file G4GDMLParser.hh.

Referenced by G4GDMLParser(), and ~G4GDMLParser().

◆ reader

G4GDMLReadStructure* G4GDMLParser::reader = nullptr
private

Definition at line 151 of file G4GDMLParser.hh.

Referenced by ExportRegions(), G4GDMLParser(), ImportRegions(), and ~G4GDMLParser().

◆ rexp

G4bool G4GDMLParser::rexp = false
private

Definition at line 155 of file G4GDMLParser.hh.

◆ rlist

G4GDMLAuxListType* G4GDMLParser::rlist = nullptr
private

Definition at line 153 of file G4GDMLParser.hh.

Referenced by ExportRegions(), and ~G4GDMLParser().

◆ strip

G4bool G4GDMLParser::strip = false
private

Definition at line 155 of file G4GDMLParser.hh.

Referenced by ImportRegions().

◆ ullist

G4GDMLAuxListType * G4GDMLParser::ullist = nullptr
private

Definition at line 153 of file G4GDMLParser.hh.

Referenced by ExportRegions(), and ~G4GDMLParser().

◆ urcode

G4bool G4GDMLParser::urcode = false
private

Definition at line 155 of file G4GDMLParser.hh.

Referenced by ~G4GDMLParser().

◆ uwcode

G4bool G4GDMLParser::uwcode = false
private

Definition at line 155 of file G4GDMLParser.hh.

Referenced by ~G4GDMLParser().

◆ writer

G4GDMLWriteStructure* G4GDMLParser::writer = nullptr
private

Definition at line 152 of file G4GDMLParser.hh.

Referenced by ExportRegions(), G4GDMLParser(), and ~G4GDMLParser().


The documentation for this class was generated from the following files: