Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
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 delete ullist;
87 delete rlist;
88 }
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 return;
102
103 G4String name = iaux->value;
104 if(strip)
105 {
106 reader->StripName(name);
107 }
108 if(name.contains("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 }
132 G4LogicalVolume* lvol =
134 aRegion->AddRootLogicalVolume(lvol);
135 if(reflFactory->IsConstituent(lvol))
136 aRegion->AddRootLogicalVolume(reflFactory->GetReflectedLV(lvol));
137 }
138 else if(tag == "pcut")
139 {
140 const G4String& cvalue = raux->value;
141 const G4String& cunit = raux->unit;
142 if(G4UnitDefinition::GetCategory(cunit) != "Length")
143 {
144 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
145 FatalException, "Invalid unit for length!");
146 }
147 G4double cut =
148 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
149 pcuts->SetProductionCut(cut, "proton");
150 }
151 else if(tag == "ecut")
152 {
153 const G4String& cvalue = raux->value;
154 const G4String& cunit = raux->unit;
155 if(G4UnitDefinition::GetCategory(cunit) != "Length")
156 {
157 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
158 FatalException, "Invalid unit for length!");
159 }
160 G4double cut =
161 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
162 pcuts->SetProductionCut(cut, "e-");
163 }
164 else if(tag == "poscut")
165 {
166 const G4String& cvalue = raux->value;
167 const G4String& cunit = raux->unit;
168 if(G4UnitDefinition::GetCategory(cunit) != "Length")
169 {
170 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
171 FatalException, "Invalid unit for length!");
172 }
173 G4double cut =
174 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
175 pcuts->SetProductionCut(cut, "e+");
176 }
177 else if(tag == "gamcut")
178 {
179 const G4String& cvalue = raux->value;
180 const G4String& cunit = raux->unit;
181 if(G4UnitDefinition::GetCategory(cunit) != "Length")
182 {
183 G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
184 FatalException, "Invalid unit for length!");
185 }
186 G4double cut =
187 eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
188 pcuts->SetProductionCut(cut, "gamma");
189 }
190 else if(tag == "ulimits")
191 {
192 G4double ustepMax = DBL_MAX, utrakMax = DBL_MAX, utimeMax = DBL_MAX;
193 G4double uekinMin = 0., urangMin = 0.;
194 const G4String& ulname = raux->value;
195 for(auto uaux = raux->auxList->cbegin();
196 uaux != raux->auxList->cend(); ++uaux)
197 {
198 const G4String& ultag = uaux->type;
199 const G4String& uvalue = uaux->value;
200 const G4String& uunit = uaux->unit;
201 G4double ulvalue = eval.Evaluate(uvalue) * eval.Evaluate(uunit);
202 if(ultag == "ustepMax")
203 {
204 ustepMax = ulvalue;
205 }
206 else if(ultag == "utrakMax")
207 {
208 utrakMax = ulvalue;
209 }
210 else if(ultag == "utimeMax")
211 {
212 utimeMax = ulvalue;
213 }
214 else if(ultag == "uekinMin")
215 {
216 uekinMin = ulvalue;
217 }
218 else if(ultag == "urangMin")
219 {
220 urangMin = ulvalue;
221 }
222 else
223 {
224 G4Exception("G4GDMLParser::ImportRegions()", "ReadError",
225 FatalException, "Invalid definition of user-limits!");
226 }
227 }
228 G4UserLimits* ulimits = new G4UserLimits(
229 ulname, ustepMax, utrakMax, utimeMax, uekinMin, urangMin);
230 aRegion->SetUserLimits(ulimits);
231 }
232 else
233 continue; // Ignore unknown tags
234 }
235 }
236 }
237}
238
239// --------------------------------------------------------------------
240void G4GDMLParser::ExportRegions(G4bool storeReferences)
241{
244 for(std::size_t i = 0; i < rstore->size(); ++i)
245 // Skip default regions associated to worlds
246 {
247 const G4String& tname = (*rstore)[i]->GetName();
248 if(tname.contains("DefaultRegionForParallelWorld"))
249 continue;
250 const G4String& rname = writer->GenerateName(tname, (*rstore)[i]);
251 rlist = new G4GDMLAuxListType();
252 G4GDMLAuxStructType raux = { "Region", rname, "", rlist };
253 auto rlvol_iter = (*rstore)[i]->GetRootLogicalVolumeIterator();
254 for(std::size_t j = 0; j < (*rstore)[i]->GetNumberOfRootVolumes(); ++j)
255 {
256 G4LogicalVolume* rlvol = *rlvol_iter;
257 if(reflFactory->IsReflected(rlvol))
258 continue;
259 G4String vname = writer->GenerateName(rlvol->GetName(), rlvol);
260 if(!storeReferences)
261 {
262 reader->StripName(vname);
263 }
264 G4GDMLAuxStructType rsubaux = { "volume", vname, "", 0 };
265 rlist->push_back(rsubaux);
266 ++rlvol_iter;
267 }
268 G4double gam_cut
269 = (*rstore)[i]->GetProductionCuts()->GetProductionCut("gamma");
271 = { "gamcut", eval.ConvertToString(gam_cut), "mm", 0 };
272 rlist->push_back(caux1);
273 G4double e_cut = (*rstore)[i]->GetProductionCuts()->GetProductionCut("e-");
275 = { "ecut", eval.ConvertToString(e_cut), "mm", 0 };
276 rlist->push_back(caux2);
277 G4double pos_cut
278 = (*rstore)[i]->GetProductionCuts()->GetProductionCut("e+");
280 = { "poscut", eval.ConvertToString(pos_cut), "mm", 0 };
281 rlist->push_back(caux3);
282 G4double p_cut
283 = (*rstore)[i]->GetProductionCuts()->GetProductionCut("proton");
285 = { "pcut", eval.ConvertToString(p_cut), "mm", 0 };
286 rlist->push_back(caux4);
287 if((*rstore)[i]->GetUserLimits())
288 {
289 const G4Track fake_trk;
290 ullist = new G4GDMLAuxListType();
291 const G4String& utype = (*rstore)[i]->GetUserLimits()->GetType();
292 G4GDMLAuxStructType uaux = { "ulimits", utype, "mm", ullist };
293 G4double max_step
294 = (*rstore)[i]->GetUserLimits()->GetMaxAllowedStep(fake_trk);
296 = { "ustepMax", eval.ConvertToString(max_step), "mm", 0 };
297 ullist->push_back(ulaux1);
298 G4double max_trk
299 = (*rstore)[i]->GetUserLimits()->GetUserMaxTrackLength(fake_trk);
301 = { "utrakMax", eval.ConvertToString(max_trk), "mm", 0 };
302 ullist->push_back(ulaux2);
303 G4double max_time
304 = (*rstore)[i]->GetUserLimits()->GetUserMaxTime(fake_trk);
306 = { "utimeMax", eval.ConvertToString(max_time), "mm", 0 };
307 ullist->push_back(ulaux3);
308 G4double min_ekin
309 = (*rstore)[i]->GetUserLimits()->GetUserMinEkine(fake_trk);
311 = { "uekinMin", eval.ConvertToString(min_ekin), "mm", 0 };
312 ullist->push_back(ulaux4);
313 G4double min_rng
314 = (*rstore)[i]->GetUserLimits()->GetUserMinRange(fake_trk);
316 = { "urangMin", eval.ConvertToString(min_rng), "mm", 0 };
317 ullist->push_back(ulaux5);
318 rlist->push_back(uaux);
319 }
320 AddAuxiliary(raux);
321 }
322}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
G4double Evaluate(const G4String &)
G4String ConvertToString(G4int ival)
const G4GDMLAuxListType * GetAuxList() const
void AddAuxiliary(G4GDMLAuxStructType myaux)
void StripName(G4String &) const
Definition: G4GDMLRead.cc:111
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:132
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4LogicalVolumeStore * GetInstance()
const G4String & GetName() const
void SetProductionCut(G4double cut, G4int index=-1)
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:283
G4bool contains(const std::string &) const
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