source: trunk/source/persistency/gdml/src/G4GDMLWriteDefine.cc@ 1316

Last change on this file since 1316 was 1228, checked in by garnier, 16 years ago

update geant4.9.3 tag

File size: 5.3 KB
Line 
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//
27// $Id: G4GDMLWriteDefine.cc,v 1.19 2009/03/24 15:47:33 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-03 $
29//
30// class G4GDMLWriteDefine Implementation
31//
32// Original author: Zoltan Torzsok, November 2007
33//
34// --------------------------------------------------------------------
35
36#include "G4GDMLWriteDefine.hh"
37
38const G4double G4GDMLWriteDefine::kRelativePrecision = DBL_EPSILON;
39const G4double G4GDMLWriteDefine::kAngularPrecision = DBL_EPSILON;
40const G4double G4GDMLWriteDefine::kLinearPrecision = DBL_EPSILON;
41
42G4GDMLWriteDefine::G4GDMLWriteDefine() : G4GDMLWrite()
43{
44}
45
46G4GDMLWriteDefine::~G4GDMLWriteDefine()
47{
48}
49
50G4ThreeVector G4GDMLWriteDefine::GetAngles(const G4RotationMatrix& mat)
51{
52 G4double x,y,z;
53
54 const G4double cosb = std::sqrt(mat.xx()*mat.xx()+mat.yx()*mat.yx());
55
56 if (cosb > kRelativePrecision)
57 {
58 x = std::atan2(mat.zy(),mat.zz());
59 y = std::atan2(-mat.zx(),cosb);
60 z = std::atan2(mat.yx(),mat.xx());
61 }
62 else
63 {
64 x = std::atan2(-mat.yz(),mat.yy());
65 y = std::atan2(-mat.zx(),cosb);
66 z = 0.0;
67 }
68
69 return G4ThreeVector(x,y,z);
70}
71
72void G4GDMLWriteDefine::
73Scale_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
74 const G4String& name, const G4ThreeVector& scl)
75{
76 const G4double x = (std::fabs(scl.x()-1.0) < kRelativePrecision)
77 ? 1.0 : scl.x();
78 const G4double y = (std::fabs(scl.y()-1.0) < kRelativePrecision)
79 ? 1.0 : scl.y();
80 const G4double z = (std::fabs(scl.z()-1.0) < kRelativePrecision)
81 ? 1.0 : scl.z();
82
83 xercesc::DOMElement* scaleElement = NewElement(tag);
84 scaleElement->setAttributeNode(NewAttribute("name",name));
85 scaleElement->setAttributeNode(NewAttribute("x",x));
86 scaleElement->setAttributeNode(NewAttribute("y",y));
87 scaleElement->setAttributeNode(NewAttribute("z",z));
88 element->appendChild(scaleElement);
89}
90
91void G4GDMLWriteDefine::
92Rotation_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
93 const G4String& name, const G4ThreeVector& rot)
94{
95 const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x();
96 const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y();
97 const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z();
98
99 xercesc::DOMElement* rotationElement = NewElement(tag);
100 rotationElement->setAttributeNode(NewAttribute("name",name));
101 rotationElement->setAttributeNode(NewAttribute("x",x/degree));
102 rotationElement->setAttributeNode(NewAttribute("y",y/degree));
103 rotationElement->setAttributeNode(NewAttribute("z",z/degree));
104 rotationElement->setAttributeNode(NewAttribute("unit","deg"));
105 element->appendChild(rotationElement);
106}
107
108void G4GDMLWriteDefine::
109Position_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
110 const G4String& name, const G4ThreeVector& pos)
111{
112 const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x();
113 const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y();
114 const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z();
115
116 xercesc::DOMElement* positionElement = NewElement(tag);
117 positionElement->setAttributeNode(NewAttribute("name",name));
118 positionElement->setAttributeNode(NewAttribute("x",x/mm));
119 positionElement->setAttributeNode(NewAttribute("y",y/mm));
120 positionElement->setAttributeNode(NewAttribute("z",z/mm));
121 positionElement->setAttributeNode(NewAttribute("unit","mm"));
122 element->appendChild(positionElement);
123}
124
125void G4GDMLWriteDefine::DefineWrite(xercesc::DOMElement* element)
126{
127 G4cout << "G4GDML: Writing definitions..." << G4endl;
128
129 defineElement = NewElement("define");
130 element->appendChild(defineElement);
131}
Note: See TracBrowser for help on using the repository browser.