Geant4 9.6.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4gsrotm.cc
Go to the documentation of this file.
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$
28//
29#include "G3toG4.hh"
30#include "G3RotTable.hh"
33#include "G4ThreeVector.hh"
34
35void PG4gsrotm(G4String *tokens)
36{
37 // fill the parameter containers
38 G3fillParams(tokens,PTgsrotm);
39
40 // interpret the parameters
41 G4int irot = Ipar[0];
42
43 // the angles in Geant are in degrees
44 G4double theta1 = Rpar[0];
45 G4double phi1 = Rpar[1];
46 G4double theta2 = Rpar[2];
47 G4double phi2 = Rpar[3];
48 G4double theta3 = Rpar[4];
49 G4double phi3 = Rpar[5];
50
51 G4gsrotm(irot, theta1,phi1, theta2,phi2, theta3,phi3);
52}
53
54void G4gsrotm(G4int irot, G4double theta1, G4double phi1,
55 G4double theta2, G4double phi2, G4double theta3, G4double phi3)
56{
57 G4double degrad = pi/180;
58
59 G4double th1r = theta1*degrad;
60 G4double th2r = theta2*degrad;
61 G4double th3r = theta3*degrad;
62
63 G4double phi1r = phi1*degrad;
64 G4double phi2r = phi2*degrad;
65 G4double phi3r = phi3*degrad;
66
67 // Construct unit vectors
68
69 G4ThreeVector x(std::sin(th1r)*std::cos(phi1r), std::sin(th1r)*std::sin(phi1r), std::cos(th1r));
70 G4ThreeVector y(std::sin(th2r)*std::cos(phi2r), std::sin(th2r)*std::sin(phi2r), std::cos(th2r));
71 G4ThreeVector z(std::sin(th3r)*std::cos(phi3r), std::sin(th3r)*std::sin(phi3r), std::cos(th3r));
72
73 // check for orthonormality and left-handedness
74
75 G4double check = (x.cross(y))*z;
76 G4double tol = 1.0e-3;
77
78 if (1-std::abs(check)>tol) {
79 G4cerr << "Coordinate axes forming rotation matrix "
80 << irot << " are not orthonormal.(" << 1-std::abs(check) << ")"
81 << G4endl;
82 G4cerr << " theta1=" << theta1;
83 G4cerr << " phi1=" << phi1;
84 G4cerr << " theta2=" << theta2;
85 G4cerr << " phi2=" << phi2;
86 G4cerr << " theta3=" << theta3;
87 G4cerr << " phi3=" << phi3;
88 G4cerr << G4endl;
89 G4Exception("G4gsrotm()", "G3toG40023", FatalException,
90 "Non orthogonal axes!");
91 return;
92 }
93 //else if (1+check<=tol) {
94 // G4cerr << "G4gsrotm warning: coordinate axes forming rotation "
95 // << "matrix " << irot << " are left-handed" << G4endl;
96 //}
97
99
100 rotp->SetRotationMatrixByRow(x, y, z);
101
102 // add it to the List
103
104 G3Rot.Put(irot, rotp);
105}
G3G4DLL_API G3RotTable G3Rot
Definition: clparse.cc:57
#define PTgsrotm
Definition: G3toG4.hh:56
G3G4DLL_API G4int Ipar[1000]
Definition: clparse.cc:66
void G3fillParams(G4String *tokens, const char *ptypes)
Definition: clparse.cc:219
G3G4DLL_API G4double Rpar[1000]
Definition: clparse.cc:67
@ FatalException
double G4double
Definition: G4Types.hh:64
int G4int
Definition: G4Types.hh:66
void PG4gsrotm(G4String *tokens)
Definition: G4gsrotm.cc:35
void G4gsrotm(G4int irot, G4double theta1, G4double phi1, G4double theta2, G4double phi2, G4double theta3, G4double phi3)
Definition: G4gsrotm.cc:54
#define G4endl
Definition: G4ios.hh:52
G4DLLIMPORT std::ostream G4cerr
Hep3Vector cross(const Hep3Vector &) const
void Put(G4int id, G4RotationMatrix *matrix)
Definition: G3RotTable.cc:53
void SetRotationMatrixByRow(const G4ThreeVector &Row1, const G4ThreeVector &Row2, const G4ThreeVector &Row3)
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *comments)
Definition: G4Exception.cc:41