Geant4-11
Public Member Functions | Protected Attributes | Private Member Functions | Private Attributes
G4CrystalUnitCell Class Reference

#include <G4CrystalUnitCell.hh>

Public Member Functions

G4double ComputeCellVolume ()
 
G4bool FillAtomicPos (G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
 
G4bool FillAtomicUnitPos (G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
 
G4bool FillElReduced (G4double Cij[6][6])
 
 G4CrystalUnitCell (G4double sizeA, G4double sizeB, G4double sizeC, G4double alpha, G4double beta, G4double gamma, G4int spacegroup)
 
G4ThreeVector GetAngle () const
 
const G4ThreeVectorGetBasis (G4int idx) const
 
theBravaisLatticeType GetBravaisLattice ()
 
G4double GetIntCosAng (G4int h1, G4int k1, G4int l1, G4int h2, G4int k2, G4int l2)
 
G4double GetIntSp2 (G4int h, G4int k, G4int l)
 
theLatticeSystemType GetLatticeSystem ()
 
G4ThreeVector GetRecAngle () const
 
const G4ThreeVectorGetRecBasis (G4int idx) const
 
G4double GetRecIntSp2 (G4int h, G4int k, G4int l)
 
G4ThreeVector GetRecSize () const
 
const G4ThreeVectorGetRecUnitBasis (G4int idx) const
 
G4double GetRecVolume () const
 
G4ThreeVector GetSize () const
 
G4int GetSpaceGroup () const
 
const G4ThreeVectorGetUnitBasis (G4int idx) const
 
G4ThreeVector GetUnitBasisTrigonal ()
 
G4double GetVolume () const
 
void SetSpaceGroup (G4int aInt)
 
virtual ~G4CrystalUnitCell ()
 

Protected Attributes

G4ThreeVector nullVec
 
G4ThreeVector theAngle
 
G4ThreeVector theBasis [3]
 
G4ThreeVector theRecAngle
 
G4ThreeVector theRecBasis [3]
 
G4ThreeVector theRecSize
 
G4ThreeVector theRecUnitBasis [3]
 
G4ThreeVector theSize
 
G4ThreeVector theUnitBasis [3]
 

Private Member Functions

G4bool FillAmorphous (G4double Cij[6][6]) const
 
G4bool FillCubic (G4double Cij[6][6]) const
 
G4bool FillHexagonal (G4double Cij[6][6]) const
 
G4bool FillMonoclinic (G4double Cij[6][6]) const
 
G4bool FillOrthorhombic (G4double Cij[6][6]) const
 
G4bool FillRhombohedral (G4double Cij[6][6]) const
 
G4bool FillTetragonal (G4double Cij[6][6]) const
 
G4bool FillTriclinic (G4double Cij[6][6]) const
 
theBravaisLatticeType GetBravaisLattice (G4int aGroup)
 
theLatticeSystemType GetLatticeSystem (G4int aGroup)
 
G4bool ReflectElReduced (G4double Cij[6][6]) const
 

Private Attributes

G4double cosa
 
G4double cosar
 
G4double cosb
 
G4double cosbr
 
G4double cosg
 
G4double cosgr
 
G4double sina
 
G4double sinb
 
G4double sing
 
G4double theRecVolume
 
G4int theSpaceGroup
 
G4double theVolume
 

Detailed Description

Definition at line 47 of file G4CrystalUnitCell.hh.

Constructor & Destructor Documentation

◆ G4CrystalUnitCell()

G4CrystalUnitCell::G4CrystalUnitCell ( G4double  sizeA,
G4double  sizeB,
G4double  sizeC,
G4double  alpha,
G4double  beta,
G4double  gamma,
G4int  spacegroup 
)

Definition at line 38 of file G4CrystalUnitCell.cc.

44 :
45theSpaceGroup(spacegroup),
46theSize(G4ThreeVector(sizeA,sizeB,sizeC)),
48{
49
50 nullVec = G4ThreeVector(0.,0.,0.);
54
58
59 cosa=std::cos(alpha), cosb=std::cos(beta), cosg=std::cos(gamma);
60 sina=std::sin(alpha), sinb=std::sin(beta), sing=std::sin(gamma);
61
65
68
69 theRecSize[0] = sizeB * sizeC * sina / theVolume;
70 theRecSize[1] = sizeC * sizeA * sinb / theVolume;
71 theRecSize[2] = sizeA * sizeB * sing / theVolume;
72
73 theRecAngle[0] = std::acos(cosar);
74 theRecAngle[1] = std::acos(cosbr);
75 theRecAngle[2] = std::acos(cosgr);
76
77 G4double x3,y3,z3;
78
80 case Amorphous:
81 break;
82 case Cubic: // Cubic, C44 set
83 break;
84 case Tetragonal:
85 break;
86 case Orthorhombic:
87 break;
88 case Rhombohedral:
89 theUnitBasis[1].rotateZ(gamma-CLHEP::halfpi); // X-Y opening angle
90 // Z' axis computed by hand to get both opening angles right
91 // X'.Z' = cos(alpha), Y'.Z' = cos(beta), solve for Z' components
92 x3=cosa, y3=(cosb-cosa*cosg)/sing, z3=std::sqrt(1.-x3*x3-y3*y3);
93 theUnitBasis[2] = G4ThreeVector(x3, y3, z3).unit();
94 break;
95 case Monoclinic:
96 theUnitBasis[2].rotateX(beta-CLHEP::halfpi); // Z-Y opening angle
97 break;
98 case Triclinic:
99 theUnitBasis[1].rotateZ(gamma-CLHEP::halfpi); // X-Y opening angle
100 // Z' axis computed by hand to get both opening angles right
101 // X'.Z' = cos(alpha), Y'.Z' = cos(beta), solve for Z' components
102 x3=cosa, y3=(cosb-cosa*cosg)/sing, z3=std::sqrt(1.-x3*x3-y3*y3);
103 theUnitBasis[2] = G4ThreeVector(x3, y3, z3).unit();
104 break;
105 case Hexagonal: // Tetragonal, C16=0
106 theUnitBasis[1].rotateZ(30.*CLHEP::deg); // X-Y opening angle
107 break;
108 default:
109 break;
110 }
111
112 for(auto i:{0,1,2}){
113 theBasis[i] = theUnitBasis[i] * theSize[i];
115 }
116
117 // Initialize sgInfo
118 /* at first some initialization for SgInfo */
119 /*
120 const T_TabSgName *tsgn = NULL;
121
122 SgInfo.MaxList = 192;
123 SgInfo.ListSeitzMx = malloc( SgInfo.MaxList * sizeof(*SgInfo.ListSeitzMx) );
124
125 // no list info needed here
126 SgInfo.ListRotMxInfo = NULL;
127 tsgn = FindTabSgNameEntry(SchoenfliesSymbols[theSpaceGroup], 'A');
128
129 // initialize SgInfo struct
130 InitSgInfo( &SgInfo );
131 SgInfo.TabSgName = tsgn;
132 if ( tsgn ){
133 SgInfo.GenOption = 1;
134 }
135
136 ParseHallSymbol( SchoenfliesSymbols[theSpaceGroup], &SgInfo );
137 CompleteSgInfo( &SgInfo );
138 Set_si( &SgInfo );
139 */
140
141}
static const G4double alpha
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
Hep3Vector & rotateX(double)
Definition: ThreeVector.cc:87
Hep3Vector unit() const
Hep3Vector & rotateZ(double)
Definition: ThreeVector.cc:107
G4ThreeVector theBasis[3]
G4ThreeVector theSize
G4ThreeVector theRecSize
G4ThreeVector theRecUnitBasis[3]
G4ThreeVector theRecBasis[3]
G4ThreeVector nullVec
G4ThreeVector theUnitBasis[3]
G4double ComputeCellVolume()
G4ThreeVector theRecAngle
theLatticeSystemType GetLatticeSystem()
G4ThreeVector theAngle
static constexpr double deg
DLL_API const Hep3Vector HepZHat
Definition: ThreeVector.h:419
DLL_API const Hep3Vector HepXHat
static constexpr double halfpi
Definition: SystemOfUnits.h:57
DLL_API const Hep3Vector HepYHat
Definition: ThreeVector.h:419

References alpha, Amorphous, anonymous_namespace{G4PionRadiativeDecayChannel.cc}::beta, ComputeCellVolume(), cosa, cosar, cosb, cosbr, cosg, cosgr, Cubic, CLHEP::deg, GetLatticeSystem(), CLHEP::halfpi, CLHEP::HepXHat, CLHEP::HepYHat, CLHEP::HepZHat, Hexagonal, Monoclinic, nullVec, Orthorhombic, Rhombohedral, CLHEP::Hep3Vector::rotateX(), CLHEP::Hep3Vector::rotateZ(), sina, sinb, sing, Tetragonal, theBasis, theRecAngle, theRecBasis, theRecSize, theRecUnitBasis, theRecVolume, theSize, theSpaceGroup, theUnitBasis, theVolume, Triclinic, and CLHEP::Hep3Vector::unit().

◆ ~G4CrystalUnitCell()

G4CrystalUnitCell::~G4CrystalUnitCell ( )
virtual

Definition at line 145 of file G4CrystalUnitCell.cc.

145{;}

Member Function Documentation

◆ ComputeCellVolume()

G4double G4CrystalUnitCell::ComputeCellVolume ( )

Definition at line 386 of file G4CrystalUnitCell.cc.

386 {
387 G4double a = theSize[0], b = theSize[1], c = theSize[2];
388
389 switch(GetLatticeSystem())
390 {
391 case Amorphous:
392 return 0.;
393 break;
394 case Cubic:
395 return a * a * a;
396 break;
397 case Tetragonal:
398 return a * a * c;
399 break;
400 case Orthorhombic:
401 return a * b * c;
402 break;
403 case Rhombohedral:
404 return a*a*a*std::sqrt(1.-3.*cosa*cosa+2.*cosa*cosa*cosa);
405 break;
406 case Monoclinic:
407 return a*b*c*sinb;
408 break;
409 case Triclinic:
410 return a*b*c*std::sqrt(1.-cosa*cosa-cosb*cosb-cosg*cosg*2.*cosa*cosb*cosg);
411 break;
412 case Hexagonal:
413 return std::sqrt(3.0)/2.*a*a*c;
414 break;
415 default:
416 break;
417 }
418
419 return 0.;
420}

References Amorphous, cosa, cosb, cosg, Cubic, GetLatticeSystem(), Hexagonal, Monoclinic, Orthorhombic, Rhombohedral, sinb, Tetragonal, theSize, and Triclinic.

Referenced by G4CrystalUnitCell().

◆ FillAmorphous()

G4bool G4CrystalUnitCell::FillAmorphous ( G4double  Cij[6][6]) const
private

Definition at line 264 of file G4CrystalUnitCell.cc.

264 {
265 Cij[3][3] = 0.5*(Cij[0][0]-Cij[0][1]);
266 return true;
267}

Referenced by FillElReduced().

◆ FillAtomicPos()

G4bool G4CrystalUnitCell::FillAtomicPos ( G4ThreeVector pos,
std::vector< G4ThreeVector > &  vecout 
)

Definition at line 217 of file G4CrystalUnitCell.cc.

217 {
218 FillAtomicUnitPos(posin,vecout);
219 for(auto &vec:vecout){
220 vec.setX(vec.x()*theSize[0]);
221 vec.setY(vec.y()*theSize[1]);
222 vec.setZ(vec.z()*theSize[2]);
223 }
224 return true;
225}
G4bool FillAtomicUnitPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)

References FillAtomicUnitPos(), and theSize.

Referenced by G4CrystalExtension::GetAtomPos().

◆ FillAtomicUnitPos()

G4bool G4CrystalUnitCell::FillAtomicUnitPos ( G4ThreeVector pos,
std::vector< G4ThreeVector > &  vecout 
)

Definition at line 207 of file G4CrystalUnitCell.cc.

207 {
208 // Just for testing the infrastructure
209 G4ThreeVector aaa = pos;
210 vecout.push_back(aaa);
211 vecout.push_back(G4ThreeVector(2.,5.,3.));
212 return true;
213}
static const G4double pos

References pos.

Referenced by FillAtomicPos().

◆ FillCubic()

G4bool G4CrystalUnitCell::FillCubic ( G4double  Cij[6][6]) const
private

Definition at line 271 of file G4CrystalUnitCell.cc.

271 {
272 G4double C11=Cij[0][0], C12=Cij[0][1], C44=Cij[3][3];
273
274 for (size_t i=0; i<6; i++) {
275 for (size_t j=i; j<6; j++) {
276 if (i<3 && j<3) Cij[i][j] = (i==j) ? C11 : C12;
277 else if (i==j && i>=3) Cij[i][i] = C44;
278 else Cij[i][j] = 0.;
279 }
280 }
281
282 ReflectElReduced(Cij);
283
284 return (C11!=0. && C12!=0. && C44!=0.);
285}
G4bool ReflectElReduced(G4double Cij[6][6]) const

References ReflectElReduced().

Referenced by FillElReduced().

◆ FillElReduced()

G4bool G4CrystalUnitCell::FillElReduced ( G4double  Cij[6][6])

Definition at line 229 of file G4CrystalUnitCell.cc.

229 {
230 switch (GetLatticeSystem()) {
231 case Amorphous:
232 return FillAmorphous(Cij);
233 break;
234 case Cubic: // Cubic, C44 set
235 return FillCubic(Cij);
236 break;
237 case Tetragonal:
238 return FillTetragonal(Cij);
239 break;
240 case Orthorhombic:
241 return FillOrthorhombic(Cij);
242 break;
243 case Rhombohedral:
244 return FillRhombohedral(Cij);
245 break;
246 case Monoclinic:
247 return FillMonoclinic(Cij);
248 break;
249 case Triclinic:
250 return FillTriclinic(Cij);
251 break;
252 case Hexagonal: // Tetragonal, C16=0
253 return FillHexagonal(Cij);
254 break;
255 default:
256 break;
257 }
258
259 return false;
260}
G4bool FillMonoclinic(G4double Cij[6][6]) const
G4bool FillRhombohedral(G4double Cij[6][6]) const
G4bool FillAmorphous(G4double Cij[6][6]) const
G4bool FillHexagonal(G4double Cij[6][6]) const
G4bool FillCubic(G4double Cij[6][6]) const
G4bool FillOrthorhombic(G4double Cij[6][6]) const
G4bool FillTriclinic(G4double Cij[6][6]) const
G4bool FillTetragonal(G4double Cij[6][6]) const

References Amorphous, Cubic, FillAmorphous(), FillCubic(), FillHexagonal(), FillMonoclinic(), FillOrthorhombic(), FillRhombohedral(), FillTetragonal(), FillTriclinic(), GetLatticeSystem(), Hexagonal, Monoclinic, Orthorhombic, Rhombohedral, Tetragonal, and Triclinic.

◆ FillHexagonal()

G4bool G4CrystalUnitCell::FillHexagonal ( G4double  Cij[6][6]) const
private

Definition at line 366 of file G4CrystalUnitCell.cc.

366 {
367 Cij[0][5] = 0.;
368 Cij[4][5] = 0.5*(Cij[0][0] - Cij[0][1]);
369 return true;
370}

Referenced by FillElReduced().

◆ FillMonoclinic()

G4bool G4CrystalUnitCell::FillMonoclinic ( G4double  Cij[6][6]) const
private

Definition at line 340 of file G4CrystalUnitCell.cc.

340 {
341 // The monoclinic matrix has 13 independent elements with no degeneracies
342 // Sanity condition is same as orthorhombic, plus C45, C(1,2,3)6
343
344 return (FillOrthorhombic(Cij) && Cij[0][5]!=0. && Cij[1][5]!=0. &&
345 Cij[2][5] != 0. && Cij[3][4]!=0.);
346}

References FillOrthorhombic().

Referenced by FillElReduced().

◆ FillOrthorhombic()

G4bool G4CrystalUnitCell::FillOrthorhombic ( G4double  Cij[6][6]) const
private

Definition at line 306 of file G4CrystalUnitCell.cc.

306 {
307 // No degenerate elements; just check for all non-zero
308 ReflectElReduced(Cij);
309
310 G4bool good = true;
311 for (size_t i=0; i<6; i++) {
312 for (size_t j=i+1; j<3; j++)
313 good &= (Cij[i][j] != 0);
314 }
315
316 return good;
317}
bool G4bool
Definition: G4Types.hh:86

References ReflectElReduced().

Referenced by FillElReduced(), and FillMonoclinic().

◆ FillRhombohedral()

G4bool G4CrystalUnitCell::FillRhombohedral ( G4double  Cij[6][6]) const
private

Definition at line 321 of file G4CrystalUnitCell.cc.

321 {
322 G4double C11=Cij[0][0], C12=Cij[0][1], C13=Cij[0][2], C14=Cij[0][3];
323 G4double C15=Cij[0][4], C33=Cij[2][2], C44=Cij[3][3], C66=0.5*(C11-C12);
324
325 Cij[1][1] = C11; // Copy small number of individual elements
326 Cij[1][2] = C13;
327 Cij[1][3] = -C14;
328 Cij[1][4] = -C15;
329 Cij[3][5] = -C15;
330 Cij[4][4] = C44;
331 Cij[4][5] = C14;
332
333 // NOTE: C15 may be zero (c.f. rhombohedral(I) vs. (II))
334 return (C11!=0 && C12!=0 && C13!=0 && C14!=0. &&
335 C33!=0. && C44!=0. && C66!=0.);
336}

Referenced by FillElReduced().

◆ FillTetragonal()

G4bool G4CrystalUnitCell::FillTetragonal ( G4double  Cij[6][6]) const
private

Definition at line 289 of file G4CrystalUnitCell.cc.

289 {
290 G4double C11=Cij[0][0], C12=Cij[0][1], C13=Cij[0][2], C16=Cij[0][5];
291 G4double C33=Cij[2][2], C44=Cij[3][3], C66=Cij[5][5];
292
293 Cij[1][1] = C11; // Copy small number of individual elements
294 Cij[1][2] = C13;
295 Cij[1][5] = -C16;
296 Cij[4][4] = C44;
297
298 ReflectElReduced(Cij);
299
300 // NOTE: Do not test for C16 != 0., to allow calling from Hexagonal
301 return (C11!=0. && C12!=0. && C13!=0. && C33!=0. && C44!=0. && C66!=0.);
302}

References ReflectElReduced().

Referenced by FillElReduced().

◆ FillTriclinic()

G4bool G4CrystalUnitCell::FillTriclinic ( G4double  Cij[6][6]) const
private

Definition at line 350 of file G4CrystalUnitCell.cc.

350 {
351 // The triclinic matrix has the entire upper half filled (21 elements)
352
353 ReflectElReduced(Cij);
354
355 G4bool good = true;
356 for (size_t i=0; i<6; i++) {
357 for (size_t j=i; j<6; j++) good &= (Cij[i][j] != 0);
358 }
359
360 return good;
361}

References ReflectElReduced().

Referenced by FillElReduced().

◆ GetAngle()

G4ThreeVector G4CrystalUnitCell::GetAngle ( ) const
inline

Definition at line 102 of file G4CrystalUnitCell.hh.

102{return theAngle;}

References theAngle.

◆ GetBasis()

const G4ThreeVector & G4CrystalUnitCell::GetBasis ( G4int  idx) const

Definition at line 180 of file G4CrystalUnitCell.cc.

180 {
181 return (idx>=0 && idx<3 ? theBasis[idx] : nullVec);
182}

References nullVec, and theBasis.

Referenced by G4LogicalCrystalVolume::GetBasis().

◆ GetBravaisLattice() [1/2]

theBravaisLatticeType G4CrystalUnitCell::GetBravaisLattice ( )
inline

Definition at line 74 of file G4CrystalUnitCell.hh.

74 {
76 }
theBravaisLatticeType GetBravaisLattice()

References GetBravaisLattice(), and theSpaceGroup.

Referenced by GetBravaisLattice().

◆ GetBravaisLattice() [2/2]

theBravaisLatticeType G4CrystalUnitCell::GetBravaisLattice ( G4int  aGroup)
private

◆ GetIntCosAng()

G4double G4CrystalUnitCell::GetIntCosAng ( G4int  h1,
G4int  k1,
G4int  l1,
G4int  h2,
G4int  k2,
G4int  l2 
)

Definition at line 550 of file G4CrystalUnitCell.cc.

555 {
556
557 /* Reference:
558 Table 2.4, pag. 65
559
560 @Inbook{Kelly2012,
561 author="Anthony A. Kelly and Kevin M. Knowles",
562 title="Appendix 3 Interplanar Spacings and Interplanar Angles",
563 bookTitle="Crystallography and Crystal Defects, 2nd Edition",
564 year="2012",
565 publisher="John Wiley & Sons, Ltd.",
566 isbn="978-0-470-75014-8",
567 doi="10.1002/9781119961468",
568 url="http://onlinelibrary.wiley.com/book/10.1002/9781119961468"
569 }
570 */
571
572 G4double a = theRecSize[0], b = theRecSize[1], c = theRecSize[2];
573 G4double a2 = a*a, b2 = b*b, c2 = c*c;
574 G4double dsp1dsp2;
575 switch(GetLatticeSystem())
576 {
577 case Amorphous:
578 return 0.;
579 break;
580 case Cubic:
581 return (h1*h2 + k1*k2 + l1+l2) / (std::sqrt(h1*h1 + k1*k1 + l1*l1) * std::sqrt(h2*h2 + k2*k2 + l2*l2));
582 break;
583 case Tetragonal:
584 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
585 return 0. ;
586 break;
587 case Orthorhombic:
588 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
589 return dsp1dsp2 * (h1*h2*a2 + k1*k2*a2 + l1*l2*c2);
590 break;
591 case Rhombohedral:
592 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
593 return dsp1dsp2 * (h1*h2*a2 + k1*k2*b2 + l1*l2*c2+
594 (k1*l2+k2*l1)*b*c*cosar+
595 (h1*l2+h2*l1)*a*c*cosbr+
596 (h1*k2+h2*k1)*a*b*cosgr);
597 break;
598 case Monoclinic:
599 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
600 return dsp1dsp2 * (h1*h2*a2 + k1*k2*b2 + l1*l2*c2+
601 (k1*l2+k2*l1)*b*c*cosar+
602 (h1*l2+h2*l1)*a*c*cosbr+
603 (h1*k2+h2*k1)*a*b*cosgr);
604 break;
605 case Triclinic:
606 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
607 return dsp1dsp2 * (h1*h2*a2 + k1*k2*b2 + l1*l2*c2+
608 (k1*l2+k2*l1)*b*c*cosar+
609 (h1*l2+h2*l1)*a*c*cosbr+
610 (h1*k2+h2*k1)*a*b*cosgr);
611 break;
612 case Hexagonal:
613 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
614 return dsp1dsp2 *( (h1*h2 + k1*k2 + 0.5*(h1*k2+k1*h2))*a2 + l1*l2*c2);
615 break;
616 default:
617 break;
618 }
619
620 return 0.;
621}
G4double GetIntSp2(G4int h, G4int k, G4int l)

References Amorphous, cosar, cosbr, cosgr, Cubic, GetIntSp2(), GetLatticeSystem(), Hexagonal, Monoclinic, Orthorhombic, Rhombohedral, Tetragonal, theRecSize, and Triclinic.

◆ GetIntSp2()

G4double G4CrystalUnitCell::GetIntSp2 ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 424 of file G4CrystalUnitCell.cc.

426 {
427
428 /* Reference:
429 Table 2.4, pag. 65
430
431 @Inbook{Ladd2003,
432 author="Ladd, Mark and Palmer, Rex",
433 title="Lattices and Space-Group Theory",
434 bookTitle="Structure Determination by X-ray Crystallography",
435 year="2003",
436 publisher="Springer US",
437 address="Boston, MA",
438 pages="51--116",
439 isbn="978-1-4615-0101-5",
440 doi="10.1007/978-1-4615-0101-5_2",
441 url="http://dx.doi.org/10.1007/978-1-4615-0101-5_2"
442 }
443 */
444
445 G4double a = theSize[0], b = theSize[1], c = theSize[2];
446 G4double a2 = a*a, b2 = b*b, c2 = c*c;
447 G4double h2 = h*h, k2 = k*k, l2 = l*l;
448
449 G4double cos2a,sin2a,sin2b;
450 G4double R,T;
451
452 switch(GetLatticeSystem())
453 {
454 case Amorphous:
455 return 0.;
456 break;
457 case Cubic:
458 return a2 / ( h2+k2+l2 );
459 break;
460 case Tetragonal:
461 return 1.0 / ( (h2 + k2)/a2 + l2/c2 );
462 break;
463 case Orthorhombic:
464 return 1.0 / ( h2/a2 + k2/b2 + l2/c2 );
465 break;
466 case Rhombohedral:
467 cos2a=cosa*cosa; sin2a=sina*sina;
468 T = h2+k2+l2+2.*(h*k+k*l+h*l) * ((cos2a-cosa)/sin2a);
469 R = sin2a / (1. - 3*cos2a + 2.*cos2a*cosa);
470 return a*a / (T*R);
471 break;
472 case Monoclinic:
473 sin2b=sinb*sinb;
474 return 1./(1./sin2b * (h2/a2+l2/c2-2*h*l*cosb/(a*c)) + k2/b2);
475 break;
476 case Triclinic:
477 return 1./GetRecIntSp2(h,k,l);
478 break;
479 case Hexagonal:
480 return 1. / ( (4.*(h2+k2+h*k) / (3.*a2)) + l2/c2 );
481 break;
482 default:
483 break;
484 }
485
486 return 0.;
487}
G4double GetRecIntSp2(G4int h, G4int k, G4int l)

References Amorphous, cosa, cosb, Cubic, GetLatticeSystem(), GetRecIntSp2(), Hexagonal, Monoclinic, Orthorhombic, Rhombohedral, sina, sinb, Tetragonal, theSize, and Triclinic.

Referenced by GetIntCosAng().

◆ GetLatticeSystem() [1/2]

theLatticeSystemType G4CrystalUnitCell::GetLatticeSystem ( )
inline

◆ GetLatticeSystem() [2/2]

theLatticeSystemType G4CrystalUnitCell::GetLatticeSystem ( G4int  aGroup)
private

Definition at line 149 of file G4CrystalUnitCell.cc.

149 {
150
151 if( aGroup >= 1 && aGroup <= 2 ) {return Triclinic;}
152 else if(aGroup >= 3 && aGroup <= 15 ) {return Monoclinic;}
153 else if(aGroup >= 16 && aGroup <= 74 ) {return Orthorhombic;}
154 else if(aGroup >= 75 && aGroup <= 142) {return Tetragonal;}
155 else if(aGroup == 146 || aGroup == 148 ||
156 aGroup == 155 || aGroup == 160 ||
157 aGroup == 161 || aGroup == 166 ||
158 aGroup == 167) {return Rhombohedral;}
159 else if(aGroup >= 143 && aGroup <= 167) {return Hexagonal;}
160 else if(aGroup >= 168 && aGroup <= 194) {return Hexagonal;}
161 else if(aGroup >= 195 && aGroup <= 230) {return Cubic;}
162
163 return Amorphous;
164}

References Amorphous, Cubic, Hexagonal, Monoclinic, Orthorhombic, Rhombohedral, Tetragonal, and Triclinic.

◆ GetRecAngle()

G4ThreeVector G4CrystalUnitCell::GetRecAngle ( ) const
inline

Definition at line 118 of file G4CrystalUnitCell.hh.

118{return theRecAngle;}

References theRecAngle.

◆ GetRecBasis()

const G4ThreeVector & G4CrystalUnitCell::GetRecBasis ( G4int  idx) const

Definition at line 192 of file G4CrystalUnitCell.cc.

192 {
193 return (idx>=0 && idx<3 ? theRecBasis[idx] : nullVec);
194}

References nullVec, and theRecBasis.

◆ GetRecIntSp2()

G4double G4CrystalUnitCell::GetRecIntSp2 ( G4int  h,
G4int  k,
G4int  l 
)

Definition at line 491 of file G4CrystalUnitCell.cc.

493 {
494 /* Reference:
495 Table 2.4, pag. 65
496
497 @Inbook{Ladd2003,
498 author="Ladd, Mark and Palmer, Rex",
499 title="Lattices and Space-Group Theory",
500 bookTitle="Structure Determination by X-ray Crystallography",
501 year="2003",
502 publisher="Springer US",
503 address="Boston, MA",
504 pages="51--116",
505 isbn="978-1-4615-0101-5",
506 doi="10.1007/978-1-4615-0101-5_2",
507 url="http://dx.doi.org/10.1007/978-1-4615-0101-5_2"
508 }
509 */
510
511 G4double a = theRecSize[0], b = theRecSize[1], c = theRecSize[2];
512 G4double a2 = a*a, b2 = b*b, c2 = c*c;
513 G4double h2 = h*h, k2 = k*k, l2 = l*l;
514
515 switch(GetLatticeSystem())
516 {
517 case Amorphous:
518 return 0.;
519 break;
520 case Cubic:
521 return a2 * (h2+k2+l2);
522 break;
523 case Tetragonal:
524 return (h2+k2)*a2 + l2*c2 ;
525 break;
526 case Orthorhombic:
527 return h2*a2 + k2+b2 + h2*c2;
528 break;
529 case Rhombohedral:
530 return (h2+k2+l2+2.*(h*k+k*l+h*l) * cosar)*a2;
531 break;
532 case Monoclinic:
533 return h2*a2+k2*b2+l2*c2+2.*h*l*a*c*cosbr;
534 break;
535 case Triclinic:
536 return h2*a2+k2*b2+l2*c2+2.*k*l*b*c*cosar+2.*l*h*c*a*cosbr+2.*h*k*a*b*cosgr;
537 break;
538 case Hexagonal:
539 return (h2+k2+h*k)*a2 + l2*c2;
540 break;
541 default:
542 break;
543 }
544
545 return 0.;
546}

References Amorphous, cosar, cosbr, cosgr, Cubic, GetLatticeSystem(), Hexagonal, Monoclinic, Orthorhombic, Rhombohedral, Tetragonal, theRecSize, and Triclinic.

Referenced by GetIntSp2().

◆ GetRecSize()

G4ThreeVector G4CrystalUnitCell::GetRecSize ( ) const
inline

Definition at line 117 of file G4CrystalUnitCell.hh.

117{return theRecSize;}

References theRecSize.

◆ GetRecUnitBasis()

const G4ThreeVector & G4CrystalUnitCell::GetRecUnitBasis ( G4int  idx) const

Definition at line 186 of file G4CrystalUnitCell.cc.

186 {
187 return (idx>=0 && idx<3 ? theRecUnitBasis[idx] : nullVec);
188}

References nullVec, and theRecUnitBasis.

◆ GetRecVolume()

G4double G4CrystalUnitCell::GetRecVolume ( ) const
inline

Definition at line 152 of file G4CrystalUnitCell.hh.

152{return theRecVolume;} //get the stored volume

References theRecVolume.

◆ GetSize()

G4ThreeVector G4CrystalUnitCell::GetSize ( ) const
inline

Definition at line 101 of file G4CrystalUnitCell.hh.

101{return theSize;}

References theSize.

◆ GetSpaceGroup()

G4int G4CrystalUnitCell::GetSpaceGroup ( ) const
inline

Definition at line 63 of file G4CrystalUnitCell.hh.

63{return theSpaceGroup;};

References theSpaceGroup.

◆ GetUnitBasis()

const G4ThreeVector & G4CrystalUnitCell::GetUnitBasis ( G4int  idx) const

Definition at line 174 of file G4CrystalUnitCell.cc.

174 {
175 return (idx>=0 && idx<3 ? theUnitBasis[idx] : nullVec);
176}

References nullVec, and theUnitBasis.

◆ GetUnitBasisTrigonal()

G4ThreeVector G4CrystalUnitCell::GetUnitBasisTrigonal ( )

Definition at line 198 of file G4CrystalUnitCell.cc.

198 {
199 // Z' axis computed by hand to get both opening angles right
200 // X'.Z' = cos(alpha), Y'.Z' = cos(beta), solve for Z' components
201 G4double x3=cosa, y3=(cosb-cosa*cosg)/sing, z3=std::sqrt(1.-x3*x3-y3*y3);
202 return G4ThreeVector(x3, y3, z3).unit();
203}

References cosa, cosb, cosg, sing, and CLHEP::Hep3Vector::unit().

◆ GetVolume()

G4double G4CrystalUnitCell::GetVolume ( ) const
inline

Definition at line 151 of file G4CrystalUnitCell.hh.

151{return theVolume;} //get the stored volume

References theVolume.

◆ ReflectElReduced()

G4bool G4CrystalUnitCell::ReflectElReduced ( G4double  Cij[6][6]) const
private

Definition at line 375 of file G4CrystalUnitCell.cc.

375 {
376 for (size_t i=1; i<6; i++) {
377 for (size_t j=i+1; j<6; j++) {
378 Cij[j][i] = Cij[i][j];
379 }
380 }
381 return true;
382}

Referenced by FillCubic(), FillOrthorhombic(), FillTetragonal(), and FillTriclinic().

◆ SetSpaceGroup()

void G4CrystalUnitCell::SetSpaceGroup ( G4int  aInt)
inline

Definition at line 64 of file G4CrystalUnitCell.hh.

64{theSpaceGroup=aInt;};

References theSpaceGroup.

Field Documentation

◆ cosa

G4double G4CrystalUnitCell::cosa
private

< struct from SgInfo library needed for further calculations see SgInfo documentation on http://cci.lbl.gov/sginfo/

Definition at line 83 of file G4CrystalUnitCell.hh.

Referenced by ComputeCellVolume(), G4CrystalUnitCell(), GetIntSp2(), and GetUnitBasisTrigonal().

◆ cosar

G4double G4CrystalUnitCell::cosar
private

Definition at line 85 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), GetIntCosAng(), and GetRecIntSp2().

◆ cosb

G4double G4CrystalUnitCell::cosb
private

◆ cosbr

G4double G4CrystalUnitCell::cosbr
private

Definition at line 85 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), GetIntCosAng(), and GetRecIntSp2().

◆ cosg

G4double G4CrystalUnitCell::cosg
private

Definition at line 83 of file G4CrystalUnitCell.hh.

Referenced by ComputeCellVolume(), G4CrystalUnitCell(), and GetUnitBasisTrigonal().

◆ cosgr

G4double G4CrystalUnitCell::cosgr
private

Definition at line 85 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), GetIntCosAng(), and GetRecIntSp2().

◆ nullVec

G4ThreeVector G4CrystalUnitCell::nullVec
protected

◆ sina

G4double G4CrystalUnitCell::sina
private

Definition at line 84 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetIntSp2().

◆ sinb

G4double G4CrystalUnitCell::sinb
private

Definition at line 84 of file G4CrystalUnitCell.hh.

Referenced by ComputeCellVolume(), G4CrystalUnitCell(), and GetIntSp2().

◆ sing

G4double G4CrystalUnitCell::sing
private

Definition at line 84 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetUnitBasisTrigonal().

◆ theAngle

G4ThreeVector G4CrystalUnitCell::theAngle
protected

Definition at line 94 of file G4CrystalUnitCell.hh.

Referenced by GetAngle().

◆ theBasis

G4ThreeVector G4CrystalUnitCell::theBasis[3]
protected

Definition at line 96 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetBasis().

◆ theRecAngle

G4ThreeVector G4CrystalUnitCell::theRecAngle
protected

Definition at line 110 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetRecAngle().

◆ theRecBasis

G4ThreeVector G4CrystalUnitCell::theRecBasis[3]
protected

Definition at line 112 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetRecBasis().

◆ theRecSize

G4ThreeVector G4CrystalUnitCell::theRecSize
protected

Definition at line 109 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), GetIntCosAng(), GetRecIntSp2(), and GetRecSize().

◆ theRecUnitBasis

G4ThreeVector G4CrystalUnitCell::theRecUnitBasis[3]
protected

Definition at line 111 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetRecUnitBasis().

◆ theRecVolume

G4double G4CrystalUnitCell::theRecVolume
private

Definition at line 156 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetRecVolume().

◆ theSize

G4ThreeVector G4CrystalUnitCell::theSize
protected

◆ theSpaceGroup

G4int G4CrystalUnitCell::theSpaceGroup
private

◆ theUnitBasis

G4ThreeVector G4CrystalUnitCell::theUnitBasis[3]
protected

Definition at line 95 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetUnitBasis().

◆ theVolume

G4double G4CrystalUnitCell::theVolume
private

Definition at line 155 of file G4CrystalUnitCell.hh.

Referenced by G4CrystalUnitCell(), and GetVolume().


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