Geant4  10.03.p02
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
XUnitCell.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 
30 #include "XUnitCell.hh"
31 #include "G4PhysicalConstants.hh"
32 #include <cmath>
33 
35  fSize = G4ThreeVector(0. * CLHEP::angstrom,
36  0. * CLHEP::angstrom,
37  0. * CLHEP::angstrom);
38  fAngle = G4ThreeVector(0.5 * CLHEP::pi * CLHEP::radian,
39  0.5 * CLHEP::pi * CLHEP::radian,
40  0.5 * CLHEP::pi * CLHEP::radian);
41  fNumberOfBases = 0;
42 }
43 
44 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
45 
47 }
48 
49 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
50 
52  return fSize;
53 }
54 
55 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
56 
58  return fAngle;
59 }
60 
61 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
62 
64  if(i<fNumberOfBases){
65  return fBase[i];
66  }
67  else{
68  G4cout << "XUnitCell::GetBase(G4int) Base "
69  << i << " does not exist" << std::endl;
70  return NULL;
71  }
72 }
73 
74 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
75 
77  fSize = vSize;
78 }
79 
80 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
81 
83  fAngle = vAngle;
84 }
85 
86 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
87 
89  if(i<fNumberOfBases){
90  fBase[i] = base;
91  }
92  else{
93  G4cout << "XUnitCell::SetBase Base(G4int,G4XLogicalBase) "
94  << i << " does not exist" << std::endl;
95  }
96 }
97 
98 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
99 
101  fNumberOfBases++;
102  //the new added basis will be in the last of the [0,fNumberOfBases-1] bases
103  fBase[fNumberOfBases-1] = base;
104 }
105 
106 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
107 
109  if(IsOrthogonal()){
110  return ( fSize.x()*fSize.y()*fSize.z() );
111  }
112  else{
113  return (fSize.x()*fSize.y()*fSize.z()*std::cos(fAngle.x())*std::sin(fAngle.z()));
114  }
115 }
116 
117 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
118 
120  G4double vAtomVolumeDensity = 0.;
121  for(G4int i=0;i< fNumberOfBases;i++){
122  vAtomVolumeDensity += (fBase[i]->GetElement()->GetZ()
123  *fBase[i]->GetLattice()->GetLatticeNumberOfAtoms());
124  }
125  vAtomVolumeDensity /= ComputeVolume();
126  return vAtomVolumeDensity;
127 }
128 
129 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
130 
132  return std::pow(h/fSize.x(),2.) + std::pow(k/fSize.y(),2.) + std::pow(l/fSize.z(),2.);
133 }
134 
135 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
136 
138  return std::pow(h*fSize.x(),2.) + std::pow(k*fSize.y(),2.) + std::pow(l*fSize.z(),2.);
139 }
140 
141 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
142 
144  G4double vReciprocalVectorSquared = 0.0;
145 
146  if(IsOrthogonal()){
147  vReciprocalVectorSquared = ComputeMillerOverSizeSquared(h,k,l);
148  }
149  else{
150  G4double vTemp[6];
151  G4double vVolume = ComputeVolume();
152 
153  vTemp[0] = fSize.y() * fSize.z() * std::sin(fAngle.y());
154  vTemp[0] = std::pow(vTemp[0] * h,2.) / vVolume;
155 
156  vTemp[1] = fSize.x() * fSize.z() * std::sin(fAngle.x());
157  vTemp[1] = std::pow(vTemp[1] * k,2.) / vVolume;
158 
159  vTemp[2] = fSize.x() * fSize.y() * std::sin(fAngle.z());
160  vTemp[2] = std::pow(vTemp[2] * l,2.) / vVolume;
161 
162  vTemp[3] = fSize.x() * fSize.y() * std::pow(fSize.z(),2.) *
163  (std::cos(fAngle.x()) * std::cos(fAngle.y()) - std::cos(fAngle.z()));
164  vTemp[3] *= (2. * h * k);
165 
166  vTemp[4] = fSize.x() * std::pow(fSize.y(),2.) * fSize.z() *
167  (std::cos(fAngle.z()) * std::cos(fAngle.x()) - std::cos(fAngle.y()));
168  vTemp[4] *= (2. * l * h);
169 
170  vTemp[5] = std::pow(fSize.x(),2.) * fSize.y() * fSize.z() *
171  (std::cos(fAngle.y()) * std::cos(fAngle.z()) - std::cos(fAngle.x()));
172  vTemp[5] *= (2. * k * l);
173 
174  vReciprocalVectorSquared = (vTemp[0]+vTemp[1]+vTemp[2]) / vVolume;
175  vReciprocalVectorSquared += (vTemp[3]+vTemp[4]+vTemp[5]);
176  vReciprocalVectorSquared /= vVolume;
177  }
178 
179  vReciprocalVectorSquared *= (4. * CLHEP::pi * CLHEP::pi);
180  return vReciprocalVectorSquared;
181 }
182 
183 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
184 
186  return std::sqrt(ComputeReciprocalVectorSquared(h,k,l));
187 }
188 
189 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
190 
192  if(IsOrthogonal()){
193  return ComputeMillerPerSizeSquared(h,k,l);
194  }
195  else{
196  double vDirectVectorSquared = 0.0;
197  vDirectVectorSquared = ComputeMillerPerSizeSquared(h,k,l);
198  vDirectVectorSquared += 2. * h * k * fSize.y() *
199  fSize.z() * std::cos(fAngle.y()) ;
200  vDirectVectorSquared += 2. * l * h * fSize.x() *
201  fSize.z() * std::cos(fAngle.x()) ;
202  vDirectVectorSquared += 2. * l * l * fSize.x() *
203  fSize.y() * std::cos(fAngle.z()) ;
204  return vDirectVectorSquared;
205  }
206 }
207 
208 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
209 
211  return std::sqrt(ComputeDirectVectorSquared(h,k,l));
212 }
213 
214 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
215 
217  return (1./ComputeMillerOverSizeSquared(h,k,l));
218 }
219 
220 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
221 
223  return std::sqrt(ComputeDirectPeriodSquared(h,k,l));
224 }
225 
226 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
227 
228 
230  if(fAngle.x() == CLHEP::pi / 2.)
231  if(fAngle.y() == CLHEP::pi / 2.)
232  if(fAngle.z() == CLHEP::pi / 2.)
233  return true;
234  return false;
235 }
236 
237 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
238 
240  if(IsOrthogonal())
241  if(fSize.x() == fSize.y())
242  if(fSize.y() == fSize.z())
243  return true;
244  return false;
245 }
246 
247 //....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
G4ThreeVector GetAngle()
Definition: XUnitCell.cc:57
G4double ComputeDirectVectorSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:191
G4Element * GetElement()
Definition: XLogicalBase.cc:56
CLHEP::Hep3Vector G4ThreeVector
double x() const
~XUnitCell()
Definition: XUnitCell.cc:46
void SetSize(G4ThreeVector)
Definition: XUnitCell.cc:76
G4bool IsCubic()
Definition: XUnitCell.cc:239
G4double GetZ() const
Definition: G4Element.hh:131
G4double ComputeMillerPerSizeSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:137
const XML_Char int const XML_Char int const XML_Char * base
Definition: expat.h:331
G4double ComputeReciprocalVectorSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:143
int G4int
Definition: G4Types.hh:78
G4double ComputeVolume()
Definition: XUnitCell.cc:108
double z() const
G4double ComputeMillerOverSizeSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:131
XLogicalAtomicLattice * GetLattice()
Definition: XLogicalBase.cc:50
G4GLOB_DLL std::ostream G4cout
bool G4bool
Definition: G4Types.hh:79
G4ThreeVector GetSize()
Definition: XUnitCell.cc:51
void AddBase(XLogicalBase *)
Definition: XUnitCell.cc:100
G4double ComputeAtomVolumeDensity()
Definition: XUnitCell.cc:119
static constexpr double radian
G4bool IsOrthogonal()
Definition: XUnitCell.cc:229
double y() const
G4double ComputeDirectPeriod(G4int, G4int, G4int)
Definition: XUnitCell.cc:222
XLogicalBase * GetBase(G4int)
Definition: XUnitCell.cc:63
void SetBase(G4int, XLogicalBase *)
Definition: XUnitCell.cc:88
G4double ComputeDirectPeriodSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:216
double G4double
Definition: G4Types.hh:76
Definition of the XUnitCell class.
void SetAngle(G4ThreeVector)
Definition: XUnitCell.cc:82
static constexpr double pi
Definition: SystemOfUnits.h:54
G4double ComputeReciprocalVector(G4int, G4int, G4int)
Definition: XUnitCell.cc:185
G4double ComputeDirectVector(G4int, G4int, G4int)
Definition: XUnitCell.cc:210
static constexpr double angstrom
Definition: SystemOfUnits.h:82