Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
XPhysicalLattice.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 
30 #include "XPhysicalLattice.hh"
31 #include "G4VPhysicalVolume.hh"
32 #include "G4PhysicalConstants.hh"
33 #include "G4SystemOfUnits.hh"
34 
35 #include "G4LogicalVolume.hh"
36 #include "G4Box.hh"
37 
39  fLattice=NULL;
40  fVolume=NULL;
41  fTheta=0.;
42  fPhi=0.;
43  fOmega=0.;
44 
45 
46  fCurvatureRadius = G4ThreeVector(0.,0.,0.); // if cr = 0 == no bending
47  fThermalVibrationAmplitude = 0.1 * angstrom; // no physical meaning
48  fMillerOrientation[0] = 2;
49  fMillerOrientation[1] = 2;
50  fMillerOrientation[2] = 0;
51 }
52 
53 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
54 
56  XLogicalLattice* Lat){
57  fLattice=Lat;
58  fVolume=Vol;
59  fA=fLattice->GetAnhDecConstant();
60  fB=fLattice->GetScatteringConstant();
61  fDosL=fLattice->GetLDOS();
62  fDosST=fLattice->GetSTDOS();
63  fDosFT=fLattice->GetFTDOS();
64  fBeta=fLattice->GetBeta();
65  fGamma=fLattice->GetGamma();
66  fLambda=fLattice->GetLambda();
67  fMu=fLattice->GetMu();
68 
69  G4RotationMatrix *rot = fVolume->GetObjectRotation();
70 
73 }
74 
75 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
76 
78 
79 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
80 
82  double Gamma,
83  double Lambda,
84  double Mu)
85 {
86  fBeta=Beta;
87  fGamma=Gamma;
88  fLambda=Lambda;
89  fMu=Mu;
90 }
91 
92 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
93 
94 
96  fA=a;
97 }
98 
99 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
100 
101 
103  fB=b;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
107 
108 
109 void XPhysicalLattice::SetLDOS(double LDOS){
110  fDosL=LDOS;
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
115 void XPhysicalLattice::SetSTDOS(double STDOS)
116 {
117  fDosST = STDOS;
118 }
119 
120 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
121 
122 void XPhysicalLattice::SetFTDOS(double FTDOS){
123  fDosFT = FTDOS;
124 }
125 
126 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
127 
128 
130  return fBeta;
131 }
132 
133 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
134 
135 
137  return fGamma;
138 }
139 
140 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
141 
143  return fLambda;
144 }
145 
146 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
147 
148 
150 {
151  return fMu;
152 }
153 
154 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
155 
156 
158 {
159  return fB;
160 }
161 
162 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
163 
164 
166 {
167  return fA;
168 }
169 
170 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
171 
172 
174  return fDosL;
175 }
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178 
179 
181  return fDosST;
182 }
183 
184 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
185 
186 
188 {
189  return fDosFT;
190 }
191 
192 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
193 
194 
195 
197 //Loads the group velocity in m/s
199 double XPhysicalLattice::MapKtoV(int polarizationState, G4ThreeVector k){
200  double groupVelocity;
201 
202  k.rotate(G4ThreeVector(0,1,0), fTheta).rotate(G4ThreeVector(0,0,1), fPhi);
203  groupVelocity = fLattice->MapKtoV(polarizationState, k);
204  k.rotate(G4ThreeVector(0,0,1), -fPhi).rotate(G4ThreeVector(0,1,0), -fTheta);
205 
206  return groupVelocity;
207 }
208 
209 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
210 
211 
212 
214 //Loads the normalized direction vector along VG
217  G4ThreeVector k){
218 
219  G4ThreeVector GroupVelocity;
220 
221  k=k.rotate(G4ThreeVector(0,1,0), fTheta).rotate(G4ThreeVector(0,0,1), fPhi);
222  GroupVelocity = fLattice->MapKtoVDir(polarizationState, k);
223 
224  return GroupVelocity.rotate(G4ThreeVector(0,0,1), -fPhi)
225  .rotate(G4ThreeVector(0,1,0), -fTheta).unit();
226 }
227 
228 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
229 
230 
232  return fVolume;
233 }
234 
235 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
236 
237 
239  fVolume=Vol;
240 }
241 
242 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
243 
244 
246  fTheta=t_rot;
247  fPhi= p_rot;
248 }
249 
250 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
251 
253  G4double o_rot,
254  G4double p_rot){
255  fTheta = t_rot;
256  fOmega = o_rot;
257  fPhi = p_rot;
258 }
259 
260 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
261 
263  //fTheta=pi/2-std::atan2(n+0.000001,l+0.000001)*rad;
264  //fPhi=pi/2-std::atan2(l+0.000001,k+0.000001)*rad;
265 
266  // // // added for channeling // // //
267  fMillerOrientation[0]=l;
268  fMillerOrientation[1]=k;
269  fMillerOrientation[2]=n;
270  // // // // // // // // // // // // //
271 }
272 
273 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
274 
275 
277  fLattice=Lat;
278 }
279 
280 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
281 // Begin Channeling specific code
282 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
283 
284 
285 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
286 
288 SetThermalVibrationAmplitude(G4double vThermalVibrationAmplitude){
289  fThermalVibrationAmplitude = vThermalVibrationAmplitude;
290 }
291 
292 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
293 
295  return fThermalVibrationAmplitude;
296 }
297 
298 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
299 
302  G4ThreeVector& vPosition){
303  vMomentum.rotate(G4ThreeVector(1.,0.,0.),fOmega)
304  .rotate(G4ThreeVector(0.,1.,0.), fTheta)
305  .rotate(G4ThreeVector(0.,0.,1.), fPhi);
306 
307  if(IsBent() ){
308  G4ThreeVector vBendingAngle = ComputeBendingAngle(vPosition);
309  vMomentum.rotate(G4ThreeVector(1.,0.,0.), vBendingAngle.z())
310  .rotate(G4ThreeVector(0.,1.,0.),vBendingAngle.x());
311  }
312 
313  return vMomentum;
314 }
315 
316 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
317 
320  G4ThreeVector& vPosition){
321  vMomentum.rotate(G4ThreeVector(0.,0.,1.), -fPhi)
322  .rotate(G4ThreeVector(0.,1.,0.), -fTheta)
323  .rotate(G4ThreeVector(1.,0.,0.), fOmega);
324 
325  if(IsBent() ){
326  G4ThreeVector vBendingAngle = ComputeBendingAngle(vPosition);
327  vMomentum.rotate(G4ThreeVector(0.,1.,0.), -vBendingAngle.x())
328  .rotate(G4ThreeVector(1.,0.,0.), -vBendingAngle.z());
329  }
330 
331  return vMomentum;
332 }
333 
334 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
335 
337  G4ThreeVector dir = G4ThreeVector(0.,0.,1.);
339  vPosition);
340 }
341 
342 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
343 
345  fUnitCell = cell;
346 }
347 
348 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
349 
351  return fUnitCell;
352 }
353 
354 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
355 
357  fCurvatureRadius = cr;
358 }
359 
360 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
361 
363  return fCurvatureRadius;
364 }
365 
366 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
367 
369  if(fCurvatureRadius.x() != 0.) {
370  return true;
371  }
372  else {
373  return false;
374  }
375 }
376 
377 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
378 
380 
381  G4double vAngleX = 0.;
382  G4double vAngleY = 0.;
383 
384  if(GetCurvatureRadius().x() != 0){
385  vAngleX = vPosition.phi();
386  }
387 
388  return G4ThreeVector(vAngleX,0.,vAngleY);
389 }
390 
391 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
392 
394  return fLattice;
395 }
396 
397 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
398 
400  if(vIndex<3 && vIndex>=0)
401  return fMillerOrientation[vIndex];
402  else return -1;
403 }
404 
405 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
406 
408  return G4ThreeVector(fPhi,fTheta,fOmega);
409 }
410 
411 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
412 
414  G4double vInterplanarPeriod =
416  GetMiller(1),
417  GetMiller(2));
418  return vInterplanarPeriod;
419 }
420 
421 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
void SetLDOS(double)
G4double GetAnhDecConstant()
void SetXLogicalLattice(XLogicalLattice *)
void SetFTDOS(double)
G4int GetMiller(G4int)
double MapKtoV(int, G4ThreeVector)
void SetThermalVibrationAmplitude(G4double)
G4ThreeVector GetLatticeAngles()
CLHEP::Hep3Vector G4ThreeVector
double x() const
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
G4ThreeVector GetCurvatureRadius()
G4ThreeVector ComputeBendingAngle(G4ThreeVector &)
tuple x
Definition: test.py:50
int G4int
Definition: G4Types.hh:78
double z() const
G4ThreeVector MapKtoVDir(int, G4ThreeVector)
double MapKtoV(int, G4ThreeVector)
void SetUnitCell(XUnitCell *)
G4ThreeVector ProjectMomentumVectorFromWorldToLattice(G4ThreeVector &, G4ThreeVector &)
G4AffineTransform fLocalToGlobal
tuple b
Definition: test.py:12
G4AffineTransform & Invert()
void SetSTDOS(double)
G4double GetAnhDecConstant()
bool G4bool
Definition: G4Types.hh:79
double phi() const
const G4int n
G4double GetThermalVibrationAmplitude()
XUnitCell * GetXUnitCell()
G4ThreeVector ProjectMomentumVectorFromLatticeToWorld(G4ThreeVector &, G4ThreeVector &)
void SetMillerOrientation(int, int, int)
void SetCurvatureRadius(G4ThreeVector)
void SetLatticeOrientation(G4double, G4double)
Hep3Vector unit() const
XLogicalLattice * GetLogicalLattice()
G4double ComputeDirectPeriod(G4int, G4int, G4int)
Definition: XUnitCell.cc:222
G4RotationMatrix * GetObjectRotation() const
void SetDynamicalConstants(double, double, double, double)
static constexpr double angstrom
Definition: G4SIunits.hh:102
Definition of the XPhysicalLattice class.
void SetAnhDecConstant(G4double)
double G4double
Definition: G4Types.hh:76
void SetPhysicalVolume(G4VPhysicalVolume *)
void SetScatteringConstant(G4double)
G4ThreeVector GetLatticeDirection(G4ThreeVector &)
G4AffineTransform fGlobalToLocal
G4double ComputeInterplanarPeriod()
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
G4ThreeVector MapKtoVDir(int, G4ThreeVector)
G4VPhysicalVolume * GetVolume()
G4double GetScatteringConstant()
G4double GetScatteringConstant()