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