Geant4 11.1.1
Toolkit for the simulation of the passage of particles through matter
Loading...
Searching...
No Matches
G4CrystalUnitCell.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//
27
28//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
29
30// 21-04-16, created by E.Bagli
31
32//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
33
34#include "G4CrystalUnitCell.hh"
36#include <cmath>
37
39 G4double sizeB,
40 G4double sizeC,
41 G4double alpha,
42 G4double beta,
43 G4double gamma,
44 G4int spacegroup):
45theSpaceGroup(spacegroup),
46theSize(G4ThreeVector(sizeA,sizeB,sizeC)),
47theAngle(G4ThreeVector(alpha,beta,gamma))
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
62 cosar = (cosb*cosg-cosa)/(sinb*sing);
63 cosbr = (cosa*cosg-cosb)/(sina*sing);
64 cosgr = (cosa*cosb-cosg)/(sina*sinb);
65
66 theVolume = ComputeCellVolume();
67 theRecVolume = 1. / theVolume;
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
79 switch (GetLatticeSystem(theSpaceGroup)) {
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}
142
143//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
144
146
147//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
148
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}
165
166//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
167/*
168theBravaisLatticeType G4CrystalUnitCell::GetBravaisLattice(G4int aGroup){
169 ;
170}
171*/
172//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
173
175 return (idx>=0 && idx<3 ? theUnitBasis[idx] : nullVec);
176}
177
178//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
179
181 return (idx>=0 && idx<3 ? theBasis[idx] : nullVec);
182}
183
184//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
185
187 return (idx>=0 && idx<3 ? theRecUnitBasis[idx] : nullVec);
188}
189
190//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
191
193 return (idx>=0 && idx<3 ? theRecBasis[idx] : nullVec);
194}
195
196//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
197
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}
204
205//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
206
207G4bool G4CrystalUnitCell::FillAtomicUnitPos(G4ThreeVector& pos, std::vector<G4ThreeVector>& vecout){
208 // Just for testing the infrastructure
209 G4ThreeVector aaa = pos;
210 vecout.push_back(aaa);
211 vecout.emplace_back(2., 5., 3.);
212 return true;
213}
214
215//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
216
217G4bool G4CrystalUnitCell::FillAtomicPos(G4ThreeVector& posin, std::vector<G4ThreeVector>& vecout){
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}
226
227//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
228
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}
261
262//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
263
264G4bool G4CrystalUnitCell::FillAmorphous(G4double Cij[6][6]) const {
265 Cij[3][3] = 0.5*(Cij[0][0]-Cij[0][1]);
266 return true;
267}
268
269//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
270
271G4bool G4CrystalUnitCell::FillCubic(G4double Cij[6][6]) const {
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)
277 {
278 Cij[i][j] = (i == j) ? C11 : C12;
279 }
280 else if(i == j && i >= 3)
281 {
282 Cij[i][i] = C44;
283 }
284 else
285 {
286 Cij[i][j] = 0.;
287 }
288 }
289 }
290
291 ReflectElReduced(Cij);
292
293 return (C11!=0. && C12!=0. && C44!=0.);
294}
295
296//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
297
298G4bool G4CrystalUnitCell::FillTetragonal(G4double Cij[6][6]) const {
299 G4double C11=Cij[0][0], C12=Cij[0][1], C13=Cij[0][2], C16=Cij[0][5];
300 G4double C33=Cij[2][2], C44=Cij[3][3], C66=Cij[5][5];
301
302 Cij[1][1] = C11; // Copy small number of individual elements
303 Cij[1][2] = C13;
304 Cij[1][5] = -C16;
305 Cij[4][4] = C44;
306
307 ReflectElReduced(Cij);
308
309 // NOTE: Do not test for C16 != 0., to allow calling from Hexagonal
310 return (C11!=0. && C12!=0. && C13!=0. && C33!=0. && C44!=0. && C66!=0.);
311}
312
313//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
314
315G4bool G4CrystalUnitCell::FillOrthorhombic(G4double Cij[6][6]) const {
316 // No degenerate elements; just check for all non-zero
317 ReflectElReduced(Cij);
318
319 G4bool good = true;
320 for (size_t i=0; i<6; i++) {
321 for(size_t j = i + 1; j < 3; j++)
322 {
323 good &= (Cij[i][j] != 0);
324 }
325 }
326
327 return good;
328}
329
330//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
331
332G4bool G4CrystalUnitCell::FillRhombohedral(G4double Cij[6][6]) const {
333 G4double C11=Cij[0][0], C12=Cij[0][1], C13=Cij[0][2], C14=Cij[0][3];
334 G4double C15=Cij[0][4], C33=Cij[2][2], C44=Cij[3][3], C66=0.5*(C11-C12);
335
336 Cij[1][1] = C11; // Copy small number of individual elements
337 Cij[1][2] = C13;
338 Cij[1][3] = -C14;
339 Cij[1][4] = -C15;
340 Cij[3][5] = -C15;
341 Cij[4][4] = C44;
342 Cij[4][5] = C14;
343
344 // NOTE: C15 may be zero (c.f. rhombohedral(I) vs. (II))
345 return (C11!=0 && C12!=0 && C13!=0 && C14!=0. &&
346 C33!=0. && C44!=0. && C66!=0.);
347}
348
349//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
350
351G4bool G4CrystalUnitCell::FillMonoclinic(G4double Cij[6][6]) const {
352 // The monoclinic matrix has 13 independent elements with no degeneracies
353 // Sanity condition is same as orthorhombic, plus C45, C(1,2,3)6
354
355 return (FillOrthorhombic(Cij) && Cij[0][5]!=0. && Cij[1][5]!=0. &&
356 Cij[2][5] != 0. && Cij[3][4]!=0.);
357}
358
359//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
360
361G4bool G4CrystalUnitCell::FillTriclinic(G4double Cij[6][6]) const {
362 // The triclinic matrix has the entire upper half filled (21 elements)
363
364 ReflectElReduced(Cij);
365
366 G4bool good = true;
367 for (size_t i=0; i<6; i++) {
368 for(size_t j = i; j < 6; j++)
369 {
370 good &= (Cij[i][j] != 0);
371 }
372 }
373
374 return good;
375}
376
377
378//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
379
380G4bool G4CrystalUnitCell::FillHexagonal(G4double Cij[6][6]) const {
381 Cij[0][5] = 0.;
382 Cij[4][5] = 0.5*(Cij[0][0] - Cij[0][1]);
383 return true;
384}
385
386
387//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
388
389G4bool G4CrystalUnitCell::ReflectElReduced(G4double Cij[6][6]) const {
390 for (size_t i=1; i<6; i++) {
391 for (size_t j=i+1; j<6; j++) {
392 Cij[j][i] = Cij[i][j];
393 }
394 }
395 return true;
396}
397
398//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
399
401 G4double a = theSize[0], b = theSize[1], c = theSize[2];
402
403 switch(GetLatticeSystem())
404 {
405 case Amorphous:
406 return 0.;
407 break;
408 case Cubic:
409 return a * a * a;
410 break;
411 case Tetragonal:
412 return a * a * c;
413 break;
414 case Orthorhombic:
415 return a * b * c;
416 break;
417 case Rhombohedral:
418 return a*a*a*std::sqrt(1.-3.*cosa*cosa+2.*cosa*cosa*cosa);
419 break;
420 case Monoclinic:
421 return a*b*c*sinb;
422 break;
423 case Triclinic:
424 return a*b*c*std::sqrt(1.-cosa*cosa-cosb*cosb-cosg*cosg*2.*cosa*cosb*cosg);
425 break;
426 case Hexagonal:
427 return std::sqrt(3.0)/2.*a*a*c;
428 break;
429 default:
430 break;
431 }
432
433 return 0.;
434}
435
436//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
437
439 G4int k,
440 G4int l){
441
442 /* Reference:
443 Table 2.4, pag. 65
444
445 @Inbook{Ladd2003,
446 author="Ladd, Mark and Palmer, Rex",
447 title="Lattices and Space-Group Theory",
448 bookTitle="Structure Determination by X-ray Crystallography",
449 year="2003",
450 publisher="Springer US",
451 address="Boston, MA",
452 pages="51--116",
453 isbn="978-1-4615-0101-5",
454 doi="10.1007/978-1-4615-0101-5_2",
455 url="http://dx.doi.org/10.1007/978-1-4615-0101-5_2"
456 }
457 */
458
459 G4double a = theSize[0], b = theSize[1], c = theSize[2];
460 G4double a2 = a*a, b2 = b*b, c2 = c*c;
461 G4double h2 = h*h, k2 = k*k, l2 = l*l;
462
463 G4double cos2a,sin2a,sin2b;
464 G4double R,T;
465
466 switch(GetLatticeSystem())
467 {
468 case Amorphous:
469 return 0.;
470 break;
471 case Cubic:
472 return a2 / ( h2+k2+l2 );
473 break;
474 case Tetragonal:
475 return 1.0 / ( (h2 + k2)/a2 + l2/c2 );
476 break;
477 case Orthorhombic:
478 return 1.0 / ( h2/a2 + k2/b2 + l2/c2 );
479 break;
480 case Rhombohedral:
481 cos2a=cosa*cosa; sin2a=sina*sina;
482 T = h2+k2+l2+2.*(h*k+k*l+h*l) * ((cos2a-cosa)/sin2a);
483 R = sin2a / (1. - 3*cos2a + 2.*cos2a*cosa);
484 return a*a / (T*R);
485 break;
486 case Monoclinic:
487 sin2b=sinb*sinb;
488 return 1./(1./sin2b * (h2/a2+l2/c2-2*h*l*cosb/(a*c)) + k2/b2);
489 break;
490 case Triclinic:
491 return 1./GetRecIntSp2(h,k,l);
492 break;
493 case Hexagonal:
494 return 1. / ( (4.*(h2+k2+h*k) / (3.*a2)) + l2/c2 );
495 break;
496 default:
497 break;
498 }
499
500 return 0.;
501}
502
503//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
504
506 G4int k,
507 G4int l){
508 /* Reference:
509 Table 2.4, pag. 65
510
511 @Inbook{Ladd2003,
512 author="Ladd, Mark and Palmer, Rex",
513 title="Lattices and Space-Group Theory",
514 bookTitle="Structure Determination by X-ray Crystallography",
515 year="2003",
516 publisher="Springer US",
517 address="Boston, MA",
518 pages="51--116",
519 isbn="978-1-4615-0101-5",
520 doi="10.1007/978-1-4615-0101-5_2",
521 url="http://dx.doi.org/10.1007/978-1-4615-0101-5_2"
522 }
523 */
524
525 G4double a = theRecSize[0], b = theRecSize[1], c = theRecSize[2];
526 G4double a2 = a*a, b2 = b*b, c2 = c*c;
527 G4double h2 = h*h, k2 = k*k, l2 = l*l;
528
529 switch(GetLatticeSystem())
530 {
531 case Amorphous:
532 return 0.;
533 break;
534 case Cubic:
535 return a2 * (h2+k2+l2);
536 break;
537 case Tetragonal:
538 return (h2+k2)*a2 + l2*c2 ;
539 break;
540 case Orthorhombic:
541 return h2*a2 + k2+b2 + h2*c2;
542 break;
543 case Rhombohedral:
544 return (h2+k2+l2+2.*(h*k+k*l+h*l) * cosar)*a2;
545 break;
546 case Monoclinic:
547 return h2*a2+k2*b2+l2*c2+2.*h*l*a*c*cosbr;
548 break;
549 case Triclinic:
550 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;
551 break;
552 case Hexagonal:
553 return (h2+k2+h*k)*a2 + l2*c2;
554 break;
555 default:
556 break;
557 }
558
559 return 0.;
560}
561
562//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo....
563
565 G4int k1,
566 G4int l1,
567 G4int h2,
568 G4int k2,
569 G4int l2){
570
571 /* Reference:
572 Table 2.4, pag. 65
573
574 @Inbook{Kelly2012,
575 author="Anthony A. Kelly and Kevin M. Knowles",
576 title="Appendix 3 Interplanar Spacings and Interplanar Angles",
577 bookTitle="Crystallography and Crystal Defects, 2nd Edition",
578 year="2012",
579 publisher="John Wiley & Sons, Ltd.",
580 isbn="978-0-470-75014-8",
581 doi="10.1002/9781119961468",
582 url="http://onlinelibrary.wiley.com/book/10.1002/9781119961468"
583 }
584 */
585
586 G4double a = theRecSize[0], b = theRecSize[1], c = theRecSize[2];
587 G4double a2 = a*a, b2 = b*b, c2 = c*c;
588 G4double dsp1dsp2;
589 switch(GetLatticeSystem())
590 {
591 case Amorphous:
592 return 0.;
593 break;
594 case Cubic:
595 return (h1*h2 + k1*k2 + l1+l2) / (std::sqrt(h1*h1 + k1*k1 + l1*l1) * std::sqrt(h2*h2 + k2*k2 + l2*l2));
596 break;
597 case Tetragonal:
598 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
599 return 0. ;
600 break;
601 case Orthorhombic:
602 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
603 return dsp1dsp2 * (h1*h2*a2 + k1*k2*a2 + l1*l2*c2);
604 break;
605 case Rhombohedral:
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 Monoclinic:
613 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
614 return dsp1dsp2 * (h1*h2*a2 + k1*k2*b2 + l1*l2*c2+
615 (k1*l2+k2*l1)*b*c*cosar+
616 (h1*l2+h2*l1)*a*c*cosbr+
617 (h1*k2+h2*k1)*a*b*cosgr);
618 break;
619 case Triclinic:
620 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
621 return dsp1dsp2 * (h1*h2*a2 + k1*k2*b2 + l1*l2*c2+
622 (k1*l2+k2*l1)*b*c*cosar+
623 (h1*l2+h2*l1)*a*c*cosbr+
624 (h1*k2+h2*k1)*a*b*cosgr);
625 break;
626 case Hexagonal:
627 dsp1dsp2 = std::sqrt(GetIntSp2(h1,k1,l1)*GetIntSp2(h2,k2,l2));
628 return dsp1dsp2 *( (h1*h2 + k1*k2 + 0.5*(h1*k2+k1*h2))*a2 + l1*l2*c2);
629 break;
630 default:
631 break;
632 }
633
634 return 0.;
635}
636
637//....oooOO0OOooo........oooOO0OOooo........oooOO0OOooo........oooOO0OOooo......
638
theLatticeSystemType
CLHEP::Hep3Vector G4ThreeVector
double G4double
Definition: G4Types.hh:83
bool G4bool
Definition: G4Types.hh:86
int G4int
Definition: G4Types.hh:85
Hep3Vector & rotateX(double)
Definition: ThreeVector.cc:87
Hep3Vector unit() const
Hep3Vector & rotateZ(double)
Definition: ThreeVector.cc:107
G4ThreeVector theBasis[3]
G4ThreeVector theSize
const G4ThreeVector & GetRecUnitBasis(G4int idx) const
G4ThreeVector GetUnitBasisTrigonal()
const G4ThreeVector & GetRecBasis(G4int idx) const
G4double GetIntSp2(G4int h, G4int k, G4int l)
G4bool FillElReduced(G4double Cij[6][6])
G4CrystalUnitCell(G4double sizeA, G4double sizeB, G4double sizeC, G4double alpha, G4double beta, G4double gamma, G4int spacegroup)
G4bool FillAtomicPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
G4ThreeVector theRecSize
const G4ThreeVector & GetUnitBasis(G4int idx) const
G4ThreeVector theRecUnitBasis[3]
G4double GetRecIntSp2(G4int h, G4int k, G4int l)
G4ThreeVector theRecBasis[3]
G4bool FillAtomicUnitPos(G4ThreeVector &pos, std::vector< G4ThreeVector > &vecout)
G4ThreeVector nullVec
G4ThreeVector theUnitBasis[3]
G4double ComputeCellVolume()
G4ThreeVector theRecAngle
theLatticeSystemType GetLatticeSystem()
const G4ThreeVector & GetBasis(G4int idx) const
G4double GetIntCosAng(G4int h1, G4int k1, G4int l1, G4int h2, G4int k2, G4int l2)
DLL_API const Hep3Vector HepZHat
Definition: ThreeVector.h:419
DLL_API const Hep3Vector HepXHat
DLL_API const Hep3Vector HepYHat
Definition: ThreeVector.h:419