Geant4  10.03.p01
 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 //
27 // $Id: G4GDMLParser.cc 96190 2016-03-29 08:07:36Z gcosmo $
28 //
29 //
30 // class G4GDMLParser Implementation
31 //
32 // -------------------------------------------------------------------------
33 
34 #include "G4GDMLParser.hh"
35 
36 #include "G4UnitsTable.hh"
37 #include "G4LogicalVolumeStore.hh"
38 #include "G4RegionStore.hh"
39 #include "G4UserLimits.hh"
40 #include "G4ProductionCuts.hh"
41 #include "G4ReflectionFactory.hh"
42 #include "G4Track.hh"
43 
45  : rlist(0), ullist(0), urcode(false), uwcode(false), strip(true), rexp(false)
46 {
47  reader = new G4GDMLReadStructure;
48  writer = new G4GDMLWriteStructure;
49  messenger = new G4GDMLMessenger(this);
50 
51  xercesc::XMLPlatformUtils::Initialize();
52 }
53 
55  : rlist(0), ullist(0), urcode(true), uwcode(false), strip(true), rexp(false)
56 {
57  reader = extr;
58  writer = new G4GDMLWriteStructure;
59  messenger = new G4GDMLMessenger(this);
60 
61  xercesc::XMLPlatformUtils::Initialize();
62 }
63 
66  : rlist(0), ullist(0), urcode(true), uwcode(true), strip(true), rexp(false)
67 {
68  reader = extr;
69  writer = extw;
70  messenger = new G4GDMLMessenger(this);
71 
72  xercesc::XMLPlatformUtils::Initialize();
73 }
74 
76 {
77  xercesc::XMLPlatformUtils::Terminate();
78  if (!urcode) { delete reader; }
79  if (!uwcode) { delete writer; delete ullist; delete rlist; }
80 
81  delete messenger;
82 }
83 
84 void G4GDMLParser::ImportRegions()
85 {
87  const G4GDMLAuxListType* auxInfoList = GetAuxList();
88  for(std::vector<G4GDMLAuxStructType>::const_iterator
89  iaux = auxInfoList->begin(); iaux != auxInfoList->end(); iaux++ )
90  {
91  if (iaux->type != "Region") return;
92 
93  G4String name = iaux->value;
94  if (strip) { reader->StripName(name); }
95  if (name.contains("DefaultRegionForTheWorld")) continue;
96 
97  if (!iaux->auxList)
98  {
99  G4Exception("G4GDMLParser::ImportRegions()",
100  "ReadError", FatalException,
101  "Invalid definition of geometrical region!");
102  }
103  else // Create region and loop over all region attributes
104  {
105  G4Region* aRegion = new G4Region(name);
106  G4ProductionCuts* pcuts = new G4ProductionCuts();
107  aRegion->SetProductionCuts(pcuts);
108  for(std::vector<G4GDMLAuxStructType>::const_iterator
109  raux = iaux->auxList->begin(); raux != iaux->auxList->end(); raux++ )
110  {
111  const G4String& tag = raux->type;
112  if (tag=="volume")
113  {
114  G4String volname = raux->value;
115  if (strip) { reader->StripName(volname); }
117  aRegion->AddRootLogicalVolume(lvol);
118  if(reflFactory->IsConstituent(lvol)) aRegion->AddRootLogicalVolume(reflFactory->GetReflectedLV(lvol));
119  } else
120  if (tag=="pcut")
121  {
122  const G4String& cvalue = raux->value;
123  const G4String& cunit = raux->unit;
124  if (G4UnitDefinition::GetCategory(cunit)!="Length") {
125  G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
126  FatalException, "Invalid unit for length!"); }
127  G4double cut = eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
128  pcuts->SetProductionCut(cut,"proton");
129  } else
130  if (tag=="ecut")
131  {
132  const G4String& cvalue = raux->value;
133  const G4String& cunit = raux->unit;
134  if (G4UnitDefinition::GetCategory(cunit)!="Length") {
135  G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
136  FatalException, "Invalid unit for length!"); }
137  G4double cut = eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
138  pcuts->SetProductionCut(cut,"e-");
139  } else
140  if (tag=="poscut")
141  {
142  const G4String& cvalue = raux->value;
143  const G4String& cunit = raux->unit;
144  if (G4UnitDefinition::GetCategory(cunit)!="Length") {
145  G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
146  FatalException, "Invalid unit for length!"); }
147  G4double cut = eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
148  pcuts->SetProductionCut(cut,"e+");
149  } else
150  if (tag=="gamcut")
151  {
152  const G4String& cvalue = raux->value;
153  const G4String& cunit = raux->unit;
154  if (G4UnitDefinition::GetCategory(cunit)!="Length") {
155  G4Exception("G4GDMLParser::ImportRegions()", "InvalidRead",
156  FatalException, "Invalid unit for length!"); }
157  G4double cut = eval.Evaluate(cvalue) * G4UnitDefinition::GetValueOf(cunit);
158  pcuts->SetProductionCut(cut,"gamma");
159  } else
160  if (tag=="ulimits")
161  {
162  G4double ustepMax=DBL_MAX, utrakMax=DBL_MAX, utimeMax=DBL_MAX;
163  G4double uekinMin=0., urangMin=0.;
164  const G4String& ulname = raux->value;
165  for(std::vector<G4GDMLAuxStructType>::const_iterator
166  uaux=raux->auxList->begin(); uaux!=raux->auxList->end(); uaux++ )
167  {
168  const G4String& ultag = uaux->type;
169  const G4String& uvalue = uaux->value;
170  const G4String& uunit = uaux->unit;
171  G4double ulvalue = eval.Evaluate(uvalue) * eval.Evaluate(uunit);
172  if (ultag=="ustepMax") { ustepMax = ulvalue; } else
173  if (ultag=="utrakMax") { utrakMax = ulvalue; } else
174  if (ultag=="utimeMax") { utimeMax = ulvalue; } else
175  if (ultag=="uekinMin") { uekinMin = ulvalue; } else
176  if (ultag=="urangMin") { urangMin = ulvalue; }
177  else
178  {
179  G4Exception("G4GDMLParser::ImportRegions()",
180  "ReadError", FatalException,
181  "Invalid definition of user-limits!");
182  }
183  }
184  G4UserLimits* ulimits = new G4UserLimits(ulname, ustepMax, utrakMax,
185  utimeMax, uekinMin, urangMin);
186  aRegion->SetUserLimits(ulimits);
187  }
188  else continue; // Ignore unknown tags
189  }
190  }
191  }
192 }
193 
194 void G4GDMLParser::ExportRegions(G4bool storeReferences)
195 {
198  for (size_t i=0; i<rstore->size(); ++i) // Skip default regions associated to worlds
199  {
200  const G4String& tname = (*rstore)[i]->GetName();
201  if (tname.contains("DefaultRegionForParallelWorld")) continue;
202  const G4String& rname = writer->GenerateName(tname,(*rstore)[i]);
203  rlist = new G4GDMLAuxListType();
204  G4GDMLAuxStructType raux = {"Region", rname, "", rlist};
205  std::vector<G4LogicalVolume*>::iterator rlvol_iter
206  = (*rstore)[i]->GetRootLogicalVolumeIterator();
207  for (size_t j=0; j<(*rstore)[i]->GetNumberOfRootVolumes(); ++j)
208  {
209  G4LogicalVolume* rlvol = *rlvol_iter;
210  if(reflFactory->IsReflected(rlvol)) continue;
211  G4String vname = writer->GenerateName(rlvol->GetName(), rlvol);
212  if (!storeReferences) { reader->StripName(vname); }
213  G4GDMLAuxStructType rsubaux = {"volume", vname, "", 0};
214  rlist->push_back(rsubaux);
215  rlvol_iter++;
216  }
217  G4double gam_cut = (*rstore)[i]->GetProductionCuts()->GetProductionCut("gamma");
218  G4GDMLAuxStructType caux1 = {"gamcut", eval.ConvertToString(gam_cut), "mm", 0};
219  rlist->push_back(caux1);
220  G4double e_cut = (*rstore)[i]->GetProductionCuts()->GetProductionCut("e-");
221  G4GDMLAuxStructType caux2 = {"ecut", eval.ConvertToString(e_cut), "mm", 0};
222  rlist->push_back(caux2);
223  G4double pos_cut = (*rstore)[i]->GetProductionCuts()->GetProductionCut("e+");
224  G4GDMLAuxStructType caux3 = {"poscut", eval.ConvertToString(pos_cut), "mm", 0};
225  rlist->push_back(caux3);
226  G4double p_cut = (*rstore)[i]->GetProductionCuts()->GetProductionCut("proton");
227  G4GDMLAuxStructType caux4 = {"pcut", eval.ConvertToString(p_cut), "mm", 0};
228  rlist->push_back(caux4);
229  if ((*rstore)[i]->GetUserLimits())
230  {
231  const G4Track fake_trk;
232  ullist = new G4GDMLAuxListType();
233  const G4String& utype = (*rstore)[i]->GetUserLimits()->GetType();
234  G4GDMLAuxStructType uaux = {"ulimits", utype, "mm", ullist};
235  G4double max_step = (*rstore)[i]->GetUserLimits()->GetMaxAllowedStep(fake_trk);
236  G4GDMLAuxStructType ulaux1 = {"ustepMax", eval.ConvertToString(max_step), "mm", 0};
237  ullist->push_back(ulaux1);
238  G4double max_trk = (*rstore)[i]->GetUserLimits()->GetUserMaxTrackLength(fake_trk);
239  G4GDMLAuxStructType ulaux2 = {"utrakMax", eval.ConvertToString(max_trk), "mm", 0};
240  ullist->push_back(ulaux2);
241  G4double max_time = (*rstore)[i]->GetUserLimits()->GetUserMaxTime(fake_trk);
242  G4GDMLAuxStructType ulaux3 = {"utimeMax", eval.ConvertToString(max_time), "mm", 0};
243  ullist->push_back(ulaux3);
244  G4double min_ekin = (*rstore)[i]->GetUserLimits()->GetUserMinEkine(fake_trk);
245  G4GDMLAuxStructType ulaux4 = {"uekinMin", eval.ConvertToString(min_ekin), "mm", 0};
246  ullist->push_back(ulaux4);
247  G4double min_rng = (*rstore)[i]->GetUserLimits()->GetUserMinRange(fake_trk);
248  G4GDMLAuxStructType ulaux5 = {"urangMin", eval.ConvertToString(min_rng), "mm", 0};
249  ullist->push_back(ulaux5);
250  rlist->push_back(uaux);
251  }
252  AddAuxiliary(raux);
253  }
254 }
const XML_Char * name
Definition: expat.h:151
void AddRootLogicalVolume(G4LogicalVolume *lv)
Definition: G4Region.cc:290
Definition: xmlparse.cc:187
void SetProductionCut(G4double cut, G4int index=-1)
G4bool IsReflected(G4LogicalVolume *lv) const
const G4GDMLAuxListType * GetAuxList() const
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4ReflectionFactory * Instance()
static G4RegionStore * GetInstance()
static G4double GetValueOf(const G4String &)
G4String ConvertToString(G4int ival)
bool G4bool
Definition: G4Types.hh:79
static G4LogicalVolumeStore * GetInstance()
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void StripName(G4String &) const
Definition: G4GDMLRead.cc:103
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41
static G4String GetCategory(const G4String &)
G4bool contains(const std::string &) const
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
G4LogicalVolume * GetReflectedLV(G4LogicalVolume *lv) const
void SetProductionCuts(G4ProductionCuts *cut)
void SetUserLimits(G4UserLimits *ul)
void AddAuxiliary(G4GDMLAuxStructType myaux)
const G4String & GetName() const
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4bool IsConstituent(G4LogicalVolume *lv) const
G4double Evaluate(const G4String &)