Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GDMLParser.cc
Go to the documentation of this file.
1//
2// ********************************************************************
3// * License and Disclaimer *
4// * *
5// * The Geant4 software is copyright of the Copyright Holders of *
6// * the Geant4 Collaboration. It is provided under the terms and *
7// * conditions of the Geant4 Software License, included in the file *
8// * LICENSE and available at http://cern.ch/geant4/license . These *
9// * include a list of copyright holders. *
10// * *
11// * Neither the authors of this software system, nor their employing *
12// * institutes,nor the agencies providing financial support for this *
13// * work make any representation or warranty, express or implied, *
14// * regarding this software system or assume any liability for its *
15// * use. Please see the license in the file LICENSE and URL above *
16// * for the full disclaimer and the limitation of liability. *
17// * *
18// * This code implementation is the result of the scientific and *
19// * technical work of the GEANT4 collaboration. *
20// * By using, copying, modifying or distributing the software (or *
21// * any work based on the software) you agree to acknowledge its *
22// * use in resulting scientific publications, and indicate your *
23// * acceptance of all terms of the Geant4 Software license. *
24// ********************************************************************
25//
26// G4GDMLParser implementation
27//
28// Author: Zoltan Torzsok, November 2007
29// --------------------------------------------------------------------
30
31#include "G4GDMLParser.hh"
32
33#include "G4UnitsTable.hh"
35#include "G4RegionStore.hh"
36#include "G4UserLimits.hh"
37#include "G4ProductionCuts.hh"
39#include "G4Track.hh"
40
41// --------------------------------------------------------------------
43 : strip(true)
44{
45 reader = new G4GDMLReadStructure;
46 writer = new G4GDMLWriteStructure;
47 messenger = new G4GDMLMessenger(this);
48
49 xercesc::XMLPlatformUtils::Initialize();
50}
51
52// --------------------------------------------------------------------
54 : urcode(true), strip(true)
55{
56 reader = extr;
57 writer = new G4GDMLWriteStructure;
58 messenger = new G4GDMLMessenger(this);
59
60 xercesc::XMLPlatformUtils::Initialize();
61}
62
63// --------------------------------------------------------------------
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}
74
75// --------------------------------------------------------------------
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}
92
93// --------------------------------------------------------------------
94void G4GDMLParser::ImportRegions()
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 {
106 reader->StripName(name);
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 =
180 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
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 =
193 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
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 =
206 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
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 =
219 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
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}
270
271// --------------------------------------------------------------------
272void G4GDMLParser::ExportRegions(G4bool storeReferences)
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;
322 ullist = new G4GDMLAuxListType();
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}
@ JustWarning
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
double G4double
Definition G4Types.hh:83
bool G4bool
Definition G4Types.hh:86
#define G4endl
Definition G4ios.hh:67
G4double Evaluate(const G4String &)
G4String ConvertToString(G4int ival)
const G4GDMLAuxListType * GetAuxList() const
void AddAuxiliary(G4GDMLAuxStructType myaux)
void StripName(G4String &) const
G4String GenerateName(const G4String &, const void *const)
const std::map< G4String, std::vector< G4LogicalVolume * > > & GetMap() const
static G4LogicalVolumeStore * GetInstance()
const G4String & GetName() const
void SetProductionCut(G4double cut, G4int index)
G4bool IsReflected(G4LogicalVolume *lv) const
G4LogicalVolume * GetReflectedLV(G4LogicalVolume *lv) const
static G4ReflectionFactory * Instance()
G4bool IsConstituent(G4LogicalVolume *lv) const
static G4RegionStore * GetInstance()
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)
#define DBL_MAX
Definition templates.hh:62