Geant4-11
Functions
G4gsrotm.cc File Reference
#include "G3toG4.hh"
#include "G3RotTable.hh"
#include "G3toG4RotationMatrix.hh"
#include "G4PhysicalConstants.hh"
#include "G4ThreeVector.hh"

Go to the source code of this file.

Functions

void G4gsrotm (G4int irot, G4double theta1, G4double phi1, G4double theta2, G4double phi2, G4double theta3, G4double phi3)
 
void PG4gsrotm (G4String *tokens)
 

Function Documentation

◆ G4gsrotm()

void G4gsrotm ( G4int  irot,
G4double  theta1,
G4double  phi1,
G4double  theta2,
G4double  phi2,
G4double  theta3,
G4double  phi3 
)

Definition at line 53 of file G4gsrotm.cc.

55{
56 G4double degrad = pi/180;
57
58 G4double th1r = theta1*degrad;
59 G4double th2r = theta2*degrad;
60 G4double th3r = theta3*degrad;
61
62 G4double phi1r = phi1*degrad;
63 G4double phi2r = phi2*degrad;
64 G4double phi3r = phi3*degrad;
65
66 // Construct unit vectors
67
68 G4ThreeVector x(std::sin(th1r)*std::cos(phi1r), std::sin(th1r)*std::sin(phi1r), std::cos(th1r));
69 G4ThreeVector y(std::sin(th2r)*std::cos(phi2r), std::sin(th2r)*std::sin(phi2r), std::cos(th2r));
70 G4ThreeVector z(std::sin(th3r)*std::cos(phi3r), std::sin(th3r)*std::sin(phi3r), std::cos(th3r));
71
72 // check for orthonormality and left-handedness
73
74 G4double check = (x.cross(y))*z;
75 G4double tol = 1.0e-3;
76
77 if (1-std::abs(check)>tol) {
78 G4cerr << "Coordinate axes forming rotation matrix "
79 << irot << " are not orthonormal.(" << 1-std::abs(check) << ")"
80 << G4endl;
81 G4cerr << " theta1=" << theta1;
82 G4cerr << " phi1=" << phi1;
83 G4cerr << " theta2=" << theta2;
84 G4cerr << " phi2=" << phi2;
85 G4cerr << " theta3=" << theta3;
86 G4cerr << " phi3=" << phi3;
87 G4cerr << G4endl;
88 G4Exception("G4gsrotm()", "G3toG40023", FatalException,
89 "Non orthogonal axes!");
90 return;
91 }
92 //else if (1+check<=tol) {
93 // G4cerr << "G4gsrotm warning: coordinate axes forming rotation "
94 // << "matrix " << irot << " are left-handed" << G4endl;
95 //}
96
98
99 rotp->SetRotationMatrixByRow(x, y, z);
100
101 // add it to the List
102
103 G3Rot.Put(irot, rotp);
104}
G3G4DLL_API G3RotTable G3Rot
Definition: clparse.cc:56
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
static constexpr double pi
Definition: G4SIunits.hh:55
double G4double
Definition: G4Types.hh:83
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
void Put(G4int id, G4RotationMatrix *matrix)
Definition: G3RotTable.cc:52
void SetRotationMatrixByRow(const G4ThreeVector &Row1, const G4ThreeVector &Row2, const G4ThreeVector &Row3)

References CLHEP::Hep3Vector::cross(), FatalException, G3Rot, G4cerr, G4endl, G4Exception(), pi, G3RotTable::Put(), and G3toG4RotationMatrix::SetRotationMatrixByRow().

Referenced by G4BuildGeom(), and PG4gsrotm().

◆ PG4gsrotm()

void PG4gsrotm ( G4String tokens)

Definition at line 34 of file G4gsrotm.cc.

35{
36 // fill the parameter containers
37 G3fillParams(tokens,PTgsrotm);
38
39 // interpret the parameters
40 G4int irot = Ipar[0];
41
42 // the angles in Geant are in degrees
43 G4double theta1 = Rpar[0];
44 G4double phi1 = Rpar[1];
45 G4double theta2 = Rpar[2];
46 G4double phi2 = Rpar[3];
47 G4double theta3 = Rpar[4];
48 G4double phi3 = Rpar[5];
49
50 G4gsrotm(irot, theta1,phi1, theta2,phi2, theta3,phi3);
51}
#define PTgsrotm
Definition: G3toG4.hh:55
G3G4DLL_API G4int Ipar[1000]
Definition: clparse.cc:65
void G3fillParams(G4String *tokens, const char *ptypes)
Definition: clparse.cc:218
G3G4DLL_API G4double Rpar[1000]
Definition: clparse.cc:66
int G4int
Definition: G4Types.hh:85
void G4gsrotm(G4int irot, G4double theta1, G4double phi1, G4double theta2, G4double phi2, G4double theta3, G4double phi3)
Definition: G4gsrotm.cc:53

References G3fillParams(), G4gsrotm(), Ipar, PTgsrotm, and Rpar.

Referenced by G3CLEval().