Back to home page

sPhenix code displayed by LXR

 
 

    


File indexing completed on 2025-12-17 09:21:54

0001 //
0002 // ********************************************************************
0003 // * License and Disclaimer                                           *
0004 // *                                                                  *
0005 // * The  Geant4 software  is  copyright of the Copyright Holders  of *
0006 // * the Geant4 Collaboration.  It is provided  under  the terms  and *
0007 // * conditions of the Geant4 Software License,  included in the file *
0008 // * LICENSE and available at  http://cern.ch/geant4/license .  These *
0009 // * include a list of copyright holders.                             *
0010 // *                                                                  *
0011 // * Neither the authors of this software system, nor their employing *
0012 // * institutes,nor the agencies providing financial support for this *
0013 // * work  make  any representation or  warranty, express or implied, *
0014 // * regarding  this  software system or assume any liability for its *
0015 // * use.  Please see the license in the file  LICENSE  and URL above *
0016 // * for the full disclaimer and the limitation of liability.         *
0017 // *                                                                  *
0018 // * This  code  implementation is the result of  the  scientific and *
0019 // * technical work of the GEANT4 collaboration.                      *
0020 // * By using,  copying,  modifying or  distributing the software (or *
0021 // * any work based  on the software)  you  agree  to acknowledge its *
0022 // * use  in  resulting  scientific  publications,  and indicate your *
0023 // * acceptance of all terms of the Geant4 Software license.          *
0024 // ********************************************************************
0025 //
0026 //
0027 // $Id: PHG4GDMLWriteDefine.cc 68053 2013-03-13 14:39:51Z gcosmo $
0028 //
0029 // class PHG4GDMLWriteDefine Implementation
0030 //
0031 // Original author: Zoltan Torzsok, November 2007
0032 //
0033 // --------------------------------------------------------------------
0034 
0035 #include "PHG4GDMLWriteDefine.hh"
0036 #include <Geant4/G4SystemOfUnits.hh>
0037 
0038 const G4double PHG4GDMLWriteDefine::kRelativePrecision = DBL_EPSILON;
0039 const G4double PHG4GDMLWriteDefine::kAngularPrecision = DBL_EPSILON;
0040 const G4double PHG4GDMLWriteDefine::kLinearPrecision = DBL_EPSILON;
0041 
0042 PHG4GDMLWriteDefine::PHG4GDMLWriteDefine()
0043   :  defineElement(nullptr)
0044 {
0045 }
0046 
0047 G4ThreeVector PHG4GDMLWriteDefine::GetAngles(const G4RotationMatrix& mtx)
0048 {
0049    G4double x;
0050    G4double y;
0051    G4double z;
0052    G4RotationMatrix mat = mtx;
0053    mat.rectify();   // Rectify matrix from possible roundoff errors
0054 
0055    // Direction of rotation given by left-hand rule; clockwise rotation
0056 
0057    static const G4double kMatrixPrecision = 10E-10;
0058    const G4double cosb = std::sqrt(mtx.xx()*mtx.xx()+mtx.yx()*mtx.yx());
0059 
0060    if (cosb > kMatrixPrecision)
0061    {
0062       x = std::atan2(mtx.zy(),mtx.zz());
0063       y = std::atan2(-mtx.zx(),cosb);
0064       z = std::atan2(mtx.yx(),mtx.xx());
0065    }
0066    else
0067    {
0068       x = std::atan2(-mtx.yz(),mtx.yy());
0069       y = std::atan2(-mtx.zx(),cosb);
0070       z = 0.0;
0071    }
0072 
0073    return G4ThreeVector(x,y,z);
0074 }
0075 
0076 void PHG4GDMLWriteDefine::
0077 Scale_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
0078                   const G4String& name, const G4ThreeVector& scl)
0079 {
0080    const G4double x = (std::fabs(scl.x()-1.0) < kRelativePrecision)
0081                     ? 1.0 : scl.x();
0082    const G4double y = (std::fabs(scl.y()-1.0) < kRelativePrecision)
0083                     ? 1.0 : scl.y();
0084    const G4double z = (std::fabs(scl.z()-1.0) < kRelativePrecision)
0085                     ? 1.0 : scl.z();
0086 
0087    xercesc::DOMElement* scaleElement = NewElement(tag);
0088    scaleElement->setAttributeNode(NewAttribute("name",name));
0089    scaleElement->setAttributeNode(NewAttribute("x",x));
0090    scaleElement->setAttributeNode(NewAttribute("y",y));
0091    scaleElement->setAttributeNode(NewAttribute("z",z));
0092    element->appendChild(scaleElement);
0093 }
0094 
0095 void PHG4GDMLWriteDefine::
0096 Rotation_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
0097                      const G4String& name, const G4ThreeVector& rot)
0098 {
0099    const G4double x = (std::fabs(rot.x()) < kAngularPrecision) ? 0.0 : rot.x();
0100    const G4double y = (std::fabs(rot.y()) < kAngularPrecision) ? 0.0 : rot.y();
0101    const G4double z = (std::fabs(rot.z()) < kAngularPrecision) ? 0.0 : rot.z();
0102 
0103    xercesc::DOMElement* rotationElement = NewElement(tag);
0104    rotationElement->setAttributeNode(NewAttribute("name",name));
0105    rotationElement->setAttributeNode(NewAttribute("x",x/degree));
0106    rotationElement->setAttributeNode(NewAttribute("y",y/degree));
0107    rotationElement->setAttributeNode(NewAttribute("z",z/degree));
0108    rotationElement->setAttributeNode(NewAttribute("unit","deg"));
0109    element->appendChild(rotationElement);
0110 }
0111 
0112 void PHG4GDMLWriteDefine::
0113 Position_vectorWrite(xercesc::DOMElement* element, const G4String& tag,
0114                      const G4String& name, const G4ThreeVector& pos)
0115 {
0116    const G4double x = (std::fabs(pos.x()) < kLinearPrecision) ? 0.0 : pos.x();
0117    const G4double y = (std::fabs(pos.y()) < kLinearPrecision) ? 0.0 : pos.y();
0118    const G4double z = (std::fabs(pos.z()) < kLinearPrecision) ? 0.0 : pos.z();
0119 
0120    xercesc::DOMElement* positionElement = NewElement(tag);
0121    positionElement->setAttributeNode(NewAttribute("name",name));
0122    positionElement->setAttributeNode(NewAttribute("x",x/mm));
0123    positionElement->setAttributeNode(NewAttribute("y",y/mm));
0124    positionElement->setAttributeNode(NewAttribute("z",z/mm));
0125    positionElement->setAttributeNode(NewAttribute("unit","mm"));
0126    element->appendChild(positionElement);
0127 }
0128 
0129 void PHG4GDMLWriteDefine::DefineWrite(xercesc::DOMElement* element)
0130 {
0131   std::cout << "PHG4GDML: Writing definitions..." << std::endl;
0132 
0133    defineElement = NewElement("define");
0134    element->appendChild(defineElement);
0135 }