Geant4 10.7.0
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4tgbRotationMatrix Class Reference

#include <G4tgbRotationMatrix.hh>

Public Member Functions

 G4tgbRotationMatrix ()
 
 ~G4tgbRotationMatrix ()
 
 G4tgbRotationMatrix (G4tgrRotationMatrix *tgr)
 
G4RotationMatrixBuildG4RotMatrix ()
 
G4RotationMatrixBuildG4RotMatrixFrom3 (std::vector< G4double > &values)
 
G4RotationMatrixBuildG4RotMatrixFrom6 (std::vector< G4double > &values)
 
G4RotationMatrixBuildG4RotMatrixFrom9 (std::vector< G4double > &values)
 
const G4StringGetName ()
 

Detailed Description

Definition at line 46 of file G4tgbRotationMatrix.hh.

Constructor & Destructor Documentation

◆ G4tgbRotationMatrix() [1/2]

G4tgbRotationMatrix::G4tgbRotationMatrix ( )

Definition at line 38 of file G4tgbRotationMatrix.cc.

39{
40}

◆ ~G4tgbRotationMatrix()

G4tgbRotationMatrix::~G4tgbRotationMatrix ( )

Definition at line 43 of file G4tgbRotationMatrix.cc.

44{
45}

◆ G4tgbRotationMatrix() [2/2]

G4tgbRotationMatrix::G4tgbRotationMatrix ( G4tgrRotationMatrix tgr)

Definition at line 48 of file G4tgbRotationMatrix.cc.

49 : theTgrRM(tgr)
50{
51}

Member Function Documentation

◆ BuildG4RotMatrix()

G4RotationMatrix * G4tgbRotationMatrix::BuildG4RotMatrix ( )

Definition at line 54 of file G4tgbRotationMatrix.cc.

55{
56 std::vector<G4double> values = theTgrRM->GetValues();
57
58 if(values.size() == 3)
59 {
60 return BuildG4RotMatrixFrom3(values);
61 }
62 else if(values.size() == 6)
63 {
64 return BuildG4RotMatrixFrom6(values);
65 }
66 else if(values.size() == 9)
67 {
68 return BuildG4RotMatrixFrom9(values);
69 }
70 else
71 {
72 G4String ErrMessage = "Number of values is: " +
73 G4UIcommand::ConvertToString(G4int(values.size())) +
74 G4String(". It should be 3, 6, or 9 !");
75 G4Exception("G4tgbRotationMatrix::BuildG4RotMatrix()", "InvalidData",
76 FatalException, ErrMessage);
77 }
78 return nullptr;
79}
@ FatalException
void G4Exception(const char *originOfException, const char *exceptionCode, G4ExceptionSeverity severity, const char *description)
Definition: G4Exception.cc:35
int G4int
Definition: G4Types.hh:85
static G4String ConvertToString(G4bool boolVal)
Definition: G4UIcommand.cc:430
G4RotationMatrix * BuildG4RotMatrixFrom9(std::vector< G4double > &values)
G4RotationMatrix * BuildG4RotMatrixFrom6(std::vector< G4double > &values)
G4RotationMatrix * BuildG4RotMatrixFrom3(std::vector< G4double > &values)
std::vector< G4double > & GetValues()

Referenced by G4tgbRotationMatrixMgr::FindOrBuildG4RotMatrix().

◆ BuildG4RotMatrixFrom3()

G4RotationMatrix * G4tgbRotationMatrix::BuildG4RotMatrixFrom3 ( std::vector< G4double > &  values)

Definition at line 83 of file G4tgbRotationMatrix.cc.

84{
85 G4RotationMatrix* rotMat = new G4RotationMatrix();
86
87 rotMat->rotateX(values[0]);
88 rotMat->rotateY(values[1]);
89 rotMat->rotateZ(values[2]);
90
91#ifdef G4VERBOSE
93 {
94 G4cout << " Constructing new G4RotationMatrix from 3 numbers " << GetName()
95 << " : " << *rotMat << G4endl;
96 }
97#endif
98
99 return rotMat;
100}
CLHEP::HepRotation G4RotationMatrix
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
HepRotation & rotateX(double delta)
Definition: Rotation.cc:61
HepRotation & rotateZ(double delta)
Definition: Rotation.cc:87
HepRotation & rotateY(double delta)
Definition: Rotation.cc:74
const G4String & GetName()
static G4int GetVerboseLevel()

Referenced by BuildG4RotMatrix().

◆ BuildG4RotMatrixFrom6()

G4RotationMatrix * G4tgbRotationMatrix::BuildG4RotMatrixFrom6 ( std::vector< G4double > &  values)

Definition at line 104 of file G4tgbRotationMatrix.cc.

105{
106 G4double thetaX = values[0];
107 G4double phiX = values[1];
108 G4double thetaY = values[2];
109 G4double phiY = values[3];
110 G4double thetaZ = values[4];
111 G4double phiZ = values[5];
112
113 // build the 3 axis from the values
114 G4ThreeVector colx(std::sin(thetaX) * std::cos(phiX),
115 std::sin(thetaX) * std::sin(phiX), std::cos(thetaX));
116 G4ThreeVector coly(std::sin(thetaY) * std::cos(phiY),
117 std::sin(thetaY) * std::sin(phiY), std::cos(thetaY));
118 G4ThreeVector colz(std::sin(thetaZ) * std::cos(phiZ),
119 std::sin(thetaZ) * std::sin(phiZ), std::cos(thetaZ));
120
121 // Now create a G4RotationMatrix (HepRotation), which can be left handed.
122 // This is not foreseen in CLHEP, but can be achieved using the
123 // constructor which does not check its input arguments!
124
125 G4Rep3x3 rottemp(colx.x(), coly.x(), colz.x(), // matrix representation
126 colx.y(), coly.y(), colz.y(), // (inverted)
127 colx.z(), coly.z(), colz.z());
128
129 G4RotationMatrix* rotMat = new G4RotationMatrix(rottemp);
130
131#ifdef G4VERBOSE
133 {
134 G4cout << " Constructing new G4RotationMatrix from 6 numbers " << GetName()
135 << " : " << *rotMat << G4endl;
136 }
137#endif
138
139 return rotMat;
140}
double G4double
Definition: G4Types.hh:83

Referenced by BuildG4RotMatrix().

◆ BuildG4RotMatrixFrom9()

G4RotationMatrix * G4tgbRotationMatrix::BuildG4RotMatrixFrom9 ( std::vector< G4double > &  values)

Definition at line 144 of file G4tgbRotationMatrix.cc.

145{
146 // build the 3 axis from the values
147 G4ThreeVector colx(values[0], values[1], values[2]);
148 G4ThreeVector coly(values[3], values[4], values[5]);
149 G4ThreeVector colz(values[6], values[7], values[8]);
150
151 // Now create a G4RotationMatrix (HepRotation), which can be left handed.
152 // This is not foreseen in CLHEP, but can be achieved using the
153 // constructor which does not check its input arguments!
154
155 G4Rep3x3 rottemp(colx.x(), coly.x(), colz.x(), // matrix representation
156 colx.y(), coly.y(), colz.y(), // (inverted)
157 colx.z(), coly.z(), colz.z());
158
159 G4RotationMatrix* rotMat = new G4RotationMatrix(rottemp);
160
161#ifdef G4VERBOSE
163 {
164 G4cout << " Constructing new G4RotationMatrix from 9 numbers " << GetName()
165 << " : " << *rotMat << G4endl;
166 }
167#endif
168
169 return rotMat;
170}

Referenced by BuildG4RotMatrix().

◆ GetName()

const G4String & G4tgbRotationMatrix::GetName ( )
inline

Definition at line 63 of file G4tgbRotationMatrix.hh.

63{ return theTgrRM->GetName(); }
const G4String & GetName()

Referenced by BuildG4RotMatrixFrom3(), BuildG4RotMatrixFrom6(), and BuildG4RotMatrixFrom9().


The documentation for this class was generated from the following files: