Geant4  10.02.p02
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 //
28 //
29 // $Id: G4LatticePhysical.cc 76883 2013-11-18 12:50:08Z gcosmo $
30 //
31 // 20131115 Save rotation results in local variable, report verbosely
32 // 20131116 Replace G4Transform3D with G4RotationMatrix
33 
34 #include "G4LatticePhysical.hh"
35 #include "G4LatticeLogical.hh"
36 #include "G4PhysicalConstants.hh"
37 #include "G4RotationMatrix.hh"
38 #include "G4SystemOfUnits.hh"
39 
40 
41 // Unit vectors defined for convenience (avoid memory churn)
42 
43 namespace {
44  G4ThreeVector xhat(1,0,0), yhat(0,1,0), zhat(0,0,1), nullVec(0,0,0);
45 }
46 
47 
48 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
49 
51  const G4RotationMatrix* Rot)
52  : verboseLevel(0), fTheta(0), fPhi(0), fLattice(Lat) {
54 }
55 
57 
58 
59 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
60 
62  if (!Rot) { // No orientation specified
63  fLocalToGlobal = fGlobalToLocal = G4RotationMatrix::IDENTITY;
64  } else {
65  fLocalToGlobal = fGlobalToLocal = *Rot; // Frame rotation
66  fGlobalToLocal.invert();
67  }
68 
69  if (verboseLevel) {
70  G4cout << "G4LatticePhysical::SetPhysicalOrientation " << *Rot
71  << "\nfLocalToGlobal: " << fLocalToGlobal
72  << "\nfGlobalToLocal: " << fGlobalToLocal
73  << G4endl;
74  }
75 }
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
80  fTheta = t_rot;
81  fPhi = p_rot;
82 
83  if (verboseLevel)
84  G4cout << "G4LatticePhysical::SetLatticeOrientation " << fTheta << " "
85  << fPhi << G4endl;
86 }
87 
88 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
89 
91  fTheta = halfpi - std::atan2(n+0.000001,l+0.000001);
92  fPhi = halfpi - std::atan2(l+0.000001,k+0.000001);
93 
94  if (verboseLevel)
95  G4cout << "G4LatticePhysical::SetMillerOrientation(" << l << k << n
96  << ") : " << fTheta << " " << fPhi << G4endl;
97 }
98 
99 
100 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
101 
103 //Loads the group velocity in m/s
106  G4ThreeVector k) const {
107  if (verboseLevel>1) G4cout << "G4LatticePhysical::MapKtoV " << k << G4endl;
108 
109  k.rotate(yhat,fTheta).rotate(zhat, fPhi);
110  return fLattice->MapKtoV(polarizationState, k);
111 }
112 
114 //Loads the normalized direction vector along VG
117  G4ThreeVector k) const {
118  if (verboseLevel>1) G4cout << "G4LatticePhysical::MapKtoVDir " << k << G4endl;
119 
120  k.rotate(yhat,fTheta).rotate(zhat,fPhi);
121 
122  G4ThreeVector VG = fLattice->MapKtoVDir(polarizationState, k);
123 
124  return VG.rotate(zhat,-fPhi).rotate(yhat,-fTheta);
125 }
126 
127 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
128 
129 // Apply orientation transforms to specified vector
130 
133  if (verboseLevel>1) {
134  G4cout << "G4LatticePhysical::RotateToGlobal " << dir
135  << "\nusing fLocalToGlobal " << fLocalToGlobal
136  << G4endl;
137  }
138 
139  G4ThreeVector result = fLocalToGlobal*dir;
140  if (verboseLevel>1) G4cout << " result " << result << G4endl;
141 
142  return result;
143 }
144 
147  if (verboseLevel>1) {
148  G4cout << "G4LatticePhysical::RotateToLocal " << dir
149  << "\nusing fGlobalToLocal " << fGlobalToLocal
150  << G4endl;
151  }
152 
153  G4ThreeVector result = fGlobalToLocal*dir;
154  if (verboseLevel>1) G4cout << " result " << result << G4endl;
155 
156  return result;
157 }
virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector &) const
virtual ~G4LatticePhysical()
G4LatticePhysical(const G4LatticeLogical *Lat=0, const G4RotationMatrix *Rot=0)
static const double halfpi
Definition: G4SIunits.hh:76
CLHEP::Hep3Vector G4ThreeVector
CLHEP::HepRotation G4RotationMatrix
const G4LatticeLogical * fLattice
virtual G4double MapKtoV(G4int, const G4ThreeVector &) const
G4double MapKtoV(G4int, G4ThreeVector) const
G4RotationMatrix fGlobalToLocal
G4ThreeVector RotateToLocal(const G4ThreeVector &dir) const
int G4int
Definition: G4Types.hh:78
void SetPhysicalOrientation(const G4RotationMatrix *Rot)
G4GLOB_DLL std::ostream G4cout
G4ThreeVector MapKtoVDir(G4int, G4ThreeVector) const
const G4int n
G4RotationMatrix fLocalToGlobal
void SetMillerOrientation(G4int, G4int, G4int)
void SetLatticeOrientation(G4double, G4double)
Definition of the G4LatticePhysical class.
Definition of the G4LatticeLogical class.
G4ThreeVector RotateToGlobal(const G4ThreeVector &dir) const
#define G4endl
Definition: G4ios.hh:61
double G4double
Definition: G4Types.hh:76