Geant4  10.02.p03
XPhysicalLattice Class Reference

#include <XPhysicalLattice.hh>

Collaboration diagram for XPhysicalLattice:

Public Member Functions

void SetDynamicalConstants (double, double, double, double)
 
void SetScatteringConstant (G4double)
 
void SetAnhDecConstant (G4double)
 
void SetLDOS (double)
 
void SetSTDOS (double)
 
void SetFTDOS (double)
 
double GetBeta ()
 
double GetGamma ()
 
double GetLambda ()
 
double GetMu ()
 
G4double GetScatteringConstant ()
 
G4double GetAnhDecConstant ()
 
double GetLDOS ()
 
double GetSTDOS ()
 
double GetFTDOS ()
 
 XPhysicalLattice ()
 
 XPhysicalLattice (G4VPhysicalVolume *, XLogicalLattice *)
 
 ~XPhysicalLattice ()
 
double MapKtoV (int, G4ThreeVector)
 
G4ThreeVector MapKtoVDir (int, G4ThreeVector)
 
G4VPhysicalVolumeGetVolume ()
 
void SetPhysicalVolume (G4VPhysicalVolume *)
 
void SetXLogicalLattice (XLogicalLattice *)
 
void SetLatticeOrientation (G4double, G4double)
 
void SetLatticeOrientation (G4double, G4double, G4double)
 
void SetMillerOrientation (int, int, int)
 
void SetUnitCell (XUnitCell *)
 
G4ThreeVector ProjectMomentumVectorFromWorldToLattice (G4ThreeVector &, G4ThreeVector &)
 
G4ThreeVector ProjectMomentumVectorFromLatticeToWorld (G4ThreeVector &, G4ThreeVector &)
 
G4ThreeVector GetLatticeDirection (G4ThreeVector &)
 
XUnitCellGetXUnitCell ()
 
XLogicalLatticeGetLogicalLattice ()
 
G4int GetMiller (G4int)
 
G4ThreeVector GetLatticeAngles ()
 
G4ThreeVector GetCurvatureRadius ()
 
void SetCurvatureRadius (G4ThreeVector)
 
G4ThreeVector ComputeBendingAngle (G4ThreeVector &)
 
G4bool IsBent ()
 
G4double ComputeInterplanarPeriod ()
 
void SetThermalVibrationAmplitude (G4double)
 
G4double GetThermalVibrationAmplitude ()
 

Public Attributes

G4AffineTransform fLocalToGlobal
 
G4AffineTransform fGlobalToLocal
 

Private Attributes

G4double fTheta
 
G4double fPhi
 
G4double fOmega
 
XLogicalLatticefLattice
 
G4VPhysicalVolumefVolume
 
double fA
 
double fB
 
double fDosL
 
double fDosST
 
double fDosFT
 
double fBeta
 
double fGamma
 
double fLambda
 
double fMu
 
G4ThreeVector fCurvatureRadius
 
G4double fThermalVibrationAmplitude
 
G4int fMillerOrientation [3]
 
XUnitCellfUnitCell
 

Detailed Description

Definition at line 37 of file XPhysicalLattice.hh.

Constructor & Destructor Documentation

◆ XPhysicalLattice() [1/2]

XPhysicalLattice::XPhysicalLattice ( )

Definition at line 35 of file XPhysicalLattice.cc.

35  {
36  fLattice=NULL;
37  fVolume=NULL;
38  fTheta=0.;
39  fPhi=0.;
40  fOmega=0.;
41 
42 
43  fCurvatureRadius = G4ThreeVector(0.,0.,0.); // if cr = 0 == no bending
44  fThermalVibrationAmplitude = 0.1 * angstrom; // no physical meaning
45  fMillerOrientation[0] = 2;
46  fMillerOrientation[1] = 2;
47  fMillerOrientation[2] = 0;
48 }
G4double fThermalVibrationAmplitude
CLHEP::Hep3Vector G4ThreeVector
XLogicalLattice * fLattice
G4VPhysicalVolume * fVolume
G4ThreeVector fCurvatureRadius
static const double angstrom
Definition: G4SIunits.hh:101

◆ XPhysicalLattice() [2/2]

XPhysicalLattice::XPhysicalLattice ( G4VPhysicalVolume Vol,
XLogicalLattice Lat 
)

Definition at line 52 of file XPhysicalLattice.cc.

53  {
54  fLattice=Lat;
55  fVolume=Vol;
64  fMu=fLattice->GetMu();
65 
67 
70 }
XLogicalLattice * fLattice
G4AffineTransform fLocalToGlobal
G4AffineTransform & Invert()
G4double GetAnhDecConstant()
G4VPhysicalVolume * fVolume
G4RotationMatrix * GetObjectRotation() const
G4AffineTransform fGlobalToLocal
G4double GetScatteringConstant()
Here is the call graph for this function:

◆ ~XPhysicalLattice()

XPhysicalLattice::~XPhysicalLattice ( )

Definition at line 74 of file XPhysicalLattice.cc.

74 {;}

Member Function Documentation

◆ ComputeBendingAngle()

G4ThreeVector XPhysicalLattice::ComputeBendingAngle ( G4ThreeVector vPosition)

Definition at line 376 of file XPhysicalLattice.cc.

376  {
377 
378  G4double vAngleX = 0.;
379  G4double vAngleY = 0.;
380 
381  if(GetCurvatureRadius().x() != 0){
382  vAngleX = vPosition.phi();
383  }
384 
385  return G4ThreeVector(vAngleX,0.,vAngleY);
386 }
double phi() const
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector GetCurvatureRadius()
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeInterplanarPeriod()

G4double XPhysicalLattice::ComputeInterplanarPeriod ( )

Definition at line 410 of file XPhysicalLattice.cc.

410  {
411  G4double vInterplanarPeriod =
413  GetMiller(1),
414  GetMiller(2));
415  return vInterplanarPeriod;
416 }
G4int GetMiller(G4int)
XUnitCell * GetXUnitCell()
G4double ComputeDirectPeriod(G4int, G4int, G4int)
Definition: XUnitCell.cc:219
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAnhDecConstant()

G4double XPhysicalLattice::GetAnhDecConstant ( )

Definition at line 162 of file XPhysicalLattice.cc.

163 {
164  return fA;
165 }

◆ GetBeta()

double XPhysicalLattice::GetBeta ( )

Definition at line 126 of file XPhysicalLattice.cc.

126  {
127  return fBeta;
128 }

◆ GetCurvatureRadius()

G4ThreeVector XPhysicalLattice::GetCurvatureRadius ( )

Definition at line 359 of file XPhysicalLattice.cc.

359  {
360  return fCurvatureRadius;
361 }
G4ThreeVector fCurvatureRadius
Here is the caller graph for this function:

◆ GetFTDOS()

double XPhysicalLattice::GetFTDOS ( )

Definition at line 184 of file XPhysicalLattice.cc.

185 {
186  return fDosFT;
187 }

◆ GetGamma()

double XPhysicalLattice::GetGamma ( )

Definition at line 133 of file XPhysicalLattice.cc.

133  {
134  return fGamma;
135 }

◆ GetLambda()

double XPhysicalLattice::GetLambda ( )

Definition at line 139 of file XPhysicalLattice.cc.

139  {
140  return fLambda;
141 }

◆ GetLatticeAngles()

G4ThreeVector XPhysicalLattice::GetLatticeAngles ( )

Definition at line 404 of file XPhysicalLattice.cc.

404  {
406 }
CLHEP::Hep3Vector G4ThreeVector

◆ GetLatticeDirection()

G4ThreeVector XPhysicalLattice::GetLatticeDirection ( G4ThreeVector vPosition)

Definition at line 333 of file XPhysicalLattice.cc.

333  {
334  G4ThreeVector dir = G4ThreeVector(0.,0.,1.);
336  vPosition);
337 }
CLHEP::Hep3Vector G4ThreeVector
TDirectory * dir
G4ThreeVector ProjectMomentumVectorFromLatticeToWorld(G4ThreeVector &, G4ThreeVector &)
Here is the call graph for this function:

◆ GetLDOS()

double XPhysicalLattice::GetLDOS ( )

Definition at line 170 of file XPhysicalLattice.cc.

170  {
171  return fDosL;
172 }

◆ GetLogicalLattice()

XLogicalLattice * XPhysicalLattice::GetLogicalLattice ( )

Definition at line 390 of file XPhysicalLattice.cc.

390  {
391  return fLattice;
392 }
XLogicalLattice * fLattice
Here is the caller graph for this function:

◆ GetMiller()

G4int XPhysicalLattice::GetMiller ( G4int  vIndex)

Definition at line 396 of file XPhysicalLattice.cc.

396  {
397  if(vIndex<3 && vIndex>=0)
398  return fMillerOrientation[vIndex];
399  else return -1;
400 }
Here is the caller graph for this function:

◆ GetMu()

double XPhysicalLattice::GetMu ( )

Definition at line 146 of file XPhysicalLattice.cc.

147 {
148  return fMu;
149 }

◆ GetScatteringConstant()

G4double XPhysicalLattice::GetScatteringConstant ( )

Definition at line 154 of file XPhysicalLattice.cc.

155 {
156  return fB;
157 }

◆ GetSTDOS()

double XPhysicalLattice::GetSTDOS ( )

Definition at line 177 of file XPhysicalLattice.cc.

177  {
178  return fDosST;
179 }

◆ GetThermalVibrationAmplitude()

G4double XPhysicalLattice::GetThermalVibrationAmplitude ( )

Definition at line 291 of file XPhysicalLattice.cc.

291  {
293 }
G4double fThermalVibrationAmplitude
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetVolume()

G4VPhysicalVolume * XPhysicalLattice::GetVolume ( )

Definition at line 228 of file XPhysicalLattice.cc.

228  {
229  return fVolume;
230 }
G4VPhysicalVolume * fVolume
Here is the caller graph for this function:

◆ GetXUnitCell()

XUnitCell * XPhysicalLattice::GetXUnitCell ( )

Definition at line 347 of file XPhysicalLattice.cc.

347  {
348  return fUnitCell;
349 }
Here is the caller graph for this function:

◆ IsBent()

G4bool XPhysicalLattice::IsBent ( )

Definition at line 365 of file XPhysicalLattice.cc.

365  {
366  if(fCurvatureRadius.x() != 0.) {
367  return true;
368  }
369  else {
370  return false;
371  }
372 }
double x() const
G4ThreeVector fCurvatureRadius
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MapKtoV()

double XPhysicalLattice::MapKtoV ( int  polarizationState,
G4ThreeVector  k 
)

Definition at line 196 of file XPhysicalLattice.cc.

196  {
197  double groupVelocity;
198 
199  k.rotate(G4ThreeVector(0,1,0), fTheta).rotate(G4ThreeVector(0,0,1), fPhi);
200  groupVelocity = fLattice->MapKtoV(polarizationState, k);
201  k.rotate(G4ThreeVector(0,0,1), -fPhi).rotate(G4ThreeVector(0,1,0), -fTheta);
202 
203  return groupVelocity;
204 }
CLHEP::Hep3Vector G4ThreeVector
XLogicalLattice * fLattice
double MapKtoV(int, G4ThreeVector)
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MapKtoVDir()

G4ThreeVector XPhysicalLattice::MapKtoVDir ( int  polarizationState,
G4ThreeVector  k 
)

Definition at line 213 of file XPhysicalLattice.cc.

214  {
215 
216  G4ThreeVector GroupVelocity;
217 
218  k=k.rotate(G4ThreeVector(0,1,0), fTheta).rotate(G4ThreeVector(0,0,1), fPhi);
219  GroupVelocity = fLattice->MapKtoVDir(polarizationState, k);
220 
221  return GroupVelocity.rotate(G4ThreeVector(0,0,1), -fPhi)
222  .rotate(G4ThreeVector(0,1,0), -fTheta).unit();
223 }
CLHEP::Hep3Vector G4ThreeVector
XLogicalLattice * fLattice
G4ThreeVector MapKtoVDir(int, G4ThreeVector)
Hep3Vector unit() const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProjectMomentumVectorFromLatticeToWorld()

G4ThreeVector XPhysicalLattice::ProjectMomentumVectorFromLatticeToWorld ( G4ThreeVector vMomentum,
G4ThreeVector vPosition 
)

Definition at line 316 of file XPhysicalLattice.cc.

317  {
318  vMomentum.rotate(G4ThreeVector(0.,0.,1.), -fPhi)
319  .rotate(G4ThreeVector(0.,1.,0.), -fTheta)
320  .rotate(G4ThreeVector(1.,0.,0.), fOmega);
321 
322  if(IsBent() ){
323  G4ThreeVector vBendingAngle = ComputeBendingAngle(vPosition);
324  vMomentum.rotate(G4ThreeVector(0.,1.,0.), -vBendingAngle.x())
325  .rotate(G4ThreeVector(1.,0.,0.), -vBendingAngle.z());
326  }
327 
328  return vMomentum;
329 }
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector ComputeBendingAngle(G4ThreeVector &)
double x() const
double z() const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ProjectMomentumVectorFromWorldToLattice()

G4ThreeVector XPhysicalLattice::ProjectMomentumVectorFromWorldToLattice ( G4ThreeVector vMomentum,
G4ThreeVector vPosition 
)

Definition at line 298 of file XPhysicalLattice.cc.

299  {
300  vMomentum.rotate(G4ThreeVector(1.,0.,0.),fOmega)
301  .rotate(G4ThreeVector(0.,1.,0.), fTheta)
302  .rotate(G4ThreeVector(0.,0.,1.), fPhi);
303 
304  if(IsBent() ){
305  G4ThreeVector vBendingAngle = ComputeBendingAngle(vPosition);
306  vMomentum.rotate(G4ThreeVector(1.,0.,0.), vBendingAngle.z())
307  .rotate(G4ThreeVector(0.,1.,0.),vBendingAngle.x());
308  }
309 
310  return vMomentum;
311 }
CLHEP::Hep3Vector G4ThreeVector
G4ThreeVector ComputeBendingAngle(G4ThreeVector &)
double x() const
double z() const
Hep3Vector & rotate(double, const Hep3Vector &)
Definition: ThreeVectorR.cc:28
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnhDecConstant()

void XPhysicalLattice::SetAnhDecConstant ( G4double  b)

Definition at line 99 of file XPhysicalLattice.cc.

99  {
100  fB=b;
101 }

◆ SetCurvatureRadius()

void XPhysicalLattice::SetCurvatureRadius ( G4ThreeVector  cr)

Definition at line 353 of file XPhysicalLattice.cc.

353  {
354  fCurvatureRadius = cr;
355 }
G4ThreeVector fCurvatureRadius
Here is the caller graph for this function:

◆ SetDynamicalConstants()

void XPhysicalLattice::SetDynamicalConstants ( double  Beta,
double  Gamma,
double  Lambda,
double  Mu 
)

Definition at line 78 of file XPhysicalLattice.cc.

82 {
83  fBeta=Beta;
84  fGamma=Gamma;
85  fLambda=Lambda;
86  fMu=Mu;
87 }

◆ SetFTDOS()

void XPhysicalLattice::SetFTDOS ( double  FTDOS)

Definition at line 119 of file XPhysicalLattice.cc.

119  {
120  fDosFT = FTDOS;
121 }

◆ SetLatticeOrientation() [1/2]

void XPhysicalLattice::SetLatticeOrientation ( G4double  t_rot,
G4double  p_rot 
)

Definition at line 242 of file XPhysicalLattice.cc.

242  {
243  fTheta=t_rot;
244  fPhi= p_rot;
245 }
Here is the caller graph for this function:

◆ SetLatticeOrientation() [2/2]

void XPhysicalLattice::SetLatticeOrientation ( G4double  t_rot,
G4double  o_rot,
G4double  p_rot 
)

Definition at line 249 of file XPhysicalLattice.cc.

251  {
252  fTheta = t_rot;
253  fOmega = o_rot;
254  fPhi = p_rot;
255 }

◆ SetLDOS()

void XPhysicalLattice::SetLDOS ( double  LDOS)

Definition at line 106 of file XPhysicalLattice.cc.

106  {
107  fDosL=LDOS;
108 }

◆ SetMillerOrientation()

void XPhysicalLattice::SetMillerOrientation ( int  l,
int  k,
int  n 
)

Definition at line 259 of file XPhysicalLattice.cc.

259  {
260  //fTheta=pi/2-std::atan2(n+0.000001,l+0.000001)*rad;
261  //fPhi=pi/2-std::atan2(l+0.000001,k+0.000001)*rad;
262 
263  // // // added for channeling // // //
264  fMillerOrientation[0]=l;
265  fMillerOrientation[1]=k;
267  // // // // // // // // // // // // //
268 }
Char_t n[5]
Here is the caller graph for this function:

◆ SetPhysicalVolume()

void XPhysicalLattice::SetPhysicalVolume ( G4VPhysicalVolume Vol)

Definition at line 235 of file XPhysicalLattice.cc.

235  {
236  fVolume=Vol;
237 }
G4VPhysicalVolume * fVolume

◆ SetScatteringConstant()

void XPhysicalLattice::SetScatteringConstant ( G4double  a)

Definition at line 92 of file XPhysicalLattice.cc.

92  {
93  fA=a;
94 }

◆ SetSTDOS()

void XPhysicalLattice::SetSTDOS ( double  STDOS)

Definition at line 112 of file XPhysicalLattice.cc.

113 {
114  fDosST = STDOS;
115 }

◆ SetThermalVibrationAmplitude()

void XPhysicalLattice::SetThermalVibrationAmplitude ( G4double  vThermalVibrationAmplitude)

Definition at line 285 of file XPhysicalLattice.cc.

285  {
286  fThermalVibrationAmplitude = vThermalVibrationAmplitude;
287 }
G4double fThermalVibrationAmplitude
Here is the caller graph for this function:

◆ SetUnitCell()

void XPhysicalLattice::SetUnitCell ( XUnitCell cell)

Definition at line 341 of file XPhysicalLattice.cc.

341  {
342  fUnitCell = cell;
343 }
Here is the caller graph for this function:

◆ SetXLogicalLattice()

void XPhysicalLattice::SetXLogicalLattice ( XLogicalLattice Lat)

Definition at line 273 of file XPhysicalLattice.cc.

273  {
274  fLattice=Lat;
275 }
XLogicalLattice * fLattice
Here is the call graph for this function:

Member Data Documentation

◆ fA

double XPhysicalLattice::fA
private

Definition at line 44 of file XPhysicalLattice.hh.

◆ fB

double XPhysicalLattice::fB
private

Definition at line 45 of file XPhysicalLattice.hh.

◆ fBeta

double XPhysicalLattice::fBeta
private

Definition at line 49 of file XPhysicalLattice.hh.

◆ fCurvatureRadius

G4ThreeVector XPhysicalLattice::fCurvatureRadius
private

Definition at line 115 of file XPhysicalLattice.hh.

◆ fDosFT

double XPhysicalLattice::fDosFT
private

Definition at line 48 of file XPhysicalLattice.hh.

◆ fDosL

double XPhysicalLattice::fDosL
private

Definition at line 46 of file XPhysicalLattice.hh.

◆ fDosST

double XPhysicalLattice::fDosST
private

Definition at line 47 of file XPhysicalLattice.hh.

◆ fGamma

double XPhysicalLattice::fGamma
private

Definition at line 49 of file XPhysicalLattice.hh.

◆ fGlobalToLocal

G4AffineTransform XPhysicalLattice::fGlobalToLocal

Definition at line 54 of file XPhysicalLattice.hh.

◆ fLambda

double XPhysicalLattice::fLambda
private

Definition at line 49 of file XPhysicalLattice.hh.

◆ fLattice

XLogicalLattice* XPhysicalLattice::fLattice
private

Definition at line 41 of file XPhysicalLattice.hh.

◆ fLocalToGlobal

G4AffineTransform XPhysicalLattice::fLocalToGlobal

Definition at line 53 of file XPhysicalLattice.hh.

◆ fMillerOrientation

G4int XPhysicalLattice::fMillerOrientation[3]
private

Definition at line 117 of file XPhysicalLattice.hh.

◆ fMu

double XPhysicalLattice::fMu
private

Definition at line 49 of file XPhysicalLattice.hh.

◆ fOmega

G4double XPhysicalLattice::fOmega
private

Definition at line 40 of file XPhysicalLattice.hh.

◆ fPhi

G4double XPhysicalLattice::fPhi
private

Definition at line 40 of file XPhysicalLattice.hh.

◆ fThermalVibrationAmplitude

G4double XPhysicalLattice::fThermalVibrationAmplitude
private

Definition at line 116 of file XPhysicalLattice.hh.

◆ fTheta

G4double XPhysicalLattice::fTheta
private

Definition at line 40 of file XPhysicalLattice.hh.

◆ fUnitCell

XUnitCell* XPhysicalLattice::fUnitCell
private

Definition at line 118 of file XPhysicalLattice.hh.

◆ fVolume

G4VPhysicalVolume* XPhysicalLattice::fVolume
private

Definition at line 42 of file XPhysicalLattice.hh.


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