Geant4  9.6.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
XLogicalLattice.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$
30 //
31 #include "XLogicalLattice.hh"
32 #include "G4PhysicalConstants.hh"
33 #include "Randomize.hh"
34 #include <cmath>
35 
37  fVresTheta=0;
38  fVresPhi=0;
39  fDresTheta=0;
40  fDresPhi=0;
41  fA=0;
42  fB=0;
43  fDosL=0;
44  fDosST=0;
45  fDosFT=0;
46  fBeta=0;
47  fLambda=0;
48  fGamma=0;
49  fMu=0;
50 }
51 
52 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
53 
54 
56 
57 
58 void XLogicalLattice::SetDynamicalConstants(double Beta,double Gamma, double Lambda, double Mu)
59 {
60  fBeta=Beta;
61  fGamma=Gamma;
62  fLambda=Lambda;
63  fMu=Mu;
64 
65 }
66 
67 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
68 
69 
71 {
72  //Constant governing the rate of isotope scattering, for use with
73  //XPhononScatteringProcess
74  fB=b;
75 }
76 
77 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
78 
79 
81 {
82  //Constant governing rate of anharmonic down conversion of L-phonons,
83  //for use with XPhononDownconversionProcess
84  fA=a;
85 }
86 
87 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
88 
89 void XLogicalLattice::SetLDOS(double LDOS)
90 {
91  //Longitudinal phonon density of states
92  fDosL=LDOS;
93 }
94 
95 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96 
97 
98 void XLogicalLattice::SetSTDOS(double STDOS)
99 {
100  //Slow transverse phonon density of states
101  fDosST=STDOS;
102 }
103 
104 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
105 
106 
107 void XLogicalLattice::SetFTDOS(double FTDOS)
108 {
109  //Fast transverse phonon density of states
110  fDosFT=FTDOS;
111 }
112 
113 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
114 
115 
117 {
118  return fBeta;
119 }
120 
121 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
122 
123 
125 {
126  return fGamma;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
130 
131 
133 {
134  return fLambda;
135 }
136 
137 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138 
139 
141 {
142  return fMu;
143 }
144 
145 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
146 
147 
149 {
150  return fB;
151 }
152 
153 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
154 
155 
157 {
158  return fA;
159 }
160 
161 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
162 
163 
165 {
166  return fDosL;
167 }
168 
169 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
170 
171 
173 {
174  return fDosST;
175 }
176 
177 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
178 
179 
181 {
182  return fDosFT;
183 }
184 
185 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
186 
187 
189 //Load map of group velocity scalars (m/s)
191 bool XLogicalLattice::LoadMap(int tRes, int pRes, int polarizationState, string m)
192 {
193  if((tRes>MAXRES)||(pRes>MAXRES)){
194  G4cout<<"\nk-V fMap exceeds maximum resolution of "<<MAXRES<<" by "<<MAXRES<<". terminating."<<endl;
195  return false; //terminate if resolution out of bounds.
196  }
197 
198 
199  fMapFile.clear();
200 
201  fThetaRes=tRes;
202  fPhiRes=pRes;
203  fMapFile.open(m.data());
204  if(!fMapFile.is_open()) return false;
205 
206  for(int theta = 0; theta<fThetaRes; theta++){
207  for(int phi = 0; phi<fPhiRes; phi++){
208  fMapFile>>fMap[polarizationState][theta][phi];
209  //G4cout<<"\nXLogicalLattice::LoadMap:loading map "<<m<<" theta: "<<theta<<" phi: "<<phi<<" map: "<<fMap[polarizationState][theta][phi];
210  }
211  }
212  fMapFile.close();
213  G4cout<<"\nXLogicalLattice::LoadMap() sucessful (Vg scalars)."<<endl;
214  fVresTheta=tRes; //store map dimensions
215  fVresPhi=pRes;
216  return true;
217 }
218 
219 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
220 
221 
223 //Load map of group velocity unit vectors
225 bool XLogicalLattice::Load_NMap(int tRes, int pRes, int polarizationState, string m)
226 {
227  if((tRes>MAXRES)||(pRes>MAXRES)){
228  G4cout<<"\nk-V map exceeds maximum resolution of "<<MAXRES<<" by "<<MAXRES<<" terminating."<<endl;
229  return false; //terminate if resolution out of bounds.
230  }
231 
232 
233  fMapFile.clear();
234 
235  fThetaRes=tRes;
236  fPhiRes=pRes;
237  fMapFile.open(m.data());
238  if(!fMapFile.is_open()) return false;
239 
240  for(int theta = 0; theta<fThetaRes; theta++){
241  for(int phi = 0; phi<fPhiRes; phi++){
242  for(int coord = 0; coord<3; coord++){
243  fMapFile>>fN_map[polarizationState][theta][phi][coord];
244  }
245  //G4cout<<"\n theta: " <<theta<<" phi: "<<phi<<" [xyz]: "<<fN_map[polarizationState][theta][phi][0]<<" "<<fN_map[polarizationState][theta][phi][1]<<" "<<fN_map[polarizationState][theta][phi][2];
246  }
247  }
248  G4cout<<"\n";
249  fMapFile.close();
250  G4cout<<"\nXLogicalLattice::Load_NMap() sucessful"<<endl;
251  fDresTheta=tRes; //store map dimesnions
252  fDresPhi=pRes;
253  return true;
254 }
255 
256 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
257 
258 
259 
260 double XLogicalLattice::MapKtoV(int polarizationState,G4ThreeVector k)
261 {
262  //Given the phonon wave vector k, phonon physical volume Vol
263  //and polarizationState(0=LON, 1=FT, 2=ST),
264  //returns phonon velocity in m/s
265 
266  double theta, phi, tRes, pRes;
267 
268  tRes=pi/(fVresTheta);//The summant "-1" is required:index=[0:array length-1]
269  pRes=2*pi/(fVresPhi);
270 
271 
272  theta=k.getTheta();
273  phi=k.getPhi();
274 
275  if(phi<0) phi = phi + 2*pi;
276  if(theta>pi) theta=theta-pi;
277  //phi=[0 to 2 pi] in accordance with DMC //if(phi>pi/2) phi=phi-pi/2;
278  if(fMap[polarizationState][int(theta/tRes)][int(phi/pRes)]==0){
279  G4cout<<"\nFound v=0 for polarization "<<polarizationState<<" theta "<<theta<<" phi "<<phi<< " translating to map coords " << "theta "<< int(theta/tRes) << " phi " << int(phi/pRes)<<endl;
280  }
281 
282  // G4cout<<"\nTheta: "<<theta<<" theta index: "<<int(theta/tRes)<<" tRes: "<<tRes;
283  return fMap[polarizationState][int(theta/tRes)][int(phi/pRes)];
284 
285 }
286 
287 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
288 
289 
290 
292 {
293  //Given the phonon wave vector k, phonon physical volume Vol
294  //and polarizationState(0=LON, 1=FT, 2=ST),
295  //returns phonon propagation direction as dimensionless unit vector
296 
297  double theta, phi, tRes, pRes;
298 
299  tRes=pi/(fDresTheta-1);//The summant "-1" is required:index=[0:array length-1]
300  pRes=2*pi/(fDresPhi-1);
301 
302  theta=k.getTheta();
303  phi=k.getPhi();
304 
305  if(theta>pi) theta=theta-pi;
306  //phi=[0 to 2 pi] in accordance with DMC //if(phi>pi/2) phi=phi-pi/2;
307  if(phi<0) phi = phi + 2*pi;
308 
309  G4int iTheta = int(theta/tRes+0.5);
310  G4int iPhi = int(phi/pRes+0.5);
311 
312 
313  G4ThreeVector v(fN_map[polarizationState][iTheta][iPhi][0],
314  fN_map[polarizationState][iTheta][iPhi][1],
315  fN_map[polarizationState][iTheta][iPhi][2]);
316 
317 
319  //v.Set(sin(theta)*cos(phi),sin(theta)*sin(phi),cos(theta));
320  //v.Set(k.getX(), k.getY(), k.getZ());
321 
322  return v.unit();
323 }
324 
325 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
326