source: trunk/GSASIIlattice.py @ 584

Last change on this file since 584 was 584, checked in by vondreele, 11 years ago

add coding line everywhere
more options in sample parm copy

  • Property svn:keywords set to Date Author Revision URL Id
File size: 55.0 KB
Line 
1# -*- coding: utf-8 -*-
2'''Perform lattice-related computations'''
3########### SVN repository information ###################
4# $Date: 2012-05-03 16:35:40 +0000 (Thu, 03 May 2012) $
5# $Author: vondreele $
6# $Revision: 584 $
7# $URL: trunk/GSASIIlattice.py $
8# $Id: GSASIIlattice.py 584 2012-05-03 16:35:40Z vondreele $
9########### SVN repository information ###################
10import math
11import numpy as np
12import numpy.linalg as nl
13
14# trig functions in degrees
15sind = lambda x: np.sin(x*np.pi/180.)
16asind = lambda x: 180.*np.arcsin(x)/np.pi
17tand = lambda x: np.tan(x*np.pi/180.)
18atand = lambda x: 180.*np.arctan(x)/np.pi
19atan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
20cosd = lambda x: np.cos(x*np.pi/180.)
21acosd = lambda x: 180.*np.arccos(x)/np.pi
22rdsq2d = lambda x,p: round(1.0/np.sqrt(x),p)
23
24def sec2HMS(sec):
25    """Convert time in sec to H:M:S string
26   
27    :param sec: time in seconds
28    return: H:M:S string (to nearest 100th second)
29   
30    """
31    H = int(sec/3600)
32    M = int(sec/60-H*60)
33    S = sec-3600*H-60*M
34    return '%d:%2d:%.2f'%(H,M,S)
35   
36def rotdMat(angle,axis=0):
37    """Prepare rotation matrix for angle in degrees about axis(=0,1,2)
38
39    :param angle: angle in degrees
40    :param axis:  axis (0,1,2 = x,y,z) about which for the rotation
41    :return: rotation matrix - 3x3 numpy array
42
43    """
44    if axis == 2:
45        return np.array([[cosd(angle),-sind(angle),0],[sind(angle),cosd(angle),0],[0,0,1]])
46    elif axis == 1:
47        return np.array([[cosd(angle),0,-sind(angle)],[0,1,0],[sind(angle),0,cosd(angle)]])
48    else:
49        return np.array([[1,0,0],[0,cosd(angle),-sind(angle)],[0,sind(angle),cosd(angle)]])
50       
51def rotdMat4(angle,axis=0):
52    """Prepare rotation matrix for angle in degrees about axis(=0,1,2) with scaling for OpenGL
53
54    :param angle: angle in degrees
55    :param axis:  axis (0,1,2 = x,y,z) about which for the rotation
56    :return: rotation matrix - 4x4 numpy array (last row/column for openGL scaling)
57
58    """
59    Mat = rotdMat(angle,axis)
60    return np.concatenate((np.concatenate((Mat,[[0],[0],[0]]),axis=1),[[0,0,0,1],]),axis=0)
61   
62def fillgmat(cell):
63    """Compute lattice metric tensor from unit cell constants
64
65    :param cell: tuple with a,b,c,alpha, beta, gamma (degrees)
66    :return: 3x3 numpy array
67
68    """
69    a,b,c,alp,bet,gam = cell
70    g = np.array([
71        [a*a,  a*b*cosd(gam),  a*c*cosd(bet)],
72        [a*b*cosd(gam),  b*b,  b*c*cosd(alp)],
73        [a*c*cosd(bet) ,b*c*cosd(alp),   c*c]])
74    return g
75           
76def cell2Gmat(cell):
77    """Compute real and reciprocal lattice metric tensor from unit cell constants
78
79    :param cell: tuple with a,b,c,alpha, beta, gamma (degrees)
80    :return: reciprocal (G) & real (g) metric tensors (list of two numpy 3x3 arrays)
81
82    """
83    g = fillgmat(cell)
84    G = nl.inv(g)       
85    return G,g
86
87def A2Gmat(A,inverse=True):
88    """Fill real & reciprocal metric tensor (G) from A
89
90    :param
91        A: reciprocal metric tensor elements as [G11,G22,G33,2*G12,2*G13,2*G23]
92        inverse: if True return bot G and g; else just G
93    :return: reciprocal (G) & real (g) metric tensors (list of two numpy 3x3 arrays)
94
95    """
96    G = np.zeros(shape=(3,3))
97    G = [
98        [A[0],  A[3]/2.,  A[4]/2.], 
99        [A[3]/2.,A[1],    A[5]/2.], 
100        [A[4]/2.,A[5]/2.,    A[2]]]
101    if inverse:
102        g = nl.inv(G)
103        return G,g
104    else:
105        return G
106
107def Gmat2A(G):
108    """Extract A from reciprocal metric tensor (G)
109
110    :param G: reciprocal maetric tensor (3x3 numpy array
111    :return: A = [G11,G22,G33,2*G12,2*G13,2*G23]
112
113    """
114    return [G[0][0],G[1][1],G[2][2],2.*G[0][1],2.*G[0][2],2.*G[1][2]]
115   
116def cell2A(cell):
117    """Obtain A = [G11,G22,G33,2*G12,2*G13,2*G23] from lattice parameters
118
119    :param cell: [a,b,c,alpha,beta,gamma] (degrees)
120    :return: G reciprocal metric tensor as 3x3 numpy array
121
122    """
123    G,g = cell2Gmat(cell)
124    return Gmat2A(G)
125
126def A2cell(A):
127    """Compute unit cell constants from A
128
129    :param A: [G11,G22,G33,2*G12,2*G13,2*G23] G - reciprocal metric tensor
130    :return: a,b,c,alpha, beta, gamma (degrees) - lattice parameters
131
132    """
133    G,g = A2Gmat(A)
134    return Gmat2cell(g)
135
136def Gmat2cell(g):
137    """Compute real/reciprocal lattice parameters from real/reciprocal metric tensor (g/G)
138    The math works the same either way.
139
140    :param g (or G): real (or reciprocal) metric tensor 3x3 array
141    :return: a,b,c,alpha, beta, gamma (degrees) (or a*,b*,c*,alpha*,beta*,gamma* degrees)
142
143    """
144    oldset = np.seterr('raise')
145    a = np.sqrt(max(0,g[0][0]))
146    b = np.sqrt(max(0,g[1][1]))
147    c = np.sqrt(max(0,g[2][2]))
148    alp = acosd(g[2][1]/(b*c))
149    bet = acosd(g[2][0]/(a*c))
150    gam = acosd(g[0][1]/(a*b))
151    np.seterr(**oldset)
152    return a,b,c,alp,bet,gam
153
154def invcell2Gmat(invcell):
155    """Compute real and reciprocal lattice metric tensor from reciprocal
156       unit cell constants
157       
158    :param invcell: [a*,b*,c*,alpha*, beta*, gamma*] (degrees)
159    :return: reciprocal (G) & real (g) metric tensors (list of two 3x3 arrays)
160
161    """
162    G = fillgmat(invcell)
163    g = nl.inv(G)
164    return G,g
165       
166def calc_rVsq(A):
167    """Compute the square of the reciprocal lattice volume (V* **2) from A'
168
169    """
170    G,g = A2Gmat(A)
171    rVsq = nl.det(G)
172    if rVsq < 0:
173        return 1
174    return rVsq
175   
176def calc_rV(A):
177    """Compute the reciprocal lattice volume (V*) from A
178    """
179    return np.sqrt(calc_rVsq(A))
180   
181def calc_V(A):
182    """Compute the real lattice volume (V) from A
183    """
184    return 1./calc_rV(A)
185
186def A2invcell(A):
187    """Compute reciprocal unit cell constants from A
188    returns tuple with a*,b*,c*,alpha*, beta*, gamma* (degrees)
189    """
190    G,g = A2Gmat(A)
191    return Gmat2cell(G)
192   
193def Gmat2AB(G):
194    """Computes orthogonalization matrix from reciprocal metric tensor G
195    returns tuple of two 3x3 numpy arrays (A,B)
196       A for crystal to Cartesian transformations A*x = np.inner(A,x) = X
197       B (= inverse of A) for Cartesian to crystal transformation B*X = np.inner(B,X) = x
198    """
199    cellstar = Gmat2cell(G)
200    g = nl.inv(G)
201    cell = Gmat2cell(g)
202    A = np.zeros(shape=(3,3))
203    # from Giacovazzo (Fundamentals 2nd Ed.) p.75
204    A[0][0] = cell[0]                # a
205    A[0][1] = cell[1]*cosd(cell[5])  # b cos(gamma)
206    A[0][2] = cell[2]*cosd(cell[4])  # c cos(beta)
207    A[1][1] = cell[1]*sind(cell[5])  # b sin(gamma)
208    A[1][2] = -cell[2]*cosd(cellstar[3])*sind(cell[4]) # - c cos(alpha*) sin(beta)
209    A[2][2] = 1/cellstar[2]         # 1/c*
210    B = nl.inv(A)
211    return A,B
212   
213
214def cell2AB(cell):
215    """Computes orthogonalization matrix from unit cell constants
216    cell is tuple with a,b,c,alpha, beta, gamma (degrees)
217    returns tuple of two 3x3 numpy arrays (A,B)
218       A for crystal to Cartesian transformations A*x = np.inner(A,x) = X
219       B (= inverse of A) for Cartesian to crystal transformation B*X = np.inner(B,X) = x
220    """
221    G,g = cell2Gmat(cell) 
222    cellstar = Gmat2cell(G)
223    A = np.zeros(shape=(3,3))
224    # from Giacovazzo (Fundamentals 2nd Ed.) p.75
225    A[0][0] = cell[0]                # a
226    A[0][1] = cell[1]*cosd(cell[5])  # b cos(gamma)
227    A[0][2] = cell[2]*cosd(cell[4])  # c cos(beta)
228    A[1][1] = cell[1]*sind(cell[5])  # b sin(gamma)
229    A[1][2] = -cell[2]*cosd(cellstar[3])*sind(cell[4]) # - c cos(alpha*) sin(beta)
230    A[2][2] = 1/cellstar[2]         # 1/c*
231    B = nl.inv(A)
232    return A,B
233   
234def U6toUij(U6):
235    """Fill matrix (Uij) from U6 = [U11,U22,U33,U12,U13,U23]
236    NB: there is a non numpy version in GSASIIspc: U2Uij
237    input:
238        U6 - 6 terms of u11,u22,...
239    returns:
240        Uij - numpy [3][3] array of uij
241    """
242    U = np.zeros(shape=(3,3))
243    U = [
244        [U6[0],  U6[3],  U6[4]], 
245        [U6[3],  U6[1],  U6[5]], 
246        [U6[4],  U6[5],  U6[2]]]
247    return U
248
249def UijtoU6(U):
250    """Fill vector [U11,U22,U33,U12,U13,U23] from Uij
251    NB: there is a non numpy version in GSASIIspc: Uij2U
252    """
253    U6 = np.zeros(6)
254    U6 = [U[0][0],U[1][1],U[2][2],U[0][1],U[0][2],U[1][2]]
255    return U6
256
257def Uij2betaij(Uij,G):
258    """
259    Convert Uij to beta-ij tensors
260    input:
261        Uij - numpy array [Uij]
262        G - reciprocal metric tensor
263    returns:
264    beta-ij - numpy array [beta-ij]
265    """
266    pass
267   
268def CosSinAngle(U,V,G):
269    """ calculate sin & cos of angle betwee U & V in generalized coordinates
270    defined by metric tensor G
271    input:
272        U & V - 3-vectors assume numpy arrays
273        G - metric tensor for U & V defined space assume numpy array
274    return:
275        cos(phi) & sin(phi)
276    """
277    u = U/np.sqrt(np.inner(U,np.inner(G,U)))
278    v = V/np.sqrt(np.inner(V,np.inner(G,V)))
279    cosP = np.inner(u,np.inner(G,v))
280    sinP = np.sqrt(max(0.0,1.0-cosP**2))
281    return cosP,sinP
282   
283def criticalEllipse(prob):
284    """
285    Calculate critical values for probability ellipsoids from probability
286    """
287    if not ( 0.01 <= prob < 1.0):
288        return 1.54 
289    coeff = np.array([6.44988E-09,4.16479E-07,1.11172E-05,1.58767E-04,0.00130554,
290        0.00604091,0.0114921,-0.040301,-0.6337203,1.311582])
291    llpr = math.log(-math.log(prob))
292    return np.polyval(coeff,llpr)
293   
294def CellBlock(nCells):
295    """
296    Generate block of unit cells n*n*n on a side; [0,0,0] centered, n = 2*nCells+1
297    currently only works for nCells = 0 or 1 (not >1)
298    """
299    if nCells:
300        N = 2*nCells+1
301        N2 = N*N
302        N3 = N*N*N
303        cellArray = []
304        A = np.array(range(N3))
305        cellGen = np.array([A/N2-1,A/N%N-1,A%N-1]).T
306        for cell in cellGen:
307            cellArray.append(cell)
308        return cellArray
309    else:
310        return [0,0,0]
311       
312def CellAbsorption(ElList,Volume):
313# ElList = dictionary of element contents including mu
314    muT = 0
315    for El in ElList:
316        muT += ElList[El]['mu']*ElList[El]['FormulaNo']
317    return muT/Volume
318   
319#Permutations and Combinations
320# Four routines: combinations,uniqueCombinations, selections & permutations
321#These taken from Python Cookbook, 2nd Edition. 19.15 p724-726
322#   
323def _combinators(_handle, items, n):
324    """ factored-out common structure of all following combinators """
325    if n==0:
326        yield [ ]
327        return
328    for i, item in enumerate(items):
329        this_one = [ item ]
330        for cc in _combinators(_handle, _handle(items, i), n-1):
331            yield this_one + cc
332def combinations(items, n):
333    """ take n distinct items, order matters """
334    def skipIthItem(items, i):
335        return items[:i] + items[i+1:]
336    return _combinators(skipIthItem, items, n)
337def uniqueCombinations(items, n):
338    """ take n distinct items, order is irrelevant """
339    def afterIthItem(items, i):
340        return items[i+1:]
341    return _combinators(afterIthItem, items, n)
342def selections(items, n):
343    """ take n (not necessarily distinct) items, order matters """
344    def keepAllItems(items, i):
345        return items
346    return _combinators(keepAllItems, items, n)
347def permutations(items):
348    """ take all items, order matters """
349    return combinations(items, len(items))
350
351#reflection generation routines
352#for these: H = [h,k,l]; A is as used in calc_rDsq; G - inv metric tensor, g - metric tensor;
353#           cell - a,b,c,alp,bet,gam in A & deg
354   
355def calc_rDsq(H,A):
356    rdsq = H[0]*H[0]*A[0]+H[1]*H[1]*A[1]+H[2]*H[2]*A[2]+H[0]*H[1]*A[3]+H[0]*H[2]*A[4]+H[1]*H[2]*A[5]
357    return rdsq
358   
359def calc_rDsq2(H,G):
360    return np.inner(H,np.inner(G,H))
361   
362def calc_rDsqZ(H,A,Z,tth,lam):
363    rpd = np.pi/180.
364    rdsq = calc_rDsq(H,A)+Z*sind(tth)*2.0*rpd/lam**2
365    return rdsq
366       
367def MaxIndex(dmin,A):
368    Hmax = [0,0,0]
369    try:
370        cell = A2cell(A)
371    except:
372        cell = [1,1,1,90,90,90]
373    for i in range(3):
374        Hmax[i] = int(round(cell[i]/dmin))
375    return Hmax
376   
377def sortHKLd(HKLd,ifreverse,ifdup):
378    #HKLd is a list of [h,k,l,d,...]; ifreverse=True for largest d first
379    #ifdup = True if duplicate d-spacings allowed
380    T = []
381    for i,H in enumerate(HKLd):
382        if ifdup:
383            T.append((H[3],i))
384        else:
385            T.append(H[3])           
386    D = dict(zip(T,HKLd))
387    T.sort()
388    if ifreverse:
389        T.reverse()
390    X = []
391    okey = ''
392    for key in T: 
393        if key != okey: X.append(D[key])    #remove duplicate d-spacings
394        okey = key
395    return X
396   
397def SwapIndx(Axis,H):
398    if Axis in [1,-1]:
399        return H
400    elif Axis in [2,-3]:
401        return [H[1],H[2],H[0]]
402    else:
403        return [H[2],H[0],H[1]]
404       
405def Rh2Hx(Rh):
406    Hx = [0,0,0]
407    Hx[0] = Rh[0]-Rh[1]
408    Hx[1] = Rh[1]-Rh[2]
409    Hx[2] = np.sum(Rh)
410    return Hx
411   
412def Hx2Rh(Hx):
413        Rh = [0,0,0]
414        itk = -Hx[0]+Hx[1]+Hx[2]
415        if itk%3 != 0:
416            return 0        #error - not rhombohedral reflection
417        else:
418            Rh[1] = itk/3
419            Rh[0] = Rh[1]+Hx[0]
420            Rh[2] = Rh[1]-Hx[1]
421            if Rh[0] < 0:
422                for i in range(3):
423                    Rh[i] = -Rh[i]
424            return Rh
425       
426def CentCheck(Cent,H):
427    h,k,l = H
428    if Cent == 'A' and (k+l)%2:
429        return False
430    elif Cent == 'B' and (h+l)%2:
431        return False
432    elif Cent == 'C' and (h+k)%2:
433        return False
434    elif Cent == 'I' and (h+k+l)%2:
435        return False
436    elif Cent == 'F' and ((h+k)%2 or (h+l)%2 or (k+l)%2):
437        return False
438    elif Cent == 'R' and (-h+k+l)%3:
439        return False
440    else:
441        return True
442                                   
443def GetBraviasNum(center,system):
444    """Determine the Bravais lattice number, as used in GenHBravais
445   
446    :param center: one of: 'P', 'C', 'I', 'F', 'R' (see SGLatt from GSASIIspc.SpcGroup)
447    :param system: one of 'cubic', 'hexagonal', 'tetragonal', 'orthorhombic', 'trigonal' (for R)
448             'monoclinic', 'triclinic' (see SGSys from GSASIIspc.SpcGroup)
449    :return: a number between 0 and 13
450          or throws a ValueError exception if the combination of center, system is not found (i.e. non-standard)
451    """
452    if center.upper() == 'F' and system.lower() == 'cubic':
453        return 0
454    elif center.upper() == 'I' and system.lower() == 'cubic':
455        return 1
456    elif center.upper() == 'P' and system.lower() == 'cubic':
457        return 2
458    elif center.upper() == 'R' and system.lower() == 'trigonal':
459        return 3
460    elif center.upper() == 'P' and system.lower() == 'hexagonal':
461        return 4
462    elif center.upper() == 'I' and system.lower() == 'tetragonal':
463        return 5
464    elif center.upper() == 'P' and system.lower() == 'tetragonal':
465        return 6
466    elif center.upper() == 'F' and system.lower() == 'orthorhombic':
467        return 7
468    elif center.upper() == 'I' and system.lower() == 'orthorhombic':
469        return 8
470    elif center.upper() == 'C' and system.lower() == 'orthorhombic':
471        return 9
472    elif center.upper() == 'P' and system.lower() == 'orthorhombic':
473        return 10
474    elif center.upper() == 'C' and system.lower() == 'monoclinic':
475        return 11
476    elif center.upper() == 'P' and system.lower() == 'monoclinic':
477        return 12
478    elif center.upper() == 'P' and system.lower() == 'triclinic':
479        return 13
480    raise ValueError,'non-standard Bravais lattice center=%s, cell=%s' % (center,system)
481
482def GenHBravais(dmin,Bravais,A):
483    """Generate the positionally unique powder diffraction reflections
484     
485    :param dmin: minimum d-spacing in A
486    :param Bravais: lattice type (see GetBraviasNum)
487        Bravais is one of::
488             0 F cubic
489             1 I cubic
490             2 P cubic
491             3 R hexagonal (trigonal not rhombohedral)
492             4 P hexagonal
493             5 I tetragonal
494             6 P tetragonal
495             7 F orthorhombic
496             8 I orthorhombic
497             9 C orthorhombic
498            10 P orthorhombic
499            11 C monoclinic
500            12 P monoclinic
501            13 P triclinic
502    :param A: reciprocal metric tensor elements as [G11,G22,G33,2*G12,2*G13,2*G23]
503    :return: HKL unique d list of [h,k,l,d,-1] sorted with largest d first
504           
505    """
506    import math
507    if Bravais in [9,11]:
508        Cent = 'C'
509    elif Bravais in [1,5,8]:
510        Cent = 'I'
511    elif Bravais in [0,7]:
512        Cent = 'F'
513    elif Bravais in [3]:
514        Cent = 'R'
515    else:
516        Cent = 'P'
517    Hmax = MaxIndex(dmin,A)
518    dminsq = 1./(dmin**2)
519    HKL = []
520    if Bravais == 13:                       #triclinic
521        for l in range(-Hmax[2],Hmax[2]+1):
522            for k in range(-Hmax[1],Hmax[1]+1):
523                hmin = 0
524                if (k < 0): hmin = 1
525                if (k ==0 and l < 0): hmin = 1
526                for h in range(hmin,Hmax[0]+1):
527                    H=[h,k,l]
528                    rdsq = calc_rDsq(H,A)
529                    if 0 < rdsq <= dminsq:
530                        HKL.append([h,k,l,rdsq2d(rdsq,6),-1])
531    elif Bravais in [11,12]:                #monoclinic - b unique
532        Hmax = SwapIndx(2,Hmax)
533        for h in range(Hmax[0]+1):
534            for k in range(-Hmax[1],Hmax[1]+1):
535                lmin = 0
536                if k < 0:lmin = 1
537                for l in range(lmin,Hmax[2]+1):
538                    [h,k,l] = SwapIndx(-2,[h,k,l])
539                    H = []
540                    if CentCheck(Cent,[h,k,l]): H=[h,k,l]
541                    if H:
542                        rdsq = calc_rDsq(H,A)
543                        if 0 < rdsq <= dminsq:
544                            HKL.append([h,k,l,rdsq2d(rdsq,6),-1])
545                    [h,k,l] = SwapIndx(2,[h,k,l])
546    elif Bravais in [7,8,9,10]:            #orthorhombic
547        for h in range(Hmax[0]+1):
548            for k in range(Hmax[1]+1):
549                for l in range(Hmax[2]+1):
550                    H = []
551                    if CentCheck(Cent,[h,k,l]): H=[h,k,l]
552                    if H:
553                        rdsq = calc_rDsq(H,A)
554                        if 0 < rdsq <= dminsq:
555                            HKL.append([h,k,l,rdsq2d(rdsq,6),-1])
556    elif Bravais in [5,6]:                  #tetragonal
557        for l in range(Hmax[2]+1):
558            for k in range(Hmax[1]+1):
559                for h in range(k,Hmax[0]+1):
560                    H = []
561                    if CentCheck(Cent,[h,k,l]): H=[h,k,l]
562                    if H:
563                        rdsq = calc_rDsq(H,A)
564                        if 0 < rdsq <= dminsq:
565                            HKL.append([h,k,l,rdsq2d(rdsq,6),-1])
566    elif Bravais in [3,4]:
567        lmin = 0
568        if Bravais == 3: lmin = -Hmax[2]                  #hexagonal/trigonal
569        for l in range(lmin,Hmax[2]+1):
570            for k in range(Hmax[1]+1):
571                hmin = k
572                if l < 0: hmin += 1
573                for h in range(hmin,Hmax[0]+1):
574                    H = []
575                    if CentCheck(Cent,[h,k,l]): H=[h,k,l]
576                    if H:
577                        rdsq = calc_rDsq(H,A)
578                        if 0 < rdsq <= dminsq:
579                            HKL.append([h,k,l,rdsq2d(rdsq,6),-1])
580
581    else:                                   #cubic
582        for l in range(Hmax[2]+1):
583            for k in range(l,Hmax[1]+1):
584                for h in range(k,Hmax[0]+1):
585                    H = []
586                    if CentCheck(Cent,[h,k,l]): H=[h,k,l]
587                    if H:
588                        rdsq = calc_rDsq(H,A)
589                        if 0 < rdsq <= dminsq:
590                            HKL.append([h,k,l,rdsq2d(rdsq,6),-1])
591    return sortHKLd(HKL,True,False)
592
593def getHKLmax(dmin,SGData,A):
594    #finds maximum allowed hkl for given A within dmin
595    SGLaue = SGData['SGLaue']
596    if SGLaue in ['3R','3mR']:        #Rhombohedral axes
597        Hmax = [0,0,0]
598        cell = A2cell(A)
599        aHx = cell[0]*math.sqrt(2.0*(1.0-cosd(cell[3])))
600        cHx = cell[0]*math.sqrt(3.0*(1.0+2.0*cosd(cell[3])))
601        Hmax[0] = Hmax[1] = int(round(aHx/dmin))
602        Hmax[2] = int(round(cHx/dmin))
603        #print Hmax,aHx,cHx
604    else:                           # all others
605        Hmax = MaxIndex(dmin,A)
606    return Hmax
607   
608def GenHLaue(dmin,SGData,A):
609    """Generate the crystallographically unique powder diffraction reflections
610    for a lattice and Bravais type
611   
612    :param dmin: minimum d-spacing
613    :param SGData: space group dictionary with at least::
614   
615        'SGLaue': Laue group symbol: one of '-1','2/m','mmm','4/m','6/m','4/mmm','6/mmm',
616                 '3m1', '31m', '3', '3R', '3mR', 'm3', 'm3m'
617        'SGLatt': lattice centering: one of 'P','A','B','C','I','F'
618        'SGUniq': code for unique monoclinic axis one of 'a','b','c' (only if 'SGLaue' is '2/m')
619            otherwise ' '
620       
621    :param A: reciprocal metric tensor elements as [G11,G22,G33,2*G12,2*G13,2*G23]
622    :return: HKL = list of [h,k,l,d] sorted with largest d first and is unique
623            part of reciprocal space ignoring anomalous dispersion
624           
625    """
626    import math
627    SGLaue = SGData['SGLaue']
628    SGLatt = SGData['SGLatt']
629    SGUniq = SGData['SGUniq']
630    #finds maximum allowed hkl for given A within dmin
631    Hmax = getHKLmax(dmin,SGData,A)
632       
633    dminsq = 1./(dmin**2)
634    HKL = []
635    if SGLaue == '-1':                       #triclinic
636        for l in range(-Hmax[2],Hmax[2]+1):
637            for k in range(-Hmax[1],Hmax[1]+1):
638                hmin = 0
639                if (k < 0) or (k ==0 and l < 0): hmin = 1
640                for h in range(hmin,Hmax[0]+1):
641                    H = []
642                    if CentCheck(SGLatt,[h,k,l]): H=[h,k,l]
643                    rdsq = calc_rDsq(H,A)
644                    if 0 < rdsq <= dminsq:
645                        HKL.append([h,k,l,1/math.sqrt(rdsq)])
646    elif SGLaue == '2/m':                #monoclinic
647        axisnum = 1 + ['a','b','c'].index(SGUniq)
648        Hmax = SwapIndx(axisnum,Hmax)
649        for h in range(Hmax[0]+1):
650            for k in range(-Hmax[1],Hmax[1]+1):
651                lmin = 0
652                if k < 0:lmin = 1
653                for l in range(lmin,Hmax[2]+1):
654                    [h,k,l] = SwapIndx(-axisnum,[h,k,l])
655                    H = []
656                    if CentCheck(SGLatt,[h,k,l]): H=[h,k,l]
657                    if H:
658                        rdsq = calc_rDsq(H,A)
659                        if 0 < rdsq <= dminsq:
660                            HKL.append([h,k,l,1/math.sqrt(rdsq)])
661                    [h,k,l] = SwapIndx(axisnum,[h,k,l])
662    elif SGLaue in ['mmm','4/m','6/m']:            #orthorhombic
663        for l in range(Hmax[2]+1):
664            for h in range(Hmax[0]+1):
665                kmin = 1
666                if SGLaue == 'mmm' or h ==0: kmin = 0
667                for k in range(kmin,Hmax[1]+1):
668                    H = []
669                    if CentCheck(SGLatt,[h,k,l]): H=[h,k,l]
670                    if H:
671                        rdsq = calc_rDsq(H,A)
672                        if 0 < rdsq <= dminsq:
673                            HKL.append([h,k,l,1/math.sqrt(rdsq)])
674    elif SGLaue in ['4/mmm','6/mmm']:                  #tetragonal & hexagonal
675        for l in range(Hmax[2]+1):
676            for h in range(Hmax[0]+1):
677                for k in range(h+1):
678                    H = []
679                    if CentCheck(SGLatt,[h,k,l]): H=[h,k,l]
680                    if H:
681                        rdsq = calc_rDsq(H,A)
682                        if 0 < rdsq <= dminsq:
683                            HKL.append([h,k,l,1/math.sqrt(rdsq)])
684    elif SGLaue in ['3m1','31m','3','3R','3mR']:                  #trigonals
685        for l in range(-Hmax[2],Hmax[2]+1):
686            hmin = 0
687            if l < 0: hmin = 1
688            for h in range(hmin,Hmax[0]+1):
689                if SGLaue in ['3R','3']:
690                    kmax = h
691                    kmin = -int((h-1.)/2.)
692                else:
693                    kmin = 0
694                    kmax = h
695                    if SGLaue in ['3m1','3mR'] and l < 0: kmax = h-1
696                    if SGLaue == '31m' and l < 0: kmin = 1
697                for k in range(kmin,kmax+1):
698                    H = []
699                    if CentCheck(SGLatt,[h,k,l]): H=[h,k,l]
700                    if SGLaue in ['3R','3mR']:
701                        H = Hx2Rh(H)
702                    if H:
703                        rdsq = calc_rDsq(H,A)
704                        if 0 < rdsq <= dminsq:
705                            HKL.append([H[0],H[1],H[2],1/math.sqrt(rdsq)])
706    else:                                   #cubic
707        for h in range(Hmax[0]+1):
708            for k in range(h+1):
709                lmin = 0
710                lmax = k
711                if SGLaue =='m3':
712                    lmax = h-1
713                    if h == k: lmax += 1
714                for l in range(lmin,lmax+1):
715                    H = []
716                    if CentCheck(SGLatt,[h,k,l]): H=[h,k,l]
717                    if H:
718                        rdsq = calc_rDsq(H,A)
719                        if 0 < rdsq <= dminsq:
720                            HKL.append([h,k,l,1/math.sqrt(rdsq)])
721    return sortHKLd(HKL,True,True)
722
723#Spherical harmonics routines
724def OdfChk(SGLaue,L,M):
725    if not L%2 and abs(M) <= L:
726        if SGLaue == '0':                      #cylindrical symmetry
727            if M == 0: return True
728        elif SGLaue == '-1':
729            return True
730        elif SGLaue == '2/m':
731            if not abs(M)%2: return True
732        elif SGLaue == 'mmm':
733            if not abs(M)%2 and M >= 0: return True
734        elif SGLaue == '4/m':
735            if not abs(M)%4: return True
736        elif SGLaue == '4/mmm':
737            if not abs(M)%4 and M >= 0: return True
738        elif SGLaue in ['3R','3']:
739            if not abs(M)%3: return True
740        elif SGLaue in ['3mR','3m1','31m']:
741            if not abs(M)%3 and M >= 0: return True
742        elif SGLaue == '6/m':
743            if not abs(M)%6: return True
744        elif SGLaue == '6/mmm':
745            if not abs(M)%6 and M >= 0: return True
746        elif SGLaue == 'm3':
747            if M > 0:
748                if L%12 == 2:
749                    if M <= L/12: return True
750                else:
751                    if M <= L/12+1: return True
752        elif SGLaue == 'm3m':
753            if M > 0:
754                if L%12 == 2:
755                    if M <= L/12: return True
756                else:
757                    if M <= L/12+1: return True
758    return False
759       
760def GenSHCoeff(SGLaue,SamSym,L,IfLMN=True):
761    coeffNames = []
762    for iord in [2*i+2 for i in range(L/2)]:
763        for m in [i-iord for i in range(2*iord+1)]:
764            if OdfChk(SamSym,iord,m):
765                for n in [i-iord for i in range(2*iord+1)]:
766                    if OdfChk(SGLaue,iord,n):
767                        if IfLMN:
768                            coeffNames.append('C(%d,%d,%d)'%(iord,m,n))
769                        else:
770                            coeffNames.append('C(%d,%d)'%(iord,n))
771    return coeffNames
772
773def CrsAng(H,cell,SGData):
774    a,b,c,al,be,ga = cell
775    SQ3 = 1.732050807569
776    H1 = np.array([1,0,0])
777    H2 = np.array([0,1,0])
778    H3 = np.array([0,0,1])
779    H4 = np.array([1,1,1])
780    G,g = cell2Gmat(cell)
781    Laue = SGData['SGLaue']
782    Naxis = SGData['SGUniq']
783    DH = np.inner(H,np.inner(G,H))
784    if Laue == '2/m':
785        if Naxis == 'a':
786            DR = np.inner(H1,np.inner(G,H1))
787            DHR = np.inner(H,np.inner(G,H1))
788        elif Naxis == 'b':
789            DR = np.inner(H2,np.inner(G,H2))
790            DHR = np.inner(H,np.inner(G,H2))
791        else:
792            DR = np.inner(H3,np.inner(G,H3))
793            DHR = np.inner(H,np.inner(G,H3))
794    elif Laue in ['R3','R3m']:
795        DR = np.inner(H4,np.inner(G,H4))
796        DHR = np.inner(H,np.inner(G,H4))
797       
798    else:
799        DR = np.inner(H3,np.inner(G,H3))
800        DHR = np.inner(H,np.inner(G,H3))
801    DHR /= np.sqrt(DR*DH)
802    phi = acosd(DHR)
803    if Laue == '-1':
804        BA = H[1]*a/(b-H[0]*cosd(ga))
805        BB = H[0]*sind(ga)**2
806    elif Laue == '2/m':
807        if Naxis == 'a':
808            BA = H[2]*b/(c-H[1]*cosd(al))
809            BB = H[1]*sind(al)**2
810        elif Naxis == 'b':
811            BA = H[0]*c/(a-H[2]*cosd(be))
812            BB = H[2]*sind(be)**2
813        else:
814            BA = H[1]*a/(b-H[0]*cosd(ga))
815            BB = H[0]*sind(ga)**2
816    elif Laue in ['mmm','4/m','4/mmm']:
817        BA = H[1]*a
818        BB = H[0]*b
819   
820    elif Laue in ['3R','3mR']:
821        BA = H[0]+H[1]-2.0*H[2]
822        BB = SQ3*(H[0]-H[1])
823    elif Laue in ['m3','m3m']:
824        BA = H[1]
825        BB = H[0]
826    else:
827        BA = H[0]+2.0*H[1]
828        BB = SQ3*H[0]
829    beta = atan2d(BA,BB)
830    return phi,beta
831   
832def SamAng(Tth,Gangls,Sangl,IFCoup):
833    """Compute sample orientation angles vs laboratory coord. system
834    input:
835        Tth:        Signed theta                                   
836        Gangls:     Sample goniometer angles phi,chi,omega,azmuth 
837        Sangl:      Sample angle zeros om-0, chi-0, phi-0         
838        IFCoup:     =.TRUE. if omega & 2-theta coupled in CW scan
839    returns: 
840        psi,gam:    Sample odf angles                             
841        dPSdA,dGMdA:    Angle zero derivatives
842    """                         
843   
844    rpd = math.pi/180.
845    if IFCoup:
846        GSomeg = sind(Gangls[2]+Tth)
847        GComeg = cosd(Gangls[2]+Tth)
848    else:
849        GSomeg = sind(Gangls[2])
850        GComeg = cosd(Gangls[2])
851    GSTth = sind(Tth)
852    GCTth = cosd(Tth)     
853    GSazm = sind(Gangls[3])
854    GCazm = cosd(Gangls[3])
855    GSchi = sind(Gangls[1])
856    GCchi = cosd(Gangls[1])
857    GSphi = sind(Gangls[0]+Sangl[2])
858    GCphi = cosd(Gangls[0]+Sangl[2])
859    SSomeg = sind(Sangl[0])
860    SComeg = cosd(Sangl[0])
861    SSchi = sind(Sangl[1])
862    SCchi = cosd(Sangl[1])
863    AT = -GSTth*GComeg+GCTth*GCazm*GSomeg
864    BT = GSTth*GSomeg+GCTth*GCazm*GComeg
865    CT = -GCTth*GSazm*GSchi
866    DT = -GCTth*GSazm*GCchi
867   
868    BC1 = -AT*GSphi+(CT+BT*GCchi)*GCphi
869    BC2 = DT-BT*GSchi
870    BC3 = AT*GCphi+(CT+BT*GCchi)*GSphi
871     
872    BC = BC1*SComeg*SCchi+BC2*SComeg*SSchi-BC3*SSomeg     
873    psi = acosd(BC)
874   
875    BD = 1.0-BC**2
876    if BD > 0.:
877        C = rpd/math.sqrt(BD)
878    else:
879        C = 0.
880    dPSdA = [-C*(-BC1*SSomeg*SCchi-BC2*SSomeg*SSchi-BC3*SComeg),
881        -C*(-BC1*SComeg*SSchi+BC2*SComeg*SCchi),
882        -C*(-BC1*SSomeg-BC3*SComeg*SCchi)]
883     
884    BA = -BC1*SSchi+BC2*SCchi
885    BB = BC1*SSomeg*SCchi+BC2*SSomeg*SSchi+BC3*SComeg
886    gam = atan2d(BB,BA)
887
888    BD = (BA**2+BB**2)/rpd
889
890    dBAdO = 0
891    dBAdC = -BC1*SCchi-BC2*SSchi
892    dBAdF = BC3*SSchi
893   
894    dBBdO = BC1*SComeg*SCchi+BC2*SComeg*SSchi-BC3*SSomeg
895    dBBdC = -BC1*SSomeg*SSchi+BC2*SSomeg*SCchi
896    dBBdF = BC1*SComeg-BC3*SSomeg*SCchi
897   
898    if BD > 0.:
899        dGMdA = [(BA*dBBdO-BB*dBAdO)/BD,(BA*dBBdC-BB*dBAdC)/BD,(BA*dBBdF-BB*dBAdF)/BD]
900    else:
901        dGMdA = [0.0,0.0,0.0]
902
903       
904    return psi,gam,dPSdA,dGMdA
905
906BOH = {
907'L=2':[[],[],[]],
908'L=4':[[0.30469720,0.36418281],[],[]],
909'L=6':[[-0.14104740,0.52775103],[],[]],
910'L=8':[[0.28646862,0.21545346,0.32826995],[],[]],
911'L=10':[[-0.16413497,0.33078546,0.39371345],[],[]],
912'L=12':[[0.26141975,0.27266871,0.03277460,0.32589402],
913    [0.09298802,-0.23773812,0.49446631,0.0],[]],
914'L=14':[[-0.17557309,0.25821932,0.27709173,0.33645360],[],[]],
915'L=16':[[0.24370673,0.29873515,0.06447688,0.00377,0.32574495],
916    [0.12039646,-0.25330128,0.23950998,0.40962508,0.0],[]],
917'L=18':[[-0.16914245,0.17017340,0.34598142,0.07433932,0.32696037],
918    [-0.06901768,0.16006562,-0.24743528,0.47110273,0.0],[]],
919'L=20':[[0.23067026,0.31151832,0.09287682,0.01089683,0.00037564,0.32573563],
920    [0.13615420,-0.25048007,0.12882081,0.28642879,0.34620433,0.0],[]],
921'L=22':[[-0.16109560,0.10244188,0.36285175,0.13377513,0.01314399,0.32585583],
922    [-0.09620055,0.20244115,-0.22389483,0.17928946,0.42017231,0.0],[]],
923'L=24':[[0.22050742,0.31770654,0.11661736,0.02049853,0.00150861,0.00003426,0.32573505],
924    [0.13651722,-0.21386648,0.00522051,0.33939435,0.10837396,0.32914497,0.0],
925    [0.05378596,-0.11945819,0.16272298,-0.26449730,0.44923956,0.0,0.0]],
926'L=26':[[-0.15435003,0.05261630,0.35524646,0.18578869,0.03259103,0.00186197,0.32574594],
927    [-0.11306511,0.22072681,-0.18706142,0.05439948,0.28122966,0.35634355,0.0],[]],
928'L=28':[[0.21225019,0.32031716,0.13604702,0.03132468,0.00362703,0.00018294,0.00000294,0.32573501],
929    [0.13219496,-0.17206256,-0.08742608,0.32671661,0.17973107,0.02567515,0.32619598,0.0],
930    [0.07989184,-0.16735346,0.18839770,-0.20705337,0.12926808,0.42715602,0.0,0.0]],
931'L=30':[[-0.14878368,0.01524973,0.33628434,0.22632587,0.05790047,0.00609812,0.00022898,0.32573594],
932    [-0.11721726,0.20915005,-0.11723436,-0.07815329,0.31318947,0.13655742,0.33241385,0.0],
933    [-0.04297703,0.09317876,-0.11831248,0.17355132,-0.28164031,0.42719361,0.0,0.0]],
934'L=32':[[0.20533892,0.32087437,0.15187897,0.04249238,0.00670516,0.00054977,0.00002018,0.00000024,0.32573501],
935    [0.12775091,-0.13523423,-0.14935701,0.28227378,0.23670434,0.05661270,0.00469819,0.32578978,0.0],
936    [0.09703829,-0.19373733,0.18610682,-0.14407046,0.00220535,0.26897090,0.36633402,0.0,0.0]],
937'L=34':[[-0.14409234,-0.01343681,0.31248977,0.25557722,0.08571889,0.01351208,0.00095792,0.00002550,0.32573508],
938    [-0.11527834,0.18472133,-0.04403280,-0.16908618,0.27227021,0.21086614,0.04041752,0.32688152,0.0],
939    [-0.06773139,0.14120811,-0.15835721,0.18357456,-0.19364673,0.08377174,0.43116318,0.0,0.0]]
940}
941
942Lnorm = lambda L: 4.*np.pi/(2.0*L+1.)
943
944def GetKcl(L,N,SGLaue,phi,beta):
945    import pytexture as ptx
946    RSQ2PI = 0.3989422804014
947    SQ2 = 1.414213562373
948    if SGLaue in ['m3','m3m']:
949        Kcl = 0.0
950        for j in range(0,L+1,4):
951            im = j/4+1
952            pcrs,dum = ptx.pyplmpsi(L,j,1,phi)
953            Kcl += BOH['L='+str(L)][N-1][im-1]*pcrs*cosd(j*beta)       
954    else:
955        pcrs,dum = ptx.pyplmpsi(L,N,1,phi)
956        pcrs *= RSQ2PI
957        if N:
958            pcrs *= SQ2
959        if SGLaue in ['mmm','4/mmm','6/mmm','R3mR','3m1','31m']:
960            if SGLaue in ['3mR','3m1','31m']: 
961                if N%6 == 3:
962                    Kcl = pcrs*sind(N*beta)
963                else:
964                    Kcl = pcrs*cosd(N*beta)
965            else:
966                Kcl = pcrs*cosd(N*beta)
967        else:
968            Kcl = pcrs*(cosd(N*beta)+sind(N*beta))
969    return Kcl
970   
971def GetKsl(L,M,SamSym,psi,gam):
972    import pytexture as ptx
973    RSQPI = 0.5641895835478
974    SQ2 = 1.414213562373
975    psrs,dpdps = ptx.pyplmpsi(L,M,1,psi)
976    psrs *= RSQPI
977    dpdps *= RSQPI
978    if M == 0:
979        psrs /= SQ2
980        dpdps /= SQ2
981    if SamSym in ['mmm',]:
982        dum = cosd(M*gam)
983        Ksl = psrs*dum
984        dKsdp = dpdps*dum
985        dKsdg = -psrs*M*sind(M*gam)
986    else:
987        dum = cosd(M*gam)+sind(M*gam)
988        Ksl = psrs*dum
989        dKsdp = dpdps*dum
990        dKsdg = psrs*M*(-sind(M*gam)+cosd(M*gam))
991    return Ksl,dKsdp,dKsdg
992   
993def GetKclKsl(L,N,SGLaue,psi,phi,beta):
994    """
995    This is used for spherical harmonics description of preferred orientation;
996        cylindrical symmetry only (M=0) and no sample angle derivatives returned
997    """
998    import pytexture as ptx
999    RSQ2PI = 0.3989422804014
1000    SQ2 = 1.414213562373
1001    Ksl,x = ptx.pyplmpsi(L,0,1,psi)
1002    Ksl *= RSQ2PI
1003    if SGLaue in ['m3','m3m']:
1004        Kcl = 0.0
1005        for j in range(0,L+1,4):
1006            im = j/4+1
1007            pcrs,dum = ptx.pyplmpsi(L,j,1,phi)
1008            Kcl += BOH['L='+str(L)][N-1][im-1]*pcrs*cosd(j*beta)       
1009    else:
1010        pcrs,dum = ptx.pyplmpsi(L,N,1,phi)
1011        pcrs *= RSQ2PI
1012        if N:
1013            pcrs *= SQ2
1014        if SGLaue in ['mmm','4/mmm','6/mmm','R3mR','3m1','31m']:
1015            if SGLaue in ['3mR','3m1','31m']: 
1016                if n%6 == 3:
1017                    Kcl = pcrs*sind(N*beta)
1018                else:
1019                    Kcl = pcrs*cosd(N*beta)
1020            else:
1021                Kcl = pcrs*cosd(N*beta)
1022        else:
1023            Kcl = pcrs*(cosd(N*beta)+sind(N*beta))
1024    return Kcl*Ksl,Lnorm(L)
1025   
1026def Glnh(Start,SHCoef,psi,gam,SamSym):
1027    import pytexture as ptx
1028    RSQPI = 0.5641895835478
1029    SQ2 = 1.414213562373
1030
1031    if Start:
1032        ptx.pyqlmninit()
1033        Start = False
1034    Fln = np.zeros(len(SHCoef))
1035    for i,term in enumerate(SHCoef):
1036        l,m,n = eval(term.strip('C'))
1037        pcrs,dum = ptx.pyplmpsi(l,m,1,psi)
1038        pcrs *= RSQPI
1039        if m == 0:
1040            pcrs /= SQ2
1041        if SamSym in ['mmm',]:
1042            Ksl = pcrs*cosd(m*gam)
1043        else:
1044            Ksl = pcrs*(cosd(m*gam)+sind(m*gam))
1045        Fln[i] = SHCoef[term]*Ksl*Lnorm(l)
1046    ODFln = dict(zip(SHCoef.keys(),list(zip(SHCoef.values(),Fln))))
1047    return ODFln
1048
1049def Flnh(Start,SHCoef,phi,beta,SGData):
1050    import pytexture as ptx
1051   
1052    FORPI = 12.5663706143592
1053    RSQPI = 0.5641895835478
1054    SQ2 = 1.414213562373
1055
1056    if Start:
1057        ptx.pyqlmninit()
1058        Start = False
1059    Fln = np.zeros(len(SHCoef))
1060    for i,term in enumerate(SHCoef):
1061        l,m,n = eval(term.strip('C'))
1062        if SGData['SGLaue'] in ['m3','m3m']:
1063            Kcl = 0.0
1064            for j in range(0,l+1,4):
1065                im = j/4+1
1066                pcrs,dum = ptx.pyplmpsi(l,j,1,phi)
1067                Kcl += BOH['L='+str(l)][n-1][im-1]*pcrs*cosd(j*beta)       
1068        else:                #all but cubic
1069            pcrs,dum = ptx.pyplmpsi(l,n,1,phi)
1070            pcrs *= RSQPI
1071            if n == 0:
1072                pcrs /= SQ2
1073            if SGData['SGLaue'] in ['mmm','4/mmm','6/mmm','R3mR','3m1','31m']:
1074               if SGData['SGLaue'] in ['3mR','3m1','31m']: 
1075                   if n%6 == 3:
1076                       Kcl = pcrs*sind(n*beta)
1077                   else:
1078                       Kcl = pcrs*cosd(n*beta)
1079               else:
1080                   Kcl = pcrs*cosd(n*beta)
1081            else:
1082                Kcl = pcrs*(cosd(n*beta)+sind(n*beta))
1083        Fln[i] = SHCoef[term]*Kcl*Lnorm(l)
1084    ODFln = dict(zip(SHCoef.keys(),list(zip(SHCoef.values(),Fln))))
1085    return ODFln
1086   
1087def polfcal(ODFln,SamSym,psi,gam):
1088    import pytexture as ptx
1089    RSQPI = 0.5641895835478
1090    SQ2 = 1.414213562373
1091    PolVal = np.ones_like(gam)
1092    for term in ODFln:
1093        if abs(ODFln[term][1]) > 1.e-3:
1094            l,m,n = eval(term.strip('C'))
1095            psrs,dum = ptx.pyplmpsi(l,m,len(psi),psi)
1096            if SamSym in ['-1','2/m']:
1097                if m != 0:
1098                    Ksl = RSQPI*psrs*(cosd(m*gam)+sind(m*gam))
1099                else:
1100                    Ksl = RSQPI*psrs/SQ2
1101            else:
1102                if m != 0:
1103                    Ksl = RSQPI*psrs*cosd(m*gam)
1104                else:
1105                    Ksl = RSQPI*psrs/SQ2
1106            PolVal += ODFln[term][1]*Ksl
1107    return PolVal
1108   
1109def invpolfcal(ODFln,SGData,phi,beta):
1110    import pytexture as ptx
1111   
1112    FORPI = 12.5663706143592
1113    RSQPI = 0.5641895835478
1114    SQ2 = 1.414213562373
1115
1116    invPolVal = np.ones_like(beta)
1117    for term in ODFln:
1118        if abs(ODFln[term][1]) > 1.e-3:
1119            l,m,n = eval(term.strip('C'))
1120            if SGData['SGLaue'] in ['m3','m3m']:
1121                Kcl = 0.0
1122                for j in range(0,l+1,4):
1123                    im = j/4+1
1124                    pcrs,dum = ptx.pyplmpsi(l,j,len(beta),phi)
1125                    Kcl += BOH['L='+str(l)][n-1][im-1]*pcrs*cosd(j*beta)       
1126            else:                #all but cubic
1127                pcrs,dum = ptx.pyplmpsi(l,n,len(beta),phi)
1128                pcrs *= RSQPI
1129                if n == 0:
1130                    pcrs /= SQ2
1131                if SGData['SGLaue'] in ['mmm','4/mmm','6/mmm','R3mR','3m1','31m']:
1132                   if SGData['SGLaue'] in ['3mR','3m1','31m']: 
1133                       if n%6 == 3:
1134                           Kcl = pcrs*sind(n*beta)
1135                       else:
1136                           Kcl = pcrs*cosd(n*beta)
1137                   else:
1138                       Kcl = pcrs*cosd(n*beta)
1139                else:
1140                    Kcl = pcrs*(cosd(n*beta)+sind(n*beta))
1141            invPolVal += ODFln[term][1]*Kcl
1142    return invPolVal
1143   
1144   
1145def textureIndex(SHCoef):
1146    Tindx = 1.0
1147    for term in SHCoef:
1148        l = eval(term.strip('C'))[0]
1149        Tindx += SHCoef[term]**2/(2.0*l+1.)
1150    return Tindx
1151   
1152# output from uctbx computed on platform darwin on 2010-05-28
1153NeedTestData = True
1154def TestData():
1155    array = np.array
1156    global NeedTestData
1157    NeedTestData = False
1158    global CellTestData
1159    CellTestData = [
1160# cell, g, G, cell*, V, V*
1161  [(4, 4, 4, 90, 90, 90), 
1162   array([[  1.60000000e+01,   9.79717439e-16,   9.79717439e-16],
1163       [  9.79717439e-16,   1.60000000e+01,   9.79717439e-16],
1164       [  9.79717439e-16,   9.79717439e-16,   1.60000000e+01]]), array([[  6.25000000e-02,   3.82702125e-18,   3.82702125e-18],
1165       [  3.82702125e-18,   6.25000000e-02,   3.82702125e-18],
1166       [  3.82702125e-18,   3.82702125e-18,   6.25000000e-02]]), (0.25, 0.25, 0.25, 90.0, 90.0, 90.0), 64.0, 0.015625],
1167# cell, g, G, cell*, V, V*
1168  [(4.0999999999999996, 5.2000000000000002, 6.2999999999999998, 100, 80, 130), 
1169   array([[ 16.81      , -13.70423184,   4.48533243],
1170       [-13.70423184,  27.04      ,  -5.6887143 ],
1171       [  4.48533243,  -5.6887143 ,  39.69      ]]), array([[ 0.10206349,  0.05083339, -0.00424823],
1172       [ 0.05083339,  0.06344997,  0.00334956],
1173       [-0.00424823,  0.00334956,  0.02615544]]), (0.31947376387537696, 0.25189277536327803, 0.16172643497798223, 85.283666420376008, 94.716333579624006, 50.825714168082683), 100.98576357983838, 0.0099023858863968445],
1174# cell, g, G, cell*, V, V*
1175  [(3.5, 3.5, 6, 90, 90, 120), 
1176   array([[  1.22500000e+01,  -6.12500000e+00,   1.28587914e-15],
1177       [ -6.12500000e+00,   1.22500000e+01,   1.28587914e-15],
1178       [  1.28587914e-15,   1.28587914e-15,   3.60000000e+01]]), array([[  1.08843537e-01,   5.44217687e-02,   3.36690552e-18],
1179       [  5.44217687e-02,   1.08843537e-01,   3.36690552e-18],
1180       [  3.36690552e-18,   3.36690552e-18,   2.77777778e-02]]), (0.32991443953692895, 0.32991443953692895, 0.16666666666666669, 90.0, 90.0, 60.000000000000021), 63.652867178156257, 0.015710211406520427],
1181  ]
1182    global CoordTestData
1183    CoordTestData = [
1184# cell, ((frac, ortho),...)
1185  ((4,4,4,90,90,90,), [
1186 ((0.10000000000000001, 0.0, 0.0),(0.40000000000000002, 0.0, 0.0)),
1187 ((0.0, 0.10000000000000001, 0.0),(2.4492935982947065e-17, 0.40000000000000002, 0.0)),
1188 ((0.0, 0.0, 0.10000000000000001),(2.4492935982947065e-17, -2.4492935982947065e-17, 0.40000000000000002)),
1189 ((0.10000000000000001, 0.20000000000000001, 0.29999999999999999),(0.40000000000000013, 0.79999999999999993, 1.2)),
1190 ((0.20000000000000001, 0.29999999999999999, 0.10000000000000001),(0.80000000000000016, 1.2, 0.40000000000000002)),
1191 ((0.29999999999999999, 0.20000000000000001, 0.10000000000000001),(1.2, 0.80000000000000004, 0.40000000000000002)),
1192 ((0.5, 0.5, 0.5),(2.0, 1.9999999999999998, 2.0)),
1193]),
1194# cell, ((frac, ortho),...)
1195  ((4.1,5.2,6.3,100,80,130,), [
1196 ((0.10000000000000001, 0.0, 0.0),(0.40999999999999998, 0.0, 0.0)),
1197 ((0.0, 0.10000000000000001, 0.0),(-0.33424955703700043, 0.39834311042186865, 0.0)),
1198 ((0.0, 0.0, 0.10000000000000001),(0.10939835193016617, -0.051013289294572106, 0.6183281045774256)),
1199 ((0.10000000000000001, 0.20000000000000001, 0.29999999999999999),(0.069695941716497567, 0.64364635296002093, 1.8549843137322766)),
1200 ((0.20000000000000001, 0.29999999999999999, 0.10000000000000001),(-0.073350319180835066, 1.1440160419710339, 0.6183281045774256)),
1201 ((0.29999999999999999, 0.20000000000000001, 0.10000000000000001),(0.67089923785616512, 0.74567293154916525, 0.6183281045774256)),
1202 ((0.5, 0.5, 0.5),(0.92574397446582857, 1.7366491056364828, 3.0916405228871278)),
1203]),
1204# cell, ((frac, ortho),...)
1205  ((3.5,3.5,6,90,90,120,), [
1206 ((0.10000000000000001, 0.0, 0.0),(0.35000000000000003, 0.0, 0.0)),
1207 ((0.0, 0.10000000000000001, 0.0),(-0.17499999999999993, 0.3031088913245536, 0.0)),
1208 ((0.0, 0.0, 0.10000000000000001),(3.6739403974420595e-17, -3.6739403974420595e-17, 0.60000000000000009)),
1209 ((0.10000000000000001, 0.20000000000000001, 0.29999999999999999),(2.7675166561703527e-16, 0.60621778264910708, 1.7999999999999998)),
1210 ((0.20000000000000001, 0.29999999999999999, 0.10000000000000001),(0.17500000000000041, 0.90932667397366063, 0.60000000000000009)),
1211 ((0.29999999999999999, 0.20000000000000001, 0.10000000000000001),(0.70000000000000018, 0.6062177826491072, 0.60000000000000009)),
1212 ((0.5, 0.5, 0.5),(0.87500000000000067, 1.5155444566227676, 3.0)),
1213]),
1214]
1215    global LaueTestData             #generated by GSAS
1216    LaueTestData = {
1217    'R 3 m':[(4.,4.,6.,90.,90.,120.),((1,0,1,6),(1,0,-2,6),(0,0,3,2),(1,1,0,6),(2,0,-1,6),(2,0,2,6),
1218        (1,1,3,12),(1,0,4,6),(2,1,1,12),(2,1,-2,12),(3,0,0,6),(1,0,-5,6),(2,0,-4,6),(3,0,-3,6),(3,0,3,6),
1219        (0,0,6,2),(2,2,0,6),(2,1,4,12),(2,0,5,6),(3,1,-1,12),(3,1,2,12),(1,1,6,12),(2,2,3,12),(2,1,-5,12))],
1220    'R 3':[(4.,4.,6.,90.,90.,120.),((1,0,1,6),(1,0,-2,6),(0,0,3,2),(1,1,0,6),(2,0,-1,6),(2,0,2,6),(1,1,3,6),
1221        (1,1,-3,6),(1,0,4,6),(3,-1,1,6),(2,1,1,6),(3,-1,-2,6),(2,1,-2,6),(3,0,0,6),(1,0,-5,6),(2,0,-4,6),
1222        (2,2,0,6),(3,0,3,6),(3,0,-3,6),(0,0,6,2),(3,-1,4,6),(2,0,5,6),(2,1,4,6),(4,-1,-1,6),(3,1,-1,6),
1223        (3,1,2,6),(4,-1,2,6),(2,2,-3,6),(1,1,-6,6),(1,1,6,6),(2,2,3,6),(2,1,-5,6),(3,-1,-5,6))],
1224    'P 3':[(4.,4.,6.,90.,90.,120.),((0,0,1,2),(1,0,0,6),(1,0,1,6),(0,0,2,2),(1,0,-1,6),(1,0,2,6),(1,0,-2,6),
1225        (1,1,0,6),(0,0,3,2),(1,1,1,6),(1,1,-1,6),(1,0,3,6),(1,0,-3,6),(2,0,0,6),(2,0,-1,6),(1,1,-2,6),
1226        (1,1,2,6),(2,0,1,6),(2,0,-2,6),(2,0,2,6),(0,0,4,2),(1,1,-3,6),(1,1,3,6),(1,0,-4,6),(1,0,4,6),
1227        (2,0,-3,6),(2,1,0,6),(2,0,3,6),(3,-1,0,6),(2,1,1,6),(3,-1,-1,6),(2,1,-1,6),(3,-1,1,6),(1,1,4,6),
1228        (3,-1,2,6),(3,-1,-2,6),(1,1,-4,6),(0,0,5,2),(2,1,2,6),(2,1,-2,6),(3,0,0,6),(3,0,1,6),(2,0,4,6),
1229        (2,0,-4,6),(3,0,-1,6),(1,0,-5,6),(1,0,5,6),(3,-1,-3,6),(2,1,-3,6),(2,1,3,6),(3,-1,3,6),(3,0,-2,6),
1230        (3,0,2,6),(1,1,5,6),(1,1,-5,6),(2,2,0,6),(3,0,3,6),(3,0,-3,6),(0,0,6,2),(2,0,-5,6),(2,1,-4,6),
1231        (2,2,-1,6),(3,-1,-4,6),(2,2,1,6),(3,-1,4,6),(2,1,4,6),(2,0,5,6),(1,0,-6,6),(1,0,6,6),(4,-1,0,6),
1232        (3,1,0,6),(3,1,-1,6),(3,1,1,6),(4,-1,-1,6),(2,2,2,6),(4,-1,1,6),(2,2,-2,6),(3,1,2,6),(3,1,-2,6),
1233        (3,0,4,6),(3,0,-4,6),(4,-1,-2,6),(4,-1,2,6),(2,2,-3,6),(1,1,6,6),(1,1,-6,6),(2,2,3,6),(3,-1,5,6),
1234        (2,1,5,6),(2,1,-5,6),(3,-1,-5,6))],
1235    'P 3 m 1':[(4.,4.,6.,90.,90.,120.),((0,0,1,2),(1,0,0,6),(1,0,-1,6),(1,0,1,6),(0,0,2,2),(1,0,-2,6),
1236        (1,0,2,6),(1,1,0,6),(0,0,3,2),(1,1,1,12),(1,0,-3,6),(1,0,3,6),(2,0,0,6),(1,1,2,12),(2,0,1,6),
1237        (2,0,-1,6),(0,0,4,2),(2,0,-2,6),(2,0,2,6),(1,1,3,12),(1,0,-4,6),(1,0,4,6),(2,0,3,6),(2,1,0,12),
1238        (2,0,-3,6),(2,1,1,12),(2,1,-1,12),(1,1,4,12),(2,1,2,12),(0,0,5,2),(2,1,-2,12),(3,0,0,6),(1,0,-5,6),
1239        (3,0,1,6),(3,0,-1,6),(1,0,5,6),(2,0,4,6),(2,0,-4,6),(2,1,3,12),(2,1,-3,12),(3,0,-2,6),(3,0,2,6),
1240        (1,1,5,12),(3,0,-3,6),(0,0,6,2),(2,2,0,6),(3,0,3,6),(2,1,4,12),(2,2,1,12),(2,0,5,6),(2,1,-4,12),
1241        (2,0,-5,6),(1,0,-6,6),(1,0,6,6),(3,1,0,12),(3,1,-1,12),(3,1,1,12),(2,2,2,12),(3,1,2,12),
1242        (3,0,4,6),(3,1,-2,12),(3,0,-4,6),(1,1,6,12),(2,2,3,12))],
1243    'P 3 1 m':[(4.,4.,6.,90.,90.,120.),((0,0,1,2),(1,0,0,6),(0,0,2,2),(1,0,1,12),(1,0,2,12),(1,1,0,6),
1244        (0,0,3,2),(1,1,-1,6),(1,1,1,6),(1,0,3,12),(2,0,0,6),(2,0,1,12),(1,1,2,6),(1,1,-2,6),(2,0,2,12),
1245        (0,0,4,2),(1,1,-3,6),(1,1,3,6),(1,0,4,12),(2,1,0,12),(2,0,3,12),(2,1,1,12),(2,1,-1,12),(1,1,-4,6),
1246        (1,1,4,6),(0,0,5,2),(2,1,-2,12),(2,1,2,12),(3,0,0,6),(1,0,5,12),(2,0,4,12),(3,0,1,12),(2,1,-3,12),
1247        (2,1,3,12),(3,0,2,12),(1,1,5,6),(1,1,-5,6),(3,0,3,12),(0,0,6,2),(2,2,0,6),(2,1,-4,12),(2,0,5,12),
1248        (2,2,-1,6),(2,2,1,6),(2,1,4,12),(3,1,0,12),(1,0,6,12),(2,2,2,6),(3,1,-1,12),(2,2,-2,6),(3,1,1,12),
1249        (3,1,-2,12),(3,0,4,12),(3,1,2,12),(1,1,-6,6),(2,2,3,6),(2,2,-3,6),(1,1,6,6))],
1250    }
1251   
1252    global FLnhTestData
1253    FLnhTestData = [{
1254    'C(4,0,0)': (0.965, 0.42760447),
1255    'C(2,0,0)': (1.0122, -0.80233610),
1256    'C(2,0,2)': (0.0061, 8.37491546E-03),
1257    'C(6,0,4)': (-0.0898, 4.37985696E-02),
1258    'C(6,0,6)': (-0.1369, -9.04081762E-02),
1259    'C(6,0,0)': (0.5935, -0.18234928),
1260    'C(4,0,4)': (0.1872, 0.16358127),
1261    'C(6,0,2)': (0.6193, 0.27573633),
1262    'C(4,0,2)': (-0.1897, 0.12530720)},[1,0,0]]
1263def test0():
1264    if NeedTestData: TestData()
1265    msg = 'test cell2Gmat, fillgmat, Gmat2cell'
1266    for (cell, tg, tG, trcell, tV, trV) in CellTestData:
1267        G, g = cell2Gmat(cell)
1268        assert np.allclose(G,tG),msg
1269        assert np.allclose(g,tg),msg
1270        tcell = Gmat2cell(g)
1271        assert np.allclose(cell,tcell),msg
1272        tcell = Gmat2cell(G)
1273        assert np.allclose(tcell,trcell),msg
1274
1275def test1():
1276    if NeedTestData: TestData()
1277    msg = 'test cell2A and A2Gmat'
1278    for (cell, tg, tG, trcell, tV, trV) in CellTestData:
1279        G, g = A2Gmat(cell2A(cell))
1280        assert np.allclose(G,tG),msg
1281        assert np.allclose(g,tg),msg
1282
1283def test2():
1284    if NeedTestData: TestData()
1285    msg = 'test Gmat2A, A2cell, A2Gmat, Gmat2cell'
1286    for (cell, tg, tG, trcell, tV, trV) in CellTestData:
1287        G, g = cell2Gmat(cell)
1288        tcell = A2cell(Gmat2A(G))
1289        assert np.allclose(cell,tcell),msg
1290
1291def test3():
1292    if NeedTestData: TestData()
1293    msg = 'test invcell2Gmat'
1294    for (cell, tg, tG, trcell, tV, trV) in CellTestData:
1295        G, g = invcell2Gmat(trcell)
1296        assert np.allclose(G,tG),msg
1297        assert np.allclose(g,tg),msg
1298
1299def test4():
1300    if NeedTestData: TestData()
1301    msg = 'test calc_rVsq, calc_rV, calc_V'
1302    for (cell, tg, tG, trcell, tV, trV) in CellTestData:
1303        assert np.allclose(calc_rV(cell2A(cell)),trV), msg
1304        assert np.allclose(calc_V(cell2A(cell)),tV), msg
1305
1306def test5():
1307    if NeedTestData: TestData()
1308    msg = 'test A2invcell'
1309    for (cell, tg, tG, trcell, tV, trV) in CellTestData:
1310        rcell = A2invcell(cell2A(cell))
1311        assert np.allclose(rcell,trcell),msg
1312
1313def test6():
1314    if NeedTestData: TestData()
1315    msg = 'test cell2AB'
1316    for (cell,coordlist) in CoordTestData:
1317        A,B = cell2AB(cell)
1318        for (frac,ortho) in coordlist:
1319            to = np.inner(A,frac)
1320            tf = np.inner(B,to)
1321            assert np.allclose(ortho,to), msg
1322            assert np.allclose(frac,tf), msg
1323            to = np.sum(A*frac,axis=1)
1324            tf = np.sum(B*to,axis=1)
1325            assert np.allclose(ortho,to), msg
1326            assert np.allclose(frac,tf), msg
1327
1328# test GetBraviasNum(...) and GenHBravais(...)
1329def test7():
1330    import os.path
1331    import sys
1332    import GSASIIspc as spc
1333    testdir = os.path.join(os.path.split(os.path.abspath( __file__ ))[0],'testinp')
1334    if os.path.exists(testdir):
1335        if testdir not in sys.path: sys.path.insert(0,testdir)
1336    import sgtbxlattinp
1337    derror = 1e-4
1338    def indexmatch(hklin, hkllist, system):
1339        for hklref in hkllist:
1340            hklref = list(hklref)
1341            # these permutations are far from complete, but are sufficient to
1342            # allow the test to complete
1343            if system == 'cubic':
1344                permlist = [(1,2,3),(1,3,2),(2,1,3),(2,3,1),(3,1,2),(3,2,1),]
1345            elif system == 'monoclinic':
1346                permlist = [(1,2,3),(-1,2,-3)]
1347            else:
1348                permlist = [(1,2,3)]
1349
1350            for perm in permlist:
1351                hkl = [abs(i) * hklin[abs(i)-1] / i for i in perm]
1352                if hkl == hklref: return True
1353                if [-i for i in hkl] == hklref: return True
1354        else:
1355            return False
1356
1357    for key in sgtbxlattinp.sgtbx7:
1358        spdict = spc.SpcGroup(key)
1359        cell = sgtbxlattinp.sgtbx7[key][0]
1360        system = spdict[1]['SGSys']
1361        center = spdict[1]['SGLatt']
1362
1363        bravcode = GetBraviasNum(center, system)
1364
1365        g2list = GenHBravais(sgtbxlattinp.dmin, bravcode, cell2A(cell))
1366
1367        assert len(sgtbxlattinp.sgtbx7[key][1]) == len(g2list), 'Reflection lists differ for %s' % key
1368        for h,k,l,d,num in g2list:
1369            for hkllist,dref in sgtbxlattinp.sgtbx7[key][1]: 
1370                if abs(d-dref) < derror:
1371                    if indexmatch((h,k,l,), hkllist, system):
1372                        break
1373            else:
1374                assert 0,'No match for %s at %s (%s)' % ((h,k,l),d,key)
1375
1376def test8():
1377    import GSASIIspc as spc
1378    import sgtbxlattinp
1379    derror = 1e-4
1380    dmin = sgtbxlattinp.dmin
1381
1382    def indexmatch(hklin, hklref, system, axis):
1383        # these permutations are far from complete, but are sufficient to
1384        # allow the test to complete
1385        if system == 'cubic':
1386            permlist = [(1,2,3),(1,3,2),(2,1,3),(2,3,1),(3,1,2),(3,2,1),]
1387        elif system == 'monoclinic' and axis=='b':
1388            permlist = [(1,2,3),(-1,2,-3)]
1389        elif system == 'monoclinic' and axis=='a':
1390            permlist = [(1,2,3),(1,-2,-3)]
1391        elif system == 'monoclinic' and axis=='c':
1392            permlist = [(1,2,3),(-1,-2,3)]
1393        elif system == 'trigonal':
1394            permlist = [(1,2,3),(2,1,3),(-1,-2,3),(-2,-1,3)]
1395        elif system == 'rhombohedral':
1396            permlist = [(1,2,3),(2,3,1),(3,1,2)]
1397        else:
1398            permlist = [(1,2,3)]
1399
1400        hklref = list(hklref)
1401        for perm in permlist:
1402            hkl = [abs(i) * hklin[abs(i)-1] / i for i in perm]
1403            if hkl == hklref: return True
1404            if [-i for i in hkl] == hklref: return True
1405        return False
1406
1407    for key in sgtbxlattinp.sgtbx8:
1408        spdict = spc.SpcGroup(key)[1]
1409        cell = sgtbxlattinp.sgtbx8[key][0]
1410        center = spdict['SGLatt']
1411        Laue = spdict['SGLaue']
1412        Axis = spdict['SGUniq']
1413        system = spdict['SGSys']
1414
1415        g2list = GenHLaue(dmin,spdict,cell2A(cell))
1416        #if len(g2list) != len(sgtbxlattinp.sgtbx8[key][1]):
1417        #    print 'failed',key,':' ,len(g2list),'vs',len(sgtbxlattinp.sgtbx8[key][1])
1418        #    print 'GSAS-II:'
1419        #    for h,k,l,d in g2list: print '  ',(h,k,l),d
1420        #    print 'SGTBX:'
1421        #    for hkllist,dref in sgtbxlattinp.sgtbx8[key][1]: print '  ',hkllist,dref
1422        assert len(g2list) == len(sgtbxlattinp.sgtbx8[key][1]), (
1423            'Reflection lists differ for %s' % key
1424            )
1425        #match = True
1426        for h,k,l,d in g2list:
1427            for hkllist,dref in sgtbxlattinp.sgtbx8[key][1]: 
1428                if abs(d-dref) < derror:
1429                    if indexmatch((h,k,l,), hkllist, system, Axis): break
1430            else:
1431                assert 0,'No match for %s at %s (%s)' % ((h,k,l),d,key)
1432                #match = False
1433        #if not match:
1434            #for hkllist,dref in sgtbxlattinp.sgtbx8[key][1]: print '  ',hkllist,dref
1435            #print center, Laue, Axis, system
1436           
1437def test9():
1438    import GSASIIspc as G2spc
1439    if NeedTestData: TestData()
1440    for spc in LaueTestData:
1441        data = LaueTestData[spc]
1442        cell = data[0]
1443        hklm = np.array(data[1])
1444        H = hklm[-1][:3]
1445        hklO = hklm.T[:3].T
1446        A = cell2A(cell)
1447        dmin = 1./np.sqrt(calc_rDsq(H,A))
1448        SGData = G2spc.SpcGroup(spc)[1]
1449        hkls = np.array(GenHLaue(dmin,SGData,A))
1450        hklN = hkls.T[:3].T
1451        print spc,hklO.shape,hklN.shape
1452        for H in hklO:
1453            if H not in hklN:
1454                print H,' missing from hkl from GSASII'
1455       
1456       
1457   
1458
1459if __name__ == '__main__':
1460    test0()
1461    test1()
1462    test2()
1463    test3()
1464    test4()
1465    test5()
1466    test6()
1467    test7()
1468    test8()
1469    print "OK"
Note: See TracBrowser for help on using the repository browser.