G4GDMLWriteDefine Class Reference

#include <G4GDMLWriteDefine.hh>

Inheritance diagram for G4GDMLWriteDefine:

G4GDMLWrite G4GDMLWriteMaterials G4GDMLWriteSolids G4GDMLWriteSetup G4GDMLWriteParamvol G4GDMLWriteStructure

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 *)

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 Attributes

xercesc::DOMElement * defineElement

Static Protected Attributes

static const G4double kRelativePrecision = DBL_EPSILON
static const G4double kAngularPrecision = DBL_EPSILON
static const G4double kLinearPrecision = DBL_EPSILON

Detailed Description

Definition at line 49 of file G4GDMLWriteDefine.hh.


Constructor & Destructor Documentation

G4GDMLWriteDefine::G4GDMLWriteDefine (  )  [protected]

Definition at line 42 of file G4GDMLWriteDefine.cc.

00043   : G4GDMLWrite(), defineElement(0)
00044 {
00045 }

G4GDMLWriteDefine::~G4GDMLWriteDefine (  )  [protected, virtual]

Definition at line 47 of file G4GDMLWriteDefine.cc.

00048 {
00049 }


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().

00071          { Position_vectorWrite(defineElement,"position",name,pos); }

void G4GDMLWriteDefine::DefineWrite ( xercesc::DOMElement *   )  [virtual]

Implements G4GDMLWrite.

Definition at line 126 of file G4GDMLWriteDefine.cc.

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

00127 {
00128    G4cout << "G4GDML: Writing definitions..." << G4endl;
00129 
00130    defineElement = NewElement("define");
00131    element->appendChild(defineElement);
00132 }

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().

00069          { Position_vectorWrite(element,"firstposition",name,pos); }

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().

00066          { Rotation_vectorWrite(element,"firstrotation",name,rot); }

G4ThreeVector G4GDMLWriteDefine::GetAngles ( const G4RotationMatrix  ) 

Definition at line 51 of file G4GDMLWriteDefine.cc.

References kRelativePrecision.

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

00052 {
00053    G4double x,y,z;
00054 
00055    const G4double cosb = std::sqrt(mat.xx()*mat.xx()+mat.yx()*mat.yx());
00056 
00057    if (cosb > kRelativePrecision)
00058    {
00059       x = std::atan2(mat.zy(),mat.zz());
00060       y = std::atan2(-mat.zx(),cosb);
00061       z = std::atan2(mat.yx(),mat.xx());
00062    }
00063    else
00064    {
00065       x = std::atan2(-mat.yz(),mat.yy());
00066       y = std::atan2(-mat.zx(),cosb);
00067       z = 0.0;
00068    }
00069 
00070    return G4ThreeVector(x,y,z);
00071 }

void G4GDMLWriteDefine::Position_vectorWrite ( xercesc::DOMElement *  ,
const G4String ,
const G4String ,
const G4ThreeVector  
) [protected]

Definition at line 110 of file G4GDMLWriteDefine.cc.

References kLinearPrecision, G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

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

00112 {
00113    const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x();
00114    const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y();
00115    const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z();
00116 
00117    xercesc::DOMElement* positionElement = NewElement(tag);
00118    positionElement->setAttributeNode(NewAttribute("name",name));
00119    positionElement->setAttributeNode(NewAttribute("x",x/mm));
00120    positionElement->setAttributeNode(NewAttribute("y",y/mm));
00121    positionElement->setAttributeNode(NewAttribute("z",z/mm));
00122    positionElement->setAttributeNode(NewAttribute("unit","mm"));
00123    element->appendChild(positionElement);
00124 }

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().

00063          { Position_vectorWrite(element,"position",name,pos); }

void G4GDMLWriteDefine::Rotation_vectorWrite ( xercesc::DOMElement *  ,
const G4String ,
const G4String ,
const G4ThreeVector  
) [protected]

Definition at line 93 of file G4GDMLWriteDefine.cc.

References kAngularPrecision, G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by FirstrotationWrite(), and RotationWrite().

00095 {
00096    const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x();
00097    const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y();
00098    const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z();
00099 
00100    xercesc::DOMElement* rotationElement = NewElement(tag);
00101    rotationElement->setAttributeNode(NewAttribute("name",name));
00102    rotationElement->setAttributeNode(NewAttribute("x",x/degree));
00103    rotationElement->setAttributeNode(NewAttribute("y",y/degree));
00104    rotationElement->setAttributeNode(NewAttribute("z",z/degree));
00105    rotationElement->setAttributeNode(NewAttribute("unit","deg"));
00106    element->appendChild(rotationElement);
00107 }

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().

00060          { Rotation_vectorWrite(element,"rotation",name,rot); }

void G4GDMLWriteDefine::Scale_vectorWrite ( xercesc::DOMElement *  ,
const G4String ,
const G4String ,
const G4ThreeVector  
) [protected]

Definition at line 74 of file G4GDMLWriteDefine.cc.

References kRelativePrecision, G4GDMLWrite::NewAttribute(), and G4GDMLWrite::NewElement().

Referenced by ScaleWrite().

00076 {
00077    const G4double x = (std::fabs(scl.x()-1.0) < kRelativePrecision)
00078                     ? 1.0 : scl.x();
00079    const G4double y = (std::fabs(scl.y()-1.0) < kRelativePrecision)
00080                     ? 1.0 : scl.y();
00081    const G4double z = (std::fabs(scl.z()-1.0) < kRelativePrecision)
00082                     ? 1.0 : scl.z();
00083 
00084    xercesc::DOMElement* scaleElement = NewElement(tag);
00085    scaleElement->setAttributeNode(NewAttribute("name",name));
00086    scaleElement->setAttributeNode(NewAttribute("x",x));
00087    scaleElement->setAttributeNode(NewAttribute("y",y));
00088    scaleElement->setAttributeNode(NewAttribute("z",z));
00089    element->appendChild(scaleElement);
00090 }

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().

00057          { Scale_vectorWrite(element,"scale",name,scl); }


Field Documentation

xercesc::DOMElement* G4GDMLWriteDefine::defineElement [protected]

Definition at line 93 of file G4GDMLWriteDefine.hh.

Referenced by AddPosition(), DefineWrite(), G4GDMLWriteMaterials::PropertyVectorWrite(), and G4GDMLWriteMaterials::PropertyWrite().

const G4double G4GDMLWriteDefine::kAngularPrecision = DBL_EPSILON [static, protected]

Definition at line 90 of file G4GDMLWriteDefine.hh.

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

const G4double G4GDMLWriteDefine::kLinearPrecision = DBL_EPSILON [static, protected]

Definition at line 91 of file G4GDMLWriteDefine.hh.

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

const G4double G4GDMLWriteDefine::kRelativePrecision = DBL_EPSILON [static, protected]

Definition at line 89 of file G4GDMLWriteDefine.hh.

Referenced by GetAngles(), G4GDMLWriteStructure::PhysvolWrite(), Scale_vectorWrite(), and G4GDMLWriteStructure::TraverseVolumeTree().


The documentation for this class was generated from the following files:
Generated on Mon May 27 17:52:03 2013 for Geant4 by  doxygen 1.4.7