source: trunk/GSASIIlattice.py @ 949

Last change on this file since 949 was 949, checked in by vondreele, 10 years ago

changes to MC/SA stuff

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