Geant4 11.2.2
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LatticePhysical.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/// \file materials/src/G4LatticePhysical.cc
27/// \brief Implementation of the G4LatticePhysical class
28//
29//
30// 20131115 Save rotation results in local variable, report verbosely
31// 20131116 Replace G4Transform3D with G4RotationMatrix
32
33#include "G4LatticePhysical.hh"
34
35#include "G4LatticeLogical.hh"
37#include "G4RotationMatrix.hh"
38#include "G4SystemOfUnits.hh"
39
40// Unit vectors defined for convenience (avoid memory churn)
41
42namespace
43{
44G4ThreeVector xhat(1, 0, 0), yhat(0, 1, 0), zhat(0, 0, 1), nullVec(0, 0, 0);
45}
46
47
49 : fLattice(Lat)
50{
52}
53
54
56{
57 if (Rot == nullptr) { // No orientation specified
58 fLocalToGlobal = fGlobalToLocal = G4RotationMatrix::IDENTITY;
59 }
60 else {
61 fLocalToGlobal = fGlobalToLocal = *Rot; // Frame rotation
62 fGlobalToLocal.invert();
63 }
64
65 if (verboseLevel != 0) {
66 G4cout << "G4LatticePhysical::SetPhysicalOrientation " << *Rot
67 << "\nfLocalToGlobal: " << fLocalToGlobal << "\nfGlobalToLocal: " << fGlobalToLocal
68 << G4endl;
69 }
70}
71
72
74{
75 fTheta = t_rot;
76 fPhi = p_rot;
77
78 if (verboseLevel != 0) {
79 G4cout << "G4LatticePhysical::SetLatticeOrientation " << fTheta << " " << fPhi << G4endl;
80 }
81}
82
83
85{
86 fTheta = halfpi - std::atan2(n + 0.000001, l + 0.000001);
87 fPhi = halfpi - std::atan2(l + 0.000001, k + 0.000001);
88
89 if (verboseLevel != 0) {
90 G4cout << "G4LatticePhysical::SetMillerOrientation(" << l << k << n << ") : " << fTheta << " "
91 << fPhi << G4endl;
92 }
93}
94
95
96///////////////////////////////
97// Loads the group velocity in m/s
98/////////////////////////////
100{
101 if (verboseLevel > 1) {
102 G4cout << "G4LatticePhysical::MapKtoV " << k << G4endl;
103 }
104
105 k.rotate(yhat, fTheta).rotate(zhat, fPhi);
106 return fLattice->MapKtoV(polarizationState, k);
107}
108
109///////////////////////////////
110// Loads the normalized direction vector along VG
111///////////////////////////////
113{
114 if (verboseLevel > 1) {
115 G4cout << "G4LatticePhysical::MapKtoVDir " << k << G4endl;
116 }
117
118 k.rotate(yhat, fTheta).rotate(zhat, fPhi);
119
120 G4ThreeVector VG = fLattice->MapKtoVDir(polarizationState, k);
121
122 return VG.rotate(zhat, -fPhi).rotate(yhat, -fTheta);
123}
124
125
126// Apply orientation transforms to specified vector
127
129{
130 if (verboseLevel > 1) {
131 G4cout << "G4LatticePhysical::RotateToGlobal " << dir << "\nusing fLocalToGlobal "
132 << fLocalToGlobal << G4endl;
133 }
134
135 G4ThreeVector result = fLocalToGlobal * dir;
136 if (verboseLevel > 1) {
137 G4cout << " result " << result << G4endl;
138 }
139
140 return result;
141}
142
144{
145 if (verboseLevel > 1) {
146 G4cout << "G4LatticePhysical::RotateToLocal " << dir << "\nusing fGlobalToLocal "
147 << fGlobalToLocal << G4endl;
148 }
149
150 G4ThreeVector result = fGlobalToLocal * dir;
151 if (verboseLevel > 1) {
152 G4cout << " result " << result << G4endl;
153 }
154
155 return result;
156}
Definition of the G4LatticeLogical class.
Definition of the G4LatticePhysical class.
double G4double
Definition G4Types.hh:83
int G4int
Definition G4Types.hh:85
#define G4endl
Definition G4ios.hh:67
G4GLOB_DLL std::ostream G4cout
Hep3Vector & rotate(double, const Hep3Vector &)
static DLL_API const HepRotation IDENTITY
Definition Rotation.h:366
HepRotation & invert()
virtual G4double MapKtoV(G4int, const G4ThreeVector &) const
virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector &) const
void SetLatticeOrientation(G4double, G4double)
void SetPhysicalOrientation(const G4RotationMatrix *Rot)
G4ThreeVector RotateToLocal(const G4ThreeVector &dir) const
G4ThreeVector RotateToGlobal(const G4ThreeVector &dir) const
G4double MapKtoV(G4int, G4ThreeVector) const
G4ThreeVector MapKtoVDir(G4int, G4ThreeVector) const
void SetMillerOrientation(G4int, G4int, G4int)
G4LatticePhysical(const G4LatticeLogical *Lat=nullptr, const G4RotationMatrix *Rot=nullptr)