Geant4  10.02.p03
G4GDMLParser Class Reference

#include <G4GDMLParser.hh>

Collaboration diagram for G4GDMLParser:

Public Member Functions

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

Private Member Functions

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

Private Attributes

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

Detailed Description

Definition at line 55 of file G4GDMLParser.hh.

Constructor & Destructor Documentation

◆ G4GDMLParser() [1/3]

G4GDMLParser::G4GDMLParser ( )

Definition at line 44 of file G4GDMLParser.cc.

45  : rlist(0), ullist(0), urcode(false), uwcode(false), strip(true), rexp(false)
46 {
49  messenger = new G4GDMLMessenger(this);
50 
52 }
void Initialize()
Definition: errprop.cc:101
G4GDMLReadStructure * reader
G4GDMLAuxListType * rlist
G4GDMLAuxListType * ullist
G4GDMLMessenger * messenger
G4GDMLWriteStructure * writer
Here is the call graph for this function:

◆ G4GDMLParser() [2/3]

G4GDMLParser::G4GDMLParser ( G4GDMLReadStructure extr)

Definition at line 54 of file G4GDMLParser.cc.

55  : rlist(0), ullist(0), urcode(true), uwcode(false), strip(true), rexp(false)
56 {
57  reader = extr;
59  messenger = new G4GDMLMessenger(this);
60 
62 }
void Initialize()
Definition: errprop.cc:101
G4GDMLReadStructure * reader
G4GDMLAuxListType * rlist
G4GDMLAuxListType * ullist
G4GDMLMessenger * messenger
G4GDMLWriteStructure * writer
Here is the call graph for this function:

◆ G4GDMLParser() [3/3]

G4GDMLParser::G4GDMLParser ( G4GDMLReadStructure extr,
G4GDMLWriteStructure extw 
)

Definition at line 64 of file G4GDMLParser.cc.

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 
73 }
void Initialize()
Definition: errprop.cc:101
G4GDMLReadStructure * reader
G4GDMLAuxListType * rlist
G4GDMLAuxListType * ullist
G4GDMLMessenger * messenger
G4GDMLWriteStructure * writer
Here is the call graph for this function:

◆ ~G4GDMLParser()

G4GDMLParser::~G4GDMLParser ( )

Definition at line 75 of file G4GDMLParser.cc.

76 {
77  xercesc::XMLPlatformUtils::Terminate();
78  if (!urcode) { delete reader; }
79  if (!uwcode) { delete writer; delete ullist; delete rlist; }
80 
81  delete messenger;
82 }
G4GDMLReadStructure * reader
G4GDMLAuxListType * rlist
G4GDMLAuxListType * ullist
G4GDMLMessenger * messenger
G4GDMLWriteStructure * writer

Member Function Documentation

◆ AddAuxiliary()

void G4GDMLParser::AddAuxiliary ( G4GDMLAuxStructType  myaux)
inline
Here is the caller graph for this function:

◆ AddModule() [1/2]

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

◆ AddModule() [2/2]

void G4GDMLParser::AddModule ( const G4int  depth)
inline

◆ AddVolumeAuxiliary()

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

◆ Clear()

void G4GDMLParser::Clear ( )
inline
Here is the caller graph for this function:

◆ ExportRegions()

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

Definition at line 194 of file G4GDMLParser.cc.

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 }
G4GDMLEvaluator eval
G4bool contains(const std::string &) const
G4bool IsReflected(G4LogicalVolume *lv) const
static G4ReflectionFactory * Instance()
static G4RegionStore * GetInstance()
G4GDMLReadStructure * reader
G4GDMLAuxListType * rlist
const G4String & GetName() const
G4String ConvertToString(G4int ival)
G4GDMLAuxListType * ullist
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:123
void StripName(G4String &) const
Definition: G4GDMLRead.cc:103
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
void AddAuxiliary(G4GDMLAuxStructType myaux)
double G4double
Definition: G4Types.hh:76
G4GDMLWriteStructure * writer
Here is the call graph for this function:

◆ GetAuxList()

const G4GDMLAuxListType* G4GDMLParser::GetAuxList ( ) const
inline
Here is the caller graph for this function:

◆ GetAuxMap()

const G4GDMLAuxMapType* G4GDMLParser::GetAuxMap ( ) const
inline
Here is the caller graph for this function:

◆ GetConstant()

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

◆ GetMatrix()

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

◆ GetPosition()

G4ThreeVector G4GDMLParser::GetPosition ( const G4String name) const
inline
Here is the caller graph for this function:

◆ GetQuantity()

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

◆ GetRotation()

G4ThreeVector G4GDMLParser::GetRotation ( const G4String name) const
inline
Here is the caller graph for this function:

◆ 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
Here is the caller graph for this function:

◆ GetWorldVolume()

G4VPhysicalVolume* G4GDMLParser::GetWorldVolume ( const G4String setupName = "Default") const
inline
Here is the caller graph for this function:

◆ ImportRegions()

void G4GDMLParser::ImportRegions ( )
private

Definition at line 84 of file G4GDMLParser.cc.

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 }
const G4GDMLAuxListType * GetAuxList() const
void AddRootLogicalVolume(G4LogicalVolume *lv)
Definition: G4Region.cc:290
G4String name
Definition: TRTMaterials.hh:40
G4GDMLEvaluator eval
Definition: xmlparse.cc:187
G4bool IsConstituent(G4LogicalVolume *lv) const
void SetProductionCut(G4double cut, G4int index=-1)
G4bool contains(const std::string &) const
static G4ReflectionFactory * Instance()
G4GDMLReadStructure * reader
static G4double GetValueOf(const G4String &)
G4LogicalVolume * GetVolume(const G4String &name, G4bool verbose=true) const
static G4LogicalVolumeStore * GetInstance()
G4LogicalVolume * GetReflectedLV(G4LogicalVolume *lv) const
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 &)
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
void SetProductionCuts(G4ProductionCuts *cut)
void SetUserLimits(G4UserLimits *ul)
double G4double
Definition: G4Types.hh:76
#define DBL_MAX
Definition: templates.hh:83
G4double Evaluate(const G4String &)
Here is the call graph for this function:

◆ IsValid()

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

◆ ParseST()

G4LogicalVolume* G4GDMLParser::ParseST ( const G4String name,
G4Material medium,
G4Material solid 
)
inline
Here is the caller graph for this function:

◆ Read()

void G4GDMLParser::Read ( const G4String filename,
G4bool  Validate = true 
)
inline
Here is the caller graph for this function:

◆ ReadModule()

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

◆ SetAddPointerToName()

void G4GDMLParser::SetAddPointerToName ( G4bool  set)
inline

◆ SetEnergyCutsExport()

void G4GDMLParser::SetEnergyCutsExport ( G4bool  )
inline
Here is the caller graph for this function:

◆ SetOverlapCheck()

void G4GDMLParser::SetOverlapCheck ( G4bool  )
inline

◆ SetRegionExport()

void G4GDMLParser::SetRegionExport ( G4bool  )
inline
Here is the caller graph for this function:

◆ SetStripFlag()

void G4GDMLParser::SetStripFlag ( G4bool  )
inline

◆ StripNamePointers()

void G4GDMLParser::StripNamePointers ( ) const
inline

◆ Write() [1/2]

void G4GDMLParser::Write ( const G4String filename,
const G4VPhysicalVolume pvol = 0,
G4bool  storeReferences = true,
const G4String SchemaLocation = G4GDML_DEFAULT_SCHEMALOCATION 
)
inline
Here is the caller graph for this function:

◆ Write() [2/2]

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

Member Data Documentation

◆ eval

G4GDMLEvaluator G4GDMLParser::eval
private

Definition at line 146 of file G4GDMLParser.hh.

◆ messenger

G4GDMLMessenger* G4GDMLParser::messenger
private

Definition at line 150 of file G4GDMLParser.hh.

◆ reader

G4GDMLReadStructure* G4GDMLParser::reader
private

Definition at line 147 of file G4GDMLParser.hh.

◆ rexp

G4bool G4GDMLParser::rexp
private

Definition at line 151 of file G4GDMLParser.hh.

◆ rlist

G4GDMLAuxListType* G4GDMLParser::rlist
private

Definition at line 149 of file G4GDMLParser.hh.

◆ strip

G4bool G4GDMLParser::strip
private

Definition at line 151 of file G4GDMLParser.hh.

◆ ullist

G4GDMLAuxListType * G4GDMLParser::ullist
private

Definition at line 149 of file G4GDMLParser.hh.

◆ urcode

G4bool G4GDMLParser::urcode
private

Definition at line 151 of file G4GDMLParser.hh.

◆ uwcode

G4bool G4GDMLParser::uwcode
private

Definition at line 151 of file G4GDMLParser.hh.

◆ writer

G4GDMLWriteStructure* G4GDMLParser::writer
private

Definition at line 148 of file G4GDMLParser.hh.


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