Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4GDMLWriteStructure Class Reference

#include <G4GDMLWriteStructure.hh>

+ Inheritance diagram for G4GDMLWriteStructure:

Public Member Functions

 G4GDMLWriteStructure ()
 
virtual ~G4GDMLWriteStructure ()
 
virtual void StructureWrite (xercesc::DOMElement *)
 
void AddVolumeAuxiliary (G4GDMLAuxStructType myaux, const G4LogicalVolume *const)
 
void SetEnergyCutsExport (G4bool)
 
void SetSDExport (G4bool)
 
G4int GetMaxExportLevel () const
 
void SetMaxExportLevel (G4int)
 
- Public Member Functions inherited from G4GDMLWriteParamvol
virtual void ParamvolWrite (xercesc::DOMElement *, const G4VPhysicalVolume *const)
 
virtual void ParamvolAlgorithmWrite (xercesc::DOMElement *paramvolElement, const G4VPhysicalVolume *const paramvol)
 
- Public Member Functions inherited from G4GDMLWriteSetup
virtual void SetupWrite (xercesc::DOMElement *, const G4LogicalVolume *const)
 
- Public Member Functions inherited from G4GDMLWriteSolids
virtual void AddSolid (const G4VSolid *const)
 
virtual void SolidsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteMaterials
void AddIsotope (const G4Isotope *const)
 
void AddElement (const G4Element *const)
 
void AddMaterial (const G4Material *const)
 
virtual void MaterialsWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWriteDefine
G4ThreeVector GetAngles (const G4RotationMatrix &)
 
void ScaleWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
 
void RotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void PositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void FirstrotationWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
 
void FirstpositionWrite (xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
 
void AddPosition (const G4String &name, const G4ThreeVector &pos)
 
virtual void DefineWrite (xercesc::DOMElement *)
 
- Public Member Functions inherited from G4GDMLWrite
G4Transform3D Write (const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
 
void AddModule (const G4VPhysicalVolume *const topVol)
 
void AddModule (const G4int depth)
 
void AddAuxiliary (G4GDMLAuxStructType myaux)
 
void SetOutputFileOverwrite (G4bool flag)
 
virtual void DefineWrite (xercesc::DOMElement *)=0
 
virtual void MaterialsWrite (xercesc::DOMElement *)=0
 
virtual void SolidsWrite (xercesc::DOMElement *)=0
 
virtual void StructureWrite (xercesc::DOMElement *)=0
 
virtual G4Transform3D TraverseVolumeTree (const G4LogicalVolume *const, const G4int)=0
 
virtual void SurfacesWrite ()=0
 
virtual void SetupWrite (xercesc::DOMElement *, const G4LogicalVolume *const)=0
 
virtual void ExtensionWrite (xercesc::DOMElement *)
 
virtual void UserinfoWrite (xercesc::DOMElement *)
 
virtual void AddExtension (xercesc::DOMElement *, const G4LogicalVolume *const)
 
G4String GenerateName (const G4String &, const void *const)
 

Protected Member Functions

void DivisionvolWrite (xercesc::DOMElement *, const G4PVDivision *const)
 
void PhysvolWrite (xercesc::DOMElement *, const G4VPhysicalVolume *const topVol, const G4Transform3D &transform, const G4String &moduleName)
 
void ReplicavolWrite (xercesc::DOMElement *, const G4VPhysicalVolume *const)
 
void AssemblyWrite (xercesc::DOMElement *, const int assemblyID)
 
G4Transform3D TraverseVolumeTree (const G4LogicalVolume *const topVol, const G4int depth)
 
void SurfacesWrite ()
 
void BorderSurfaceCache (const G4LogicalBorderSurface *const)
 
void SkinSurfaceCache (const G4LogicalSkinSurface *const)
 
const G4LogicalBorderSurfaceGetBorderSurface (const G4VPhysicalVolume *const)
 
const G4LogicalSkinSurfaceGetSkinSurface (const G4LogicalVolume *const)
 
G4bool FindOpticalSurface (const G4SurfaceProperty *)
 
void ExportEnergyCuts (const G4LogicalVolume *const)
 
void ExportSD (const G4LogicalVolume *const)
 
- Protected Member Functions inherited from G4GDMLWriteParamvol
 G4GDMLWriteParamvol ()
 
virtual ~G4GDMLWriteParamvol ()
 
void Box_dimensionsWrite (xercesc::DOMElement *, const G4Box *const)
 
void Trd_dimensionsWrite (xercesc::DOMElement *, const G4Trd *const)
 
void Trap_dimensionsWrite (xercesc::DOMElement *, const G4Trap *const)
 
void Tube_dimensionsWrite (xercesc::DOMElement *, const G4Tubs *const)
 
void Cone_dimensionsWrite (xercesc::DOMElement *, const G4Cons *const)
 
void Sphere_dimensionsWrite (xercesc::DOMElement *, const G4Sphere *const)
 
void Orb_dimensionsWrite (xercesc::DOMElement *, const G4Orb *const)
 
void Torus_dimensionsWrite (xercesc::DOMElement *, const G4Torus *const)
 
void Ellipsoid_dimensionsWrite (xercesc::DOMElement *, const G4Ellipsoid *const)
 
void Para_dimensionsWrite (xercesc::DOMElement *, const G4Para *const)
 
void Hype_dimensionsWrite (xercesc::DOMElement *, const G4Hype *const)
 
void Polycone_dimensionsWrite (xercesc::DOMElement *, const G4Polycone *const)
 
void Polyhedra_dimensionsWrite (xercesc::DOMElement *, const G4Polyhedra *const)
 
void ParametersWrite (xercesc::DOMElement *, const G4VPhysicalVolume *const, const G4int &)
 
- Protected Member Functions inherited from G4GDMLWriteSetup
 G4GDMLWriteSetup ()
 
virtual ~G4GDMLWriteSetup ()
 
- Protected Member Functions inherited from G4GDMLWriteSolids
 G4GDMLWriteSolids ()
 
virtual ~G4GDMLWriteSolids ()
 
void MultiUnionWrite (xercesc::DOMElement *solElement, const G4MultiUnion *const)
 
void BooleanWrite (xercesc::DOMElement *, const G4BooleanSolid *const)
 
void ScaledWrite (xercesc::DOMElement *, const G4ScaledSolid *const)
 
void BoxWrite (xercesc::DOMElement *, const G4Box *const)
 
void ConeWrite (xercesc::DOMElement *, const G4Cons *const)
 
void ElconeWrite (xercesc::DOMElement *, const G4EllipticalCone *const)
 
void EllipsoidWrite (xercesc::DOMElement *, const G4Ellipsoid *const)
 
void EltubeWrite (xercesc::DOMElement *, const G4EllipticalTube *const)
 
void XtruWrite (xercesc::DOMElement *, const G4ExtrudedSolid *const)
 
void HypeWrite (xercesc::DOMElement *, const G4Hype *const)
 
void OrbWrite (xercesc::DOMElement *, const G4Orb *const)
 
void ParaWrite (xercesc::DOMElement *, const G4Para *const)
 
void ParaboloidWrite (xercesc::DOMElement *, const G4Paraboloid *const)
 
void PolyconeWrite (xercesc::DOMElement *, const G4Polycone *const)
 
void GenericPolyconeWrite (xercesc::DOMElement *, const G4GenericPolycone *const)
 
void PolyhedraWrite (xercesc::DOMElement *, const G4Polyhedra *const)
 
void SphereWrite (xercesc::DOMElement *, const G4Sphere *const)
 
void TessellatedWrite (xercesc::DOMElement *, const G4TessellatedSolid *const)
 
void TetWrite (xercesc::DOMElement *, const G4Tet *const)
 
void TorusWrite (xercesc::DOMElement *, const G4Torus *const)
 
void GenTrapWrite (xercesc::DOMElement *, const G4GenericTrap *const)
 
void TrapWrite (xercesc::DOMElement *, const G4Trap *const)
 
void TrdWrite (xercesc::DOMElement *, const G4Trd *const)
 
void TubeWrite (xercesc::DOMElement *, const G4Tubs *const)
 
void CutTubeWrite (xercesc::DOMElement *, const G4CutTubs *const)
 
void TwistedboxWrite (xercesc::DOMElement *, const G4TwistedBox *const)
 
void TwistedtrapWrite (xercesc::DOMElement *, const G4TwistedTrap *const)
 
void TwistedtrdWrite (xercesc::DOMElement *, const G4TwistedTrd *const)
 
void TwistedtubsWrite (xercesc::DOMElement *, const G4TwistedTubs *const)
 
void ZplaneWrite (xercesc::DOMElement *, const G4double &, const G4double &, const G4double &)
 
void RZPointWrite (xercesc::DOMElement *, const G4double &, const G4double &)
 
void OpticalSurfaceWrite (xercesc::DOMElement *, const G4OpticalSurface *const)
 
void PropertyWrite (xercesc::DOMElement *, const G4OpticalSurface *const)
 
- Protected Member Functions inherited from G4GDMLWriteMaterials
 G4GDMLWriteMaterials ()
 
virtual ~G4GDMLWriteMaterials ()
 
void AtomWrite (xercesc::DOMElement *, const G4double &)
 
void DWrite (xercesc::DOMElement *, const G4double &)
 
void PWrite (xercesc::DOMElement *, const G4double &)
 
void TWrite (xercesc::DOMElement *, const G4double &)
 
void MEEWrite (xercesc::DOMElement *, const G4double &)
 
void IsotopeWrite (const G4Isotope *const)
 
void ElementWrite (const G4Element *const)
 
void MaterialWrite (const G4Material *const)
 
void PropertyWrite (xercesc::DOMElement *, const G4Material *const)
 
void PropertyVectorWrite (const G4String &, const G4PhysicsOrderedFreeVector *const)
 
void PropertyConstWrite (const G4String &, const G4double, const G4MaterialPropertiesTable *)
 
- Protected Member Functions inherited from G4GDMLWriteDefine
 G4GDMLWriteDefine ()
 
virtual ~G4GDMLWriteDefine ()
 
void Scale_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Rotation_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
void Position_vectorWrite (xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
 
- Protected Member Functions inherited from G4GDMLWrite
 G4GDMLWrite ()
 
virtual ~G4GDMLWrite ()
 
VolumeMapType & VolumeMap ()
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4String &)
 
xercesc::DOMAttr * NewAttribute (const G4String &, const G4double &)
 
xercesc::DOMElement * NewElement (const G4String &)
 
G4String Modularize (const G4VPhysicalVolume *const topvol, const G4int depth)
 
void AddAuxInfo (G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
 
G4bool FileExists (const G4String &) const
 
PhysVolumeMapType & PvolumeMap ()
 
DepthMapType & DepthMap ()
 

Protected Attributes

xercesc::DOMElement * structureElement = nullptr
 
std::vector< xercesc::DOMElement * > borderElementVec
 
std::vector< xercesc::DOMElement * > skinElementVec
 
std::map< const G4LogicalVolume *, G4GDMLAuxListTypeauxmap
 
- Protected Attributes inherited from G4GDMLWriteSolids
std::vector< const G4VSolid * > solidList
 
xercesc::DOMElement * solidsElement = nullptr
 
- Protected Attributes inherited from G4GDMLWriteMaterials
std::vector< const G4Isotope * > isotopeList
 
std::vector< const G4Element * > elementList
 
std::vector< const G4Material * > materialList
 
std::vector< const G4PhysicsOrderedFreeVector * > propertyList
 
xercesc::DOMElement * materialsElement = nullptr
 
- Protected Attributes inherited from G4GDMLWriteDefine
xercesc::DOMElement * defineElement = nullptr
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc = nullptr
 
xercesc::DOMElement * extElement = nullptr
 
xercesc::DOMElement * userinfoElement = nullptr
 
XMLCh tempStr [10000]
 
G4GDMLAuxListType auxList
 
G4bool overwriteOutputFile = false
 

Additional Inherited Members

- Static Public Member Functions inherited from G4GDMLWrite
static void SetAddPointerToName (G4bool)
 
- Static Protected Attributes inherited from G4GDMLWriteSolids
static const G4int maxTransforms = 8
 
- Static Protected Attributes inherited from G4GDMLWriteDefine
static const G4double kRelativePrecision = DBL_EPSILON
 
static const G4double kAngularPrecision = DBL_EPSILON
 
static const G4double kLinearPrecision = DBL_EPSILON
 
- Static Protected Attributes inherited from G4GDMLWrite
static G4bool addPointerToName = true
 

Detailed Description

Definition at line 51 of file G4GDMLWriteStructure.hh.

Constructor & Destructor Documentation

◆ G4GDMLWriteStructure()

G4GDMLWriteStructure::G4GDMLWriteStructure ( )

Definition at line 61 of file G4GDMLWriteStructure.cc.

63 , maxLevel(INT_MAX)
64{
65 reflFactory = G4ReflectionFactory::Instance();
66}
static G4ReflectionFactory * Instance()
#define INT_MAX
Definition: templates.hh:90

◆ ~G4GDMLWriteStructure()

G4GDMLWriteStructure::~G4GDMLWriteStructure ( )
virtual

Definition at line 69 of file G4GDMLWriteStructure.cc.

70{
71}

Member Function Documentation

◆ AddVolumeAuxiliary()

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

Definition at line 822 of file G4GDMLWriteStructure.cc.

824{
825 auto pos = auxmap.find(lvol);
826
827 if(pos == auxmap.cend())
828 {
829 auxmap[lvol] = G4GDMLAuxListType();
830 }
831
832 auxmap[lvol].push_back(myaux);
833}
std::vector< G4GDMLAuxStructType > G4GDMLAuxListType
std::map< const G4LogicalVolume *, G4GDMLAuxListType > auxmap

Referenced by ExportEnergyCuts(), and ExportSD().

◆ AssemblyWrite()

void G4GDMLWriteStructure::AssemblyWrite ( xercesc::DOMElement *  volumeElement,
const int  assemblyID 
)
protected

Definition at line 264 of file G4GDMLWriteStructure.cc.

266{
268 G4AssemblyVolume* myassembly = assemblies->GetAssembly(assemblyID);
269
270 xercesc::DOMElement* assemblyElement = NewElement("assembly");
271 G4String name = "Assembly_" + std::to_string(assemblyID);
272
273 assemblyElement->setAttributeNode(NewAttribute("name", name));
274
275 auto vit = myassembly->GetTripletsIterator();
276
277 G4int depth = 0;
278 const G4String ModuleName;
279
280 for(std::size_t i5 = 0; i5 < myassembly->TotalTriplets(); ++i5)
281 {
282 TraverseVolumeTree((*vit).GetVolume(), depth + 1);
283
284 const G4ThreeVector rot = GetAngles((*vit).GetRotation()->inverse());
285 const G4ThreeVector pos = (*vit).GetTranslation();
286
287 const G4String pname =
288 GenerateName((*vit).GetVolume()->GetName() + "_pv", &(*vit));
289
290 xercesc::DOMElement* physvolElement = NewElement("physvol");
291 physvolElement->setAttributeNode(NewAttribute("name", pname));
292
293 assemblyElement->appendChild(physvolElement);
294
295 const G4String volumeref =
296 GenerateName((*vit).GetVolume()->GetName(), (*vit).GetVolume());
297
298 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
299 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
300 physvolElement->appendChild(volumerefElement);
301
302 if(std::fabs(pos.x()) > kLinearPrecision ||
303 std::fabs(pos.y()) > kLinearPrecision ||
304 std::fabs(pos.z()) > kLinearPrecision)
305 {
306 PositionWrite(physvolElement, "Position_" + std::to_string(i5), pos);
307 }
308
309 if(std::fabs(rot.x()) > kAngularPrecision ||
310 std::fabs(rot.y()) > kAngularPrecision ||
311 std::fabs(rot.z()) > kAngularPrecision)
312 {
313 RotationWrite(physvolElement, "Rotation_" + std::to_string(i5), rot);
314 }
315 ++vit;
316 }
317
318 volumeElement->appendChild(assemblyElement);
319}
int G4int
Definition: G4Types.hh:85
double z() const
double x() const
double y() const
static G4AssemblyStore * GetInstance()
G4AssemblyVolume * GetAssembly(unsigned int id, G4bool verbose=true) const
std::vector< G4AssemblyTriplet >::iterator GetTripletsIterator()
std::size_t TotalTriplets() const
static const G4double kLinearPrecision
void RotationWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &rot)
static const G4double kAngularPrecision
void PositionWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &pos)
G4ThreeVector GetAngles(const G4RotationMatrix &)
G4Transform3D TraverseVolumeTree(const G4LogicalVolume *const topVol, const G4int depth)
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:181
G4String GenerateName(const G4String &, const void *const)
Definition: G4GDMLWrite.cc:132
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:155
const char * name(G4int ptype)

Referenced by TraverseVolumeTree().

◆ BorderSurfaceCache()

void G4GDMLWriteStructure::BorderSurfaceCache ( const G4LogicalBorderSurface * const  bsurf)
protected

Definition at line 322 of file G4GDMLWriteStructure.cc.

324{
325 if(bsurf == nullptr)
326 {
327 return;
328 }
329
330 const G4SurfaceProperty* psurf = bsurf->GetSurfaceProperty();
331
332 // Generate the new element for border-surface
333 //
334 const G4String& bsname = GenerateName(bsurf->GetName(), bsurf);
335 const G4String& psname = GenerateName(psurf->GetName(), psurf);
336 xercesc::DOMElement* borderElement = NewElement("bordersurface");
337 borderElement->setAttributeNode(NewAttribute("name", bsname));
338 borderElement->setAttributeNode(NewAttribute("surfaceproperty", psname));
339
340 const G4String volumeref1 =
341 GenerateName(bsurf->GetVolume1()->GetName(), bsurf->GetVolume1());
342 const G4String volumeref2 =
343 GenerateName(bsurf->GetVolume2()->GetName(), bsurf->GetVolume2());
344 xercesc::DOMElement* volumerefElement1 = NewElement("physvolref");
345 xercesc::DOMElement* volumerefElement2 = NewElement("physvolref");
346 volumerefElement1->setAttributeNode(NewAttribute("ref", volumeref1));
347 volumerefElement2->setAttributeNode(NewAttribute("ref", volumeref2));
348 borderElement->appendChild(volumerefElement1);
349 borderElement->appendChild(volumerefElement2);
350
351 if(FindOpticalSurface(psurf))
352 {
353 const G4OpticalSurface* opsurf =
354 dynamic_cast<const G4OpticalSurface*>(psurf);
355 if(opsurf == nullptr)
356 {
357 G4Exception("G4GDMLWriteStructure::BorderSurfaceCache()", "InvalidSetup",
358 FatalException, "No optical surface found!");
359 return;
360 }
362 }
363
364 borderElementVec.push_back(borderElement);
365}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
xercesc::DOMElement * solidsElement
void OpticalSurfaceWrite(xercesc::DOMElement *, const G4OpticalSurface *const)
std::vector< xercesc::DOMElement * > borderElementVec
G4bool FindOpticalSurface(const G4SurfaceProperty *)
const G4VPhysicalVolume * GetVolume2() const
const G4VPhysicalVolume * GetVolume1() const
const G4String & GetName() const
G4SurfaceProperty * GetSurfaceProperty() const
const G4String & GetName() const
const G4String & GetName() const

Referenced by GetBorderSurface().

◆ DivisionvolWrite()

void G4GDMLWriteStructure::DivisionvolWrite ( xercesc::DOMElement *  volumeElement,
const G4PVDivision * const  divisionvol 
)
protected

Definition at line 74 of file G4GDMLWriteStructure.cc.

76{
77 EAxis axis = kUndefined;
78 G4int number = 0;
79 G4double width = 0.0;
80 G4double offset = 0.0;
81 G4bool consuming = false;
82
83 divisionvol->GetReplicationData(axis, number, width, offset, consuming);
84 axis = divisionvol->GetDivisionAxis();
85
86 G4String unitString("mm");
87 G4String axisString("kUndefined");
88 if(axis == kXAxis)
89 {
90 axisString = "kXAxis";
91 }
92 else if(axis == kYAxis)
93 {
94 axisString = "kYAxis";
95 }
96 else if(axis == kZAxis)
97 {
98 axisString = "kZAxis";
99 }
100 else if(axis == kRho)
101 {
102 axisString = "kRho";
103 }
104 else if(axis == kPhi)
105 {
106 axisString = "kPhi";
107 unitString = "rad";
108 }
109
110 const G4String name = GenerateName(divisionvol->GetName(), divisionvol);
111 const G4String volumeref =
112 GenerateName(divisionvol->GetLogicalVolume()->GetName(),
113 divisionvol->GetLogicalVolume());
114
115 xercesc::DOMElement* divisionvolElement = NewElement("divisionvol");
116 divisionvolElement->setAttributeNode(NewAttribute("axis", axisString));
117 divisionvolElement->setAttributeNode(NewAttribute("number", number));
118 divisionvolElement->setAttributeNode(NewAttribute("width", width));
119 divisionvolElement->setAttributeNode(NewAttribute("offset", offset));
120 divisionvolElement->setAttributeNode(NewAttribute("unit", unitString));
121 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
122 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
123 divisionvolElement->appendChild(volumerefElement);
124 volumeElement->appendChild(divisionvolElement);
125}
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
const G4String & GetName() const
EAxis GetDivisionAxis() const
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const
G4LogicalVolume * GetLogicalVolume() const
EAxis
Definition: geomdefs.hh:54
@ kPhi
Definition: geomdefs.hh:60
@ kYAxis
Definition: geomdefs.hh:56
@ kXAxis
Definition: geomdefs.hh:55
@ kZAxis
Definition: geomdefs.hh:57
@ kUndefined
Definition: geomdefs.hh:61
@ kRho
Definition: geomdefs.hh:58

Referenced by TraverseVolumeTree().

◆ ExportEnergyCuts()

void G4GDMLWriteStructure::ExportEnergyCuts ( const G4LogicalVolume * const  lvol)
protected

Definition at line 842 of file G4GDMLWriteStructure.cc.

843{
844 G4GDMLEvaluator eval;
845 G4ProductionCuts* pcuts = lvol->GetRegion()->GetProductionCuts();
847 G4Gamma* gamma = G4Gamma::Gamma();
851
852 G4double gamma_cut = ctab->ConvertRangeToEnergy(
853 gamma, lvol->GetMaterial(), pcuts->GetProductionCut("gamma"));
854 G4double eminus_cut = ctab->ConvertRangeToEnergy(
855 eminus, lvol->GetMaterial(), pcuts->GetProductionCut("e-"));
856 G4double eplus_cut = ctab->ConvertRangeToEnergy(
857 eplus, lvol->GetMaterial(), pcuts->GetProductionCut("e+"));
858 G4double proton_cut = ctab->ConvertRangeToEnergy(
859 proton, lvol->GetMaterial(), pcuts->GetProductionCut("proton"));
860
861 G4GDMLAuxStructType gammainfo = { "gammaECut",
862 eval.ConvertToString(gamma_cut), "MeV", 0 };
863 G4GDMLAuxStructType eminusinfo = { "electronECut",
864 eval.ConvertToString(eminus_cut), "MeV",
865 0 };
866 G4GDMLAuxStructType eplusinfo = { "positronECut",
867 eval.ConvertToString(eplus_cut), "MeV", 0 };
868 G4GDMLAuxStructType protinfo = { "protonECut",
869 eval.ConvertToString(proton_cut), "MeV", 0 };
870
871 AddVolumeAuxiliary(gammainfo, lvol);
872 AddVolumeAuxiliary(eminusinfo, lvol);
873 AddVolumeAuxiliary(eplusinfo, lvol);
874 AddVolumeAuxiliary(protinfo, lvol);
875}
static G4Electron * Electron()
Definition: G4Electron.cc:93
G4String ConvertToString(G4int ival)
void AddVolumeAuxiliary(G4GDMLAuxStructType myaux, const G4LogicalVolume *const)
static G4Gamma * Gamma()
Definition: G4Gamma.cc:85
G4Region * GetRegion() const
G4Material * GetMaterial() const
static G4Positron * Positron()
Definition: G4Positron.cc:93
static G4ProductionCutsTable * GetProductionCutsTable()
G4double ConvertRangeToEnergy(const G4ParticleDefinition *particle, const G4Material *material, G4double range)
G4double GetProductionCut(G4int index) const
static G4Proton * Proton()
Definition: G4Proton.cc:92
G4ProductionCuts * GetProductionCuts() const

Referenced by TraverseVolumeTree().

◆ ExportSD()

void G4GDMLWriteStructure::ExportSD ( const G4LogicalVolume * const  lvol)
protected

Definition at line 884 of file G4GDMLWriteStructure.cc.

885{
887
888 if(sd != nullptr)
889 {
890 G4String SDname = sd->GetName();
891
892 G4GDMLAuxStructType SDinfo = { "SensDet", SDname, "", 0 };
893 AddVolumeAuxiliary(SDinfo, lvol);
894 }
895}
G4VSensitiveDetector * GetSensitiveDetector() const

Referenced by TraverseVolumeTree().

◆ FindOpticalSurface()

G4bool G4GDMLWriteStructure::FindOpticalSurface ( const G4SurfaceProperty psurf)
protected

Definition at line 409 of file G4GDMLWriteStructure.cc.

410{
411 const G4OpticalSurface* osurf = dynamic_cast<const G4OpticalSurface*>(psurf);
412 auto pos = std::find(opt_vec.cbegin(), opt_vec.cend(), osurf);
413 if(pos != opt_vec.cend())
414 {
415 return false;
416 } // item already created!
417
418 opt_vec.push_back(osurf); // cache it for future reference
419 return true;
420}

Referenced by BorderSurfaceCache(), and SkinSurfaceCache().

◆ GetBorderSurface()

const G4LogicalBorderSurface * G4GDMLWriteStructure::GetBorderSurface ( const G4VPhysicalVolume * const  pvol)
protected

Definition at line 445 of file G4GDMLWriteStructure.cc.

447{
448 G4LogicalBorderSurface* surf = nullptr;
450 if(nsurf)
451 {
452 const G4LogicalBorderSurfaceTable* btable =
454 for(auto pos = btable->cbegin(); pos != btable->cend(); ++pos)
455 {
456 if(pvol == pos->first.first) // just the first in the couple
457 { // could be enough?
458 surf = pos->second; // break;
459 BorderSurfaceCache(surf);
460 }
461 }
462 }
463 return surf;
464}
std::map< std::pair< const G4VPhysicalVolume *, const G4VPhysicalVolume * >, G4LogicalBorderSurface * > G4LogicalBorderSurfaceTable
void BorderSurfaceCache(const G4LogicalBorderSurface *const)
static size_t GetNumberOfBorderSurfaces()
static const G4LogicalBorderSurfaceTable * GetSurfaceTable()

Referenced by TraverseVolumeTree().

◆ GetMaxExportLevel()

G4int G4GDMLWriteStructure::GetMaxExportLevel ( ) const

Definition at line 898 of file G4GDMLWriteStructure.cc.

899{
900 return maxLevel;
901}

◆ GetSkinSurface()

const G4LogicalSkinSurface * G4GDMLWriteStructure::GetSkinSurface ( const G4LogicalVolume * const  lvol)
protected

Definition at line 423 of file G4GDMLWriteStructure.cc.

425{
426 G4LogicalSkinSurface* surf = 0;
428 if(nsurf)
429 {
430 const G4LogicalSkinSurfaceTable* stable =
432 for(auto pos = stable->cbegin(); pos != stable->cend(); ++pos)
433 {
434 if(lvol == (*pos)->GetLogicalVolume())
435 {
436 surf = *pos;
437 break;
438 }
439 }
440 }
441 return surf;
442}
std::vector< G4LogicalSkinSurface * > G4LogicalSkinSurfaceTable
static const G4LogicalSkinSurfaceTable * GetSurfaceTable()
static std::size_t GetNumberOfSkinSurfaces()

Referenced by TraverseVolumeTree().

◆ PhysvolWrite()

void G4GDMLWriteStructure::PhysvolWrite ( xercesc::DOMElement *  volumeElement,
const G4VPhysicalVolume *const  topVol,
const G4Transform3D transform,
const G4String moduleName 
)
protected

Definition at line 128 of file G4GDMLWriteStructure.cc.

132{
133 HepGeom::Scale3D scale;
134 HepGeom::Rotate3D rotate;
135 HepGeom::Translate3D translate;
136
137 T.getDecomposition(scale, rotate, translate);
138
139 const G4ThreeVector scl(scale(0, 0), scale(1, 1), scale(2, 2));
140 const G4ThreeVector rot = GetAngles(rotate.getRotation());
141 const G4ThreeVector pos = T.getTranslation();
142
143 const G4String name = GenerateName(physvol->GetName(), physvol);
144 const G4int copynumber = physvol->GetCopyNo();
145
146 xercesc::DOMElement* physvolElement = NewElement("physvol");
147 physvolElement->setAttributeNode(NewAttribute("name", name));
148 if(copynumber)
149 {
150 physvolElement->setAttributeNode(NewAttribute("copynumber", copynumber));
151 }
152
153 volumeElement->appendChild(physvolElement);
154
155 G4LogicalVolume* lv;
156 // Is it reflected?
157 if(reflFactory->IsReflected(physvol->GetLogicalVolume()))
158 {
159 lv = reflFactory->GetConstituentLV(physvol->GetLogicalVolume());
160 }
161 else
162 {
163 lv = physvol->GetLogicalVolume();
164 }
165
166 const G4String volumeref = GenerateName(lv->GetName(), lv);
167
168 if(ModuleName.empty())
169 {
170 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
171 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
172 physvolElement->appendChild(volumerefElement);
173 }
174 else
175 {
176 xercesc::DOMElement* fileElement = NewElement("file");
177 fileElement->setAttributeNode(NewAttribute("name", ModuleName));
178 fileElement->setAttributeNode(NewAttribute("volname", volumeref));
179 physvolElement->appendChild(fileElement);
180 }
181
182 if(std::fabs(pos.x()) > kLinearPrecision ||
183 std::fabs(pos.y()) > kLinearPrecision ||
184 std::fabs(pos.z()) > kLinearPrecision)
185 {
186 PositionWrite(physvolElement, name + "_pos", pos);
187 }
188 if(std::fabs(rot.x()) > kAngularPrecision ||
189 std::fabs(rot.y()) > kAngularPrecision ||
190 std::fabs(rot.z()) > kAngularPrecision)
191 {
192 RotationWrite(physvolElement, name + "_rot", rot);
193 }
194 if(std::fabs(scl.x() - 1.0) > kRelativePrecision ||
195 std::fabs(scl.y() - 1.0) > kRelativePrecision ||
196 std::fabs(scl.z() - 1.0) > kRelativePrecision)
197 {
198 ScaleWrite(physvolElement, name + "_scl", scl);
199 }
200}
void ScaleWrite(xercesc::DOMElement *element, const G4String &name, const G4ThreeVector &scl)
static const G4double kRelativePrecision
G4bool IsReflected(G4LogicalVolume *lv) const
G4LogicalVolume * GetConstituentLV(G4LogicalVolume *reflLV) const
CLHEP::HepRotation getRotation() const
void getDecomposition(Scale3D &scale, Rotate3D &rotation, Translate3D &translation) const
Definition: Transform3D.cc:173

Referenced by TraverseVolumeTree().

◆ ReplicavolWrite()

void G4GDMLWriteStructure::ReplicavolWrite ( xercesc::DOMElement *  volumeElement,
const G4VPhysicalVolume * const  replicavol 
)
protected

Definition at line 203 of file G4GDMLWriteStructure.cc.

205{
206 EAxis axis = kUndefined;
207 G4int number = 0;
208 G4double width = 0.0;
209 G4double offset = 0.0;
210 G4bool consuming = false;
211 G4String unitString("mm");
212
213 replicavol->GetReplicationData(axis, number, width, offset, consuming);
214
215 const G4String volumeref = GenerateName(
216 replicavol->GetLogicalVolume()->GetName(), replicavol->GetLogicalVolume());
217
218 xercesc::DOMElement* replicavolElement = NewElement("replicavol");
219 replicavolElement->setAttributeNode(NewAttribute("number", number));
220 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
221 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
222 replicavolElement->appendChild(volumerefElement);
223 xercesc::DOMElement* replicateElement = NewElement("replicate_along_axis");
224 replicavolElement->appendChild(replicateElement);
225
226 xercesc::DOMElement* dirElement = NewElement("direction");
227 if(axis == kXAxis)
228 {
229 dirElement->setAttributeNode(NewAttribute("x", "1"));
230 }
231 else if(axis == kYAxis)
232 {
233 dirElement->setAttributeNode(NewAttribute("y", "1"));
234 }
235 else if(axis == kZAxis)
236 {
237 dirElement->setAttributeNode(NewAttribute("z", "1"));
238 }
239 else if(axis == kRho)
240 {
241 dirElement->setAttributeNode(NewAttribute("rho", "1"));
242 }
243 else if(axis == kPhi)
244 {
245 dirElement->setAttributeNode(NewAttribute("phi", "1"));
246 unitString = "rad";
247 }
248 replicateElement->appendChild(dirElement);
249
250 xercesc::DOMElement* widthElement = NewElement("width");
251 widthElement->setAttributeNode(NewAttribute("value", width));
252 widthElement->setAttributeNode(NewAttribute("unit", unitString));
253 replicateElement->appendChild(widthElement);
254
255 xercesc::DOMElement* offsetElement = NewElement("offset");
256 offsetElement->setAttributeNode(NewAttribute("value", offset));
257 offsetElement->setAttributeNode(NewAttribute("unit", unitString));
258 replicateElement->appendChild(offsetElement);
259
260 volumeElement->appendChild(replicavolElement);
261}
virtual void GetReplicationData(EAxis &axis, G4int &nReplicas, G4double &width, G4double &offset, G4bool &consuming) const =0

Referenced by TraverseVolumeTree().

◆ SetEnergyCutsExport()

void G4GDMLWriteStructure::SetEnergyCutsExport ( G4bool  fcuts)

Definition at line 836 of file G4GDMLWriteStructure.cc.

837{
838 cexport = fcuts;
839}

◆ SetMaxExportLevel()

void G4GDMLWriteStructure::SetMaxExportLevel ( G4int  level)

Definition at line 904 of file G4GDMLWriteStructure.cc.

905{
906 if(level <= 0)
907 {
908 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()", "InvalidSetup",
909 FatalException, "Levels to export must be greater than zero!");
910 return;
911 }
912 maxLevel = level;
913 levelNo = 0;
914}

◆ SetSDExport()

void G4GDMLWriteStructure::SetSDExport ( G4bool  fsd)

Definition at line 878 of file G4GDMLWriteStructure.cc.

879{
880 sdexport = fsd;
881}

◆ SkinSurfaceCache()

void G4GDMLWriteStructure::SkinSurfaceCache ( const G4LogicalSkinSurface * const  ssurf)
protected

Definition at line 368 of file G4GDMLWriteStructure.cc.

370{
371 if(ssurf == nullptr)
372 {
373 return;
374 }
375
376 const G4SurfaceProperty* psurf = ssurf->GetSurfaceProperty();
377
378 // Generate the new element for border-surface
379 //
380 const G4String& ssname = GenerateName(ssurf->GetName(), ssurf);
381 const G4String& psname = GenerateName(psurf->GetName(), psurf);
382 xercesc::DOMElement* skinElement = NewElement("skinsurface");
383 skinElement->setAttributeNode(NewAttribute("name", ssname));
384 skinElement->setAttributeNode(NewAttribute("surfaceproperty", psname));
385
386 const G4String volumeref = GenerateName(ssurf->GetLogicalVolume()->GetName(),
387 ssurf->GetLogicalVolume());
388 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
389 volumerefElement->setAttributeNode(NewAttribute("ref", volumeref));
390 skinElement->appendChild(volumerefElement);
391
392 if(FindOpticalSurface(psurf))
393 {
394 const G4OpticalSurface* opsurf =
395 dynamic_cast<const G4OpticalSurface*>(psurf);
396 if(opsurf == nullptr)
397 {
398 G4Exception("G4GDMLWriteStructure::SkinSurfaceCache()", "InvalidSetup",
399 FatalException, "No optical surface found!");
400 return;
401 }
403 }
404
405 skinElementVec.push_back(skinElement);
406}
std::vector< xercesc::DOMElement * > skinElementVec
const G4LogicalVolume * GetLogicalVolume() const

Referenced by TraverseVolumeTree().

◆ StructureWrite()

void G4GDMLWriteStructure::StructureWrite ( xercesc::DOMElement *  gdmlElement)
virtual

Implements G4GDMLWrite.

Definition at line 485 of file G4GDMLWriteStructure.cc.

486{
487#ifdef G4VERBOSE
488 G4cout << "G4GDML: Writing structure..." << G4endl;
489#endif
490
491 // filling the list of phys volumes that are parts of assemblies
492
494
495 for(auto it = assemblies->cbegin(); it != assemblies->cend(); ++it)
496 {
497 auto vit = (*it)->GetVolumesIterator();
498
499 for(std::size_t i5 = 0; i5 < (*it)->TotalImprintedVolumes(); ++i5)
500 {
501 G4String pvname = (*vit)->GetName();
502 std::size_t pos = pvname.find("_impr_") + 6;
503 G4String impID = pvname.substr(pos);
504
505 pos = impID.find("_");
506 impID = impID.substr(0, pos);
507
508 assemblyVolMap[*vit] = (*it)->GetAssemblyID();
509 imprintsMap[*vit] = std::atoi(impID.c_str());
510 ++vit;
511 }
512 }
513
514 structureElement = NewElement("structure");
515 gdmlElement->appendChild(structureElement);
516}
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
xercesc::DOMElement * structureElement

◆ SurfacesWrite()

void G4GDMLWriteStructure::SurfacesWrite ( )
protectedvirtual

Implements G4GDMLWrite.

Definition at line 467 of file G4GDMLWriteStructure.cc.

468{
469#ifdef G4VERBOSE
470 G4cout << "G4GDML: Writing surfaces..." << G4endl;
471#endif
472 for(auto pos = skinElementVec.cbegin();
473 pos != skinElementVec.cend(); ++pos)
474 {
475 structureElement->appendChild(*pos);
476 }
477 for(auto pos = borderElementVec.cbegin();
478 pos != borderElementVec.cend(); ++pos)
479 {
480 structureElement->appendChild(*pos);
481 }
482}

◆ TraverseVolumeTree()

G4Transform3D G4GDMLWriteStructure::TraverseVolumeTree ( const G4LogicalVolume *const  topVol,
const G4int  depth 
)
protectedvirtual

Implements G4GDMLWrite.

Definition at line 519 of file G4GDMLWriteStructure.cc.

521{
522 if(VolumeMap().find(volumePtr) != VolumeMap().cend())
523 {
524 return VolumeMap()[volumePtr]; // Volume is already processed
525 }
526
527 G4VSolid* solidPtr = volumePtr->GetSolid();
528 G4Transform3D R, invR;
529 G4int trans = 0;
530
531 std::map<const G4LogicalVolume*, G4GDMLAuxListType>::iterator auxiter;
532
533 ++levelNo;
534
535 while(true) // Solve possible displacement/reflection
536 { // of the referenced solid!
537 if(trans > maxTransforms)
538 {
539 G4String ErrorMessage = "Referenced solid in volume '" +
540 volumePtr->GetName() +
541 "' was displaced/reflected too many times!";
542 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()", "InvalidSetup",
543 FatalException, ErrorMessage);
544 }
545
546 if(G4ReflectedSolid* refl = dynamic_cast<G4ReflectedSolid*>(solidPtr))
547 {
548 R = R * refl->GetTransform3D();
549 solidPtr = refl->GetConstituentMovedSolid();
550 ++trans;
551 continue;
552 }
553
554 if(G4DisplacedSolid* disp = dynamic_cast<G4DisplacedSolid*>(solidPtr))
555 {
556 R = R * G4Transform3D(disp->GetObjectRotation(),
557 disp->GetObjectTranslation());
558 solidPtr = disp->GetConstituentMovedSolid();
559 ++trans;
560 continue;
561 }
562
563 break;
564 }
565
566 // check if it is a reflected volume
567 G4LogicalVolume* tmplv = const_cast<G4LogicalVolume*>(volumePtr);
568
569 if(reflFactory->IsReflected(tmplv))
570 {
571 tmplv = reflFactory->GetConstituentLV(tmplv);
572 if(VolumeMap().find(tmplv) != VolumeMap().cend())
573 {
574 return R; // Volume is already processed
575 }
576 }
577
578 // Only compute the inverse when necessary!
579 //
580 if(trans > 0)
581 {
582 invR = R.inverse();
583 }
584
585 const G4String name = GenerateName(tmplv->GetName(), tmplv);
586
587 G4String materialref = "NULL";
588
589 if(volumePtr->GetMaterial())
590 {
591 materialref = GenerateName(volumePtr->GetMaterial()->GetName(),
592 volumePtr->GetMaterial());
593 }
594
595 const G4String solidref = GenerateName(solidPtr->GetName(), solidPtr);
596
597 xercesc::DOMElement* volumeElement = NewElement("volume");
598 volumeElement->setAttributeNode(NewAttribute("name", name));
599 xercesc::DOMElement* materialrefElement = NewElement("materialref");
600 materialrefElement->setAttributeNode(NewAttribute("ref", materialref));
601 volumeElement->appendChild(materialrefElement);
602 xercesc::DOMElement* solidrefElement = NewElement("solidref");
603 solidrefElement->setAttributeNode(NewAttribute("ref", solidref));
604 volumeElement->appendChild(solidrefElement);
605
606 G4int daughterCount = volumePtr->GetNoDaughters();
607
608 if(levelNo == maxLevel) // Stop exporting if reached levels limit
609 {
610 daughterCount = 0;
611 }
612
613 std::vector<G4int> addedImprints;
614
615 for(G4int i = 0; i < daughterCount; ++i) // Traverse all the children!
616 {
617 const G4VPhysicalVolume* const physvol = volumePtr->GetDaughter(i);
618 const G4String ModuleName = Modularize(physvol, depth);
619
620 G4Transform3D daughterR;
621
622 if(ModuleName.empty()) // Check if subtree requested to be
623 { // a separate module!
624 daughterR = TraverseVolumeTree(physvol->GetLogicalVolume(), depth + 1);
625 }
626 else
627 {
629 daughterR = writer.Write(ModuleName, physvol->GetLogicalVolume(),
630 SchemaLocation, depth + 1);
631 }
632
633 if(const G4PVDivision* const divisionvol =
634 dynamic_cast<const G4PVDivision*>(physvol)) // Is it division?
635 {
636 if(!G4Transform3D::Identity.isNear(invR * daughterR, kRelativePrecision))
637 {
638 G4String ErrorMessage = "Division volume in '" + name +
639 "' can not be related to reflected solid!";
640 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
641 "InvalidSetup", FatalException, ErrorMessage);
642 }
643 DivisionvolWrite(volumeElement, divisionvol);
644 }
645 else
646 {
647 if(physvol->IsParameterised()) // Is it a paramvol?
648 {
649 if(!G4Transform3D::Identity.isNear(invR * daughterR,
651 {
652 G4String ErrorMessage = "Parameterised volume in '" + name +
653 "' can not be related to reflected solid!";
654 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
655 "InvalidSetup", FatalException, ErrorMessage);
656 }
657 ParamvolWrite(volumeElement, physvol);
658 }
659 else
660 {
661 if(physvol->IsReplicated()) // Is it a replicavol?
662 {
663 if(!G4Transform3D::Identity.isNear(invR * daughterR,
665 {
666 G4String ErrorMessage = "Replica volume in '" + name +
667 "' can not be related to reflected solid!";
668 G4Exception("G4GDMLWriteStructure::TraverseVolumeTree()",
669 "InvalidSetup", FatalException, ErrorMessage);
670 }
671 ReplicavolWrite(volumeElement, physvol);
672 }
673 else // Is it a physvol or an assembly?
674 {
675 if(assemblyVolMap.find(physvol) != assemblyVolMap.cend())
676 {
677 G4int assemblyID = assemblyVolMap[physvol];
678
679 G4String assemblyref = "Assembly_" + std::to_string(assemblyID);
680
681 // here I need to retrieve the imprint ID
682
683 G4int imprintID = imprintsMap[physvol];
684
685 // there are 2 steps:
686 //
687 // 1) add assembly to the structure if that has not yet been done
688 // (but after the constituents volumes have been added)
689 //
690
691 if(std::find(addedAssemblies.cbegin(), addedAssemblies.cend(),
692 assemblyID) == addedAssemblies.cend())
693 {
694 AssemblyWrite(structureElement, assemblyID);
695 addedAssemblies.push_back(assemblyID);
696 }
697
698 // 2) add the assembly (as physical volume) to the mother volume
699 // (but only once), using it's original position and rotation.
700 //
701
702 // here I need a check if assembly has been already added to the
703 // mother volume
704 if(std::find(addedImprints.cbegin(), addedImprints.cend(),
705 imprintID) == addedImprints.cend())
706 {
707 G4String imprintname = "Imprint_" + std::to_string(imprintID);
708 imprintname = GenerateName(imprintname, physvol);
709
710 // I need to get those two from the container of imprints from
711 // the assembly I have the imprint ID, I need to get pos and rot
712 //
714 ->GetAssembly(assemblyID)
715 ->GetImprintTransformation(imprintID);
716
717 HepGeom::Scale3D scale;
718 HepGeom::Rotate3D rotate;
719 HepGeom::Translate3D translate;
720
721 transf.getDecomposition(scale, rotate, translate);
722
723 const G4ThreeVector scl(scale(0, 0), scale(1, 1), scale(2, 2));
724 const G4ThreeVector rot =
725 GetAngles(rotate.getRotation().inverse());
726 const G4ThreeVector pos = transf.getTranslation();
727
728 // here I need a normal physvol referencing to my assemblyref
729
730 xercesc::DOMElement* physvolElement = NewElement("physvol");
731 physvolElement->setAttributeNode(
732 NewAttribute("name", imprintname));
733
734 xercesc::DOMElement* volumerefElement = NewElement("volumeref");
735 volumerefElement->setAttributeNode(
736 NewAttribute("ref", assemblyref));
737 physvolElement->appendChild(volumerefElement);
738
739 if(std::fabs(pos.x()) > kLinearPrecision ||
740 std::fabs(pos.y()) > kLinearPrecision ||
741 std::fabs(pos.z()) > kLinearPrecision)
742 {
743 PositionWrite(physvolElement, imprintname + "_pos", pos);
744 }
745 if(std::fabs(rot.x()) > kAngularPrecision ||
746 std::fabs(rot.y()) > kAngularPrecision ||
747 std::fabs(rot.z()) > kAngularPrecision)
748 {
749 RotationWrite(physvolElement, imprintname + "_rot", rot);
750 }
751 if(std::fabs(scl.x() - 1.0) > kRelativePrecision ||
752 std::fabs(scl.y() - 1.0) > kRelativePrecision ||
753 std::fabs(scl.z() - 1.0) > kRelativePrecision)
754 {
755 ScaleWrite(physvolElement, name + "_scl", scl);
756 }
757
758 volumeElement->appendChild(physvolElement);
759 //
760 addedImprints.push_back(imprintID);
761 }
762 }
763 else // not part of assembly, so a normal physical volume
764 {
766 if(physvol->GetFrameRotation() != nullptr)
767 {
768 rot = *(physvol->GetFrameRotation());
769 }
770 G4Transform3D P(rot, physvol->GetObjectTranslation());
771
772 PhysvolWrite(volumeElement, physvol, invR * P * daughterR,
773 ModuleName);
774 }
775 }
776 }
777 }
778 // BorderSurfaceCache(GetBorderSurface(physvol));
779 GetBorderSurface(physvol);
780 }
781
782 if(cexport)
783 {
784 ExportEnergyCuts(volumePtr);
785 }
786 // Add optional energy cuts
787
788 if(sdexport)
789 {
790 ExportSD(volumePtr);
791 }
792 // Add optional SDs
793
794 // Here write the auxiliary info
795 //
796 auxiter = auxmap.find(volumePtr);
797 if(auxiter != auxmap.cend())
798 {
799 AddAuxInfo(&(auxiter->second), volumeElement);
800 }
801
802 structureElement->appendChild(volumeElement);
803 // Append the volume AFTER traversing the children so that
804 // the order of volumes will be correct!
805
806 VolumeMap()[tmplv] = R;
807
808 AddExtension(volumeElement, volumePtr);
809 // Add any possible user defined extension attached to a volume
810
811 AddMaterial(volumePtr->GetMaterial());
812 // Add the involved materials and solids!
813
814 AddSolid(solidPtr);
815
817
818 return R;
819}
HepGeom::Transform3D G4Transform3D
HepRotation inverse() const
G4Transform3D & GetImprintTransformation(unsigned int imprintID)
void AddMaterial(const G4Material *const)
virtual void ParamvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
virtual void AddSolid(const G4VSolid *const)
static const G4int maxTransforms
const G4LogicalBorderSurface * GetBorderSurface(const G4VPhysicalVolume *const)
void AssemblyWrite(xercesc::DOMElement *, const int assemblyID)
void SkinSurfaceCache(const G4LogicalSkinSurface *const)
const G4LogicalSkinSurface * GetSkinSurface(const G4LogicalVolume *const)
void PhysvolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const topVol, const G4Transform3D &transform, const G4String &moduleName)
void ExportEnergyCuts(const G4LogicalVolume *const)
void ExportSD(const G4LogicalVolume *const)
void DivisionvolWrite(xercesc::DOMElement *, const G4PVDivision *const)
void ReplicavolWrite(xercesc::DOMElement *, const G4VPhysicalVolume *const)
G4String SchemaLocation
Definition: G4GDMLWrite.hh:130
G4String Modularize(const G4VPhysicalVolume *const topvol, const G4int depth)
Definition: G4GDMLWrite.cc:358
virtual void AddExtension(xercesc::DOMElement *, const G4LogicalVolume *const)
Definition: G4GDMLWrite.cc:80
G4Transform3D Write(const G4String &filename, const G4LogicalVolume *const topLog, const G4String &schemaPath, const G4int depth, G4bool storeReferences=true)
Definition: G4GDMLWrite.cc:188
void AddAuxInfo(G4GDMLAuxListType *auxInfoList, xercesc::DOMElement *element)
Definition: G4GDMLWrite.cc:94
VolumeMapType & VolumeMap()
Definition: G4GDMLWrite.cc:60
virtual G4bool IsReplicated() const =0
const G4RotationMatrix * GetFrameRotation() const
G4ThreeVector GetObjectTranslation() const
virtual G4bool IsParameterised() const =0
G4String GetName() const
static DLL_API const Transform3D Identity
Definition: Transform3D.h:196
CLHEP::Hep3Vector getTranslation() const
Transform3D inverse() const
Definition: Transform3D.cc:141

Referenced by AssemblyWrite(), and TraverseVolumeTree().

Member Data Documentation

◆ auxmap

std::map<const G4LogicalVolume*, G4GDMLAuxListType> G4GDMLWriteStructure::auxmap
protected

Definition at line 94 of file G4GDMLWriteStructure.hh.

Referenced by AddVolumeAuxiliary(), and TraverseVolumeTree().

◆ borderElementVec

std::vector<xercesc::DOMElement*> G4GDMLWriteStructure::borderElementVec
protected

Definition at line 92 of file G4GDMLWriteStructure.hh.

Referenced by BorderSurfaceCache(), and SurfacesWrite().

◆ skinElementVec

std::vector<xercesc::DOMElement*> G4GDMLWriteStructure::skinElementVec
protected

Definition at line 93 of file G4GDMLWriteStructure.hh.

Referenced by SkinSurfaceCache(), and SurfacesWrite().

◆ structureElement

xercesc::DOMElement* G4GDMLWriteStructure::structureElement = nullptr
protected

Definition at line 91 of file G4GDMLWriteStructure.hh.

Referenced by StructureWrite(), SurfacesWrite(), and TraverseVolumeTree().


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