Geant4.10
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | Static Protected Attributes
G4GDMLWriteDefine Class Reference

#include <G4GDMLWriteDefine.hh>

Inheritance diagram for G4GDMLWriteDefine:
G4GDMLWrite G4GDMLWriteMaterials G4GDMLWriteSolids G4GDMLWriteSetup G4GDMLWriteParamvol G4GDMLWriteStructure G03ColorWriter

Public Member Functions

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)
 
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 AddExtension (xercesc::DOMElement *, const G4LogicalVolume *const)
 

Protected Member Functions

 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 ()
 
G4String GenerateName (const G4String &, const void *const)
 
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)
 
G4bool FileExists (const G4String &) const
 
PhysVolumeMapType & PvolumeMap ()
 
DepthMapType & DepthMap ()
 

Protected Attributes

xercesc::DOMElement * defineElement
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc
 
xercesc::DOMElement * extElement
 
XMLCh tempStr [10000]
 

Static Protected Attributes

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
 

Additional Inherited Members

- Static Public Member Functions inherited from G4GDMLWrite
static void SetAddPointerToName (G4bool)
 

Detailed Description

Definition at line 49 of file G4GDMLWriteDefine.hh.

Constructor & Destructor Documentation

G4GDMLWriteDefine::G4GDMLWriteDefine ( )
protected

Definition at line 42 of file G4GDMLWriteDefine.cc.

44 {
45 }
xercesc::DOMElement * defineElement
G4GDMLWriteDefine::~G4GDMLWriteDefine ( )
protectedvirtual

Definition at line 47 of file G4GDMLWriteDefine.cc.

48 {
49 }

Member Function Documentation

void G4GDMLWriteDefine::AddPosition ( const G4String name,
const G4ThreeVector pos 
)
inline

Definition at line 70 of file G4GDMLWriteDefine.hh.

References defineElement, and Position_vectorWrite().

Referenced by G4GDMLWriteSolids::TessellatedWrite(), and G4GDMLWriteSolids::TetWrite().

71  { Position_vectorWrite(defineElement,"position",name,pos); }
xercesc::DOMElement * defineElement
void Position_vectorWrite(xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
void G4GDMLWriteDefine::DefineWrite ( xercesc::DOMElement *  element)
virtual

Implements G4GDMLWrite.

Definition at line 126 of file G4GDMLWriteDefine.cc.

References defineElement, G4cout, G4endl, and G4GDMLWrite::NewElement().

127 {
128  G4cout << "G4GDML: Writing definitions..." << G4endl;
129 
130  defineElement = NewElement("define");
131  element->appendChild(defineElement);
132 }
xercesc::DOMElement * defineElement
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:127
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61
void G4GDMLWriteDefine::FirstpositionWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector pos 
)
inline

Definition at line 67 of file G4GDMLWriteDefine.hh.

References Position_vectorWrite().

Referenced by G4GDMLWriteSolids::BooleanWrite().

69  { Position_vectorWrite(element,"firstposition",name,pos); }
void Position_vectorWrite(xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
void G4GDMLWriteDefine::FirstrotationWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector rot 
)
inline

Definition at line 64 of file G4GDMLWriteDefine.hh.

References Rotation_vectorWrite().

Referenced by G4GDMLWriteSolids::BooleanWrite().

66  { Rotation_vectorWrite(element,"firstrotation",name,rot); }
void Rotation_vectorWrite(xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
G4ThreeVector G4GDMLWriteDefine::GetAngles ( const G4RotationMatrix mat)

Definition at line 51 of file G4GDMLWriteDefine.cc.

References kRelativePrecision, test::x, CLHEP::HepRotation::xx(), CLHEP::HepRotation::yx(), CLHEP::HepRotation::yy(), CLHEP::HepRotation::yz(), z, CLHEP::HepRotation::zx(), CLHEP::HepRotation::zy(), and CLHEP::HepRotation::zz().

Referenced by G4GDMLWriteSolids::BooleanWrite(), G4GDMLWriteParamvol::ParametersWrite(), and G4GDMLWriteStructure::PhysvolWrite().

52 {
53  G4double x,y,z;
54 
55  const G4double cosb = std::sqrt(mat.xx()*mat.xx()+mat.yx()*mat.yx());
56 
57  if (cosb > kRelativePrecision)
58  {
59  x = std::atan2(mat.zy(),mat.zz());
60  y = std::atan2(-mat.zx(),cosb);
61  z = std::atan2(mat.yx(),mat.xx());
62  }
63  else
64  {
65  x = std::atan2(-mat.yz(),mat.yy());
66  y = std::atan2(-mat.zx(),cosb);
67  z = 0.0;
68  }
69 
70  return G4ThreeVector(x,y,z);
71 }
double xx() const
CLHEP::Hep3Vector G4ThreeVector
G4double z
Definition: TRTMaterials.hh:39
double yy() const
double zx() const
double yz() const
static const G4double kRelativePrecision
double zy() const
double yx() const
double G4double
Definition: G4Types.hh:76
double zz() const
void G4GDMLWriteDefine::Position_vectorWrite ( xercesc::DOMElement *  element,
const G4String tag,
const G4String name,
const G4ThreeVector pos 
)
protected

Definition at line 110 of file G4GDMLWriteDefine.cc.

References kLinearPrecision, python.hepunit::mm, G4GDMLWrite::NewAttribute(), G4GDMLWrite::NewElement(), test::x, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), z, and CLHEP::Hep3Vector::z().

Referenced by AddPosition(), FirstpositionWrite(), and PositionWrite().

112 {
113  const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x();
114  const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y();
115  const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z();
116 
117  xercesc::DOMElement* positionElement = NewElement(tag);
118  positionElement->setAttributeNode(NewAttribute("name",name));
119  positionElement->setAttributeNode(NewAttribute("x",x/mm));
120  positionElement->setAttributeNode(NewAttribute("y",y/mm));
121  positionElement->setAttributeNode(NewAttribute("z",z/mm));
122  positionElement->setAttributeNode(NewAttribute("unit","mm"));
123  element->appendChild(positionElement);
124 }
double x() const
G4double z
Definition: TRTMaterials.hh:39
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:127
double z() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:103
static const G4double kLinearPrecision
double y() const
double G4double
Definition: G4Types.hh:76
void G4GDMLWriteDefine::PositionWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector pos 
)
inline

Definition at line 61 of file G4GDMLWriteDefine.hh.

References Position_vectorWrite().

Referenced by G4GDMLWriteSolids::BooleanWrite(), G4GDMLWriteParamvol::ParametersWrite(), and G4GDMLWriteStructure::PhysvolWrite().

63  { Position_vectorWrite(element,"position",name,pos); }
void Position_vectorWrite(xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
void G4GDMLWriteDefine::Rotation_vectorWrite ( xercesc::DOMElement *  element,
const G4String tag,
const G4String name,
const G4ThreeVector rot 
)
protected

Definition at line 93 of file G4GDMLWriteDefine.cc.

References python.hepunit::degree, kAngularPrecision, G4GDMLWrite::NewAttribute(), G4GDMLWrite::NewElement(), test::x, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), z, and CLHEP::Hep3Vector::z().

Referenced by FirstrotationWrite(), and RotationWrite().

95 {
96  const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x();
97  const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y();
98  const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z();
99 
100  xercesc::DOMElement* rotationElement = NewElement(tag);
101  rotationElement->setAttributeNode(NewAttribute("name",name));
102  rotationElement->setAttributeNode(NewAttribute("x",x/degree));
103  rotationElement->setAttributeNode(NewAttribute("y",y/degree));
104  rotationElement->setAttributeNode(NewAttribute("z",z/degree));
105  rotationElement->setAttributeNode(NewAttribute("unit","deg"));
106  element->appendChild(rotationElement);
107 }
static const G4double kAngularPrecision
double x() const
G4double z
Definition: TRTMaterials.hh:39
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:127
double z() const
tuple degree
Definition: hepunit.py:69
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:103
double y() const
double G4double
Definition: G4Types.hh:76
void G4GDMLWriteDefine::RotationWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector rot 
)
inline

Definition at line 58 of file G4GDMLWriteDefine.hh.

References Rotation_vectorWrite().

Referenced by G4GDMLWriteSolids::BooleanWrite(), G4GDMLWriteParamvol::ParametersWrite(), and G4GDMLWriteStructure::PhysvolWrite().

60  { Rotation_vectorWrite(element,"rotation",name,rot); }
void Rotation_vectorWrite(xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)
void G4GDMLWriteDefine::Scale_vectorWrite ( xercesc::DOMElement *  element,
const G4String tag,
const G4String name,
const G4ThreeVector scl 
)
protected

Definition at line 74 of file G4GDMLWriteDefine.cc.

References kRelativePrecision, G4GDMLWrite::NewAttribute(), G4GDMLWrite::NewElement(), test::x, CLHEP::Hep3Vector::x(), CLHEP::Hep3Vector::y(), z, and CLHEP::Hep3Vector::z().

Referenced by ScaleWrite().

76 {
77  const G4double x = (std::fabs(scl.x()-1.0) < kRelativePrecision)
78  ? 1.0 : scl.x();
79  const G4double y = (std::fabs(scl.y()-1.0) < kRelativePrecision)
80  ? 1.0 : scl.y();
81  const G4double z = (std::fabs(scl.z()-1.0) < kRelativePrecision)
82  ? 1.0 : scl.z();
83 
84  xercesc::DOMElement* scaleElement = NewElement(tag);
85  scaleElement->setAttributeNode(NewAttribute("name",name));
86  scaleElement->setAttributeNode(NewAttribute("x",x));
87  scaleElement->setAttributeNode(NewAttribute("y",y));
88  scaleElement->setAttributeNode(NewAttribute("z",z));
89  element->appendChild(scaleElement);
90 }
double x() const
G4double z
Definition: TRTMaterials.hh:39
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:127
double z() const
static const G4double kRelativePrecision
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:103
double y() const
double G4double
Definition: G4Types.hh:76
void G4GDMLWriteDefine::ScaleWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector scl 
)
inline

Definition at line 55 of file G4GDMLWriteDefine.hh.

References Scale_vectorWrite().

Referenced by G4GDMLWriteStructure::PhysvolWrite().

57  { Scale_vectorWrite(element,"scale",name,scl); }
void Scale_vectorWrite(xercesc::DOMElement *, const G4String &, const G4String &, const G4ThreeVector &)

Field Documentation

xercesc::DOMElement* G4GDMLWriteDefine::defineElement
protected
const G4double G4GDMLWriteDefine::kAngularPrecision = DBL_EPSILON
staticprotected
const G4double G4GDMLWriteDefine::kLinearPrecision = DBL_EPSILON
staticprotected
const G4double G4GDMLWriteDefine::kRelativePrecision = DBL_EPSILON
staticprotected

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