-
Notifications
You must be signed in to change notification settings - Fork 7
/
Copy pathtools.py
1653 lines (1352 loc) · 54.8 KB
/
tools.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
"""
xfab.tools module is a collection of functions
for doing calculation useful in crystallography
"""
from __future__ import absolute_import
from __future__ import print_function
import numpy as n
from math import degrees
from six.moves import range
def find_omega_general(g_w, twoth, w_x, w_y):
"""
For gw find the omega rotation (in radians) around an axis
tilted by w_x radians around x (chi) and w_y radians around y (wedge)
Furthermore find eta (in radians)
Soren Schmidt, implemented by Jette Oddershede
"""
assert abs(n.dot(g_w, g_w)-n.sin(twoth/2)**2) < 1e-9, \
'g-vector must have length sin(theta)'
w_mat_x = n.array([[1, 0 , 0 ],
[0, n.cos(w_x), -n.sin(w_x)],
[0, n.sin(w_x), n.cos(w_x)]])
w_mat_y = n.array([[ n.cos(w_y), 0, n.sin(w_y)],
[0 , 1, 0 ],
[-n.sin(w_y), 0, n.cos(w_y)]])
r_mat = n.dot(w_mat_x,w_mat_y)
a = g_w[0]*r_mat[0][0] + g_w[1]*r_mat[0][1]
b = g_w[0]*r_mat[1][0] - g_w[1]*r_mat[0][0]
c = - n.dot(g_w, g_w) - g_w[2]*r_mat[0][2]
d = a*a + b*b - c*c
omega = []
eta = []
if d < 0:
pass
else:
sq_d = n.sqrt(d)
for i in range(2):
cosomega = (a*c + b*sq_d*(-1)**i)/(a*a+b*b)
sinomega = (b*c + a*sq_d*(-1)**(i+1))/(a*a+b*b)
omega.append(n.arctan2(sinomega, cosomega))
if omega[i] > n.pi:
omega[i] = omega[i] - 2*n.pi
omega_mat = form_omega_mat_general(omega[i], w_x, w_y)
g_t = n.dot(omega_mat, g_w)
sineta = -2*g_t[1]/n.sin(twoth)
coseta = 2*g_t[2]/n.sin(twoth)
eta.append(n.arctan2(sineta, coseta))
return n.array(omega), n.array(eta)
def find_omega_quart(g_w, twoth, w_x, w_y):
"""
For gw find the omega rotation (in radians) around an axis
tilted by w_x radians around x (chi) and w_y radians around y (wedge)
Furthermore find eta (in radians)
Soren Schmidt, implemented by Jette Oddershede
"""
assert abs(n.dot(g_w, g_w)-n.sin(twoth/2)**2) < 1e-9, \
'g-vector must have length sin(theta)'
w_mat_x = n.array([[1, 0 , 0 ],
[0, n.cos(w_x), -n.sin(w_x)],
[0, n.sin(w_x), n.cos(w_x)]])
w_mat_y = n.array([[ n.cos(w_y), 0, n.sin(w_y)],
[0 , 1, 0 ],
[-n.sin(w_y), 0, n.cos(w_y)]])
normal = n.dot(w_mat_x, n.dot(w_mat_y, n.array([0, 0, 1])))
a = g_w[0]*(1-normal[0]**2) - \
g_w[1]*normal[0]*normal[1] - \
g_w[2]*normal[0]*normal[2]
b = g_w[2]*normal[1] - g_w[1]*normal[2]
c = - n.dot(g_w, g_w) - g_w[0]*normal[0]**2 - \
g_w[1]*normal[0]*normal[1] - \
g_w[2]*normal[0]*normal[2]
d = a*a + b*b - c*c
omega = []
eta = []
if d < 0:
pass
else:
sq_d = n.sqrt(d)
for i in range(2):
cosomega = (a*c + b*sq_d*(-1)**i)/(a*a+b*b)
sinomega = (b*c + a*sq_d*(-1)**(i+1))/(a*a+b*b)
omega.append(n.arctan2(sinomega, cosomega))
if omega[i] > n.pi:
omega[i] = omega[i] - 2*n.pi
omega_mat = quart_to_omega(omega[i]*180./n.pi, w_x, w_y)
g_t = n.dot(omega_mat, g_w)
sineta = -2*g_t[1]/n.sin(twoth)
coseta = 2*g_t[2]/n.sin(twoth)
eta.append(n.arctan2(sineta, coseta))
return n.array(omega), n.array(eta)
def find_omega_wedge(g_w, twoth, wedge):
"""
calculate the omega and eta angles for a g-vector g_w given twotheta of
g_w is provided. The calculation takes a possible wedge angle intp account.
This code is a translation of c-code from GrainSpotter by Soren Schmidt
"""
#Normalize G-vector
g_w = g_w/n.sqrt(n.dot(g_w, g_w))
costth = n.cos(twoth)
sintth = n.sin(twoth)
cosfactor = costth -1
length = n.sqrt(-2*cosfactor)
sinwedge = n.sin(wedge)
coswedge = n.cos(wedge)
# Determine cos(eta)
coseta = ( g_w[2]*length + sinwedge * cosfactor ) / coswedge / sintth
omega = []
eta = []
if (abs(coseta) > 1.):
return omega, eta
# else calc the two eta values
eta = n.array([n.arccos(coseta), -n.arccos(coseta)])
# replaced this by the above to make find_omega_wedge
# and find_omega_quart give same eta
# eta = n.array([n.arccos(coseta), 2*n.pi-n.arccos(coseta)])
#print eta*180.0/n.pi
# Now find the Omega value(s)
# A slight change in the code from GrainSpotter: the lenght
# Here the original a and b is divided by length since
# these can be only scale somega and comega equally
a = (coswedge * cosfactor + sinwedge * sintth * coseta)
for i in range(2):
b = -sintth * n.sin(eta[i])
somega = (b*g_w[0] - a*g_w[1])/(a*a + b*b)
comega = (g_w[0] - b*somega)/a
omega.append(n.arctan2(somega, comega))
if omega[i] > n.pi:
omega[i] = omega[i] - 2*n.pi
return n.array(omega), eta
def find_omega(g_w, twoth):
"""
Calculate the omega angles for a g-vector gw given twotheta using Soeren
Schmidts algorithm.
Solves an equation of type a*cos(w)+b*sin(w) = c by the fixpoint method.
"""
g_g = n.sqrt(n.dot(g_w, g_w))
costth = n.cos(twoth)
a = g_w[0]/g_g
b = -g_w[1]/g_g
c = (costth - 1)/n.sqrt(2*(1 - costth))
d = a**2 + b**2
sq_d = d - c**2
omega = []
if sq_d > 0:
sq_d = n.sqrt(sq_d)
comega = (a*c + b*sq_d)/d
somega = (b*c - a*sq_d)/d
omega.append(n.arccos(comega))
# if omega[0] > n.pi:
# omega[0] = omega[0] - 2*n.pi
if somega < 0:
omega[0] = -omega[0]
comega = comega - 2*b*sq_d/d
somega = somega + 2*a*sq_d/d
omega.append(n.arccos(comega))
# if omega[1] > n.pi:
# omega[1] = omega[1] - 2*n.pi
if somega < 0:
omega[1] = -omega[1]
return n.array(omega)
def cell_invert(unit_cell):
"""
cell_invert calculates the reciprocal unit cell parameters
from the direct space unit cell
INPUT: unit_cell = [a, b, c, alpha, beta, gamma]
OUTPUT unit_cell_star = [astar, bstar, cstar, alpastar, betastar, gammastar]
"""
a = unit_cell[0]
b = unit_cell[1]
c = unit_cell[2]
calp = n.cos(unit_cell[3]*n.pi/180.)
cbet = n.cos(unit_cell[4]*n.pi/180.)
cgam = n.cos(unit_cell[5]*n.pi/180.)
salp = n.sin(unit_cell[3]*n.pi/180.)
sbet = n.sin(unit_cell[4]*n.pi/180.)
sgam = n.sin(unit_cell[5]*n.pi/180.)
V = cell_volume(unit_cell)
astar = b*c*salp/V
bstar = a*c*sbet/V
cstar = a*b*sgam/V
#salpstar = V/(a*b*c*sbet*sgam)
#sbetstar = V/(a*b*c*salp*sgam)
#sgamstar = V/(a*b*c*salp*sbet)
calpstar = (cbet*cgam-calp)/(sbet*sgam)
cbetstar = (calp*cgam-cbet)/(salp*sgam)
cgamstar = (calp*cbet-cgam)/(salp*sbet)
alpstar = n.arccos(calpstar)*180./n.pi
betstar = n.arccos(cbetstar)*180./n.pi
gamstar = n.arccos(cgamstar)*180./n.pi
return [astar, bstar, cstar, alpstar, betstar, gamstar]
def form_omega_mat(omega):
"""
Calc Omega rotation matrix having an omega angle of "omega"
INPUT: omega (in radians)
OUTPUT: Omega rotation matrix
"""
Om = n.array([[n.cos(omega), -n.sin(omega), 0],
[n.sin(omega), n.cos(omega), 0],
[ 0 , 0 , 1]])
return Om
def form_omega_mat_general(omega,chi,wedge):
"""
Calc Omega rotation matrix having an omega angle of "omega"
INPUT: omega,chi and wedge (in radians)
OUTPUT: Omega rotation matrix
"""
phi_x = n.array([[1, 0 , 0 ],
[0, n.cos(chi), -n.sin(chi)],
[0, n.sin(chi), n.cos(chi)]])
phi_y = n.array([[ n.cos(wedge), 0, n.sin(wedge)],
[0 , 1, 0 ],
[-n.sin(wedge), 0, n.cos(wedge)]])
Om = form_omega_mat(omega)
Om = n.dot(phi_x,n.dot(phi_y,Om))
return Om
def cell_volume(unit_cell):
"""
cell_volume calculates the volume of the unit cell in AA^3
from the direct space unit cell parameters
INPUT: unit_cell = [a, b, c, alpha, beta, gamma]
OUTPUT: volume
"""
a = unit_cell[0]
b = unit_cell[1]
c = unit_cell[2]
calp = n.cos(unit_cell[3]*n.pi/180.)
cbet = n.cos(unit_cell[4]*n.pi/180.)
cgam = n.cos(unit_cell[5]*n.pi/180.)
angular = n.sqrt(1 - calp*calp - cbet*cbet - cgam*cgam + 2*calp*cbet*cgam)
#Volume of unit cell
V = a*b*c*angular
return V
def form_b_mat(unit_cell):
"""
calculate B matrix of (Gcart = B Ghkl) following eq. 3.4 in
H.F. Poulsen.
Three-dimensional X-ray diffraction microscopy.
Mapping polycrystals and their dynamics.
Springer Tracts in Modern Physics, v. 205), (Springer, Berlin, 2004).
INPUT: unit_cell - unit_cell = [a, b, c, alpha, beta, gamma]
OUTPUT: B - a 3x3 matrix
Henning Osholm Sorensen, Risoe-DTU, June 11, 2007.
"""
a = unit_cell[0]
b = unit_cell[1]
c = unit_cell[2]
calp = n.cos(unit_cell[3]*n.pi/180.)
cbet = n.cos(unit_cell[4]*n.pi/180.)
cgam = n.cos(unit_cell[5]*n.pi/180.)
salp = n.sin(unit_cell[3]*n.pi/180.)
sbet = n.sin(unit_cell[4]*n.pi/180.)
sgam = n.sin(unit_cell[5]*n.pi/180.)
#Volume of unit cell
V = cell_volume(unit_cell)
# Calculate reciprocal lattice parameters:
# NOTICE PHYSICIST DEFINITION of recip axes with 2*pi
astar = 2*n.pi*b*c*salp/V
bstar = 2*n.pi*a*c*sbet/V
cstar = 2*n.pi*a*b*sgam/V
#salpstar = V/(a*b*c*sbet*sgam)
sbetstar = V/(a*b*c*salp*sgam)
sgamstar = V/(a*b*c*salp*sbet)
#calpstar = (cbet*cgam-calp)/(sbet*sgam)
cbetstar = (calp*cgam-cbet)/(salp*sgam)
cgamstar = (calp*cbet-cgam)/(salp*sbet)
# Form B matrix following eq. 3.4 in H.F Poulsen
B = n.array([[astar, bstar*cgamstar, cstar*cbetstar ],
[0, bstar*sgamstar, -cstar*sbetstar*calp ],
[0, 0, cstar*sbetstar*salp ]])
return B
def form_a_mat(unit_cell):
"""
calculate the A matrix given in eq. 3.23 of H.F. Poulsen.
Three-dimensional X-ray diffraction microscopy.
Mapping polycrystals and their dynamics.
(Springer Tracts in Modern Physics, v. 205), (Springer, Berlin, 2004).
INPUT: unit_cell - unit_cell = [a, b, c, alpha, beta, gamma]
OUTPUT A - a 3x3 matrix
Jette Oddershede, March 7, 2008.
"""
a = unit_cell[0]
b = unit_cell[1]
c = unit_cell[2]
calp = n.cos(unit_cell[3]*n.pi/180.)
cbet = n.cos(unit_cell[4]*n.pi/180.)
cgam = n.cos(unit_cell[5]*n.pi/180.)
#salp = n.sin(unit_cell[3]*n.pi/180.)
sbet = n.sin(unit_cell[4]*n.pi/180.)
sgam = n.sin(unit_cell[5]*n.pi/180.)
#Volume of unit cell
V = cell_volume(unit_cell)
# Calculate reciprocal lattice parameters
salpstar = V/(a*b*c*sbet*sgam)
calpstar = (cbet*cgam-calp)/(sbet*sgam)
# Form A matrix following eq. 3.23 in H.F Poulsen
A = n.array([[a, b*cgam, c*cbet ],
[0, b*sgam, -c*sbet*calpstar ],
[0, 0, c*sbet*salpstar ]])
return A
def form_a_mat_inv(unit_cell):
"""
calculate the inverse of the A matrix given in eq. 3.23 of H.F. Poulsen.
Three-dimensional X-ray diffraction microscopy.
Mapping polycrystals and their dynamics.
(Springer Tracts in Modern Physics, v. 205), (Springer, Berlin, 2004).
INPUT: unit_cell - unit_cell = [a, b, c, alpha, beta, gamma]
OUTPUT A^-1 - a 3x3 matrix
Jette Oddershede, March 7, 2008.
"""
A = form_a_mat(unit_cell)
Ainv = n.linalg.inv(A)
return Ainv
def ubi_to_cell(ubi):
"""
calculate lattice constants from the UBI-matrix as
defined in H.F.Poulsen 2004 eqn.3.23
ubi_to_cell(ubi)
ubi [3x3] matrix of (U*B)^-1
in this case B = B /2pi
returns unit_cell = [a, b, c, alpha, beta, gamma]
"""
return n.array(a_to_cell(n.transpose(ubi)))
def ubi_to_u(ubi):
"""
calculate lattice constants from the UBI-matrix
defined(U*B)^-1 and is B from form_b_mat devided by 2pi
ubi_to_u(ubi)
ubi [3x3] matrix
returns U matrix
"""
unit_cell = ubi_to_cell(ubi)
B = form_b_mat(unit_cell)
U = n.transpose(n.dot(B, ubi))/(2*n.pi)
return U
def ubi_to_u_and_eps(ubi,unit_cell):
"""
calculate lattice lattice rotation and strain from the UBI-matrix
(U,eps) = ubi_to_u_and_eps(ubi,unit_cell)
ubi [3x3] matrix, (UB)^-1, where B=B/2*pi
unit_cell = [a,b,c,alpha,beta,gamma]
returns U matrix and strain tensor components
eps = [e11, e12, e13, e22, e23, e33]
"""
unit_cell_ubi = ubi_to_cell(ubi)
B_ubi = form_b_mat(unit_cell_ubi)
U = n.transpose(n.dot(B_ubi, ubi))/(2*n.pi)
eps = b_to_epsilon(B_ubi, unit_cell)
return (U,eps)
def a_to_cell(A):
"""
calculate lattice constants from the A-matix as
defined in H.F.Poulsen 2004 eqn.3.23
a_to_cell(A)
A [3x3] upper triangular matrix
returns unit_cell = [a, b, c, alpha, beta, gamma]
Jette Oddershede, March 10, 2008.
"""
g = n.dot(n.transpose(A), A)
a = n.sqrt(g[0, 0])
b = n.sqrt(g[1, 1])
c = n.sqrt(g[2, 2])
alpha = degrees(n.arccos(g[1, 2]/b/c))
beta = degrees(n.arccos(g[0, 2]/a/c))
gamma = degrees(n.arccos(g[0, 1]/a/b))
unit_cell = [a, b, c, alpha, beta, gamma]
return unit_cell
def b_to_cell(B):
"""
calculate lattice constants from the B-matix as
defined in H.F.Poulsen 2004 eqn.3.4
B [3x3] upper triangular matrix
returns unit_cell = [a, b, c, alpha, beta, gamma]
Jette Oddershede, April 21, 2008.
"""
B = B/(2*n.pi)
g = n.dot(n.transpose(B), B)
astar = n.sqrt(g[0, 0])
bstar = n.sqrt(g[1, 1])
cstar = n.sqrt(g[2, 2])
alphastar = degrees(n.arccos(g[1, 2]/bstar/cstar))
betastar = degrees(n.arccos(g[0, 2]/astar/cstar))
gammastar = degrees(n.arccos(g[0, 1]/astar/bstar))
unit_cell = cell_invert([astar, bstar, cstar,
alphastar, betastar, gammastar])
return unit_cell
def epsilon_to_b_old(epsilon, unit_cell):
"""
calculate B matrix of (Gcart = B Ghkl) from epsilon and
unstrained cell as in H.F. Poulsen (2004) page 33.
INPUT: epsilon - strain tensor [e11, e12, e13, e22, e23, e33]
unit_cell - unit cell = [a, b, c, alpha, beta, gamma]
OUTPUT: B - [3x3] for strained lattice constants
Jette Oddershede, March 10, 2008.
"""
A0inv = form_a_mat_inv(unit_cell)
A = n.zeros((3, 3))
A[0, 0] = (epsilon[0]+1)/A0inv[0, 0]
A[1, 1] = (epsilon[3]+1)/A0inv[1, 1]
A[2, 2] = (epsilon[5]+1)/A0inv[2, 2]
A[0, 1] = (2*epsilon[1]-A[0, 0]*A0inv[0, 1])/A0inv[1, 1]
A[1, 2] = (2*epsilon[4]-A[1, 1]*A0inv[1, 2])/A0inv[2, 2]
A[0, 2] = (2*epsilon[2]-A[0, 0]*A0inv[0, 2]-A[0, 1]*A0inv[1, 2])/A0inv[2, 2]
strainedcell = a_to_cell(A)
B = form_b_mat(strainedcell)
return B
def b_to_epsilon_old(B, unit_cell):
"""
calculate epsilon from the the unstrained cell and
the B matrix of (Gcart = B Ghkl) as in H.F. Poulsen (2004) page 33.
INPUT: B - upper triangular 3x3 matrix of strained lattice constants
unit_cell - unit cell = [a, b, c, alpha, beta, gamma]
of unstrained lattice
OUTPUT: epsilon = [e11, e12, e13, e22, e23, e33]
Jette Oddershede, April 21, 2008.
"""
A0inv = form_a_mat_inv(unit_cell)
A = form_a_mat(b_to_cell(B))
T = n.dot(A, A0inv)
I = n.eye(3)
eps = 0.5*(T+n.transpose(T))-I
epsilon = [eps[0, 0], eps[0, 1], eps[0, 2],
eps[1, 1], eps[1, 2], eps[2, 2]]
return epsilon
def b_to_epsilon(B, unit_cell):
"""
calculate epsilon from the the unstrained cell and
the B matrix of (Gcart = B Ghkl).
The algorithm is similar to H.F. Poulsen (2004) page 33,
except A is defined such that A'B = I
T' = (A*A0inv)' = A0inv'*A' = B0*Binv
This definition of A ensures that U relating Gcart to Gsam
also relates epsilon_cart to epsilon_sam
INPUT: B - upper triangular 3x3 matrix of strained lattice constants
unit_cell - unit cell = [a, b, c, alpha, beta, gamma]
of unstrained lattice
OUTPUT: epsilon = [e11, e12, e13, e22, e23, e33]
Jette Oddershede, [email protected], January 2012.
"""
B0 = form_b_mat(unit_cell)
T = n.dot(B0,n.linalg.inv(B))
I = n.eye(3)
eps = 0.5*(T+n.transpose(T))-I
epsilon = [eps[0, 0], eps[0, 1], eps[0, 2],
eps[1, 1], eps[1, 2], eps[2, 2]]
return epsilon
def epsilon_to_b(epsilon, unit_cell):
"""
calculate B matrix of (Gcart = B Ghkl) from epsilon and
unstrained cell as in H.F. Poulsen (2004) page 33 with the
exception that A is defined such that A'B = I
2*(epsilon+I) = B0*Binv + (B0*Binv)'
INPUT: epsilon - strain tensor [e11, e12, e13, e22, e23, e33]
unit_cell - unit cell = [a, b, c, alpha, beta, gamma]
of unstrained lattice
OUTPUT: B - [3x3] for strained lattice constants
Jette Oddershede, [email protected], January 2012.
"""
B0 = form_b_mat(unit_cell)
Binv = n.zeros((3, 3))
Binv[0, 0] = (epsilon[0]+1)/B0[0, 0]
Binv[1, 1] = (epsilon[3]+1)/B0[1, 1]
Binv[2, 2] = (epsilon[5]+1)/B0[2, 2]
Binv[0, 1] = (2*epsilon[1]-B0[0, 1]*Binv[1, 1])/B0[0, 0]
Binv[1, 2] = (2*epsilon[4]-B0[1, 2]*Binv[2, 2])/B0[1, 1]
Binv[0, 2] = (2*epsilon[2]-B0[0, 1]*Binv[1, 2]-B0[0, 2]*Binv[2, 2])/B0[0, 0]
return n.linalg.inv(Binv)
def euler_to_u(phi1, PHI, phi2):
"""
U matrix from Euler angles phi1, PHI, phi2.
The formalism follows the ID11-3DXRD specs
U = euler_to_u(phi1, PHI, phi2)
INPUT: phi1, PHI, and phi2 angles in radians
OUTPUT: U = array([[U11, U12, U13],[ U21, U22, U23],[U31, U32, U33]])
Changed input angles to be in radians instead of degrees
Henning Osholm Sorensen, Riso National Laboratory, June 23, 2006.
Translated from MATLAB to python by Jette Oddershede, March 26 2008
Origingal MATLAB code from: Henning Poulsen, Risoe 15/6 2002.
"""
U = n.zeros((3, 3))
U[0, 0] = n.cos(phi1)*n.cos(phi2)-n.sin(phi1)*n.sin(phi2)*n.cos(PHI)
U[1, 0] = n.sin(phi1)*n.cos(phi2)+n.cos(phi1)*n.sin(phi2)*n.cos(PHI)
U[2, 0] = n.sin(phi2)*n.sin(PHI)
U[0, 1] = -n.cos(phi1)*n.sin(phi2)-n.sin(phi1)*n.cos(phi2)*n.cos(PHI)
U[1, 1] = -n.sin(phi1)*n.sin(phi2)+n.cos(phi1)*n.cos(phi2)*n.cos(PHI)
U[2, 1] = n.cos(phi2)*n.sin(PHI)
U[0, 2] = n.sin(phi1)*n.sin(PHI)
U[1, 2] = -n.cos(phi1)*n.sin(PHI)
U[2, 2] = n.cos(PHI)
return U
def _arctan2(y, x):
"""Modified arctan function used locally in u_to_euler().
"""
tol = 1e-8
if n.abs(x)<tol: x = 0
if n.abs(y)<tol: y = 0
if x>0:
return n.arctan(y/x)
elif x<0 and y>=0:
return n.arctan(y/x) + n.pi
elif x<0 and y<0:
return n.arctan(y/x) - n.pi
elif x==0 and y>0:
return n.pi/2
elif x==0 and y<0:
return -n.pi/2
elif x==0 and y==0:
raise ValueError('Local function _arctan2() does not accept arguments (0,0)')
def u_to_euler(U):
"""Convert unitary 3x3 rotation matrix into Euler angles in Bunge notation.
The returned Euler angles are all in the range [0, 2*pi]. If Gimbal lock occurs
(PHI=0) infinite number of solutions exists. The solution returned is phi2=0.
Implementation is based on the notes by
Depriester, Dorian. (2018). Computing Euler angles with Bunge convention from rotation matrix.
https://www.researchgate.net/publication/324088567_Computing_Euler_angles_with_Bunge_convention_from_rotation_matrix
notationwise U = g.T in in these notes.
INPUT:
U : unitary 3x3 rotation matrix as a numpy array (shape=(3,3))
RETURNS
angles : Euler angles in radians. numpy array as, [phi1, PHI, phi2].
Last Modified: Axel Henningsson, January 2021
"""
tol = 1e-8
PHI = n.arccos(U[2, 2])
if n.abs(PHI)<tol:
phi1 = _arctan2(-U[0, 1], U[0, 0])
phi2 = 0
elif n.abs(PHI-n.pi)<tol:
phi1 = _arctan2(U[0, 1], U[0, 0])
phi2 = 0
else:
phi1 = _arctan2(U[0, 2], -U[1, 2])
phi2 = _arctan2(U[2, 0], U[2, 1])
if phi1<0:
phi1 = phi1 + 2*n.pi
if phi2<0:
phi2 = phi2 + 2*n.pi
return n.array([ phi1, PHI, phi2 ])
# def U2rod_old(U):
# """
# Get Rodrigues vector from U matrix (Busing Levy)
# Taken from ImageD11.indexing.ubitoRod of Jon Wright
# added argsort of w
# """
# w, v = n.linalg.eig(U)
# #print w
# #print v
# order = n.argsort(w.real)
# #print w, order
# ehat = v[:, order[-1]]
# print order
# # if order.tolist() != range(3):
# # print 'HHFH'
# # angle = -1*n.arccos(w[order[1]].real)
# # else:
# # angle = n.arccos(w[order[1]].real)
# angle = n.arccos(w[order[1]].real)
# Rod = ehat * n.tan(angle/2.)
# return Rod.real
def u_to_rod(U):
"""
Get Rodrigues vector from U matrix (Busing Levy)
INPUT: U 3x3 matrix
OUTPUT: Rodrigues vector
Function taken from GrainsSpotter by Soeren Schmidt
"""
ttt = 1+U[0, 0]+U[1, 1]+U[2, 2]
if abs(ttt) < 1e-16:
raise ValueError('Wrong trace of U')
a = 1/ttt
r1 = (U[1, 2]-U[2, 1])*a
r2 = (U[2, 0]-U[0, 2])*a
r3 = (U[0, 1]-U[1, 0])*a
return n.array([r1, r2, r3])
def u_to_ubi(u_mat,unit_cell):
"""
Get UBI matrix from U matrix and unit cell
INPUT: U orientaion matrix and unit cell
OUTPUT: UBI 3x3 matrix
"""
b_mat = form_b_mat(unit_cell)
return n.linalg.inv(n.dot(u_mat,b_mat))*(2*n.pi)
def ubi_to_rod(ubi):
"""
Get Rodrigues vector from UBI matrix
INPUT: UBI 3x3 matrix
OUTPUT: Rodrigues vector
"""
return u_to_rod(ubi_to_u(ubi))
def ubi_to_u_b(ubi):
"""
Get Rodrigues vector from UBI matrix
INPUT: UBI 3x3 matrix
OUTPUT: U orientaion matrix and B metric matrix
"""
return ub_to_u_b(n.linalg.inv(ubi)*(2*n.pi))
def rod_to_u(r):
"""
rod_to_u calculates the U orientation matrix given an oriention
represented in Rodrigues space. r = [r1, r2, r3]
"""
g = n.zeros((3, 3))
r2 = n.dot(r , r)
for i in range(3):
for j in range(3):
if i == j:
fac = 1
else:
fac = 0
term = 0
for k in range(3):
if [i, j, k] == [0, 1, 2] or \
[i, j, k] == [1, 2, 0] or \
[i, j, k] == [2, 0, 1]:
sign = 1
elif [i, j, k] == [2, 1, 0] or \
[i, j, k] == [0, 2, 1] or \
[i, j, k] == [1, 0, 2]:
sign = -1
else:
sign = 0
term = term + 2*sign*r[k]
g[i, j] = 1/(1+r2) * ((1-r2)*fac + 2*r[i]*r[j] - term)
return n.transpose(g)
def ub_to_u_b(UB):
"""
qr decomposition to get U unitary and B upper triangular with positive
diagonal from UB
"""
(U, B) = n.linalg.qr(UB)
if B[0, 0] < 0:
B[0, 0] = -B[0, 0]
B[0, 1] = -B[0, 1]
B[0, 2] = -B[0, 2]
U[0, 0] = -U[0, 0]
U[1, 0] = -U[1, 0]
U[2, 0] = -U[2, 0]
if B[1, 1] < 0:
B[1, 1] = -B[1, 1]
B[1, 2] = -B[1, 2]
U[0, 1] = -U[0, 1]
U[1, 1] = -U[1, 1]
U[2, 1] = -U[2, 1]
if B[2, 2] < 0:
B[2, 2] = -B[2, 2]
U[0, 2] = -U[0, 2]
U[1, 2] = -U[1, 2]
U[2, 2] = -U[2, 2]
return (U, B)
def reduce_cell(unit_cell,uvw = 3):
"""
reduce unit cell
INPUT: unit_cell - unit_cell = [a, b, c, alpha, beta, gamma]
OUTPUT unit_cell reduced - array([a, b, c, alpha, beta, gamma])
"""
res = n.zeros((0, 4))
red_a_mat = n.zeros((3, 3))
a_mat = form_a_mat(unit_cell)
for i in n.arange(-uvw, uvw):
for j in n.arange(-uvw, uvw):
for k in n.arange(-uvw, uvw):
tmp = n.dot(a_mat, n.array([i, j, k]))
res = n.concatenate((res, [[i, j, k, n.linalg.norm(tmp)]]))
res = res[n.argsort(res[:, 3]), :]
red_a_mat[0] = n.dot(a_mat, res[1, :3])
for i in range(2, len(res)):
tmp = n.dot(a_mat, res[i, :3])
kryds = n.cross(tmp, red_a_mat[0])
if n.sum(n.abs(kryds)) > 0.00001:
red_a_mat[1] = tmp
break
for j in range(i, len(res)):
tmp = n.dot(a_mat, res[j,:3])
dist = n.dot(kryds, tmp)/n.linalg.norm(kryds)
if dist > 0.00001:
red_a_mat[2] = tmp
break
return a_to_cell(red_a_mat)
def detect_tilt(tilt_x, tilt_y, tilt_z):
"""
Calculate the tilt matrix
tiltR(tilt_x,tilt_y,tilt_z)
input tilt_x, tilt_y, tilt_z
Henning Osholm Sorensen 2006
"""
Rx = n.array([[ 1, 0, 0],
[ 0, n.cos(tilt_x), -n.sin(tilt_x)],
[ 0, n.sin(tilt_x), n.cos(tilt_x)]])
Ry = n.array([[ n.cos(tilt_y), 0, n.sin(tilt_y)],
[ 0, 1, 0],
[ -n.sin(tilt_y), 0, n.cos(tilt_y)]])
Rz = n.array([[ n.cos(tilt_z), -n.sin(tilt_z), 0],
[ n.sin(tilt_z), n.cos(tilt_z), 0],
[ 0, 0, 1]])
R = n.dot(Rx, n.dot(Ry, Rz))
return R
def quart_to_omega(w, w_x, w_y):
"""
Calculate the Omega rotation matrix given w (the motorised rotation in
degrees, usually around the z-axis).
wx and wy (the rotations around x and y bringing the z-axis to the true
rotation axis, in radians).
Quarternions are used for the calculations to avoid singularities in
subsequent refinements.
"""
whalf = w*n.pi/360.
w_mat_x = n.array([[1, 0 , 0 ],
[0, n.cos(w_x), -n.sin(w_x)],
[0, n.sin(w_x), n.cos(w_x)]])
w_mat_y = n.array([[ n.cos(w_y), 0, n.sin(w_y)],
[0 , 1, 0 ],
[-n.sin(w_y), 0, n.cos(w_y)]])
qua = n.dot(w_mat_x, n.dot(w_mat_y, n.array([0, 0, n.sin(whalf)])))
q = [n.cos(whalf), qua[0], qua[1], qua[2]]
omega_mat = n.array([[1-2*q[2]**2-2*q[3]**2 ,
2*q[1]*q[2]-2*q[3]*q[0],
2*q[1]*q[3]+2*q[2]*q[0]],
[2*q[1]*q[2]+2*q[3]*q[0],
1-2*q[1]**2-2*q[3]**2 ,
2*q[2]*q[3]-2*q[1]*q[0]],
[2*q[1]*q[3]-2*q[2]*q[0],
2*q[2]*q[3]+2*q[1]*q[0],
1-2*q[1]**2-2*q[2]**2]])
return omega_mat
def sintl(unit_cell, hkl):
"""
sintl calculate sin(theta)/lambda of the reflection "hkl" given
the unit cell "unit_cell"
sintl(unit_cell,hkl)
INPUT: unit_cell = [a, b, c, alpha, beta, gamma]
hkl = [h, k, l]
OUTPUT: sin(theta)/lambda
Henning Osholm Sorensen, Risoe National Laboratory, June 23, 2006.
"""
a = float(unit_cell[0])
b = float(unit_cell[1])
c = float(unit_cell[2])
calp = n.cos(unit_cell[3]*n.pi/180.)
cbet = n.cos(unit_cell[4]*n.pi/180.)
cgam = n.cos(unit_cell[5]*n.pi/180.)
(h, k, l) = hkl
part1 = (h*h/a**2) * (1-calp**2) + (k*k/b**2) *\
(1-cbet**2) + (l*l/c**2) * (1-cgam**2) +\
2*h*k*(calp*cbet-cgam)/(a*b) + 2*h*l*(calp*cgam-cbet)/(a*c) +\
2*k*l*(cbet*cgam-calp)/(b*c)
part2 = 1 - (calp**2 + cbet**2 + cgam**2) + 2*calp*cbet*cgam
stl = n.sqrt(part1) / (2*n.sqrt(part2))
return stl
def tth(unit_cell, hkl, wavelength):
"""
tth calculate two theta of the reflection given
the unit cell and wavelenght
INPUT: unit_cell = [a, b, c, alpha, beta, gamma] (in Angstroem and degrees)
hkl = [h, k, l]
wavelenth (in Angstroem)
OUTPUT: twotheta (in radians)
Henning Osholm Sorensen, Risoe-DTU, July 16, 2008.
"""
stl = sintl(unit_cell, hkl) # calls sintl function in tools
twotheta = 2*n.arcsin(wavelength*stl)
return twotheta
def tth2(gve, wavelength):
"""
calculates two theta for a scattering vector given the wavelenght
INPUT: gve: scattering vector
(defined in reciprocal space (2*pi/lambda))
wavelenth (in Angstroem)
OUTPUT: twotheta (in radians)
Henning Osholm Sorensen, Risoe DTU, July 17, 2008.
"""
length = n.sqrt(n.dot(gve, gve))
twotheta = 2.0*n.arcsin(length*wavelength/(4*n.pi))
return twotheta
def genhkl_all(unit_cell, sintlmin, sintlmax, sgname=None, sgno=None, cell_choice='standard', output_stl=False):
"""
Generate the full set of reflections given a unit cell and space group up to maximum sin(theta)/lambda (sintlmax)
The function is using the function genhkl_base for the actual generation
INPUT: unit cell : [a , b, c, alpha, beta, gamma]
sintlmin : minimum sin(theta)/lambda for generated reflections
sintlmax : maximum sin(theta)/lambda for generated reflections
sgno/sgname : provide either the space group number or its name
e.g. sgno=225 or equivalently
sgname='Fm-3m'
output_stl : Should sin(theta)/lambda be output (True/False)
default=False
OUTPUT: list of reflections (n by 3) or (n by 4)
if sin(theta)/lambda is chosen to be output
The algorithm follows the method described in:
Le Page and Gabe (1979) J. Appl. Cryst., 12, 464-466
Henning Osholm Sorensen, University of Copenhagen, July 22, 2010.
"""
from xfab import sg
if sgname != None:
spg = sg.sg(sgname=sgname,cell_choice=cell_choice)