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

Last change on this file since 1177 was 987, checked in by garnier, 15 years ago

fichiers manquants

File size: 5.2 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.18 2008/07/16 15:46:34 gcosmo Exp $
28// GEANT4 tag $Name: geant4-09-02-ref-02 $
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
42G4ThreeVector G4GDMLWriteDefine::GetAngles(const G4RotationMatrix& mat)
43{
44   G4double x,y,z;
45
46   const G4double cosb = std::sqrt(mat.xx()*mat.xx()+mat.yx()*mat.yx());
47
48   if (cosb > kRelativePrecision)
49   {
50      x = std::atan2(mat.zy(),mat.zz());
51      y = std::atan2(-mat.zx(),cosb);
52      z = std::atan2(mat.yx(),mat.xx());
53   }
54   else
55   {
56      x = std::atan2(-mat.yz(),mat.yy());
57      y = std::atan2(-mat.zx(),cosb);
58      z = 0.0;
59   }
60
61   return G4ThreeVector(x,y,z);
62}
63
64void G4GDMLWriteDefine::
65Scale_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
66                  const G4String& name, const G4ThreeVector& scl)
67{
68   const G4double x = (std::fabs(scl.x()-1.0) < kRelativePrecision)
69                    ? 1.0 : scl.x();
70   const G4double y = (std::fabs(scl.y()-1.0) < kRelativePrecision)
71                    ? 1.0 : scl.y();
72   const G4double z = (std::fabs(scl.z()-1.0) < kRelativePrecision)
73                    ? 1.0 : scl.z();
74
75   xercesc::DOMElement* scaleElement = NewElement(tag);
76   scaleElement->setAttributeNode(NewAttribute("name",name));
77   scaleElement->setAttributeNode(NewAttribute("x",x));
78   scaleElement->setAttributeNode(NewAttribute("y",y));
79   scaleElement->setAttributeNode(NewAttribute("z",z));
80   element->appendChild(scaleElement);
81}
82
83void G4GDMLWriteDefine::
84Rotation_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
85                     const G4String& name, const G4ThreeVector& rot)
86{
87   const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x();
88   const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y();
89   const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z();
90
91   xercesc::DOMElement* rotationElement = NewElement(tag);
92   rotationElement->setAttributeNode(NewAttribute("name",name));
93   rotationElement->setAttributeNode(NewAttribute("x",x/degree));
94   rotationElement->setAttributeNode(NewAttribute("y",y/degree));
95   rotationElement->setAttributeNode(NewAttribute("z",z/degree));
96   rotationElement->setAttributeNode(NewAttribute("unit","deg"));
97   element->appendChild(rotationElement);
98}
99
100void G4GDMLWriteDefine::
101Position_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
102                     const G4String& name, const G4ThreeVector& pos)
103{
104   const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x();
105   const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y();
106   const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z();
107
108   xercesc::DOMElement* positionElement = NewElement(tag);
109   positionElement->setAttributeNode(NewAttribute("name",name));
110   positionElement->setAttributeNode(NewAttribute("x",x/mm));
111   positionElement->setAttributeNode(NewAttribute("y",y/mm));
112   positionElement->setAttributeNode(NewAttribute("z",z/mm));
113   positionElement->setAttributeNode(NewAttribute("unit","mm"));
114   element->appendChild(positionElement);
115}
116
117void G4GDMLWriteDefine::DefineWrite(xercesc::DOMElement* element)
118{
119   G4cout << "G4GDML: Writing definitions..." << G4endl;
120
121   defineElement = NewElement("define");
122   element->appendChild(defineElement);
123}
Note: See TracBrowser for help on using the repository browser.