Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4LatticeLogical.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//
26/// \file materials/src/G4LatticeLogical.cc
27/// \brief Implementation of the G4LatticeLogical class
28//
29//
30#include "G4LatticeLogical.hh"
31#include "G4SystemOfUnits.hh"
33#include <cmath>
34#include <fstream>
35
36
37//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
38
40 : verboseLevel(0), fVresTheta(0), fVresPhi(0), fDresTheta(0), fDresPhi(0),
41 fA(0), fB(0), fLDOS(0), fSTDOS(0), fFTDOS(0),
42 fBeta(0), fGamma(0), fLambda(0), fMu(0) {
43 for (G4int i=0; i<3; i++) {
44 for (G4int j=0; j<MAXRES; j++) {
45 for (G4int k=0; k<MAXRES; k++) {
46 fMap[i][j][k] = 0.;
47 fN_map[i][j][k].set(0.,0.,0.);
48 }
49 }
50 }
51}
52
54
55//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56
57
58///////////////////////////////////////////
59//Load map of group velocity scalars (m/s)
60////////////////////////////////////////////
62 G4int polarizationState, G4String map) {
63 if (tRes>MAXRES || pRes>MAXRES) {
64 G4cerr << "G4LatticeLogical::LoadMap exceeds maximum resolution of "
65 << MAXRES << " by " << MAXRES << ". terminating." << G4endl;
66 return false; //terminate if resolution out of bounds.
67 }
68
69 std::ifstream fMapFile(map.data());
70 if(!fMapFile.is_open())
71 {
72 return false;
73 }
74
75 G4double vgrp = 0.;
76 for (G4int theta = 0; theta<tRes; theta++) {
77 for (G4int phi = 0; phi<pRes; phi++) {
78 fMapFile >> vgrp;
79 fMap[polarizationState][theta][phi] = vgrp * (m/s);
80 }
81 }
82
83 if(verboseLevel != 0)
84 {
85 G4cout << "\nG4LatticeLogical::LoadMap(" << map << ") successful"
86 << " (Vg scalars " << tRes << " x " << pRes << " for polarization "
87 << polarizationState << ")." << G4endl;
88 }
89
90 fVresTheta=tRes; //store map dimensions
91 fVresPhi=pRes;
92 return true;
93}
94
95//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
96
97
98////////////////////////////////////
99//Load map of group velocity unit vectors
100///////////////////////////////////
102 G4int polarizationState, G4String map) {
103 if (tRes>MAXRES || pRes>MAXRES) {
104 G4cerr << "G4LatticeLogical::LoadMap exceeds maximum resolution of "
105 << MAXRES << " by " << MAXRES << ". terminating." << G4endl;
106 return false; //terminate if resolution out of bounds.
107 }
108
109 std::ifstream fMapFile(map.data());
110 if(!fMapFile.is_open())
111 {
112 return false;
113 }
114
115 G4double x,y,z; // Buffers to read coordinates from file
116 G4ThreeVector dir;
117 for (G4int theta = 0; theta<tRes; theta++) {
118 for (G4int phi = 0; phi<pRes; phi++) {
119 fMapFile >> x >> y >> z;
120 dir.set(x,y,z);
121 fN_map[polarizationState][theta][phi] = dir.unit(); // Enforce unity
122 }
123 }
124
125 if(verboseLevel != 0)
126 {
127 G4cout << "\nG4LatticeLogical::Load_NMap(" << map << ") successful"
128 << " (Vdir " << tRes << " x " << pRes << " for polarization "
129 << polarizationState << ")." << G4endl;
130 }
131
132 fDresTheta=tRes; //store map dimensions
133 fDresPhi=pRes;
134 return true;
135}
136
137//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
138
139//Given the phonon wave vector k, phonon physical volume Vol
140//and polarizationState(0=LON, 1=FT, 2=ST),
141//returns phonon velocity in m/s
142
144 const G4ThreeVector& k) const {
145 G4double theta, phi, tRes, pRes;
146
147 tRes=pi/fVresTheta;
148 pRes=twopi/fVresPhi;
149
150 theta=k.getTheta();
151 phi=k.getPhi();
152
153 if(phi < 0)
154 {
155 phi = phi + twopi;
156 }
157 if(theta > pi)
158 {
159 theta = theta - pi;
160 }
161
162 G4double Vg = fMap[polarizationState][int(theta/tRes)][int(phi/pRes)];
163
164 if(Vg == 0){
165 G4cout<<"\nFound v=0 for polarization "<<polarizationState
166 <<" theta "<<theta<<" phi "<<phi<< " translating to map coords "
167 <<"theta "<< int(theta/tRes) << " phi " << int(phi/pRes)<<G4endl;
168 }
169
170 if (verboseLevel>1) {
171 G4cout << "G4LatticeLogical::MapKtoV theta,phi=" << theta << " " << phi
172 << " : ith,iph " << int(theta/tRes) << " " << int(phi/pRes)
173 << " : V " << Vg << G4endl;
174 }
175
176 return Vg;
177}
178
179//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
180
181//Given the phonon wave vector k, phonon physical volume Vol
182//and polarizationState(0=LON, 1=FT, 2=ST),
183//returns phonon propagation direction as dimensionless unit vector
184
186 const G4ThreeVector& k) const {
187 G4double theta, phi, tRes, pRes;
188
189 tRes=pi/(fDresTheta-1);//The summant "-1" is required:index=[0:array length-1]
190 pRes=2*pi/(fDresPhi-1);
191
192 theta=k.getTheta();
193 phi=k.getPhi();
194
195 if(theta > pi)
196 {
197 theta = theta - pi;
198 }
199 //phi=[0 to 2 pi] in accordance with DMC //if(phi>pi/2) phi=phi-pi/2;
200 if(phi < 0)
201 {
202 phi = phi + 2 * pi;
203 }
204
205 G4int iTheta = int(theta/tRes+0.5);
206 G4int iPhi = int(phi/pRes+0.5);
207
208 if (verboseLevel>1) {
209 G4cout << "G4LatticeLogical::MapKtoVDir theta,phi=" << theta << " " << phi
210 << " : ith,iph " << iTheta << " " << iPhi
211 << " : dir " << fN_map[polarizationState][iTheta][iPhi] << G4endl;
212 }
213
214 return fN_map[polarizationState][iTheta][iPhi];
215}
216
217//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
218
219// Dump structure in format compatible with reading back
220
221void G4LatticeLogical::Dump(std::ostream& os) const {
222 os << "dyn " << fBeta << " " << fGamma << " " << fLambda << " " << fMu
223 << "\nscat " << fB << " decay " << fA
224 << "\nLDOS " << fLDOS << " STDOS " << fSTDOS << " FTDOS " << fFTDOS
225 << std::endl;
226
227 Dump_NMap(os, 0, "LVec.ssv");
228 Dump_NMap(os, 1, "FTVec.ssv");
229 Dump_NMap(os, 2, "STVec.ssv");
230
231 DumpMap(os, 0, "L.ssv");
232 DumpMap(os, 1, "FT.ssv");
233 DumpMap(os, 2, "ST.ssv");
234}
235
236void G4LatticeLogical::DumpMap(std::ostream& os, G4int pol,
237 const G4String& name) const {
238 os << "VG " << name << " " << (pol==0?"L":pol==1?"FT":pol==2?"ST":"??")
239 << " " << fVresTheta << " " << fVresPhi << std::endl;
240
241 for (G4int iTheta=0; iTheta<fVresTheta; iTheta++) {
242 for (G4int iPhi=0; iPhi<fVresPhi; iPhi++) {
243 os << fMap[pol][iTheta][iPhi] << std::endl;
244 }
245 }
246}
247
248void G4LatticeLogical::Dump_NMap(std::ostream& os, G4int pol,
249 const G4String& name) const {
250 os << "VDir " << name << " " << (pol==0?"L":pol==1?"FT":pol==2?"ST":"??")
251 << " " << fDresTheta << " " << fDresPhi << std::endl;
252
253 for (G4int iTheta=0; iTheta<fDresTheta; iTheta++) {
254 for (G4int iPhi=0; iPhi<fDresPhi; iPhi++) {
255 os << fN_map[pol][iTheta][iPhi].x()
256 << " " << fN_map[pol][iTheta][iPhi].y()
257 << " " << fN_map[pol][iTheta][iPhi].z()
258 << std::endl;
259 }
260 }
261}
262
Definition of the G4LatticeLogical class.
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
G4GLOB_DLL std::ostream G4cerr
#define G4endl
Definition: G4ios.hh:57
G4GLOB_DLL std::ostream G4cout
double z() const
Hep3Vector unit() const
double x() const
double y() const
double getTheta() const
void set(double x, double y, double z)
double getPhi() const
virtual G4double MapKtoV(G4int, const G4ThreeVector &) const
void DumpMap(std::ostream &os, G4int pol, const G4String &name) const
G4bool LoadMap(G4int, G4int, G4int, G4String)
void Dump_NMap(std::ostream &os, G4int pol, const G4String &name) const
virtual G4ThreeVector MapKtoVDir(G4int, const G4ThreeVector &) const
void Dump(std::ostream &os) const
G4bool Load_NMap(G4int, G4int, G4int, G4String)
virtual ~G4LatticeLogical()