Geant4  10.02.p03
XUnitCell Class Reference

#include <XUnitCell.hh>

Collaboration diagram for XUnitCell:

Public Member Functions

G4ThreeVector GetSize ()
 
G4ThreeVector GetAngle ()
 
XLogicalBaseGetBase (G4int)
 
void SetSize (G4ThreeVector)
 
void SetAngle (G4ThreeVector)
 
void SetBase (G4int, XLogicalBase *)
 
void AddBase (XLogicalBase *)
 
G4double ComputeVolume ()
 
G4double ComputeMillerOverSizeSquared (G4int, G4int, G4int)
 
G4double ComputeMillerPerSizeSquared (G4int, G4int, G4int)
 
G4double ComputeReciprocalVectorSquared (G4int, G4int, G4int)
 
G4double ComputeReciprocalVector (G4int, G4int, G4int)
 
G4double ComputeDirectVectorSquared (G4int, G4int, G4int)
 
G4double ComputeDirectVector (G4int, G4int, G4int)
 
G4double ComputeDirectPeriodSquared (G4int, G4int, G4int)
 
G4double ComputeDirectPeriod (G4int, G4int, G4int)
 
G4double ComputeAtomVolumeDensity ()
 
G4complex ComputeStructureFactor (G4int, G4int, G4int)
 
G4bool IsOrthogonal ()
 
G4bool IsCubic ()
 
 XUnitCell ()
 
 ~XUnitCell ()
 

Private Attributes

G4int fNumberOfBases
 
XLogicalBasefBase [MAXBASENUMBER]
 
G4ThreeVector fSize
 
G4ThreeVector fAngle
 

Detailed Description

Definition at line 43 of file XUnitCell.hh.

Constructor & Destructor Documentation

◆ XUnitCell()

XUnitCell::XUnitCell ( )

Definition at line 31 of file XUnitCell.cc.

31  {
33  0. * CLHEP::angstrom,
34  0. * CLHEP::angstrom);
36  0.5 * CLHEP::pi * CLHEP::radian,
37  0.5 * CLHEP::pi * CLHEP::radian);
38  fNumberOfBases = 0;
39 }
static const double radian
CLHEP::Hep3Vector G4ThreeVector
static const double angstrom
Definition: SystemOfUnits.h:81
static const double pi
Definition: SystemOfUnits.h:53
G4ThreeVector fSize
Definition: XUnitCell.hh:49
G4ThreeVector fAngle
Definition: XUnitCell.hh:50
G4int fNumberOfBases
Definition: XUnitCell.hh:46

◆ ~XUnitCell()

XUnitCell::~XUnitCell ( )

Definition at line 43 of file XUnitCell.cc.

43  {
44 }

Member Function Documentation

◆ AddBase()

void XUnitCell::AddBase ( XLogicalBase base)

Definition at line 97 of file XUnitCell.cc.

97  {
99  //the new added basis will be in the last of the [0,fNumberOfBases-1] bases
100  fBase[fNumberOfBases-1] = base;
101 }
XLogicalBase * fBase[MAXBASENUMBER]
Definition: XUnitCell.hh:47
G4int fNumberOfBases
Definition: XUnitCell.hh:46
Here is the caller graph for this function:

◆ ComputeAtomVolumeDensity()

G4double XUnitCell::ComputeAtomVolumeDensity ( )

Definition at line 116 of file XUnitCell.cc.

116  {
117  G4double vAtomVolumeDensity = 0.;
118  for(G4int i=0;i< fNumberOfBases;i++){
119  vAtomVolumeDensity += (fBase[i]->GetElement()->GetZ()
121  }
122  vAtomVolumeDensity /= ComputeVolume();
123  return vAtomVolumeDensity;
124 }
XLogicalBase * fBase[MAXBASENUMBER]
Definition: XUnitCell.hh:47
G4Element * GetElement()
Definition: XLogicalBase.cc:53
int G4int
Definition: G4Types.hh:78
G4double ComputeVolume()
Definition: XUnitCell.cc:105
XLogicalAtomicLattice * GetLattice()
Definition: XLogicalBase.cc:47
double G4double
Definition: G4Types.hh:76
G4int fNumberOfBases
Definition: XUnitCell.hh:46
G4double GetZ() const
Definition: G4Element.hh:131
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeDirectPeriod()

G4double XUnitCell::ComputeDirectPeriod ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 219 of file XUnitCell.cc.

219  {
220  return std::sqrt(ComputeDirectPeriodSquared(h,k,l));
221 }
G4double ComputeDirectPeriodSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:213
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeDirectPeriodSquared()

G4double XUnitCell::ComputeDirectPeriodSquared ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 213 of file XUnitCell.cc.

213  {
214  return (1./ComputeMillerOverSizeSquared(h,k,l));
215 }
G4double ComputeMillerOverSizeSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:128
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeDirectVector()

G4double XUnitCell::ComputeDirectVector ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 207 of file XUnitCell.cc.

207  {
208  return std::sqrt(ComputeDirectVectorSquared(h,k,l));
209 }
G4double ComputeDirectVectorSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:188
Here is the call graph for this function:

◆ ComputeDirectVectorSquared()

G4double XUnitCell::ComputeDirectVectorSquared ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 188 of file XUnitCell.cc.

188  {
189  if(IsOrthogonal()){
190  return ComputeMillerPerSizeSquared(h,k,l);
191  }
192  else{
193  double vDirectVectorSquared = 0.0;
194  vDirectVectorSquared = ComputeMillerPerSizeSquared(h,k,l);
195  vDirectVectorSquared += 2. * h * k * fSize.y() *
196  fSize.z() * std::cos(fAngle.y()) ;
197  vDirectVectorSquared += 2. * l * h * fSize.x() *
198  fSize.z() * std::cos(fAngle.x()) ;
199  vDirectVectorSquared += 2. * l * l * fSize.x() *
200  fSize.y() * std::cos(fAngle.z()) ;
201  return vDirectVectorSquared;
202  }
203 }
G4double ComputeMillerPerSizeSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:134
double x() const
G4bool IsOrthogonal()
Definition: XUnitCell.cc:226
double y() const
double z() const
G4ThreeVector fSize
Definition: XUnitCell.hh:49
G4ThreeVector fAngle
Definition: XUnitCell.hh:50
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeMillerOverSizeSquared()

G4double XUnitCell::ComputeMillerOverSizeSquared ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 128 of file XUnitCell.cc.

128  {
129  return std::pow(h/fSize.x(),2.) + std::pow(k/fSize.y(),2.) + std::pow(l/fSize.z(),2.);
130 }
double x() const
double y() const
double z() const
G4ThreeVector fSize
Definition: XUnitCell.hh:49
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeMillerPerSizeSquared()

G4double XUnitCell::ComputeMillerPerSizeSquared ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 134 of file XUnitCell.cc.

134  {
135  return std::pow(h*fSize.x(),2.) + std::pow(k*fSize.y(),2.) + std::pow(l*fSize.z(),2.);
136 }
double x() const
double y() const
double z() const
G4ThreeVector fSize
Definition: XUnitCell.hh:49
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeReciprocalVector()

G4double XUnitCell::ComputeReciprocalVector ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 182 of file XUnitCell.cc.

182  {
183  return std::sqrt(ComputeReciprocalVectorSquared(h,k,l));
184 }
G4double ComputeReciprocalVectorSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:140
Here is the call graph for this function:

◆ ComputeReciprocalVectorSquared()

G4double XUnitCell::ComputeReciprocalVectorSquared ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 140 of file XUnitCell.cc.

140  {
141  G4double vReciprocalVectorSquared = 0.0;
142 
143  if(IsOrthogonal()){
144  vReciprocalVectorSquared = ComputeMillerOverSizeSquared(h,k,l);
145  }
146  else{
147  G4double vTemp[6];
148  G4double vVolume = ComputeVolume();
149 
150  vTemp[0] = fSize.y() * fSize.z() * std::sin(fAngle.y());
151  vTemp[0] = std::pow(vTemp[0] * h,2.) / vVolume;
152 
153  vTemp[1] = fSize.x() * fSize.z() * std::sin(fAngle.x());
154  vTemp[1] = std::pow(vTemp[1] * k,2.) / vVolume;
155 
156  vTemp[2] = fSize.x() * fSize.y() * std::sin(fAngle.z());
157  vTemp[2] = std::pow(vTemp[2] * l,2.) / vVolume;
158 
159  vTemp[3] = fSize.x() * fSize.y() * std::pow(fSize.z(),2.) *
160  (std::cos(fAngle.x()) * std::cos(fAngle.y()) - std::cos(fAngle.z()));
161  vTemp[3] *= (2. * h * k);
162 
163  vTemp[4] = fSize.x() * std::pow(fSize.y(),2.) * fSize.z() *
164  (std::cos(fAngle.z()) * std::cos(fAngle.x()) - std::cos(fAngle.y()));
165  vTemp[4] *= (2. * l * h);
166 
167  vTemp[5] = std::pow(fSize.x(),2.) * fSize.y() * fSize.z() *
168  (std::cos(fAngle.y()) * std::cos(fAngle.z()) - std::cos(fAngle.x()));
169  vTemp[5] *= (2. * k * l);
170 
171  vReciprocalVectorSquared = (vTemp[0]+vTemp[1]+vTemp[2]) / vVolume;
172  vReciprocalVectorSquared += (vTemp[3]+vTemp[4]+vTemp[5]);
173  vReciprocalVectorSquared /= vVolume;
174  }
175 
176  vReciprocalVectorSquared *= (4. * CLHEP::pi * CLHEP::pi);
177  return vReciprocalVectorSquared;
178 }
G4double ComputeVolume()
Definition: XUnitCell.cc:105
G4double ComputeMillerOverSizeSquared(G4int, G4int, G4int)
Definition: XUnitCell.cc:128
static const double pi
Definition: SystemOfUnits.h:53
double x() const
G4bool IsOrthogonal()
Definition: XUnitCell.cc:226
double y() const
double z() const
G4ThreeVector fSize
Definition: XUnitCell.hh:49
G4ThreeVector fAngle
Definition: XUnitCell.hh:50
double G4double
Definition: G4Types.hh:76
Here is the call graph for this function:
Here is the caller graph for this function:

◆ ComputeStructureFactor()

G4complex XUnitCell::ComputeStructureFactor ( G4int  ,
G4int  ,
G4int   
)

◆ ComputeVolume()

G4double XUnitCell::ComputeVolume ( )

Definition at line 105 of file XUnitCell.cc.

105  {
106  if(IsOrthogonal()){
107  return ( fSize.x()*fSize.y()*fSize.z() );
108  }
109  else{
110  return (fSize.x()*fSize.y()*fSize.z()*std::cos(fAngle.x())*std::sin(fAngle.z()));
111  }
112 }
double x() const
G4bool IsOrthogonal()
Definition: XUnitCell.cc:226
double y() const
double z() const
G4ThreeVector fSize
Definition: XUnitCell.hh:49
G4ThreeVector fAngle
Definition: XUnitCell.hh:50
Here is the call graph for this function:
Here is the caller graph for this function:

◆ GetAngle()

G4ThreeVector XUnitCell::GetAngle ( )

Definition at line 54 of file XUnitCell.cc.

54  {
55  return fAngle;
56 }
G4ThreeVector fAngle
Definition: XUnitCell.hh:50

◆ GetBase()

XLogicalBase * XUnitCell::GetBase ( G4int  i)

Definition at line 60 of file XUnitCell.cc.

60  {
61  if(i<fNumberOfBases){
62  return fBase[i];
63  }
64  else{
65  G4cout << "XUnitCell::GetBase(G4int) Base "
66  << i << " does not exist" << std::endl;
67  return NULL;
68  }
69 }
XLogicalBase * fBase[MAXBASENUMBER]
Definition: XUnitCell.hh:47
G4GLOB_DLL std::ostream G4cout
G4int fNumberOfBases
Definition: XUnitCell.hh:46
Here is the caller graph for this function:

◆ GetSize()

G4ThreeVector XUnitCell::GetSize ( )

Definition at line 48 of file XUnitCell.cc.

48  {
49  return fSize;
50 }
G4ThreeVector fSize
Definition: XUnitCell.hh:49
Here is the caller graph for this function:

◆ IsCubic()

G4bool XUnitCell::IsCubic ( )

Definition at line 236 of file XUnitCell.cc.

236  {
237  if(IsOrthogonal())
238  if(fSize.x() == fSize.y())
239  if(fSize.y() == fSize.z())
240  return true;
241  return false;
242 }
double x() const
G4bool IsOrthogonal()
Definition: XUnitCell.cc:226
double y() const
double z() const
G4ThreeVector fSize
Definition: XUnitCell.hh:49
Here is the call graph for this function:

◆ IsOrthogonal()

G4bool XUnitCell::IsOrthogonal ( )

Definition at line 226 of file XUnitCell.cc.

226  {
227  if(fAngle.x() == CLHEP::pi / 2.)
228  if(fAngle.y() == CLHEP::pi / 2.)
229  if(fAngle.z() == CLHEP::pi / 2.)
230  return true;
231  return false;
232 }
static const double pi
Definition: SystemOfUnits.h:53
double x() const
double y() const
double z() const
G4ThreeVector fAngle
Definition: XUnitCell.hh:50
Here is the call graph for this function:
Here is the caller graph for this function:

◆ SetAngle()

void XUnitCell::SetAngle ( G4ThreeVector  vAngle)

Definition at line 79 of file XUnitCell.cc.

79  {
80  fAngle = vAngle;
81 }
G4ThreeVector fAngle
Definition: XUnitCell.hh:50

◆ SetBase()

void XUnitCell::SetBase ( G4int  i,
XLogicalBase base 
)

Definition at line 85 of file XUnitCell.cc.

85  {
86  if(i<fNumberOfBases){
87  fBase[i] = base;
88  }
89  else{
90  G4cout << "XUnitCell::SetBase Base(G4int,G4XLogicalBase) "
91  << i << " does not exist" << std::endl;
92  }
93 }
XLogicalBase * fBase[MAXBASENUMBER]
Definition: XUnitCell.hh:47
G4GLOB_DLL std::ostream G4cout
G4int fNumberOfBases
Definition: XUnitCell.hh:46

◆ SetSize()

void XUnitCell::SetSize ( G4ThreeVector  vSize)

Definition at line 73 of file XUnitCell.cc.

73  {
74  fSize = vSize;
75 }
G4ThreeVector fSize
Definition: XUnitCell.hh:49
Here is the caller graph for this function:

Member Data Documentation

◆ fAngle

G4ThreeVector XUnitCell::fAngle
private

Definition at line 50 of file XUnitCell.hh.

◆ fBase

XLogicalBase* XUnitCell::fBase[MAXBASENUMBER]
private

Definition at line 47 of file XUnitCell.hh.

◆ fNumberOfBases

G4int XUnitCell::fNumberOfBases
private

Definition at line 46 of file XUnitCell.hh.

◆ fSize

G4ThreeVector XUnitCell::fSize
private

Definition at line 49 of file XUnitCell.hh.


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