Geant4  10.02.p03
XLogicalLattice Class Reference

#include <XLogicalLattice.hh>

Collaboration diagram for XLogicalLattice:

Public Member Functions

 XLogicalLattice ()
 
 ~XLogicalLattice ()
 
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 ()
 
bool LoadMap (int, int, int, string)
 
bool Load_NMap (int, int, int, string)
 
double MapKtoV (int, G4ThreeVector)
 
G4ThreeVector MapKtoVDir (int, G4ThreeVector)
 

Private Attributes

int fVresTheta
 
int fVresPhi
 
int fDresTheta
 
int fDresPhi
 
double fMap [3][MAXRES][MAXRES]
 
double fN_map [3][MAXRES][MAXRES][3]
 
double fA
 
double fB
 
double fDosL
 
double fDosST
 
double fDosFT
 
double fBeta
 
double fGamma
 
double fLambda
 
double fMu
 
ifstream fMapFile
 
int fThetaRes
 
int fPhiRes
 

Detailed Description

Definition at line 40 of file XLogicalLattice.hh.

Constructor & Destructor Documentation

◆ XLogicalLattice()

XLogicalLattice::XLogicalLattice ( )

Definition at line 32 of file XLogicalLattice.cc.

32  {
33  fVresTheta=0;
34  fVresPhi=0;
35  fDresTheta=0;
36  fDresPhi=0;
37  fA=0;
38  fB=0;
39  fDosL=0;
40  fDosST=0;
41  fDosFT=0;
42  fBeta=0;
43  fLambda=0;
44  fGamma=0;
45  fMu=0;
46 }

◆ ~XLogicalLattice()

XLogicalLattice::~XLogicalLattice ( )

Definition at line 51 of file XLogicalLattice.cc.

51 {;}

Member Function Documentation

◆ GetAnhDecConstant()

G4double XLogicalLattice::GetAnhDecConstant ( )

Definition at line 155 of file XLogicalLattice.cc.

156 {
157  return fA;
158 }
Here is the caller graph for this function:

◆ GetBeta()

double XLogicalLattice::GetBeta ( )

Definition at line 115 of file XLogicalLattice.cc.

116 {
117  return fBeta;
118 }
Here is the caller graph for this function:

◆ GetFTDOS()

double XLogicalLattice::GetFTDOS ( )

Definition at line 179 of file XLogicalLattice.cc.

180 {
181  return fDosFT;
182 }
Here is the caller graph for this function:

◆ GetGamma()

double XLogicalLattice::GetGamma ( )

Definition at line 123 of file XLogicalLattice.cc.

124 {
125  return fGamma;
126 }
Here is the caller graph for this function:

◆ GetLambda()

double XLogicalLattice::GetLambda ( )

Definition at line 131 of file XLogicalLattice.cc.

132 {
133  return fLambda;
134 }
Here is the caller graph for this function:

◆ GetLDOS()

double XLogicalLattice::GetLDOS ( )

Definition at line 163 of file XLogicalLattice.cc.

164 {
165  return fDosL;
166 }
Here is the caller graph for this function:

◆ GetMu()

double XLogicalLattice::GetMu ( )

Definition at line 139 of file XLogicalLattice.cc.

140 {
141  return fMu;
142 }
Here is the caller graph for this function:

◆ GetScatteringConstant()

G4double XLogicalLattice::GetScatteringConstant ( )

Definition at line 147 of file XLogicalLattice.cc.

148 {
149  return fB;
150 }
Here is the caller graph for this function:

◆ GetSTDOS()

double XLogicalLattice::GetSTDOS ( )

Definition at line 171 of file XLogicalLattice.cc.

172 {
173  return fDosST;
174 }
Here is the caller graph for this function:

◆ Load_NMap()

bool XLogicalLattice::Load_NMap ( int  tRes,
int  pRes,
int  polarizationState,
string  m 
)

Definition at line 227 of file XLogicalLattice.cc.

231 {
232  if((tRes>MAXRES)||(pRes>MAXRES)){
233  G4cout<<"\nk-V map exceeds maximum resolution of "<<MAXRES<<
234  " by "<<MAXRES<<" terminating."<<endl;
235  return false; //terminate if resolution out of bounds.
236  }
237 
238 
239  fMapFile.clear();
240 
241  fThetaRes=tRes;
242  fPhiRes=pRes;
243  fMapFile.open(m.data());
244  if(!fMapFile.is_open()) return false;
245 
246  for(int theta = 0; theta<fThetaRes; theta++){
247  for(int phi = 0; phi<fPhiRes; phi++){
248  for(int coord = 0; coord<3; coord++){
249  fMapFile>>fN_map[polarizationState][theta][phi][coord];
250  }
251  }
252  }
253  G4cout<<"\n";
254  fMapFile.close();
255  G4cout<<"\nXLogicalLattice::Load_NMap() sucessful"<<endl;
256  fDresTheta=tRes; //store map dimesnions
257  fDresPhi=pRes;
258  return true;
259 }
#define MAXRES
G4GLOB_DLL std::ostream G4cout
double fN_map[3][MAXRES][MAXRES][3]
static const double m
Definition: G4SIunits.hh:128

◆ LoadMap()

bool XLogicalLattice::LoadMap ( int  tRes,
int  pRes,
int  polarizationState,
string  m 
)

Definition at line 190 of file XLogicalLattice.cc.

194 {
195  if((tRes>MAXRES)||(pRes>MAXRES)){
196  G4cout<<"\nk-V fMap exceeds maximum resolution of "<<MAXRES<<
197  " by "<<MAXRES<<". terminating."<<endl;
198  return false; //terminate if resolution out of bounds.
199  }
200 
201 
202  fMapFile.clear();
203 
204  fThetaRes=tRes;
205  fPhiRes=pRes;
206  fMapFile.open(m.data());
207  if(!fMapFile.is_open()) return false;
208 
209  for(int theta = 0; theta<fThetaRes; theta++){
210  for(int phi = 0; phi<fPhiRes; phi++){
211  fMapFile>>fMap[polarizationState][theta][phi];
212  }
213  }
214  fMapFile.close();
215  G4cout<<"\nXLogicalLattice::LoadMap() sucessful (Vg scalars)."<<endl;
216  fVresTheta=tRes; //store map dimensions
217  fVresPhi=pRes;
218  return true;
219 }
#define MAXRES
double fMap[3][MAXRES][MAXRES]
G4GLOB_DLL std::ostream G4cout
static const double m
Definition: G4SIunits.hh:128

◆ MapKtoV()

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

Definition at line 265 of file XLogicalLattice.cc.

266 {
267  //Given the phonon wave vector k, phonon physical volume Vol
268  //and polarizationState(0=LON, 1=FT, 2=ST),
269  //returns phonon velocity in m/s
270 
271  double theta, phi, tRes, pRes;
272 
273  tRes=pi/(fVresTheta);//The summant "-1" is required:index=[0:array length-1]
274  pRes=2*pi/(fVresPhi);
275 
276 
277  theta=k.getTheta();
278  phi=k.getPhi();
279 
280  if(phi<0) phi = phi + 2*pi;
281  if(theta>pi) theta=theta-pi;
282  //phi=[0 to 2 pi] in accordance with DMC //if(phi>pi/2) phi=phi-pi/2;
283  if(fMap[polarizationState][int(theta/tRes)][int(phi/pRes)]==0){
284  G4cout<<"\nFound v=0 for polarization "<<polarizationState
285  <<" theta "<<theta<<" phi "<<phi<< " translating to map coords "
286  << "theta "<< int(theta/tRes) << " phi " << int(phi/pRes)<<endl;
287  }
288 
289  return fMap[polarizationState][int(theta/tRes)][int(phi/pRes)];
290 
291 }
double getTheta() const
double getPhi() const
double fMap[3][MAXRES][MAXRES]
G4GLOB_DLL std::ostream G4cout
static const double pi
Definition: G4SIunits.hh:74
Here is the call graph for this function:
Here is the caller graph for this function:

◆ MapKtoVDir()

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

Definition at line 296 of file XLogicalLattice.cc.

297 {
298  //Given the phonon wave vector k, phonon physical volume Vol
299  //and polarizationState(0=LON, 1=FT, 2=ST),
300  //returns phonon propagation direction as dimensionless unit vector
301 
302  double theta, phi, tRes, pRes;
303 
304  tRes=pi/(fDresTheta-1);//The summant "-1" is required:index=[0:array length-1]
305  pRes=2*pi/(fDresPhi-1);
306 
307  theta=k.getTheta();
308  phi=k.getPhi();
309 
310  if(theta>pi) theta=theta-pi;
311  //phi=[0 to 2 pi] in accordance with DMC //if(phi>pi/2) phi=phi-pi/2;
312  if(phi<0) phi = phi + 2*pi;
313 
314  G4int iTheta = int(theta/tRes+0.5);
315  G4int iPhi = int(phi/pRes+0.5);
316 
317 
318  G4ThreeVector v(fN_map[polarizationState][iTheta][iPhi][0],
319  fN_map[polarizationState][iTheta][iPhi][1],
320  fN_map[polarizationState][iTheta][iPhi][2]);
321 
322 
324  //v.Set(sin(theta)*cos(phi),sin(theta)*sin(phi),cos(theta));
325  //v.Set(k.getX(), k.getY(), k.getZ());
326 
327  return v.unit();
328 }
double getTheta() const
double getPhi() const
int G4int
Definition: G4Types.hh:78
double fN_map[3][MAXRES][MAXRES][3]
static const double pi
Definition: G4SIunits.hh:74
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAnhDecConstant()

void XLogicalLattice::SetAnhDecConstant ( G4double  a)

Definition at line 79 of file XLogicalLattice.cc.

80 {
81  //Constant governing rate of anharmonic down conversion of L-phonons,
82  //for use with XPhononDownconversionProcess
83  fA=a;
84 }

◆ SetDynamicalConstants()

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

Definition at line 54 of file XLogicalLattice.cc.

58 {
59  fBeta=Beta;
60  fGamma=Gamma;
61  fLambda=Lambda;
62  fMu=Mu;
63 
64 }

◆ SetFTDOS()

void XLogicalLattice::SetFTDOS ( double  FTDOS)

Definition at line 106 of file XLogicalLattice.cc.

107 {
108  //Fast transverse phonon density of states
109  fDosFT=FTDOS;
110 }

◆ SetLDOS()

void XLogicalLattice::SetLDOS ( double  LDOS)

Definition at line 88 of file XLogicalLattice.cc.

89 {
90  //Longitudinal phonon density of states
91  fDosL=LDOS;
92 }

◆ SetScatteringConstant()

void XLogicalLattice::SetScatteringConstant ( G4double  b)

Definition at line 69 of file XLogicalLattice.cc.

70 {
71  //Constant governing the rate of isotope scattering, for use with
72  //XPhononScatteringProcess
73  fB=b;
74 }
Here is the caller graph for this function:

◆ SetSTDOS()

void XLogicalLattice::SetSTDOS ( double  STDOS)

Definition at line 97 of file XLogicalLattice.cc.

98 {
99  //Slow transverse phonon density of states
100  fDosST=STDOS;
101 }

Member Data Documentation

◆ fA

double XLogicalLattice::fA
private

Definition at line 53 of file XLogicalLattice.hh.

◆ fB

double XLogicalLattice::fB
private

Definition at line 54 of file XLogicalLattice.hh.

◆ fBeta

double XLogicalLattice::fBeta
private

Definition at line 58 of file XLogicalLattice.hh.

◆ fDosFT

double XLogicalLattice::fDosFT
private

Definition at line 57 of file XLogicalLattice.hh.

◆ fDosL

double XLogicalLattice::fDosL
private

Definition at line 55 of file XLogicalLattice.hh.

◆ fDosST

double XLogicalLattice::fDosST
private

Definition at line 56 of file XLogicalLattice.hh.

◆ fDresPhi

int XLogicalLattice::fDresPhi
private

Definition at line 47 of file XLogicalLattice.hh.

◆ fDresTheta

int XLogicalLattice::fDresTheta
private

Definition at line 46 of file XLogicalLattice.hh.

◆ fGamma

double XLogicalLattice::fGamma
private

Definition at line 58 of file XLogicalLattice.hh.

◆ fLambda

double XLogicalLattice::fLambda
private

Definition at line 58 of file XLogicalLattice.hh.

◆ fMap

double XLogicalLattice::fMap[3][MAXRES][MAXRES]
private

Definition at line 49 of file XLogicalLattice.hh.

◆ fMapFile

ifstream XLogicalLattice::fMapFile
private

Definition at line 61 of file XLogicalLattice.hh.

◆ fMu

double XLogicalLattice::fMu
private

Definition at line 58 of file XLogicalLattice.hh.

◆ fN_map

double XLogicalLattice::fN_map[3][MAXRES][MAXRES][3]
private

Definition at line 50 of file XLogicalLattice.hh.

◆ fPhiRes

int XLogicalLattice::fPhiRes
private

Definition at line 62 of file XLogicalLattice.hh.

◆ fThetaRes

int XLogicalLattice::fThetaRes
private

Definition at line 62 of file XLogicalLattice.hh.

◆ fVresPhi

int XLogicalLattice::fVresPhi
private

Definition at line 45 of file XLogicalLattice.hh.

◆ fVresTheta

int XLogicalLattice::fVresTheta
private

Definition at line 44 of file XLogicalLattice.hh.


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