Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
XLogicalLattice Class Reference

#include <XLogicalLattice.hh>

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)
 

Detailed Description

Definition at line 43 of file XLogicalLattice.hh.

Constructor & Destructor Documentation

XLogicalLattice::XLogicalLattice ( )

Definition at line 35 of file XLogicalLattice.cc.

35  {
36  fVresTheta=0;
37  fVresPhi=0;
38  fDresTheta=0;
39  fDresPhi=0;
40  fA=0;
41  fB=0;
42  fDosL=0;
43  fDosST=0;
44  fDosFT=0;
45  fBeta=0;
46  fLambda=0;
47  fGamma=0;
48  fMu=0;
49 }
XLogicalLattice::~XLogicalLattice ( )

Definition at line 54 of file XLogicalLattice.cc.

54 {;}

Member Function Documentation

G4double XLogicalLattice::GetAnhDecConstant ( )

Definition at line 158 of file XLogicalLattice.cc.

159 {
160  return fA;
161 }

Here is the caller graph for this function:

double XLogicalLattice::GetBeta ( )

Definition at line 118 of file XLogicalLattice.cc.

119 {
120  return fBeta;
121 }

Here is the caller graph for this function:

double XLogicalLattice::GetFTDOS ( )

Definition at line 182 of file XLogicalLattice.cc.

183 {
184  return fDosFT;
185 }

Here is the caller graph for this function:

double XLogicalLattice::GetGamma ( )

Definition at line 126 of file XLogicalLattice.cc.

127 {
128  return fGamma;
129 }

Here is the caller graph for this function:

double XLogicalLattice::GetLambda ( )

Definition at line 134 of file XLogicalLattice.cc.

135 {
136  return fLambda;
137 }

Here is the caller graph for this function:

double XLogicalLattice::GetLDOS ( )

Definition at line 166 of file XLogicalLattice.cc.

167 {
168  return fDosL;
169 }

Here is the caller graph for this function:

double XLogicalLattice::GetMu ( )

Definition at line 142 of file XLogicalLattice.cc.

143 {
144  return fMu;
145 }

Here is the caller graph for this function:

G4double XLogicalLattice::GetScatteringConstant ( )

Definition at line 150 of file XLogicalLattice.cc.

151 {
152  return fB;
153 }

Here is the caller graph for this function:

double XLogicalLattice::GetSTDOS ( )

Definition at line 174 of file XLogicalLattice.cc.

175 {
176  return fDosST;
177 }

Here is the caller graph for this function:

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

Definition at line 230 of file XLogicalLattice.cc.

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

Definition at line 193 of file XLogicalLattice.cc.

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

Definition at line 268 of file XLogicalLattice.cc.

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

Here is the call graph for this function:

Here is the caller graph for this function:

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

Definition at line 299 of file XLogicalLattice.cc.

300 {
301  //Given the phonon wave vector k, phonon physical volume Vol
302  //and polarizationState(0=LON, 1=FT, 2=ST),
303  //returns phonon propagation direction as dimensionless unit vector
304 
305  double theta, phi, tRes, pRes;
306 
307  tRes=pi/(fDresTheta-1);//The summant "-1" is required:index=[0:array length-1]
308  pRes=2*pi/(fDresPhi-1);
309 
310  theta=k.getTheta();
311  phi=k.getPhi();
312 
313  if(theta>pi) theta=theta-pi;
314  //phi=[0 to 2 pi] in accordance with DMC //if(phi>pi/2) phi=phi-pi/2;
315  if(phi<0) phi = phi + 2*pi;
316 
317  G4int iTheta = int(theta/tRes+0.5);
318  G4int iPhi = int(phi/pRes+0.5);
319 
320 
321  G4ThreeVector v(fN_map[polarizationState][iTheta][iPhi][0],
322  fN_map[polarizationState][iTheta][iPhi][1],
323  fN_map[polarizationState][iTheta][iPhi][2]);
324 
325 
327  //v.Set(sin(theta)*cos(phi),sin(theta)*sin(phi),cos(theta));
328  //v.Set(k.getX(), k.getY(), k.getZ());
329 
330  return v.unit();
331 }
double getTheta() const
int G4int
Definition: G4Types.hh:78
typedef int(XMLCALL *XML_NotStandaloneHandler)(void *userData)
tuple v
Definition: test.py:18
double getPhi() const
static constexpr double pi
Definition: G4SIunits.hh:75

Here is the call graph for this function:

Here is the caller graph for this function:

void XLogicalLattice::SetAnhDecConstant ( G4double  a)

Definition at line 82 of file XLogicalLattice.cc.

83 {
84  //Constant governing rate of anharmonic down conversion of L-phonons,
85  //for use with XPhononDownconversionProcess
86  fA=a;
87 }
std::vector< ExP01TrackerHit * > a
Definition: ExP01Classes.hh:33
void XLogicalLattice::SetDynamicalConstants ( double  Beta,
double  Gamma,
double  Lambda,
double  Mu 
)

Definition at line 57 of file XLogicalLattice.cc.

61 {
62  fBeta=Beta;
63  fGamma=Gamma;
64  fLambda=Lambda;
65  fMu=Mu;
66 
67 }
void XLogicalLattice::SetFTDOS ( double  FTDOS)

Definition at line 109 of file XLogicalLattice.cc.

110 {
111  //Fast transverse phonon density of states
112  fDosFT=FTDOS;
113 }
void XLogicalLattice::SetLDOS ( double  LDOS)

Definition at line 91 of file XLogicalLattice.cc.

92 {
93  //Longitudinal phonon density of states
94  fDosL=LDOS;
95 }
void XLogicalLattice::SetScatteringConstant ( G4double  b)

Definition at line 72 of file XLogicalLattice.cc.

73 {
74  //Constant governing the rate of isotope scattering, for use with
75  //XPhononScatteringProcess
76  fB=b;
77 }
tuple b
Definition: test.py:12
void XLogicalLattice::SetSTDOS ( double  STDOS)

Definition at line 100 of file XLogicalLattice.cc.

101 {
102  //Slow transverse phonon density of states
103  fDosST=STDOS;
104 }

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