1 | "GSASII - Space group interpretion routines" |
---|
2 | |
---|
3 | ########### SVN repository information ################### |
---|
4 | # $Date: 2012-03-01 16:28:56 +0000 (Thu, 01 Mar 2012) $ |
---|
5 | # $Author: vondreele $ |
---|
6 | # $Revision: 504 $ |
---|
7 | # $URL: trunk/GSASIIspc.py $ |
---|
8 | # $Id: GSASIIspc.py 504 2012-03-01 16:28:56Z vondreele $ |
---|
9 | ########### SVN repository information ################### |
---|
10 | import numpy as np |
---|
11 | import numpy.ma as ma |
---|
12 | import numpy.linalg as nl |
---|
13 | import math |
---|
14 | import sys |
---|
15 | import os.path as ospath |
---|
16 | |
---|
17 | import GSASIIpath |
---|
18 | import pyspg |
---|
19 | |
---|
20 | def SpcGroup(SGSymbol): |
---|
21 | ''' |
---|
22 | Determines cell and symmetry information from a short H-M space group name |
---|
23 | input: |
---|
24 | SGSymbol - space group symbol (string) with spaces between axial fields |
---|
25 | returns: |
---|
26 | SGError = 0 for no errors; >0 for errors (see SGErrors below for details) |
---|
27 | SGData - dictionary with entries: |
---|
28 | 'SpGrp': space group symbol slightly cleaned up |
---|
29 | 'Laue': one of '-1','2/m','mmm','4/m','4/mmm','3R','3mR','3', |
---|
30 | '3m1','31m','6/m','6/mmm','m3','m3m' |
---|
31 | 'SGInv': boolean; True if centrosymmetric, False if not |
---|
32 | 'SGLatt': one of 'P','A','B','C','I','F','R' |
---|
33 | 'SGUniq': one of 'a','b','c' if monoclinic, '' otherwise |
---|
34 | 'SGCen': cell centering vectors [0,0,0] at least |
---|
35 | 'SGOps': symmetry operations as [M,T] so that M*x+T = x' |
---|
36 | 'SGSys': one of 'triclinic','monoclinic','orthorhombic','tetragonal','rhombohedral','trigonal','hexagonal','cubic' |
---|
37 | ''' |
---|
38 | LaueSym = ('-1','2/m','mmm','4/m','4/mmm','3R','3mR','3','3m1','31m','6/m','6/mmm','m3','m3m') |
---|
39 | LattSym = ('P','A','B','C','I','F','R') |
---|
40 | UniqSym = ('','','a','b','c','',) |
---|
41 | SysSym = ('triclinic','monoclinic','orthorhombic','tetragonal','rhombohedral','trigonal','hexagonal','cubic') |
---|
42 | SGData = {} |
---|
43 | SGData['SpGrp'] = SGSymbol.strip().lower().capitalize() |
---|
44 | SGInfo = pyspg.sgforpy(SGSymbol) |
---|
45 | SGData['SGLaue'] = LaueSym[SGInfo[0]-1] |
---|
46 | SGData['SGInv'] = bool(SGInfo[1]) |
---|
47 | SGData['SGLatt'] = LattSym[SGInfo[2]-1] |
---|
48 | SGData['SGUniq'] = UniqSym[SGInfo[3]+1] |
---|
49 | if SGData['SGLatt'] == 'P': |
---|
50 | SGData['SGCen'] = np.array(([0,0,0],)) |
---|
51 | elif SGData['SGLatt'] == 'A': |
---|
52 | SGData['SGCen'] = np.array(([0,0,0],[0,.5,.5])) |
---|
53 | elif SGData['SGLatt'] == 'B': |
---|
54 | SGData['SGCen'] = np.array(([0,0,0],[.5,0,.5])) |
---|
55 | elif SGData['SGLatt'] == 'C': |
---|
56 | SGData['SGCen'] = np.array(([0,0,0],[.5,.5,0,])) |
---|
57 | elif SGData['SGLatt'] == 'I': |
---|
58 | SGData['SGCen'] = np.array(([0,0,0],[.5,.5,.5])) |
---|
59 | elif SGData['SGLatt'] == 'F': |
---|
60 | SGData['SGCen'] = np.array(([0,0,0],[0,.5,.5],[.5,0,.5],[.5,.5,0,])) |
---|
61 | elif SGData['SGLatt'] == 'R': |
---|
62 | SGData['SGCen'] = np.array(([0,0,0],[1./3.,2./3.,2./3.],[2./3.,1./3.,1./3.])) |
---|
63 | SGData['SGOps'] = [] |
---|
64 | for i in range(SGInfo[5]): |
---|
65 | Mat = np.array(SGInfo[6][i]) |
---|
66 | Trns = np.array(SGInfo[7][i]) |
---|
67 | SGData['SGOps'].append([Mat,Trns]) |
---|
68 | if SGData['SGLaue'] in '-1': |
---|
69 | SGData['SGSys'] = SysSym[0] |
---|
70 | elif SGData['SGLaue'] in '2/m': |
---|
71 | SGData['SGSys'] = SysSym[1] |
---|
72 | elif SGData['SGLaue'] in 'mmm': |
---|
73 | SGData['SGSys'] = SysSym[2] |
---|
74 | elif SGData['SGLaue'] in ['4/m','4/mmm']: |
---|
75 | SGData['SGSys'] = SysSym[3] |
---|
76 | elif SGData['SGLaue'] in ['3R','3mR']: |
---|
77 | SGData['SGSys'] = SysSym[4] |
---|
78 | elif SGData['SGLaue'] in ['3','3m1','31m']: |
---|
79 | SGData['SGSys'] = SysSym[5] |
---|
80 | elif SGData['SGLaue'] in ['6/m','6/mmm']: |
---|
81 | SGData['SGSys'] = SysSym[6] |
---|
82 | elif SGData['SGLaue'] in ['m3','m3m']: |
---|
83 | SGData['SGSys'] = SysSym[7] |
---|
84 | SGData['SGPolax'] = SGpolar(SGData) |
---|
85 | return SGInfo[8],SGData |
---|
86 | |
---|
87 | def SGErrors(IErr): |
---|
88 | ''' |
---|
89 | Interprets the error message code from SpcGroup. Used in SpaceGroup. |
---|
90 | input: |
---|
91 | SGError - from SpcGroup |
---|
92 | returns: |
---|
93 | ErrString - a string with the error message or "Unknown error" |
---|
94 | ''' |
---|
95 | |
---|
96 | ErrString = [' ', |
---|
97 | 'Less than 2 operator fields were found', |
---|
98 | 'Illegal Lattice type, not P, A, B, C, I, F or R', |
---|
99 | 'Rhombohedral lattice requires a 3-axis', |
---|
100 | 'Minus sign does not preceed 1, 2, 3, 4 or 6', |
---|
101 | 'Either a 5-axis anywhere or a 3-axis in field not allowed', |
---|
102 | ' ', |
---|
103 | 'I for COMPUTED GO TO out of range.', |
---|
104 | 'An a-glide mirror normal to A not allowed', |
---|
105 | 'A b-glide mirror normal to B not allowed', |
---|
106 | 'A c-glide mirror normal to C not allowed', |
---|
107 | 'D-glide in a primitive lattice not allowed', |
---|
108 | 'A 4-axis not allowed in the 2nd operator field', |
---|
109 | 'A 6-axis not allowed in the 2nd operator field', |
---|
110 | 'More than 24 matrices needed to define group', |
---|
111 | ' ', |
---|
112 | 'Improper construction of a rotation operator', |
---|
113 | 'Mirror following a / not allowed', |
---|
114 | 'A translation conflict between operators', |
---|
115 | 'The 2bar operator is not allowed', |
---|
116 | '3 fields are legal only in R & m3 cubic groups', |
---|
117 | 'Syntax error. Expected I -4 3 d at this point', |
---|
118 | ' ', |
---|
119 | 'A or B centered tetragonal not allowed', |
---|
120 | ' ','unknown error in sgroup',' ',' ',' ', |
---|
121 | 'Illegal character in the space group symbol', |
---|
122 | ] |
---|
123 | try: |
---|
124 | return ErrString[IErr] |
---|
125 | except: |
---|
126 | return "Unknown error" |
---|
127 | def SGpolar(SGData): |
---|
128 | ''' |
---|
129 | Determine identity of polar axes if any |
---|
130 | ''' |
---|
131 | POL = ('','x','y','x y','z','x z','y z','xyz','111') |
---|
132 | NP = [1,2,4] |
---|
133 | NPZ = [0,1] |
---|
134 | for M,T in SGData['SGOps']: |
---|
135 | for i in range(3): |
---|
136 | if M[i][i] <= 0.: NP[i] = 0 |
---|
137 | if M[0][2] > 0: NPZ[0] = 8 |
---|
138 | if M[1][2] > 0: NPZ[1] = 0 |
---|
139 | NPol = (NP[0]+NP[1]+NP[2]+NPZ[0]*NPZ[1])*(1-int(SGData['SGInv'])) |
---|
140 | return POL[NPol] |
---|
141 | |
---|
142 | def SGPrint(SGData): |
---|
143 | ''' |
---|
144 | Print the output of SpcGroup in a nicely formatted way. Used in SpaceGroup |
---|
145 | input: |
---|
146 | SGData - from SpcGroup |
---|
147 | returns: |
---|
148 | SGText - list of strings with the space group details |
---|
149 | ''' |
---|
150 | Mult = len(SGData['SGCen'])*len(SGData['SGOps'])*(int(SGData['SGInv'])+1) |
---|
151 | SGText = [] |
---|
152 | SGText.append(' Space Group: '+SGData['SpGrp']) |
---|
153 | CentStr = 'centrosymmetric' |
---|
154 | if not SGData['SGInv']: |
---|
155 | CentStr = 'non'+CentStr |
---|
156 | if SGData['SGLatt'] in 'ABCIFR': |
---|
157 | SGText.append(' The lattice is '+CentStr+' '+SGData['SGLatt']+'-centered '+SGData['SGSys'].lower()) |
---|
158 | else: |
---|
159 | SGText.append(' The lattice is '+CentStr+' '+'primitive '+SGData['SGSys'].lower()) |
---|
160 | SGText.append(' Multiplicity of a general site is '+str(Mult)) |
---|
161 | SGText.append(' The Laue symmetry is '+SGData['SGLaue']) |
---|
162 | if SGData['SGUniq'] in ['a','b','c']: |
---|
163 | SGText.append(' The unique monoclinic axis is '+SGData['SGUniq']) |
---|
164 | if SGData['SGInv']: |
---|
165 | SGText.append(' The inversion center is located at 0,0,0') |
---|
166 | if SGData['SGPolax']: |
---|
167 | SGText.append(' The location of the origin is arbitrary in '+SGData['SGPolax']) |
---|
168 | SGText.append('\n'+' The equivalent positions are:') |
---|
169 | if SGData['SGLatt'] != 'P': |
---|
170 | SGText.append('\n ('+Latt2text(SGData['SGLatt'])+')+\n') |
---|
171 | Ncol = 2 |
---|
172 | line = ' ' |
---|
173 | col = 0 |
---|
174 | for iop,[M,T] in enumerate(SGData['SGOps']): |
---|
175 | OPtxt = MT2text(M,T) |
---|
176 | Fld = '(%2i) '%(iop+1)+OPtxt+'\t' |
---|
177 | line += Fld |
---|
178 | if '/' not in Fld: |
---|
179 | line += '\t' |
---|
180 | col += 1 |
---|
181 | if col == Ncol: |
---|
182 | SGText.append(line) |
---|
183 | line = ' ' |
---|
184 | col = 0 |
---|
185 | SGText.append(line) |
---|
186 | return SGText |
---|
187 | |
---|
188 | def MT2text(M,T): |
---|
189 | #From space group matrix/translation operator returns text version |
---|
190 | XYZ = ('-Z ','-Y ','-X ','X-Y','ERR','Y-X',' X ',' Y ',' Z ','+X ','+Y ','+Z ') |
---|
191 | TRA = (' ','ERR','1/6','1/4','1/3','ERR','1/2','ERR','2/3','3/4','5/6','ERR') |
---|
192 | Fld = '' |
---|
193 | for j in range(3): |
---|
194 | IJ = int(round(2*M[j][0]+3*M[j][1]+4*M[j][2]+4))%12 |
---|
195 | IK = int(round(T[j]*12))%12 |
---|
196 | if IK > 0 and IJ > 4: IJ += 3 |
---|
197 | Fld += TRA[IK]+XYZ[IJ] |
---|
198 | if j != 2: Fld += ',' |
---|
199 | return Fld |
---|
200 | |
---|
201 | def Latt2text(Latt): |
---|
202 | #From lattice type ('P',A', etc.) returns ';' delimited cell centering vectors |
---|
203 | lattTxt = {'A':'0,0,0; 0,1/2,1/2','B':'0,0,0; 1/2,0,1/2', |
---|
204 | 'C':'0,0,0; 1/2,1/2,0','I':'0,0,0; 1/2,1/2,1/2', |
---|
205 | 'F':'0,0,0; 0,1/2,1/2; 1/2,0,1/2; 1/2,1/2,0', |
---|
206 | 'R':'0,0,0; 1/3,2/3,2/3; 2/3,1/3,1/3','P':'0,0,0'} |
---|
207 | return lattTxt[Latt] |
---|
208 | |
---|
209 | def SpaceGroup(SGSymbol): |
---|
210 | ''' |
---|
211 | Print the output of SpcGroup in a nicely formatted way. |
---|
212 | input: |
---|
213 | SGSymbol - space group symbol (string) with spaces between axial fields |
---|
214 | returns: |
---|
215 | nothing |
---|
216 | ''' |
---|
217 | E,A = SpcGroup(SGSymbol) |
---|
218 | if E > 0: |
---|
219 | print SGErrors(E) |
---|
220 | return |
---|
221 | for l in SGPrint(A): |
---|
222 | print l |
---|
223 | |
---|
224 | def MoveToUnitCell(xyz): |
---|
225 | ''' |
---|
226 | Translates a set of coordinates so that all values are >=0 and < 1 |
---|
227 | input: |
---|
228 | xyz - a list or numpy array of fractional coordinates |
---|
229 | returns: |
---|
230 | XYZ - numpy array of new coordinates inside 0-1 |
---|
231 | ''' |
---|
232 | XYZ = np.zeros(3) |
---|
233 | for i,x in enumerate(xyz): |
---|
234 | XYZ[i] = (x-int(x))%1.0 |
---|
235 | return XYZ |
---|
236 | |
---|
237 | def Opposite(XYZ,toler=0.0002): |
---|
238 | ''' |
---|
239 | Gives opposite corner, edge or face of unit cell for position within tolerance. |
---|
240 | Result may be just outside the cell within tolerance |
---|
241 | input: |
---|
242 | XYZ: 0 >= np.array[x,y,z] > 1 as by MoveToUnitCell |
---|
243 | toler: unit cell fraction tolerance making opposite |
---|
244 | returns: |
---|
245 | XYZ: array of opposite positions; always contains XYZ |
---|
246 | ''' |
---|
247 | perm3 = [[1,1,1],[0,1,1],[1,0,1],[1,1,0],[1,0,0],[0,1,0],[0,0,1],[0,0,0]] |
---|
248 | TB = np.where(abs(XYZ-1)<toler,-1,0)+np.where(abs(XYZ)<toler,1,0) |
---|
249 | perm = TB*perm3 |
---|
250 | cperm = ['%d%d%d'%(i,j,k) for i,j,k in perm] |
---|
251 | D = dict(zip(cperm,perm)) |
---|
252 | new = [] |
---|
253 | for key in D: |
---|
254 | new.append(np.array(D[key])+np.array(XYZ)) |
---|
255 | return new |
---|
256 | |
---|
257 | def GenAtom(XYZ,SGData,All=False,Uij=[],Move=True): |
---|
258 | ''' |
---|
259 | Generates the equivalent positions for a specified coordinate and space group |
---|
260 | input: |
---|
261 | XYZ an array, tuple or list containing 3 elements: x, y & z |
---|
262 | SGData, from SpcGroup |
---|
263 | All = True return all equivalent positions including duplicates |
---|
264 | = False return only unique positions |
---|
265 | Uij = [U11,U22,U33,U12,U13,U23] or [] if no Uij |
---|
266 | return: [[XYZEquiv],Idup,[UijEquiv]] |
---|
267 | [XYZEquiv] is list of equivalent positions (XYZ is first entry) |
---|
268 | Idup = [-][C]SS where SS is the symmetry operator number (1-24), C (if not 0,0,0) |
---|
269 | is centering operator number (1-4) and - is for inversion |
---|
270 | Cell = unit cell translations needed to put new positions inside cell |
---|
271 | [UijEquiv] - equivalent Uij; absent if no Uij given |
---|
272 | ''' |
---|
273 | XYZEquiv = [] |
---|
274 | UijEquiv = [] |
---|
275 | Idup = [] |
---|
276 | Cell = [] |
---|
277 | X = np.array(XYZ) |
---|
278 | if Move: |
---|
279 | X = MoveToUnitCell(X) |
---|
280 | for ic,cen in enumerate(SGData['SGCen']): |
---|
281 | C = np.array(cen) |
---|
282 | for invers in range(int(SGData['SGInv']+1)): |
---|
283 | for io,[M,T] in enumerate(SGData['SGOps']): |
---|
284 | idup = ((io+1)+100*ic)*(1-2*invers) |
---|
285 | XT = np.inner(M,X)+T |
---|
286 | if len(Uij): |
---|
287 | U = Uij2U(Uij) |
---|
288 | U = np.inner(M,np.inner(U,M).T) |
---|
289 | newUij = U2Uij(U) |
---|
290 | if invers: |
---|
291 | XT = -XT |
---|
292 | XT += C |
---|
293 | newX = MoveToUnitCell(XT) |
---|
294 | cell = np.asarray(np.rint(newX-XT),dtype=np.int32) |
---|
295 | if All: |
---|
296 | if np.allclose(newX,X,atol=0.0002): |
---|
297 | idup = False |
---|
298 | else: |
---|
299 | if True in [np.allclose(newX,oldX,atol=0.0002) for oldX in XYZEquiv]: |
---|
300 | idup = False |
---|
301 | if All or idup: |
---|
302 | XYZEquiv.append(newX) |
---|
303 | Idup.append(idup) |
---|
304 | Cell.append(cell) |
---|
305 | if len(Uij): |
---|
306 | UijEquiv.append(newUij) |
---|
307 | if len(Uij): |
---|
308 | return zip(XYZEquiv,UijEquiv,Idup,Cell) |
---|
309 | else: |
---|
310 | return zip(XYZEquiv,Idup,Cell) |
---|
311 | |
---|
312 | def GenHKLf(HKL,SGData,Friedel=False): |
---|
313 | ''' |
---|
314 | Uses old GSAS Fortran routine genhkl.for |
---|
315 | input: |
---|
316 | HKL - [h,k,l] |
---|
317 | SGData - space group data obtained from SpcGroup |
---|
318 | Friedel = True to retain Friedel pairs for centrosymmetric case |
---|
319 | returns: |
---|
320 | iabsnt = True is reflection is forbidden by symmetry |
---|
321 | mulp = reflection multiplicity including Fridel pairs |
---|
322 | Uniq = numpy array of equivalent hkl in descending order of h,k,l |
---|
323 | ''' |
---|
324 | hklf = HKL+[0,] |
---|
325 | Ops = SGData['SGOps'] |
---|
326 | OpM = np.array([op[0] for op in Ops]) |
---|
327 | OpT = np.array([op[1] for op in Ops]) |
---|
328 | Inv = SGData['SGInv'] |
---|
329 | Cen = np.array([cen for cen in SGData['SGCen']]) |
---|
330 | |
---|
331 | Nuniq,Uniq,iabsnt,mulp = pyspg.genhklpy(hklf,len(Ops),OpM,OpT,SGData['SGInv'],len(Cen),Cen) |
---|
332 | h,k,l,f = Uniq |
---|
333 | Uniq=np.array(zip(h[:Nuniq],k[:Nuniq],l[:Nuniq])) |
---|
334 | phi = f[:Nuniq] |
---|
335 | Uniq = np.array(Uniq) |
---|
336 | |
---|
337 | return iabsnt,2*mulp,Uniq,phi #include Friedel pairs in powder mulp |
---|
338 | |
---|
339 | def GetOprPtrName(key): |
---|
340 | OprPtrName = { |
---|
341 | '-6643':[ 2,' 1bar ', 1],'6479' :[ 10,' 2z ', 2],'-6479':[ 9,' mz ', 3], |
---|
342 | '6481' :[ 7,' my ', 4],'-6481':[ 6,' 2y ', 5],'6641' :[ 4,' mx ', 6], |
---|
343 | '-6641':[ 3,' 2x ', 7],'6591' :[ 28,' m+-0 ', 8],'-6591':[ 27,' 2+-0 ', 9], |
---|
344 | '6531' :[ 25,' m110 ',10],'-6531':[ 24,' 2110 ',11],'6537' :[ 61,' 4z ',12], |
---|
345 | '-6537':[ 62,' -4z ',13],'975' :[ 68,' 3+++1',14],'6456' :[ 114,' 3z1 ',15], |
---|
346 | '-489' :[ 73,' 3+-- ',16],'483' :[ 78,' 3-+- ',17],'-969' :[ 83,' 3--+ ',18], |
---|
347 | '819' :[ 22,' m+0- ',19],'-819' :[ 21,' 2+0- ',20],'2431' :[ 16,' m0+- ',21], |
---|
348 | '-2431':[ 15,' 20+- ',22],'-657' :[ 19,' m101 ',23],'657' :[ 18,' 2101 ',24], |
---|
349 | '1943' :[ 48,' -4x ',25],'-1943':[ 47,' 4x ',26],'-2429':[ 13,' m011 ',27], |
---|
350 | '2429' :[ 12,' 2011 ',28],'639' :[ 55,' -4y ',29],'-639' :[ 54,' 4y ',30], |
---|
351 | '-6484':[ 146,' 2010 ', 4],'6484' :[ 139,' m010 ', 5],'-6668':[ 145,' 2100 ', 6], |
---|
352 | '6668' :[ 138,' m100 ', 7],'-6454':[ 148,' 2120 ',18],'6454' :[ 141,' m120 ',19], |
---|
353 | '-6638':[ 149,' 2210 ',20],'6638' :[ 142,' m210 ',21], #search ends here |
---|
354 | '2223' :[ 68,' 3+++2',39], |
---|
355 | '6538' :[ 106,' 6z1 ',40],'-2169':[ 83,' 3--+2',41],'2151' :[ 73,' 3+--2',42], |
---|
356 | '2205' :[ 79,'-3-+-2',43],'-2205':[ 78,' 3-+-2',44],'489' :[ 74,'-3+--1',45], |
---|
357 | '801' :[ 53,' 4y1 ',46],'1945' :[ 47,' 4x3 ',47],'-6585':[ 62,' -4z3 ',48], |
---|
358 | '6585' :[ 61,' 4z3 ',49],'6584' :[ 114,' 3z2 ',50],'6666' :[ 106,' 6z5 ',51], |
---|
359 | '6643' :[ 1,' Iden ',52],'-801' :[ 55,' -4y1 ',53],'-1945':[ 48,' -4x3 ',54], |
---|
360 | '-6666':[ 105,' -6z5 ',55],'-6538':[ 105,' -6z1 ',56],'-2223':[ 69,'-3+++2',57], |
---|
361 | '-975' :[ 69,'-3+++1',58],'-6456':[ 113,' -3z1 ',59],'-483' :[ 79,'-3-+-1',60], |
---|
362 | '969' :[ 84,'-3--+1',61],'-6584':[ 113,' -3z2 ',62],'2169' :[ 84,'-3--+2',63], |
---|
363 | '-2151':[ 74,'-3+--2',64],'0':[0,' ????',0] |
---|
364 | } |
---|
365 | return OprPtrName[key] |
---|
366 | |
---|
367 | def GetKNsym(key): |
---|
368 | KNsym = { |
---|
369 | '0' :' 1 ','1' :' -1 ','64' :' 2(100)','32' :' m(100)', |
---|
370 | '97' :'2/m(100)','16' :' 2(010)','8' :' m(010)','25' :'2/m(010)', |
---|
371 | '2' :' 2(001)','4' :' m(001)','7' :'2/m(001)','134217728' :' 2(011)', |
---|
372 | '67108864' :' m(011)','201326593' :'2/m(011)','2097152' :' 2(0+-)','1048576' :' m(0+-)', |
---|
373 | '3145729' :'2/m(0+-)','8388608' :' 2(101)','4194304' :' m(101)','12582913' :'2/m(101)', |
---|
374 | '524288' :' 2(+0-)','262144' :' m(+0-)','796433' :'2/m(+0-)','1024' :' 2(110)', |
---|
375 | '512' :' m(110)','1537' :'2/m(110)','256' :' 2(+-0)','128' :' m(+-0)', |
---|
376 | '385' :'2/m(+-0)','76' :'mm2(100)','52' :'mm2(010)','42' :'mm2(001)', |
---|
377 | '135266336' :'mm2(011)','69206048' :'mm2(0+-)','8650760' :'mm2(101)','4718600' :'mm2(+0-)', |
---|
378 | '1156' :'mm2(110)','772' :'mm2(+-0)','82' :' 222 ','136314944' :'222(100)', |
---|
379 | '8912912' :'222(010)','1282' :'222(001)','127' :' mmm ','204472417' :'mmm(100)', |
---|
380 | '13369369' :'mmm(010)','1927' :'mmm(001)','33554496' :' 4(100)','16777280' :' -4(100)', |
---|
381 | '50331745' :'4/m(100)','169869394' :'422(100)','84934738' :'-42m 100','101711948' :'4mm(100)', |
---|
382 | '254804095' :'4/mmm100','536870928 ':' 4(010)','268435472' :' -4(010)','805306393' :'4/m (10)', |
---|
383 | '545783890' :'422(010)','272891986' :'-42m 010','541327412' :'4mm(010)','818675839' :'4/mmm010', |
---|
384 | '2050' :' 4(001)','4098' :' -4(001)','6151' :'4/m(001)','3410' :'422(001)', |
---|
385 | '4818' :'-42m 001','2730' :'4mm(001)','8191' :'4/mmm001','8192' :' 3(111)', |
---|
386 | '8193' :' -3(111)','2629888' :' 32(111)','1319040' :' 3m(111)','3940737' :'-3m(111)', |
---|
387 | '32768' :' 3(+--)','32769' :' -3(+--)','10519552' :' 32(+--)','5276160' :' 3m(+--)', |
---|
388 | '15762945' :'-3m(+--)','65536' :' 3(-+-)','65537' :' -3(-+-)','134808576' :' 32(-+-)', |
---|
389 | '67437056' :' 3m(-+-)','202180097' :'-3m(-+-)','131072' :' 3(--+)','131073' :' -3(--+)', |
---|
390 | '142737664' :' 32(--+)','71434368' :' 3m(--+)','214040961' :'-3m(--+)','237650' :' 23 ', |
---|
391 | '237695' :' m3 ','715894098' :' 432 ','358068946' :' -43m ','1073725439':' m3m ', |
---|
392 | '68157504' :' mm2d100','4456464' :' mm2d010','642' :' mm2d001','153092172' :'-4m2 100', |
---|
393 | '277348404' :'-4m2 010','5418' :'-4m2 001','1075726335':' 6/mmm ','1074414420':'-6m2 100', |
---|
394 | '1075070124':'-6m2 120','1075069650':' 6mm ','1074414890':' 622 ','1073758215':' 6/m ', |
---|
395 | '1073758212':' -6 ','1073758210':' 6 ','1073759865':'-3m(100)','1075724673':'-3m(120)', |
---|
396 | '1073758800':' 3m(100)','1075069056':' 3m(120)','1073759272':' 32(100)','1074413824':' 32(120)', |
---|
397 | '1073758209':' -3 ','1073758208':' 3 ','1074135143':'mmm(100)','1075314719':'mmm(010)', |
---|
398 | '1073743751':'mmm(110)','1074004034':' mm2z100','1074790418':' mm2z010','1073742466':' mm2z110', |
---|
399 | '1074004004':'mm2(100)','1074790412':'mm2(010)','1073742980':'mm2(110)','1073872964':'mm2(120)', |
---|
400 | '1074266132':'mm2(210)','1073742596':'mm2(+-0)','1073872930':'222(100)','1074266122':'222(010)', |
---|
401 | '1073743106':'222(110)','1073741831':'2/m(001)','1073741921':'2/m(100)','1073741849':'2/m(010)', |
---|
402 | '1073743361':'2/m(110)','1074135041':'2/m(120)','1075314689':'2/m(210)','1073742209':'2/m(+-0)', |
---|
403 | '1073741828':' m(001) ','1073741888':' m(100) ','1073741840':' m(010) ','1073742336':' m(110) ', |
---|
404 | '1074003968':' m(120) ','1074790400':' m(210) ','1073741952':' m(+-0) ','1073741826':' 2(001) ', |
---|
405 | '1073741856':' 2(100) ','1073741832':' 2(010) ','1073742848':' 2(110) ','1073872896':' 2(120) ', |
---|
406 | '1074266112':' 2(210) ','1073742080':' 2(+-0) ','1073741825':' -1 ' |
---|
407 | } |
---|
408 | return KNsym[key] |
---|
409 | |
---|
410 | def GetNXUPQsym(siteSym): |
---|
411 | NXUPQsym = { |
---|
412 | ' 1 ':(28,29,28,28),' -1 ':( 1,29,28, 0),' 2(100)':(12,18,12,25),' m(100)':(25,18,12,25), |
---|
413 | '2/m(100)':( 1,18, 0,-1),' 2(010)':(13,17,13,24),' m(010)':(24,17,13,24),'2/m(010)':( 1,17, 0,-1), |
---|
414 | ' 2(001)':(14,16,14,23),' m(001)':(23,16,14,23),'2/m(001)':( 1,16, 0,-1),' 2(011)':(10,23,10,22), |
---|
415 | ' m(011)':(22,23,10,22),'2/m(011)':( 1,23, 0,-1),' 2(0+-)':(11,24,11,21),' m(0+-)':(21,24,11,21), |
---|
416 | '2/m(0+-)':( 1,24, 0,-1),' 2(101)':( 8,21, 8,20),' m(101)':(20,21, 8,20),'2/m(101)':( 1,21, 0,-1), |
---|
417 | ' 2(+0-)':( 9,22, 9,19),' m(+0-)':(19,22, 9,19),'2/m(+0-)':( 1,22, 0,-1),' 2(110)':( 6,19, 6,18), |
---|
418 | ' m(110)':(18,19, 6,18),'2/m(110)':( 1,19, 0,-1),' 2(+-0)':( 7,20, 7,17),' m(+-0)':(17,20, 7,17), |
---|
419 | '2/m(+-0)':( 1,20, 0,-1),'mm2(100)':(12,10, 0,-1),'mm2(010)':(13,10, 0,-1),'mm2(001)':(14,10, 0,-1), |
---|
420 | 'mm2(011)':(10,13, 0,-1),'mm2(0+-)':(11,13, 0,-1),'mm2(101)':( 8,12, 0,-1),'mm2(+0-)':( 9,12, 0,-1), |
---|
421 | 'mm2(110)':( 6,11, 0,-1),'mm2(+-0)':( 7,11, 0,-1),' 222 ':( 1,10, 0,-1),'222(100)':( 1,13, 0,-1), |
---|
422 | '222(010)':( 1,12, 0,-1),'222(001)':( 1,11, 0,-1),' mmm ':( 1,10, 0,-1),'mmm(100)':( 1,13, 0,-1), |
---|
423 | 'mmm(010)':( 1,12, 0,-1),'mmm(001)':( 1,11, 0,-1),' 4(100)':(12, 4,12, 0),' -4(100)':( 1, 4,12, 0), |
---|
424 | '4/m(100)':( 1, 4,12,-1),'422(100)':( 1, 4, 0,-1),'-42m 100':( 1, 4, 0,-1),'4mm(100)':(12, 4, 0,-1), |
---|
425 | '4/mmm100':( 1, 4, 0,-1),' 4(010)':(13, 3,13, 0),' -4(010)':( 1, 3,13, 0),'4/m (10)':( 1, 3,13,-1), |
---|
426 | '422(010)':( 1, 3, 0,-1),'-42m 010':( 1, 3, 0,-1),'4mm(010)':(13, 3, 0,-1),'4/mmm010':(1, 3, 0,-1,), |
---|
427 | ' 4(001)':(14, 2,14, 0),' -4(001)':( 1, 2,14, 0),'4/m(001)':( 1, 2,14,-1),'422(001)':( 1, 2, 0,-1), |
---|
428 | '-42m 001':( 1, 2, 0,-1),'4mm(001)':(14, 2, 0,-1),'4/mmm001':( 1, 2, 0,-1),' 3(111)':( 2, 5, 2, 0), |
---|
429 | ' -3(111)':( 1, 5, 2, 0),' 32(111)':( 1, 5, 0, 2),' 3m(111)':( 2, 5, 0, 2),'-3m(111)':( 1, 5, 0,-1), |
---|
430 | ' 3(+--)':( 5, 8, 5, 0),' -3(+--)':( 1, 8, 5, 0),' 32(+--)':( 1, 8, 0, 5),' 3m(+--)':( 5, 8, 0, 5), |
---|
431 | '-3m(+--)':( 1, 8, 0,-1),' 3(-+-)':( 4, 7, 4, 0),' -3(-+-)':( 1, 7, 4, 0),' 32(-+-)':( 1, 7, 0, 4), |
---|
432 | ' 3m(-+-)':( 4, 7, 0, 4),'-3m(-+-)':( 1, 7, 0,-1),' 3(--+)':( 3, 6, 3, 0),' -3(--+)':( 1, 6, 3, 0), |
---|
433 | ' 32(--+)':( 1, 6, 0, 3),' 3m(--+)':( 3, 6, 0, 3),'-3m(--+)':( 1, 6, 0,-1),' 23 ':( 1, 1, 0, 0), |
---|
434 | ' m3 ':( 1, 1, 0, 0),' 432 ':( 1, 1, 0, 0),' -43m ':( 1, 1, 0, 0),' m3m ':( 1, 1, 0, 0), |
---|
435 | ' mm2d100':(12,13, 0,-1),' mm2d010':(13,12, 0,-1),' mm2d001':(14,11, 0,-1),'-4m2 100':( 1, 4, 0,-1), |
---|
436 | '-4m2 010':( 1, 3, 0,-1),'-4m2 001':( 1, 2, 0,-1),' 6/mmm ':( 1, 9, 0,-1),'-6m2 100':( 1, 9, 0,-1), |
---|
437 | '-6m2 120':( 1, 9, 0,-1),' 6mm ':(14, 9, 0,-1),' 622 ':( 1, 9, 0,-1),' 6/m ':( 1, 9,14,-1), |
---|
438 | ' -6 ':( 1, 9,14, 0),' 6 ':(14, 9,14, 0),'-3m(100)':( 1, 9, 0,-1),'-3m(120)':( 1, 9, 0,-1), |
---|
439 | ' 3m(100)':(14, 9, 0,14),' 3m(120)':(14, 9, 0,14),' 32(100)':( 1, 9, 0,14),' 32(120)':( 1, 9, 0,14), |
---|
440 | ' -3 ':( 1, 9,14, 0),' 3 ':(14, 9,14, 0),'mmm(100)':( 1,14, 0,-1),'mmm(010)':( 1,15, 0,-1), |
---|
441 | 'mmm(110)':( 1,11, 0,-1),' mm2z100':(14,14, 0,-1),' mm2z010':(14,15, 0,-1),' mm2z110':(14,11, 0,-1), |
---|
442 | 'mm2(100)':(12,14, 0,-1),'mm2(010)':(13,15, 0,-1),'mm2(110)':( 6,11, 0,-1),'mm2(120)':(15,14, 0,-1), |
---|
443 | 'mm2(210)':(16,15, 0,-1),'mm2(+-0)':( 7,11, 0,-1),'222(100)':( 1,14, 0,-1),'222(010)':( 1,15, 0,-1), |
---|
444 | '222(110)':( 1,11, 0,-1),'2/m(001)':( 1,16,14,-1),'2/m(100)':( 1,25,12,-1),'2/m(010)':( 1,28,13,-1), |
---|
445 | '2/m(110)':( 1,19, 6,-1),'2/m(120)':( 1,27,15,-1),'2/m(210)':( 1,26,16,-1),'2/m(+-0)':( 1,20,17,-1), |
---|
446 | ' m(001) ':(23,16,14,23),' m(100) ':(26,25,12,26),' m(010) ':(27,28,13,27),' m(110) ':(18,19, 6,18), |
---|
447 | ' m(120) ':(24,27,15,24),' m(210) ':(25,26,16,25),' m(+-0) ':(17,20, 7,17),' 2(001) ':(14,16,14,23), |
---|
448 | ' 2(100) ':(12,25,12,26),' 2(010) ':(13,28,13,27),' 2(110) ':( 6,19, 6,18),' 2(120) ':(15,27,15,24), |
---|
449 | ' 2(210) ':(16,26,16,25),' 2(+-0) ':( 7,20, 7,17),' -1 ':( 1,29,28, 0) |
---|
450 | } |
---|
451 | return NXUPQsym[siteSym] |
---|
452 | |
---|
453 | def GetCSxinel(siteSym): |
---|
454 | CSxinel = [[], # 0th empty - indices are Fortran style |
---|
455 | [[0,0,0],[ 0.0, 0.0, 0.0]], # 0 0 0 |
---|
456 | [[1,1,1],[ 1.0, 1.0, 1.0]], # X X X |
---|
457 | [[1,1,1],[ 1.0, 1.0,-1.0]], # X X -X |
---|
458 | [[1,1,1],[ 1.0,-1.0, 1.0]], # X -X X |
---|
459 | [[1,1,1],[ 1.0,-1.0,-1.0]], # -X X X |
---|
460 | [[1,1,0],[ 1.0, 1.0, 0.0]], # X X 0 |
---|
461 | [[1,1,0],[ 1.0,-1.0, 0.0]], # X -X 0 |
---|
462 | [[1,0,1],[ 1.0, 0.0, 1.0]], # X 0 X |
---|
463 | [[1,0,1],[ 1.0, 0.0,-1.0]], # X 0 -X |
---|
464 | [[0,1,1],[ 0.0, 1.0, 1.0]], # 0 Y Y |
---|
465 | [[0,1,1],[ 0.0, 1.0,-1.0]], # 0 Y -Y |
---|
466 | [[1,0,0],[ 1.0, 0.0, 0.0]], # X 0 0 |
---|
467 | [[0,1,0],[ 0.0, 1.0, 0.0]], # 0 Y 0 |
---|
468 | [[0,0,1],[ 0.0, 0.0, 1.0]], # 0 0 Z |
---|
469 | [[1,1,0],[ 1.0, 2.0, 0.0]], # X 2X 0 |
---|
470 | [[1,1,0],[ 2.0, 1.0, 0.0]], # 2X X 0 |
---|
471 | [[1,1,2],[ 1.0, 1.0, 1.0]], # X X Z |
---|
472 | [[1,1,2],[ 1.0,-1.0, 1.0]], # X -X Z |
---|
473 | [[1,2,1],[ 1.0, 1.0, 1.0]], # X Y X |
---|
474 | [[1,2,1],[ 1.0, 1.0,-1.0]], # X Y -X |
---|
475 | [[1,2,2],[ 1.0, 1.0, 1.0]], # X Y Y |
---|
476 | [[1,2,2],[ 1.0, 1.0,-1.0]], # X Y -Y |
---|
477 | [[1,2,0],[ 1.0, 1.0, 0.0]], # X Y 0 |
---|
478 | [[1,0,2],[ 1.0, 0.0, 1.0]], # X 0 Z |
---|
479 | [[0,1,2],[ 0.0, 1.0, 1.0]], # 0 Y Z |
---|
480 | [[1,1,2],[ 1.0, 2.0, 1.0]], # X 2X Z |
---|
481 | [[1,1,2],[ 2.0, 1.0, 1.0]], # 2X X Z |
---|
482 | [[1,2,3],[ 1.0, 1.0, 1.0]], # X Y Z |
---|
483 | ] |
---|
484 | indx = GetNXUPQsym(siteSym) |
---|
485 | return CSxinel[indx[0]] |
---|
486 | |
---|
487 | def GetCSuinel(siteSym): |
---|
488 | # returns Uij terms, multipliers, GUI flags & Uiso2Uij multipliers |
---|
489 | CSuinel = [[], # 0th empty - indices are Fortran style |
---|
490 | [[1,1,1,0,0,0],[ 1.0, 1.0, 1.0, 0.0, 0.0, 0.0],[1,0,0,0,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A A 0 0 0 |
---|
491 | [[1,1,2,0,0,0],[ 1.0, 1.0, 1.0, 0.0, 0.0, 0.0],[1,0,1,0,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A C 0 0 0 |
---|
492 | [[1,2,1,0,0,0],[ 1.0, 1.0, 1.0, 0.0, 0.0, 0.0],[1,1,0,0,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B A 0 0 0 |
---|
493 | [[1,2,2,0,0,0],[ 1.0, 1.0, 1.0, 0.0, 0.0, 0.0],[1,1,0,0,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B B 0 0 0 |
---|
494 | [[1,1,1,2,2,2],[ 1.0, 1.0, 1.0, 1.0, 1.0, 1.0],[1,0,0,1,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A A D D D |
---|
495 | [[1,1,1,2,2,2],[ 1.0, 1.0, 1.0, 1.0,-1.0,-1.0],[1,0,0,1,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A A D -D -D |
---|
496 | [[1,1,1,2,2,2],[ 1.0, 1.0, 1.0, 1.0,-1.0, 1.0],[1,0,0,1,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A A D -D D |
---|
497 | [[1,1,1,2,2,2],[ 1.0, 1.0, 1.0, 1.0, 1.0,-1.0],[1,0,0,1,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A A D D -D |
---|
498 | [[1,1,2,1,0,0],[ 1.0, 1.0, 1.0, 0.5, 0.0, 0.0],[1,0,1,0,0,0],[1.0,1.0,1.0,0.5,0.0,0.0]], # A A C A/2 0 0 |
---|
499 | [[1,2,3,0,0,0],[ 1.0, 1.0, 1.0, 0.0, 0.0, 0.0],[1,1,1,0,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B C 0 0 0 |
---|
500 | [[1,1,2,3,0,0],[ 1.0, 1.0, 1.0, 1.0, 0.0, 0.0],[1,0,1,1,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A C D 0 0 |
---|
501 | [[1,2,1,0,3,0],[ 1.0, 1.0, 1.0, 0.0, 1.0, 0.0],[1,1,0,0,1,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B A 0 E 0 |
---|
502 | [[1,2,2,0,0,3],[ 1.0, 1.0, 1.0, 0.0, 0.0, 1.0],[1,1,0,0,0,1],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B B 0 0 F |
---|
503 | [[1,2,3,2,0,0],[ 1.0, 1.0, 1.0, 0.5, 0.0, 0.0],[1,1,1,0,0,0],[1.0,1.0,1.0,0.0,0.5,0.0]], # A B C B/2 0 0 |
---|
504 | [[1,2,3,1,0,0],[ 1.0, 1.0, 1.0, 0.5, 0.0, 0.0],[1,1,1,0,0,0],[1.0,1.0,1.0,0.0,0.5,0.0]], # A B C A/2 0 0 |
---|
505 | [[1,2,3,4,0,0],[ 1.0, 1.0, 1.0, 1.0, 0.0, 0.0],[1,1,1,1,0,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B C D 0 0 |
---|
506 | [[1,2,3,0,4,0],[ 1.0, 1.0, 1.0, 0.0, 1.0, 0.0],[1,1,1,0,1,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B C 0 E 0 |
---|
507 | [[1,2,3,0,0,4],[ 1.0, 1.0, 1.0, 0.0, 0.0, 1.0],[1,1,1,0,0,1],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B C 0 0 F |
---|
508 | [[1,1,2,3,4,4],[ 1.0, 1.0, 1.0, 1.0, 1.0,-1.0],[1,0,1,1,1,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A C D E -E |
---|
509 | [[1,1,2,3,4,4],[ 1.0, 1.0, 1.0, 1.0, 1.0, 1.0],[1,0,1,1,1,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A A C D E E |
---|
510 | [[1,2,1,3,4,3],[ 1.0, 1.0, 1.0, 1.0, 1.0,-1.0],[1,1,0,1,1,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B A D E -D |
---|
511 | [[1,2,1,3,4,3],[ 1.0, 1.0, 1.0, 1.0, 1.0, 1.0],[1,1,0,1,1,0],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B A D E D |
---|
512 | [[1,2,2,3,3,4],[ 1.0, 1.0, 1.0, 1.0,-1.0, 1.0],[1,1,0,1,0,1],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B B D -D F |
---|
513 | [[1,2,2,3,3,4],[ 1.0, 1.0, 1.0, 1.0, 1.0, 1.0],[1,1,0,1,0,1],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B B D D F |
---|
514 | [[1,2,3,2,4,4],[ 1.0, 1.0, 1.0, 0.5, 0.5, 1.0],[1,1,1,0,0,1],[1.0,1.0,1.0,0.5,0.0,0.0]], # A B C B/2 F/2 F |
---|
515 | [[1,2,3,1,0,4],[ 1.0, 1.0, 1.0, 0.5, 0.0, 1.0],[1,1,1,0,0,1],[1.0,1.0,1.0,0.5,0.0,0.0]], # A B C A/2 0 F |
---|
516 | [[1,2,3,2,4,0],[ 1.0, 1.0, 1.0, 0.5, 1.0, 0.0],[1,1,1,0,1,0],[1.0,1.0,1.0,0.5,0.0,0.0]], # A B C B/2 E 0 |
---|
517 | [[1,2,3,1,4,4],[ 1.0, 1.0, 1.0, 0.5, 1.0, 0.5],[1,1,1,0,1,0],[1.0,1.0,1.0,0.5,0.0,0.0]], # A B C A/2 E E/2 |
---|
518 | [[1,2,3,4,5,6],[ 1.0, 1.0, 1.0, 1.0, 1.0, 1.0],[1,1,1,1,1,1],[1.0,1.0,1.0,0.0,0.0,0.0]], # A B C D E F |
---|
519 | ] |
---|
520 | indx = GetNXUPQsym(siteSym) |
---|
521 | return CSuinel[indx[1]] |
---|
522 | |
---|
523 | def MustrainNames(SGData): |
---|
524 | laue = SGData['SGLaue'] |
---|
525 | uniq = SGData['SGUniq'] |
---|
526 | if laue in ['m3','m3m']: |
---|
527 | return ['S400','S220'] |
---|
528 | elif laue in ['6/m','6/mmm','3m1']: |
---|
529 | return ['S400','S004','S202'] |
---|
530 | elif laue in ['31m','3']: |
---|
531 | return ['S400','S004','S202','S211'] |
---|
532 | elif laue in ['3R','3mR']: |
---|
533 | return ['S400','S220','S310','S211'] |
---|
534 | elif laue in ['4/m','4/mmm']: |
---|
535 | return ['S400','S004','S220','S022'] |
---|
536 | elif laue in ['mmm']: |
---|
537 | return ['S400','S040','S004','S220','S202','S022'] |
---|
538 | elif laue in ['2/m']: |
---|
539 | SHKL = ['S400','S040','S004','S220','S202','S022'] |
---|
540 | if uniq == 'a': |
---|
541 | SHKL += ['S013','S031','S211'] |
---|
542 | elif uniq == 'b': |
---|
543 | SHKL += ['S301','S103','S121'] |
---|
544 | elif uniq == 'c': |
---|
545 | SHKL += ['S130','S310','S112'] |
---|
546 | return SHKL |
---|
547 | else: |
---|
548 | SHKL = ['S400','S040','S004','S220','S202','S022'] |
---|
549 | SHKL += ['S310','S103','S031','S130','S301','S013'] |
---|
550 | SHKL += ['S211','S121','S112'] |
---|
551 | return SHKL |
---|
552 | |
---|
553 | def HStrainNames(SGData): |
---|
554 | laue = SGData['SGLaue'] |
---|
555 | uniq = SGData['SGUniq'] |
---|
556 | if laue in ['m3','m3m']: |
---|
557 | return ['D11','eA'] #add cubic strain term |
---|
558 | elif laue in ['6/m','6/mmm','3m1','31m','3']: |
---|
559 | return ['D11','D33'] |
---|
560 | elif laue in ['3R','3mR']: |
---|
561 | return ['D11','D12'] |
---|
562 | elif laue in ['4/m','4/mmm']: |
---|
563 | return ['D11','D33'] |
---|
564 | elif laue in ['mmm']: |
---|
565 | return ['D11','D22','D33'] |
---|
566 | elif laue in ['2/m']: |
---|
567 | Dij = ['D11','D22','D33'] |
---|
568 | if uniq == 'a': |
---|
569 | Dij += ['D23'] |
---|
570 | elif uniq == 'b': |
---|
571 | Dij += ['D13'] |
---|
572 | elif uniq == 'c': |
---|
573 | Dij += ['D12'] |
---|
574 | return Dij |
---|
575 | else: |
---|
576 | Dij = ['D11','D22','D33','D12','D13','D23'] |
---|
577 | return Dij |
---|
578 | |
---|
579 | def MustrainCoeff(HKL,SGData): |
---|
580 | #NB: order of terms is the same as returned by MustrainNames |
---|
581 | laue = SGData['SGLaue'] |
---|
582 | uniq = SGData['SGUniq'] |
---|
583 | h,k,l = HKL |
---|
584 | Strm = [] |
---|
585 | if laue in ['m3','m3m']: |
---|
586 | Strm.append(h**4+k**4+l**4) |
---|
587 | Strm.append(3.0*((h*k)**2+(h*l)**2+(k*l)**2)) |
---|
588 | elif laue in ['6/m','6/mmm','3m1']: |
---|
589 | Strm.append(h**4+k**4+2.0*k*h**3+2.0*h*k**3+3.0*(h*k)**2) |
---|
590 | Strm.append(l**4) |
---|
591 | Strm.append(3.0*((h*l)**2+(k*l)**2+h*k*l**2)) |
---|
592 | elif laue in ['31m','3']: |
---|
593 | Strm.append(h**4+k**4+2.0*k*h**3+2.0*h*k**3+3.0*(h*k)**2) |
---|
594 | Strm.append(l**4) |
---|
595 | Strm.append(3.0*((h*l)**2+(k*l)**2+h*k*l**2)) |
---|
596 | Strm.append(4.0*h*k*l*(h+k)) |
---|
597 | elif laue in ['3R','3mR']: |
---|
598 | Strm.append(h**4+k**4+l**4) |
---|
599 | Strm.append(3.0*((h*k)**2+(h*l)**2+(k*l)**2)) |
---|
600 | Strm.append(2.0*(h*l**3+l*k**3+k*h**3)+2.0*(l*h**3+k*l**3+l*k**3)) |
---|
601 | Strm.append(4.0*(k*l*h**2+h*l*k**2+h*k*l**2)) |
---|
602 | elif laue in ['4/m','4/mmm']: |
---|
603 | Strm.append(h**4+k**4) |
---|
604 | Strm.append(l**4) |
---|
605 | Strm.append(3.0*(h*k)**2) |
---|
606 | Strm.append(3.0*((h*l)**2+(k*l)**2)) |
---|
607 | elif laue in ['mmm']: |
---|
608 | Strm.append(h**4) |
---|
609 | Strm.append(k**4) |
---|
610 | Strm.append(l**4) |
---|
611 | Strm.append(3.0*(h*k)**2) |
---|
612 | Strm.append(3.0*(h*l)**2) |
---|
613 | Strm.append(3.0*(k*l)**2) |
---|
614 | elif laue in ['2/m']: |
---|
615 | Strm.append(h**4) |
---|
616 | Strm.append(k**4) |
---|
617 | Strm.append(l**4) |
---|
618 | Strm.append(3.0*(h*k)**2) |
---|
619 | Strm.append(3.0*(h*l)**2) |
---|
620 | Strm.append(3.0*(k*l)**2) |
---|
621 | if uniq == 'a': |
---|
622 | Strm.append(2.0*k*l**3) |
---|
623 | Strm.append(2.0*l*k**3) |
---|
624 | Strm.append(4.0*k*l*h**2) |
---|
625 | elif uniq == 'b': |
---|
626 | Strm.append(2.0*l*h**3) |
---|
627 | Strm.append(2.0*h*l**3) |
---|
628 | Strm.append(4.0*h*l*k**2) |
---|
629 | elif uniq == 'c': |
---|
630 | Strm.append(2.0*h*k**3) |
---|
631 | Strm.append(2.0*k*h**3) |
---|
632 | Strm.append(4.0*h*k*l**2) |
---|
633 | else: |
---|
634 | Strm.append(h**4) |
---|
635 | Strm.append(k**4) |
---|
636 | Strm.append(l**4) |
---|
637 | Strm.append(3.0*(h*k)**2) |
---|
638 | Strm.append(3.0*(h*l)**2) |
---|
639 | Strm.append(3.0*(k*l)**2) |
---|
640 | Strm.append(2.0*k*h**3) |
---|
641 | Strm.append(2.0*h*l**3) |
---|
642 | Strm.append(2.0*l*k**3) |
---|
643 | Strm.append(2.0*h*k**3) |
---|
644 | Strm.append(2.0*l*h**3) |
---|
645 | Strm.append(2.0*k*l**3) |
---|
646 | Strm.append(4.0*k*l*h**2) |
---|
647 | Strm.append(4.0*h*l*k**2) |
---|
648 | Strm.append(4.0*k*h*l**2) |
---|
649 | return Strm |
---|
650 | |
---|
651 | def Muiso2Shkl(muiso,SGData,cell): |
---|
652 | #this is to convert isotropic mustrain to generalized Shkls - doesn't work just now |
---|
653 | import GSASIIlattice as G2lat |
---|
654 | from scipy.optimize import fmin |
---|
655 | A = G2lat.cell2AB(cell)[0] |
---|
656 | def minMus(Shkl,H,muiso,SGData,A): |
---|
657 | U = np.inner(A.T,H) |
---|
658 | S = np.array(MustrainCoeff(H.T,SGData)) |
---|
659 | sum = np.sqrt(np.sum(np.multiply(S,Shkl))) |
---|
660 | return abs(muiso-sum*H) |
---|
661 | laue = SGData['SGLaue'] |
---|
662 | if laue in ['m3','m3m']: |
---|
663 | H = [[1,0,0],[1,1,0]] |
---|
664 | S0 = [0.01,0.01] |
---|
665 | elif laue in ['6/m','6/mmm','3m1']: |
---|
666 | H = [[1,0,0],[0,0,1],[1,0,1]] |
---|
667 | S0 = [0.01,0.01,0.01] |
---|
668 | elif laue in ['31m','3']: |
---|
669 | H = [[1,0,0],[0,0,1],[1,0,1],[1,1,1]] |
---|
670 | S0 = [0.01,0.01,0.01,0.01] |
---|
671 | elif laue in ['3R','3mR']: |
---|
672 | H = [[1,0,0],[1,1,0],[1,0,1],[1,1,1]] |
---|
673 | S0 = [0.01,0.01,0.01,0.01] |
---|
674 | elif laue in ['4/m','4/mmm']: |
---|
675 | H = [[1,0,0],[0,0,1],[1,1,0],[1,0,1]] |
---|
676 | S0 = [0.01,0.01,0.01,0.01] |
---|
677 | elif laue in ['mmm']: |
---|
678 | H = [[1,0,0],[0,1,0],[0,0,1],[1,1,0],[1,0,1],[0,1,1]] |
---|
679 | S0 = [0.01,0.01,0.01,0.01,0.01,0.01] |
---|
680 | elif laue in ['2/m']: |
---|
681 | H = [[1,0,0],[0,1,0],[0,0,1],[1,1,0],[1,0,1],[0,1,1]] |
---|
682 | if uniq == 'a': |
---|
683 | H.append([0,1,-1]) |
---|
684 | H.append([0,-2,1]) |
---|
685 | elif uniq == 'b': |
---|
686 | H.append([1,0,-1]) |
---|
687 | H.append([-2,0,1]) |
---|
688 | elif uniq == 'c': |
---|
689 | H.append([1,-1,0]) |
---|
690 | H.append([-2,1,0]) |
---|
691 | H.append([1,1,1]) |
---|
692 | S0 = [9*[0.01,]] |
---|
693 | else: |
---|
694 | H = [[1,0,0],[0,1,0],[0,0,1],[1,1,0],[1,0,1],[0,1,1], |
---|
695 | [-1,1,0],[1,0,-1],[0,-1,1],[1,-2,0],[-2,0,1],[0,1,-2], |
---|
696 | [1,-1,1],[-1, 1, 1],[1,-1,1]] |
---|
697 | S0 = [15*[0.01,]] |
---|
698 | H = np.array(H) |
---|
699 | S0 = np.array(S0) |
---|
700 | return fmin(minMus,S0,(H,muiso,SGData,A)) |
---|
701 | |
---|
702 | def SytSym(XYZ,SGData): |
---|
703 | ''' |
---|
704 | Generates the number of equivalent positions and a site symmetry code for a specified coordinate and space group |
---|
705 | input: |
---|
706 | XYZ: an array, tuple or list containing 3 elements: x, y & z |
---|
707 | SGData: from SpcGroup |
---|
708 | Returns a two element tuple: |
---|
709 | The 1st element is a code for the site symmetry (see GetKNsym) |
---|
710 | The 2nd element is the site multiplicity |
---|
711 | ''' |
---|
712 | def PackRot(SGOps): |
---|
713 | IRT = [] |
---|
714 | for ops in SGOps: |
---|
715 | M = ops[0] |
---|
716 | irt = 0 |
---|
717 | for j in range(2,-1,-1): |
---|
718 | for k in range(2,-1,-1): |
---|
719 | irt *= 3 |
---|
720 | irt += M[k][j] |
---|
721 | IRT.append(int(irt)) |
---|
722 | return IRT |
---|
723 | |
---|
724 | SymName = '' |
---|
725 | Mult = 1 |
---|
726 | Isym = 0 |
---|
727 | if SGData['SGLaue'] in ['3','3m1','31m','6/m','6/mmm']: |
---|
728 | Isym = 1073741824 |
---|
729 | Jdup = 0 |
---|
730 | Xeqv = GenAtom(XYZ,SGData,True) |
---|
731 | IRT = PackRot(SGData['SGOps']) |
---|
732 | L = -1 |
---|
733 | for ic,cen in enumerate(SGData['SGCen']): |
---|
734 | for invers in range(int(SGData['SGInv']+1)): |
---|
735 | for io,ops in enumerate(SGData['SGOps']): |
---|
736 | irtx = (1-2*invers)*IRT[io] |
---|
737 | L += 1 |
---|
738 | if not Xeqv[L][1]: |
---|
739 | Jdup += 1 |
---|
740 | jx = GetOprPtrName(str(irtx)) |
---|
741 | if jx[2] < 39: |
---|
742 | Isym += 2**(jx[2]-1) |
---|
743 | if Isym == 1073741824: Isym = 0 |
---|
744 | Mult = len(SGData['SGOps'])*len(SGData['SGCen'])*(int(SGData['SGInv'])+1)/Jdup |
---|
745 | |
---|
746 | return GetKNsym(str(Isym)),Mult |
---|
747 | |
---|
748 | def ElemPosition(SGData): |
---|
749 | ''' Under development |
---|
750 | Object here is to return a list of symmetry element types and locations suitable |
---|
751 | for say drawing them. |
---|
752 | So far I have the element type... getting all possible locations without lookup may be impossible! |
---|
753 | ''' |
---|
754 | SymElements = [] |
---|
755 | Inv = SGData['SGInv'] |
---|
756 | Cen = SGData['SGCen'] |
---|
757 | eleSym = {-3:['','-1'],-2:['',-6],-1:['2','-4'],0:['3','-3'],1:['4','m'],2:['6',''],3:['1','']} |
---|
758 | # get operators & expand if centrosymmetric |
---|
759 | Ops = SGData['SGOps'] |
---|
760 | opM = np.array([op[0].T for op in Ops]) |
---|
761 | opT = np.array([op[1] for op in Ops]) |
---|
762 | if Inv: |
---|
763 | opM = np.concatenate((opM,-opM)) |
---|
764 | opT = np.concatenate((opT,-opT)) |
---|
765 | opMT = zip(opM,opT) |
---|
766 | for M,T in opMT[1:]: #skip I |
---|
767 | Dt = int(nl.det(M)) |
---|
768 | Tr = int(np.trace(M)) |
---|
769 | Dt = -(Dt-1)/2 |
---|
770 | Es = eleSym[Tr][Dt] |
---|
771 | if Dt: #rotation-inversion |
---|
772 | I = np.eye(3) |
---|
773 | if Tr == 1: #mirrors/glides |
---|
774 | if np.any(T): #glide |
---|
775 | M2 = np.inner(M,M) |
---|
776 | MT = np.inner(M,T)+T |
---|
777 | print 'glide',Es,MT |
---|
778 | print M2 |
---|
779 | else: #mirror |
---|
780 | print 'mirror',Es,T |
---|
781 | print I-M |
---|
782 | X = [-1,-1,-1] |
---|
783 | elif Tr == -3: # pure inversion |
---|
784 | X = np.inner(nl.inv(I-M),T) |
---|
785 | print 'inversion',Es,X |
---|
786 | else: #other rotation-inversion |
---|
787 | M2 = np.inner(M,M) |
---|
788 | MT = np.inner(M,T)+T |
---|
789 | print 'rot-inv',Es,MT |
---|
790 | print M2 |
---|
791 | X = [-1,-1,-1] |
---|
792 | |
---|
793 | |
---|
794 | |
---|
795 | else: #rotations |
---|
796 | print 'rotation',Es |
---|
797 | X = [-1,-1,-1] |
---|
798 | #SymElements.append([Es,X]) |
---|
799 | |
---|
800 | return #SymElements |
---|
801 | |
---|
802 | def ApplyStringOps(A,SGData,X,Uij=[]): |
---|
803 | SGOps = SGData['SGOps'] |
---|
804 | SGCen = SGData['SGCen'] |
---|
805 | Ax = A.split('+') |
---|
806 | Ax[0] = int(Ax[0]) |
---|
807 | iC = 0 |
---|
808 | if Ax[0] < 0: |
---|
809 | iC = 1 |
---|
810 | Ax[0] = abs(Ax[0]) |
---|
811 | nA = Ax[0]%100-1 |
---|
812 | cA = Ax[0]/100 |
---|
813 | Cen = SGCen[cA] |
---|
814 | M,T = SGOps[nA] |
---|
815 | if len(Ax)>1: |
---|
816 | cellA = Ax[1].split(',') |
---|
817 | cellA = np.array([int(a) for a in cellA]) |
---|
818 | else: |
---|
819 | cellA = np.zeros(3) |
---|
820 | newX = (1-2*iC)*(Cen+np.inner(M,X)+T)+cellA |
---|
821 | if len(Uij): |
---|
822 | U = Uij2U(Uij) |
---|
823 | U = np.inner(M,np.inner(U,M).T) |
---|
824 | newUij = U2Uij(U) |
---|
825 | return [newX,newUij] |
---|
826 | else: |
---|
827 | return newX |
---|
828 | |
---|
829 | def StringOpsProd(A,B,SGData): |
---|
830 | ''' Find A*B where A & B are in strings '-' + '100*c+n' + '+ijk' |
---|
831 | where '-' indicates inversion, c(>0) is the cell centering operator, |
---|
832 | n is operator number from SgOps and ijk are unit cell translations (each may be <0). |
---|
833 | Should return resultant string - C. |
---|
834 | SGData - dictionary using entries: |
---|
835 | 'SGCen': cell centering vectors [0,0,0] at least |
---|
836 | 'SGOps': symmetry operations as [M,T] so that M*x+T = x' |
---|
837 | ''' |
---|
838 | SGOps = SGData['SGOps'] |
---|
839 | SGCen = SGData['SGCen'] |
---|
840 | #1st split out the cell translation part & work on the operator parts |
---|
841 | Ax = A.split('+'); Bx = B.split('+') |
---|
842 | Ax[0] = int(Ax[0]); Bx[0] = int(Bx[0]) |
---|
843 | iC = 0 |
---|
844 | if Ax[0]*Bx[0] < 0: |
---|
845 | iC = 1 |
---|
846 | Ax[0] = abs(Ax[0]); Bx[0] = abs(Bx[0]) |
---|
847 | nA = Ax[0]%100-1; nB = Bx[0]%100-1 |
---|
848 | cA = Ax[0]/100; cB = Bx[0]/100 |
---|
849 | Cen = (SGCen[cA]+SGCen[cB])%1.0 |
---|
850 | cC = np.nonzero([np.allclose(C,Cen) for C in SGCen])[0][0] |
---|
851 | Ma,Ta = SGOps[nA]; Mb,Tb = SGOps[nB] |
---|
852 | Mc = np.inner(Ma,Mb.T) |
---|
853 | # print Ma,Mb,Mc |
---|
854 | Tc = (np.add(np.inner(Mb,Ta)+1.,Tb))%1.0 |
---|
855 | # print Ta,Tb,Tc |
---|
856 | # print [np.allclose(M,Mc)&np.allclose(T,Tc) for M,T in SGOps] |
---|
857 | nC = np.nonzero([np.allclose(M,Mc)&np.allclose(T,Tc) for M,T in SGOps])[0][0] |
---|
858 | #now the cell translation part |
---|
859 | if len(Ax)>1: |
---|
860 | cellA = Ax[1].split(',') |
---|
861 | cellA = [int(a) for a in cellA] |
---|
862 | else: |
---|
863 | cellA = [0,0,0] |
---|
864 | if len(Bx)>1: |
---|
865 | cellB = Bx[1].split(',') |
---|
866 | cellB = [int(b) for b in cellB] |
---|
867 | else: |
---|
868 | cellB = [0,0,0] |
---|
869 | cellC = np.add(cellA,cellB) |
---|
870 | C = str(((nC+1)+(100*cC))*(1-2*iC))+'+'+ \ |
---|
871 | str(int(cellC[0]))+','+str(int(cellC[1]))+','+str(int(cellC[2])) |
---|
872 | return C |
---|
873 | |
---|
874 | def U2Uij(U): |
---|
875 | #returns the UIJ vector U11,U22,U33,U12,U13,U23 from tensor U |
---|
876 | return [U[0][0],U[1][1],U[2][2],U[0][1],U[0][2],U[1][2]] |
---|
877 | |
---|
878 | def Uij2U(Uij): |
---|
879 | #returns the thermal motion tensor U from Uij as numpy array |
---|
880 | return np.array([[Uij[0],Uij[3],Uij[4]],[Uij[3],Uij[1],Uij[5]],[Uij[4],Uij[5],Uij[2]]]) |
---|
881 | |
---|
882 | '''A list of space groups as ordered and named in the pre-2002 International |
---|
883 | Tables Volume A, except that spaces are used following the GSAS convention to |
---|
884 | separate the different crystallographic directions. |
---|
885 | Note that the symmetry codes here will recognize many non-standard space group |
---|
886 | symbols with different settings. |
---|
887 | ''' |
---|
888 | spglist = { |
---|
889 | 'triclinic' : ('P 1','P -1',), # 1-2 |
---|
890 | 'monoclinic': ('P 2','P 21','C 2','P m','P c','C m','C c','P 2/m','P 21/m', |
---|
891 | 'C 2/m','P 2/c','P 21/c','C 2/c',), #3-15 |
---|
892 | 'orthorhombic': ('P 2 2 2','P 2 2 21','P 21 21 2','P 21 21 21','C 2 2 21', |
---|
893 | 'C 2 2 2','F 2 2 2','I 2 2 2','I 21 21 21', |
---|
894 | 'P m m 2','P m c 21','P c c 2','P m a 2','P c a 21', |
---|
895 | 'P n c 2','P m n 21','P b a 2','P n a 21','P n n 2', |
---|
896 | 'C m m 2','C m c 21','C c c 2','A m m 2','A b m 2', |
---|
897 | 'A m a 2','A b a 2','F m m 2','F d d 2','I m m 2', |
---|
898 | 'I b a 2','I m a 2','P m m m','P n n n','P c c m', |
---|
899 | 'P b a n','P m m a','P n n a','P m n a','P c c a', |
---|
900 | 'P b a m','P c c n','P b c m','P n n m','P m m n', |
---|
901 | 'P b c n','P b c a','P n m a','C m c m','C m c a', |
---|
902 | 'C m m m','C c c m','C m m a','C c c a','F m m m', |
---|
903 | 'F d d d','I m m m','I b a m','I b c a','I m m a',), #16-74 |
---|
904 | 'tetragonal': ('P 4','P 41','P 42','P 43','I 4','I 41','P -4','I -4', |
---|
905 | 'P 4/m','P 42/m','P 4/n','P 42/n','I 4/m','I 41/a', |
---|
906 | 'P 4 2 2','P 4 21 2','P 41 2 2','P 41 21 2','P 42 2 2', |
---|
907 | 'P 42 21 2','P 43 2 2','P 43 21 2','I 4 2 2','I 41 2 2', |
---|
908 | 'P 4 m m','P 4 b m','P 42 c m','P 42 n m','P 4 c c', |
---|
909 | 'P 4 n c','P 42 m c','P 42 b c','I 4 m m','I 4 c m', |
---|
910 | 'I 41 m d','I 41 c d','P -4 2 m','P -4 2 c','P -4 21 m', |
---|
911 | 'P -4 21 c','P -4 m 2','P -4 c 2','P -4 b 2','P -4 n 2', |
---|
912 | 'I -4 m 2','I -4 c 2','I -4 2 m','I -4 2 d','P 4/m m m', |
---|
913 | 'P 4/m c c','P 4/n b m','P 4/n n c','P 4/m b m','P 4/m n c', |
---|
914 | 'P 4/n m m','P 4/n c c','P 42/m m c','P 42/m c m', |
---|
915 | 'P 42/n b c','P 42/n n m','P 42/m b c','P 42/m n m', |
---|
916 | 'P 42/n m c','P 42/n c m','I 4/m m m','I 4/m c m', |
---|
917 | 'I 41/a m d','I 41/a c d',), # 75-142 |
---|
918 | 'trigonal': ('P 3','P 31','P 32','R 3','P -3','R -3','P 3 1 2','P 3 2 1', |
---|
919 | 'P 31 1 2','P 31 2 1','P 32 1 2','P 32 2 1','R 3 2', 'P 3 m 1', |
---|
920 | 'P 3 1 m','P 3 c 1','P 3 1 c','R 3 m','R 3 c','P -3 1 m', |
---|
921 | 'P -3 1 c','P -3 m 1','P -3 c 1','R -3 m','R -3 c',), #143-167 |
---|
922 | 'hexagonal': ('P 6','P 61','P 65','P 62','P 64','P 63','P -6','P 6/m', |
---|
923 | 'P 63/m','P 6 2 2','P 61 2 2','P 65 2 2','P 62 2 2', |
---|
924 | 'P 64 2 2','P 63 2 2','P 6 m m','P 6 c c','P 63 c m', |
---|
925 | 'P 63 m c','P -6 m 2','P -6 c 2','P -6 2 m','P -6 2 c', |
---|
926 | 'P 6/m m m','P 6/m c c','P 63/m c m','P 63/m m c',), #144-194 |
---|
927 | 'cubic': ('P 2 3','F 2 3','I 2 3','P 21 3','I 21 3','P m -3','P n -3', |
---|
928 | 'F m -3','F d -3','I m -3','P a -3','I a -3','P 4 3 2','P 42 3 2', |
---|
929 | 'F 4 3 2','F 41 3 2','I 4 3 2','P 43 3 2','P 41 3 2','I 41 3 2', |
---|
930 | 'P -4 3 m','F -4 3 m','I -4 3 m','P -4 3 n','F -4 3 c','I -4 3 d', |
---|
931 | 'P m -3 m','P n -3 n','P m -3 n','P n -3 m','F m -3 m','F m -3 c', |
---|
932 | 'F d -3 m','F d -3 c','I m -3 m','I a -3 d',), #195-230 |
---|
933 | } |
---|
934 | 'A few non-standard space groups for test use' |
---|
935 | nonstandard_sglist = ('P 21 1 1','P 1 21 1','P 1 1 21','R 3 r','R 3 2 h', |
---|
936 | 'R -3 r', 'R 3 2 r','R 3 m h', 'R 3 m r', |
---|
937 | 'R 3 c r','R -3 c r','R -3 m r',), |
---|
938 | '''A list of orthorhombic space groups that were renamed in the 2002 Volume A, |
---|
939 | along with the pre-2002 name. The e designates a double glide-plane''' |
---|
940 | sgequiv_2002_orthorhombic= (('A e m 2', 'A b m 2',), |
---|
941 | ('A e a 2', 'A b a 2',), |
---|
942 | ('C m c e', 'C m c a',), |
---|
943 | ('C m m e', 'C m m a',), |
---|
944 | ('C c c e', 'C c c a'),) |
---|
945 | '''Use the space groups types in this order to list the symbols in the |
---|
946 | order they are listed in the International Tables, vol. A''' |
---|
947 | symtypelist = ('triclinic', 'monoclinic', 'orthorhombic', 'tetragonal', |
---|
948 | 'trigonal', 'hexagonal', 'cubic') |
---|
949 | |
---|
950 | # self-test materials follow. Requires files in directory testinp |
---|
951 | def test0(): |
---|
952 | '''test #0: exercise MoveToUnitCell''' |
---|
953 | msg = "MoveToUnitCell failed" |
---|
954 | assert (MoveToUnitCell([1,2,3]) == [0,0,0]).all, msg |
---|
955 | assert (MoveToUnitCell([2,-1,-2]) == [0,0,0]).all, msg |
---|
956 | assert abs(MoveToUnitCell(np.array([-.1]))[0]-0.9) < 1e-6, msg |
---|
957 | assert abs(MoveToUnitCell(np.array([.1]))[0]-0.1) < 1e-6, msg |
---|
958 | |
---|
959 | def test1(): |
---|
960 | ''' test #1: SpcGroup and SGPrint against previous results''' |
---|
961 | testdir = ospath.join(ospath.split(ospath.abspath( __file__ ))[0],'testinp') |
---|
962 | if ospath.exists(testdir): |
---|
963 | if testdir not in sys.path: sys.path.insert(0,testdir) |
---|
964 | import spctestinp |
---|
965 | def CompareSpcGroup(spc, referr, refdict, reflist): |
---|
966 | 'Compare output from GSASIIspc.SpcGroup with results from a previous run' |
---|
967 | # if an error is reported, the dictionary can be ignored |
---|
968 | msg0 = "CompareSpcGroup failed on space group %s" % spc |
---|
969 | result = SpcGroup(spc) |
---|
970 | if result[0] == referr and referr > 0: return True |
---|
971 | keys = result[1].keys() |
---|
972 | #print result[1]['SpGrp'] |
---|
973 | msg = msg0 + " in list lengths" |
---|
974 | assert len(keys) == len(refdict.keys()), msg |
---|
975 | for key in keys: |
---|
976 | if key == 'SGOps' or key == 'SGCen': |
---|
977 | msg = msg0 + (" in key %s length" % key) |
---|
978 | assert len(refdict[key]) == len(result[1][key]), msg |
---|
979 | for i in range(len(refdict[key])): |
---|
980 | msg = msg0 + (" in key %s level 0" % key) |
---|
981 | assert np.allclose(result[1][key][i][0],refdict[key][i][0]), msg |
---|
982 | msg = msg0 + (" in key %s level 1" % key) |
---|
983 | assert np.allclose(result[1][key][i][1],refdict[key][i][1]), msg |
---|
984 | else: |
---|
985 | msg = msg0 + (" in key %s" % key) |
---|
986 | assert result[1][key] == refdict[key], msg |
---|
987 | msg = msg0 + (" in key %s reflist" % key) |
---|
988 | #for (l1,l2) in zip(reflist, SGPrint(result[1])): |
---|
989 | # assert l2.replace('\t','').replace(' ','') == l1.replace(' ',''), 'SGPrint ' +msg |
---|
990 | assert reflist == SGPrint(result[1]), 'SGPrint ' +msg |
---|
991 | for spc in spctestinp.SGdat: |
---|
992 | CompareSpcGroup(spc, 0, spctestinp.SGdat[spc], spctestinp.SGlist[spc] ) |
---|
993 | |
---|
994 | def test2(): |
---|
995 | ''' test #2: SpcGroup against cctbx (sgtbx) computations''' |
---|
996 | testdir = ospath.join(ospath.split(ospath.abspath( __file__ ))[0],'testinp') |
---|
997 | if ospath.exists(testdir): |
---|
998 | if testdir not in sys.path: sys.path.insert(0,testdir) |
---|
999 | import sgtbxtestinp |
---|
1000 | def CompareWcctbx(spcname, cctbx_in, debug=0): |
---|
1001 | 'Compare output from GSASIIspc.SpcGroup with results from cctbx.sgtbx' |
---|
1002 | cctbx = cctbx_in[:] # make copy so we don't delete from the original |
---|
1003 | spc = (SpcGroup(spcname))[1] |
---|
1004 | if debug: print spc['SpGrp'] |
---|
1005 | if debug: print spc['SGCen'] |
---|
1006 | latticetype = spcname.strip().upper()[0] |
---|
1007 | # lattice type of R implies Hexagonal centering", fix the rhombohedral settings |
---|
1008 | if latticetype == "R" and len(spc['SGCen']) == 1: latticetype = 'P' |
---|
1009 | assert latticetype == spc['SGLatt'], "Failed: %s does not match Lattice: %s" % (spcname, spc['SGLatt']) |
---|
1010 | onebar = [1] |
---|
1011 | if spc['SGInv']: onebar.append(-1) |
---|
1012 | for (op,off) in spc['SGOps']: |
---|
1013 | for inv in onebar: |
---|
1014 | for cen in spc['SGCen']: |
---|
1015 | noff = off + cen |
---|
1016 | noff = MoveToUnitCell(noff) |
---|
1017 | mult = tuple((op*inv).ravel().tolist()) |
---|
1018 | if debug: print "\n%s: %s + %s" % (spcname,mult,noff) |
---|
1019 | for refop in cctbx: |
---|
1020 | if debug: print refop |
---|
1021 | # check the transform |
---|
1022 | if refop[:9] != mult: continue |
---|
1023 | if debug: print "mult match" |
---|
1024 | # check the translation |
---|
1025 | reftrans = list(refop[-3:]) |
---|
1026 | reftrans = MoveToUnitCell(reftrans) |
---|
1027 | if all(abs(noff - reftrans) < 1.e-5): |
---|
1028 | cctbx.remove(refop) |
---|
1029 | break |
---|
1030 | else: |
---|
1031 | assert False, "failed on %s:\n\t %s + %s" % (spcname,mult,noff) |
---|
1032 | for key in sgtbxtestinp.sgtbx: |
---|
1033 | CompareWcctbx(key, sgtbxtestinp.sgtbx[key]) |
---|
1034 | |
---|
1035 | def test3(): |
---|
1036 | ''' test #3: exercise SytSym (includes GetOprPtrName, GenAtom, GetKNsym) |
---|
1037 | for selected space groups against info in IT Volume A ''' |
---|
1038 | def ExerciseSiteSym (spc, crdlist): |
---|
1039 | 'compare site symmetries and multiplicities for a specified space group' |
---|
1040 | msg = "failed on site sym test for %s" % spc |
---|
1041 | (E,S) = SpcGroup(spc) |
---|
1042 | assert not E, msg |
---|
1043 | for t in crdlist: |
---|
1044 | symb, m = SytSym(t[0],S) |
---|
1045 | if symb.strip() != t[2].strip() or m != t[1]: |
---|
1046 | print spc,t[0],m,symb |
---|
1047 | assert m == t[1] |
---|
1048 | #assert symb.strip() == t[2].strip() |
---|
1049 | |
---|
1050 | ExerciseSiteSym('p 1',[ |
---|
1051 | ((0.13,0.22,0.31),1,'1'), |
---|
1052 | ((0,0,0),1,'1'), |
---|
1053 | ]) |
---|
1054 | ExerciseSiteSym('p -1',[ |
---|
1055 | ((0.13,0.22,0.31),2,'1'), |
---|
1056 | ((0,0.5,0),1,'-1'), |
---|
1057 | ]) |
---|
1058 | ExerciseSiteSym('C 2/c',[ |
---|
1059 | ((0.13,0.22,0.31),8,'1'), |
---|
1060 | ((0.0,.31,0.25),4,'2(010)'), |
---|
1061 | ((0.25,.25,0.5),4,'-1'), |
---|
1062 | ((0,0.5,0),4,'-1'), |
---|
1063 | ]) |
---|
1064 | ExerciseSiteSym('p 2 2 2',[ |
---|
1065 | ((0.13,0.22,0.31),4,'1'), |
---|
1066 | ((0,0.5,.31),2,'2(001)'), |
---|
1067 | ((0.5,.31,0.5),2,'2(010)'), |
---|
1068 | ((.11,0,0),2,'2(100)'), |
---|
1069 | ((0,0.5,0),1,'222'), |
---|
1070 | ]) |
---|
1071 | ExerciseSiteSym('p 4/n',[ |
---|
1072 | ((0.13,0.22,0.31),8,'1'), |
---|
1073 | ((0.25,0.75,.31),4,'2(001)'), |
---|
1074 | ((0.5,0.5,0.5),4,'-1'), |
---|
1075 | ((0,0.5,0),4,'-1'), |
---|
1076 | ((0.25,0.25,.31),2,'4(001)'), |
---|
1077 | ((0.25,.75,0.5),2,'-4(001)'), |
---|
1078 | ((0.25,.75,0.0),2,'-4(001)'), |
---|
1079 | ]) |
---|
1080 | ExerciseSiteSym('p 31 2 1',[ |
---|
1081 | ((0.13,0.22,0.31),6,'1'), |
---|
1082 | ((0.13,0.0,0.833333333),3,'2(100)'), |
---|
1083 | ((0.13,0.13,0.),3,'2(110)'), |
---|
1084 | ]) |
---|
1085 | ExerciseSiteSym('R 3 c',[ |
---|
1086 | ((0.13,0.22,0.31),18,'1'), |
---|
1087 | ((0.0,0.0,0.31),6,'3'), |
---|
1088 | ]) |
---|
1089 | ExerciseSiteSym('R 3 c R',[ |
---|
1090 | ((0.13,0.22,0.31),6,'1'), |
---|
1091 | ((0.31,0.31,0.31),2,'3(111)'), |
---|
1092 | ]) |
---|
1093 | ExerciseSiteSym('P 63 m c',[ |
---|
1094 | ((0.13,0.22,0.31),12,'1'), |
---|
1095 | ((0.11,0.22,0.31),6,'m(100)'), |
---|
1096 | ((0.333333,0.6666667,0.31),2,'3m(100)'), |
---|
1097 | ((0,0,0.31),2,'3m(100)'), |
---|
1098 | ]) |
---|
1099 | ExerciseSiteSym('I a -3',[ |
---|
1100 | ((0.13,0.22,0.31),48,'1'), |
---|
1101 | ((0.11,0,0.25),24,'2(100)'), |
---|
1102 | ((0.11,0.11,0.11),16,'3(111)'), |
---|
1103 | ((0,0,0),8,'-3(111)'), |
---|
1104 | ]) |
---|
1105 | |
---|
1106 | if __name__ == '__main__': |
---|
1107 | test0() |
---|
1108 | test1() |
---|
1109 | test2() |
---|
1110 | test3() |
---|
1111 | print "OK" |
---|