Geant4 Cross Reference |
1 // 2 // ******************************************************************** 3 // * License and Disclaimer * 4 // * * 5 // * The Geant4 software is copyright of the Copyright Holders of * 6 // * the Geant4 Collaboration. It is provided under the terms and * 7 // * conditions of the Geant4 Software License, included in the file * 8 // * LICENSE and available at http://cern.ch/geant4/license . These * 9 // * include a list of copyright holders. * 10 // * * 11 // * Neither the authors of this software system, nor their employing * 12 // * institutes,nor the agencies providing financial support for this * 13 // * work make any representation or warranty, express or implied, * 14 // * regarding this software system or assume any liability for its * 15 // * use. Please see the license in the file LICENSE and URL above * 16 // * for the full disclaimer and the limitation of liability. * 17 // * * 18 // * This code implementation is the result of the scientific and * 19 // * technical work of the GEANT4 collaboration. * 20 // * By using, copying, modifying or distributing the software (or * 21 // * any work based on the software) you agree to acknowledge its * 22 // * use in resulting scientific publications, and indicate your * 23 // * acceptance of all terms of the Geant4 Software license. * 24 // ******************************************************************** 25 // 26 // G4GDMLWriteDefine implementation 27 // 28 // Author: Zoltan Torzsok, November 2007 29 // -------------------------------------------------------------------- 30 31 #include "G4GDMLWriteDefine.hh" 32 #include "G4SystemOfUnits.hh" 33 34 const G4double G4GDMLWriteDefine::kRelativePrecision = DBL_EPSILON; 35 const G4double G4GDMLWriteDefine::kAngularPrecision = DBL_EPSILON; 36 const G4double G4GDMLWriteDefine::kLinearPrecision = DBL_EPSILON; 37 38 // -------------------------------------------------------------------- 39 G4GDMLWriteDefine::G4GDMLWriteDefine() 40 : G4GDMLWrite() 41 { 42 } 43 44 // -------------------------------------------------------------------- 45 G4GDMLWriteDefine::~G4GDMLWriteDefine() 46 { 47 } 48 49 // -------------------------------------------------------------------- 50 G4ThreeVector G4GDMLWriteDefine::GetAngles(const G4RotationMatrix& mtx) 51 { 52 G4double x, y, z; 53 G4RotationMatrix mat = mtx; 54 mat.rectify(); // Rectify matrix from possible roundoff errors 55 56 // Direction of rotation given by left-hand rule; clockwise rotation 57 58 static const G4double kMatrixPrecision = 10E-10; 59 const G4double cosb = std::sqrt(mtx.xx() * mtx.xx() + mtx.yx() * mtx.yx()); 60 61 if(cosb > kMatrixPrecision) 62 { 63 x = std::atan2(mtx.zy(), mtx.zz()); 64 y = std::atan2(-mtx.zx(), cosb); 65 z = std::atan2(mtx.yx(), mtx.xx()); 66 } 67 else 68 { 69 x = std::atan2(-mtx.yz(), mtx.yy()); 70 y = std::atan2(-mtx.zx(), cosb); 71 z = 0.0; 72 } 73 74 return G4ThreeVector(x, y, z); 75 } 76 77 // -------------------------------------------------------------------- 78 void G4GDMLWriteDefine::Scale_vectorWrite(xercesc::DOMElement* element, 79 const G4String& tag, 80 const G4String& name, 81 const G4ThreeVector& scl) 82 { 83 const G4double x = 84 (std::fabs(scl.x() - 1.0) < kRelativePrecision) ? 1.0 : scl.x(); 85 const G4double y = 86 (std::fabs(scl.y() - 1.0) < kRelativePrecision) ? 1.0 : scl.y(); 87 const G4double z = 88 (std::fabs(scl.z() - 1.0) < kRelativePrecision) ? 1.0 : scl.z(); 89 90 xercesc::DOMElement* scaleElement = NewElement(tag); 91 scaleElement->setAttributeNode(NewAttribute("name", name)); 92 scaleElement->setAttributeNode(NewAttribute("x", x)); 93 scaleElement->setAttributeNode(NewAttribute("y", y)); 94 scaleElement->setAttributeNode(NewAttribute("z", z)); 95 element->appendChild(scaleElement); 96 } 97 98 // -------------------------------------------------------------------- 99 void G4GDMLWriteDefine::Rotation_vectorWrite(xercesc::DOMElement* element, 100 const G4String& tag, 101 const G4String& name, 102 const G4ThreeVector& rot) 103 { 104 const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x(); 105 const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y(); 106 const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z(); 107 108 xercesc::DOMElement* rotationElement = NewElement(tag); 109 rotationElement->setAttributeNode(NewAttribute("name", name)); 110 rotationElement->setAttributeNode(NewAttribute("x", x / degree)); 111 rotationElement->setAttributeNode(NewAttribute("y", y / degree)); 112 rotationElement->setAttributeNode(NewAttribute("z", z / degree)); 113 rotationElement->setAttributeNode(NewAttribute("unit", "deg")); 114 element->appendChild(rotationElement); 115 } 116 117 // -------------------------------------------------------------------- 118 void G4GDMLWriteDefine::Position_vectorWrite(xercesc::DOMElement* element, 119 const G4String& tag, 120 const G4String& name, 121 const G4ThreeVector& pos) 122 { 123 const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x(); 124 const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y(); 125 const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z(); 126 127 xercesc::DOMElement* positionElement = NewElement(tag); 128 positionElement->setAttributeNode(NewAttribute("name", name)); 129 positionElement->setAttributeNode(NewAttribute("x", x / mm)); 130 positionElement->setAttributeNode(NewAttribute("y", y / mm)); 131 positionElement->setAttributeNode(NewAttribute("z", z / mm)); 132 positionElement->setAttributeNode(NewAttribute("unit", "mm")); 133 element->appendChild(positionElement); 134 } 135 136 // -------------------------------------------------------------------- 137 void G4GDMLWriteDefine::DefineWrite(xercesc::DOMElement* element) 138 { 139 #ifdef G4VERBOSE 140 G4cout << "G4GDML: Writing definitions..." << G4endl; 141 #endif 142 defineElement = NewElement("define"); 143 element->appendChild(defineElement); 144 } 145