Geant4  10.03.p01
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Pages
G4GDMLWriteDefine Class Reference

#include <G4GDMLWriteDefine.hh>

Inheritance diagram for G4GDMLWriteDefine:
Collaboration diagram for G4GDMLWriteDefine:

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)
 
void AddAuxiliary (G4GDMLAuxStructType myaux)
 
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

 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 * defineElement
 
- Protected Attributes inherited from G4GDMLWrite
G4String SchemaLocation
 
xercesc::DOMDocument * doc
 
xercesc::DOMElement * extElement
 
xercesc::DOMElement * userinfoElement
 
XMLCh tempStr [10000]
 
G4GDMLAuxListType auxList
 

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.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::DefineWrite ( xercesc::DOMElement *  element)
virtual

Implements G4GDMLWrite.

Definition at line 131 of file G4GDMLWriteDefine.cc.

132 {
133  G4cout << "G4GDML: Writing definitions..." << G4endl;
134 
135  defineElement = NewElement("define");
136  element->appendChild(defineElement);
137 }
xercesc::DOMElement * defineElement
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
G4GLOB_DLL std::ostream G4cout
#define G4endl
Definition: G4ios.hh:61

Here is the call graph for this function:

void G4GDMLWriteDefine::FirstpositionWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector pos 
)
inline

Definition at line 67 of file G4GDMLWriteDefine.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::FirstrotationWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector rot 
)
inline

Definition at line 64 of file G4GDMLWriteDefine.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

G4ThreeVector G4GDMLWriteDefine::GetAngles ( const G4RotationMatrix mtx)

Definition at line 51 of file G4GDMLWriteDefine.cc.

52 {
53  G4double x,y,z;
54  G4RotationMatrix mat = mtx;
55  mat.rectify(); // Rectify matrix from possible roundoff errors
56 
57  // Direction of rotation given by left-hand rule; clockwise rotation
58 
59  static const G4double kMatrixPrecision = 10E-10;
60  const G4double cosb = std::sqrt(mtx.xx()*mtx.xx()+mtx.yx()*mtx.yx());
61 
62  if (cosb > kMatrixPrecision)
63  {
64  x = std::atan2(mtx.zy(),mtx.zz());
65  y = std::atan2(-mtx.zx(),cosb);
66  z = std::atan2(mtx.yx(),mtx.xx());
67  }
68  else
69  {
70  x = std::atan2(-mtx.yz(),mtx.yy());
71  y = std::atan2(-mtx.zx(),cosb);
72  z = 0.0;
73  }
74 
75  return G4ThreeVector(x,y,z);
76 }
double xx() const
CLHEP::Hep3Vector G4ThreeVector
double yy() const
double zx() const
double yz() const
double zy() const
double yx() const
double G4double
Definition: G4Types.hh:76
double zz() const

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::Position_vectorWrite ( xercesc::DOMElement *  element,
const G4String tag,
const G4String name,
const G4ThreeVector pos 
)
protected

Definition at line 115 of file G4GDMLWriteDefine.cc.

117 {
118  const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x();
119  const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y();
120  const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z();
121 
122  xercesc::DOMElement* positionElement = NewElement(tag);
123  positionElement->setAttributeNode(NewAttribute("name",name));
124  positionElement->setAttributeNode(NewAttribute("x",x/mm));
125  positionElement->setAttributeNode(NewAttribute("y",y/mm));
126  positionElement->setAttributeNode(NewAttribute("z",z/mm));
127  positionElement->setAttributeNode(NewAttribute("unit","mm"));
128  element->appendChild(positionElement);
129 }
static constexpr double mm
Definition: G4SIunits.hh:115
double x() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
double z() const
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
static const G4double kLinearPrecision
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::PositionWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector pos 
)
inline

Definition at line 61 of file G4GDMLWriteDefine.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::Rotation_vectorWrite ( xercesc::DOMElement *  element,
const G4String tag,
const G4String name,
const G4ThreeVector rot 
)
protected

Definition at line 98 of file G4GDMLWriteDefine.cc.

100 {
101  const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x();
102  const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y();
103  const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z();
104 
105  xercesc::DOMElement* rotationElement = NewElement(tag);
106  rotationElement->setAttributeNode(NewAttribute("name",name));
107  rotationElement->setAttributeNode(NewAttribute("x",x/degree));
108  rotationElement->setAttributeNode(NewAttribute("y",y/degree));
109  rotationElement->setAttributeNode(NewAttribute("z",z/degree));
110  rotationElement->setAttributeNode(NewAttribute("unit","deg"));
111  element->appendChild(rotationElement);
112 }
static const G4double kAngularPrecision
double x() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
double z() const
static constexpr double degree
Definition: G4SIunits.hh:144
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::RotationWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector rot 
)
inline

Definition at line 58 of file G4GDMLWriteDefine.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::Scale_vectorWrite ( xercesc::DOMElement *  element,
const G4String tag,
const G4String name,
const G4ThreeVector scl 
)
protected

Definition at line 79 of file G4GDMLWriteDefine.cc.

81 {
82  const G4double x = (std::fabs(scl.x()-1.0) < kRelativePrecision)
83  ? 1.0 : scl.x();
84  const G4double y = (std::fabs(scl.y()-1.0) < kRelativePrecision)
85  ? 1.0 : scl.y();
86  const G4double z = (std::fabs(scl.z()-1.0) < kRelativePrecision)
87  ? 1.0 : scl.z();
88 
89  xercesc::DOMElement* scaleElement = NewElement(tag);
90  scaleElement->setAttributeNode(NewAttribute("name",name));
91  scaleElement->setAttributeNode(NewAttribute("x",x));
92  scaleElement->setAttributeNode(NewAttribute("y",y));
93  scaleElement->setAttributeNode(NewAttribute("z",z));
94  element->appendChild(scaleElement);
95 }
double x() const
xercesc::DOMElement * NewElement(const G4String &)
Definition: G4GDMLWrite.cc:160
double z() const
static const G4double kRelativePrecision
xercesc::DOMAttr * NewAttribute(const G4String &, const G4String &)
Definition: G4GDMLWrite.cc:136
double y() const
double G4double
Definition: G4Types.hh:76

Here is the call graph for this function:

Here is the caller graph for this function:

void G4GDMLWriteDefine::ScaleWrite ( xercesc::DOMElement *  element,
const G4String name,
const G4ThreeVector scl 
)
inline

Definition at line 55 of file G4GDMLWriteDefine.hh.

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

Here is the call graph for this function:

Here is the caller graph for this function:

Member Data Documentation

xercesc::DOMElement* G4GDMLWriteDefine::defineElement
protected

Definition at line 93 of file G4GDMLWriteDefine.hh.

const G4double G4GDMLWriteDefine::kAngularPrecision = DBL_EPSILON
staticprotected

Definition at line 90 of file G4GDMLWriteDefine.hh.

const G4double G4GDMLWriteDefine::kLinearPrecision = DBL_EPSILON
staticprotected

Definition at line 91 of file G4GDMLWriteDefine.hh.

const G4double G4GDMLWriteDefine::kRelativePrecision = DBL_EPSILON
staticprotected

Definition at line 89 of file G4GDMLWriteDefine.hh.


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