1 | # -*- coding: utf-8 -*- |
---|
2 | ''' |
---|
3 | *GSASIIstrMath - structure math routines* |
---|
4 | ----------------------------------------- |
---|
5 | ''' |
---|
6 | ########### SVN repository information ################### |
---|
7 | # $Date: 2013-10-09 15:27:52 +0000 (Wed, 09 Oct 2013) $ |
---|
8 | # $Author: vondreele $ |
---|
9 | # $Revision: 1092 $ |
---|
10 | # $URL: trunk/GSASIIstrMath.py $ |
---|
11 | # $Id: GSASIIstrMath.py 1092 2013-10-09 15:27:52Z vondreele $ |
---|
12 | ########### SVN repository information ################### |
---|
13 | import time |
---|
14 | import math |
---|
15 | import copy |
---|
16 | import numpy as np |
---|
17 | import numpy.ma as ma |
---|
18 | import numpy.linalg as nl |
---|
19 | import scipy.optimize as so |
---|
20 | import scipy.stats as st |
---|
21 | import GSASIIpath |
---|
22 | GSASIIpath.SetVersionNumber("$Revision: 1092 $") |
---|
23 | import GSASIIElem as G2el |
---|
24 | import GSASIIlattice as G2lat |
---|
25 | import GSASIIspc as G2spc |
---|
26 | import GSASIIpwd as G2pwd |
---|
27 | import GSASIImapvars as G2mv |
---|
28 | import GSASIImath as G2mth |
---|
29 | |
---|
30 | sind = lambda x: np.sin(x*np.pi/180.) |
---|
31 | cosd = lambda x: np.cos(x*np.pi/180.) |
---|
32 | tand = lambda x: np.tan(x*np.pi/180.) |
---|
33 | asind = lambda x: 180.*np.arcsin(x)/np.pi |
---|
34 | acosd = lambda x: 180.*np.arccos(x)/np.pi |
---|
35 | atan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi |
---|
36 | |
---|
37 | ateln2 = 8.0*math.log(2.0) |
---|
38 | |
---|
39 | ################################################################################ |
---|
40 | ##### Rigid Body Models |
---|
41 | ################################################################################ |
---|
42 | |
---|
43 | def ApplyRBModels(parmDict,Phases,rigidbodyDict,Update=False): |
---|
44 | ''' Takes RB info from RBModels in Phase and RB data in rigidbodyDict along with |
---|
45 | current RB values in parmDict & modifies atom contents (xyz & Uij) of parmDict |
---|
46 | ''' |
---|
47 | atxIds = ['Ax:','Ay:','Az:'] |
---|
48 | atuIds = ['AU11:','AU22:','AU33:','AU12:','AU13:','AU23:'] |
---|
49 | RBIds = rigidbodyDict.get('RBIds',{'Vector':[],'Residue':[]}) #these are lists of rbIds |
---|
50 | if not RBIds['Vector'] and not RBIds['Residue']: |
---|
51 | return |
---|
52 | VRBIds = RBIds['Vector'] |
---|
53 | RRBIds = RBIds['Residue'] |
---|
54 | if Update: |
---|
55 | RBData = rigidbodyDict |
---|
56 | else: |
---|
57 | RBData = copy.deepcopy(rigidbodyDict) # don't mess with original! |
---|
58 | if RBIds['Vector']: # first update the vector magnitudes |
---|
59 | VRBData = RBData['Vector'] |
---|
60 | for i,rbId in enumerate(VRBIds): |
---|
61 | if VRBData[rbId]['useCount']: |
---|
62 | for j in range(len(VRBData[rbId]['VectMag'])): |
---|
63 | name = '::RBV;'+str(j)+':'+str(i) |
---|
64 | VRBData[rbId]['VectMag'][j] = parmDict[name] |
---|
65 | for phase in Phases: |
---|
66 | Phase = Phases[phase] |
---|
67 | General = Phase['General'] |
---|
68 | cell = General['Cell'][1:7] |
---|
69 | Amat,Bmat = G2lat.cell2AB(cell) |
---|
70 | AtLookup = G2mth.FillAtomLookUp(Phase['Atoms']) |
---|
71 | pfx = str(Phase['pId'])+'::' |
---|
72 | if Update: |
---|
73 | RBModels = Phase['RBModels'] |
---|
74 | else: |
---|
75 | RBModels = copy.deepcopy(Phase['RBModels']) # again don't mess with original! |
---|
76 | for irb,RBObj in enumerate(RBModels.get('Vector',[])): |
---|
77 | jrb = VRBIds.index(RBObj['RBId']) |
---|
78 | rbsx = str(irb)+':'+str(jrb) |
---|
79 | for i,px in enumerate(['RBVPx:','RBVPy:','RBVPz:']): |
---|
80 | RBObj['Orig'][0][i] = parmDict[pfx+px+rbsx] |
---|
81 | for i,po in enumerate(['RBVOa:','RBVOi:','RBVOj:','RBVOk:']): |
---|
82 | RBObj['Orient'][0][i] = parmDict[pfx+po+rbsx] |
---|
83 | RBObj['Orient'][0] = G2mth.normQ(RBObj['Orient'][0]) |
---|
84 | TLS = RBObj['ThermalMotion'] |
---|
85 | if 'T' in TLS[0]: |
---|
86 | for i,pt in enumerate(['RBVT11:','RBVT22:','RBVT33:','RBVT12:','RBVT13:','RBVT23:']): |
---|
87 | TLS[1][i] = parmDict[pfx+pt+rbsx] |
---|
88 | if 'L' in TLS[0]: |
---|
89 | for i,pt in enumerate(['RBVL11:','RBVL22:','RBVL33:','RBVL12:','RBVL13:','RBVL23:']): |
---|
90 | TLS[1][i+6] = parmDict[pfx+pt+rbsx] |
---|
91 | if 'S' in TLS[0]: |
---|
92 | for i,pt in enumerate(['RBVS12:','RBVS13:','RBVS21:','RBVS23:','RBVS31:','RBVS32:','RBVSAA:','RBVSBB:']): |
---|
93 | TLS[1][i+12] = parmDict[pfx+pt+rbsx] |
---|
94 | if 'U' in TLS[0]: |
---|
95 | TLS[1][0] = parmDict[pfx+'RBVU:'+rbsx] |
---|
96 | XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Vector') |
---|
97 | UIJ = G2mth.UpdateRBUIJ(Bmat,Cart,RBObj) |
---|
98 | for i,x in enumerate(XYZ): |
---|
99 | atId = RBObj['Ids'][i] |
---|
100 | for j in [0,1,2]: |
---|
101 | parmDict[pfx+atxIds[j]+str(AtLookup[atId])] = x[j] |
---|
102 | if UIJ[i][0] == 'A': |
---|
103 | for j in range(6): |
---|
104 | parmDict[pfx+atuIds[j]+str(AtLookup[atId])] = UIJ[i][j+2] |
---|
105 | elif UIJ[i][0] == 'I': |
---|
106 | parmDict[pfx+'AUiso:'+str(AtLookup[atId])] = UIJ[i][1] |
---|
107 | |
---|
108 | for irb,RBObj in enumerate(RBModels.get('Residue',[])): |
---|
109 | jrb = RRBIds.index(RBObj['RBId']) |
---|
110 | rbsx = str(irb)+':'+str(jrb) |
---|
111 | for i,px in enumerate(['RBRPx:','RBRPy:','RBRPz:']): |
---|
112 | RBObj['Orig'][0][i] = parmDict[pfx+px+rbsx] |
---|
113 | for i,po in enumerate(['RBROa:','RBROi:','RBROj:','RBROk:']): |
---|
114 | RBObj['Orient'][0][i] = parmDict[pfx+po+rbsx] |
---|
115 | RBObj['Orient'][0] = G2mth.normQ(RBObj['Orient'][0]) |
---|
116 | TLS = RBObj['ThermalMotion'] |
---|
117 | if 'T' in TLS[0]: |
---|
118 | for i,pt in enumerate(['RBRT11:','RBRT22:','RBRT33:','RBRT12:','RBRT13:','RBRT23:']): |
---|
119 | RBObj['ThermalMotion'][1][i] = parmDict[pfx+pt+rbsx] |
---|
120 | if 'L' in TLS[0]: |
---|
121 | for i,pt in enumerate(['RBRL11:','RBRL22:','RBRL33:','RBRL12:','RBRL13:','RBRL23:']): |
---|
122 | RBObj['ThermalMotion'][1][i+6] = parmDict[pfx+pt+rbsx] |
---|
123 | if 'S' in TLS[0]: |
---|
124 | for i,pt in enumerate(['RBRS12:','RBRS13:','RBRS21:','RBRS23:','RBRS31:','RBRS32:','RBRSAA:','RBRSBB:']): |
---|
125 | RBObj['ThermalMotion'][1][i+12] = parmDict[pfx+pt+rbsx] |
---|
126 | if 'U' in TLS[0]: |
---|
127 | RBObj['ThermalMotion'][1][0] = parmDict[pfx+'RBRU:'+rbsx] |
---|
128 | for itors,tors in enumerate(RBObj['Torsions']): |
---|
129 | tors[0] = parmDict[pfx+'RBRTr;'+str(itors)+':'+rbsx] |
---|
130 | XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Residue') |
---|
131 | UIJ = G2mth.UpdateRBUIJ(Bmat,Cart,RBObj) |
---|
132 | for i,x in enumerate(XYZ): |
---|
133 | atId = RBObj['Ids'][i] |
---|
134 | for j in [0,1,2]: |
---|
135 | parmDict[pfx+atxIds[j]+str(AtLookup[atId])] = x[j] |
---|
136 | if UIJ[i][0] == 'A': |
---|
137 | for j in range(6): |
---|
138 | parmDict[pfx+atuIds[j]+str(AtLookup[atId])] = UIJ[i][j+2] |
---|
139 | elif UIJ[i][0] == 'I': |
---|
140 | parmDict[pfx+'AUiso:'+str(AtLookup[atId])] = UIJ[i][1] |
---|
141 | |
---|
142 | def ApplyRBModelDervs(dFdvDict,parmDict,rigidbodyDict,Phase): |
---|
143 | 'Needs a doc string' |
---|
144 | atxIds = ['dAx:','dAy:','dAz:'] |
---|
145 | atuIds = ['AU11:','AU22:','AU33:','AU12:','AU13:','AU23:'] |
---|
146 | TIds = ['T11:','T22:','T33:','T12:','T13:','T23:'] |
---|
147 | LIds = ['L11:','L22:','L33:','L12:','L13:','L23:'] |
---|
148 | SIds = ['S12:','S13:','S21:','S23:','S31:','S32:','SAA:','SBB:'] |
---|
149 | PIds = ['Px:','Py:','Pz:'] |
---|
150 | OIds = ['Oa:','Oi:','Oj:','Ok:'] |
---|
151 | RBIds = rigidbodyDict.get('RBIds',{'Vector':[],'Residue':[]}) #these are lists of rbIds |
---|
152 | if not RBIds['Vector'] and not RBIds['Residue']: |
---|
153 | return |
---|
154 | VRBIds = RBIds['Vector'] |
---|
155 | RRBIds = RBIds['Residue'] |
---|
156 | RBData = rigidbodyDict |
---|
157 | for item in parmDict: |
---|
158 | if 'RB' in item: |
---|
159 | dFdvDict[item] = 0. #NB: this is a vector which is no. refl. long & must be filled! |
---|
160 | General = Phase['General'] |
---|
161 | cell = General['Cell'][1:7] |
---|
162 | Amat,Bmat = G2lat.cell2AB(cell) |
---|
163 | rpd = np.pi/180. |
---|
164 | rpd2 = rpd**2 |
---|
165 | g = nl.inv(np.inner(Bmat,Bmat)) |
---|
166 | gvec = np.sqrt(np.array([g[0][0]**2,g[1][1]**2,g[2][2]**2, |
---|
167 | g[0][0]*g[1][1],g[0][0]*g[2][2],g[1][1]*g[2][2]])) |
---|
168 | AtLookup = G2mth.FillAtomLookUp(Phase['Atoms']) |
---|
169 | pfx = str(Phase['pId'])+'::' |
---|
170 | RBModels = Phase['RBModels'] |
---|
171 | for irb,RBObj in enumerate(RBModels.get('Vector',[])): |
---|
172 | VModel = RBData['Vector'][RBObj['RBId']] |
---|
173 | Q = RBObj['Orient'][0] |
---|
174 | Pos = RBObj['Orig'][0] |
---|
175 | jrb = VRBIds.index(RBObj['RBId']) |
---|
176 | rbsx = str(irb)+':'+str(jrb) |
---|
177 | dXdv = [] |
---|
178 | for iv in range(len(VModel['VectMag'])): |
---|
179 | dCdv = [] |
---|
180 | for vec in VModel['rbVect'][iv]: |
---|
181 | dCdv.append(G2mth.prodQVQ(Q,vec)) |
---|
182 | dXdv.append(np.inner(Bmat,np.array(dCdv)).T) |
---|
183 | XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Vector') |
---|
184 | for ia,atId in enumerate(RBObj['Ids']): |
---|
185 | atNum = AtLookup[atId] |
---|
186 | dx = 0.00001 |
---|
187 | for iv in range(len(VModel['VectMag'])): |
---|
188 | for ix in [0,1,2]: |
---|
189 | dFdvDict['::RBV;'+str(iv)+':'+str(jrb)] += dXdv[iv][ia][ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)] |
---|
190 | for i,name in enumerate(['RBVPx:','RBVPy:','RBVPz:']): |
---|
191 | dFdvDict[pfx+name+rbsx] += dFdvDict[pfx+atxIds[i]+str(atNum)] |
---|
192 | for iv in range(4): |
---|
193 | Q[iv] -= dx |
---|
194 | XYZ1 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q)) |
---|
195 | Q[iv] += 2.*dx |
---|
196 | XYZ2 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q)) |
---|
197 | Q[iv] -= dx |
---|
198 | dXdO = (XYZ2[ia]-XYZ1[ia])/(2.*dx) |
---|
199 | for ix in [0,1,2]: |
---|
200 | dFdvDict[pfx+'RBV'+OIds[iv]+rbsx] += dXdO[ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)] |
---|
201 | X = G2mth.prodQVQ(Q,Cart[ia]) |
---|
202 | dFdu = np.array([dFdvDict[pfx+Uid+str(AtLookup[atId])] for Uid in atuIds]).T/gvec |
---|
203 | dFdu = G2lat.U6toUij(dFdu.T) |
---|
204 | dFdu = np.tensordot(Amat,np.tensordot(Amat,dFdu,([1,0])),([0,1])) |
---|
205 | dFdu = G2lat.UijtoU6(dFdu) |
---|
206 | atNum = AtLookup[atId] |
---|
207 | if 'T' in RBObj['ThermalMotion'][0]: |
---|
208 | for i,name in enumerate(['RBVT11:','RBVT22:','RBVT33:','RBVT12:','RBVT13:','RBVT23:']): |
---|
209 | dFdvDict[pfx+name+rbsx] += dFdu[i] |
---|
210 | if 'L' in RBObj['ThermalMotion'][0]: |
---|
211 | dFdvDict[pfx+'RBVL11:'+rbsx] += rpd2*(dFdu[1]*X[2]**2+dFdu[2]*X[1]**2-dFdu[5]*X[1]*X[2]) |
---|
212 | dFdvDict[pfx+'RBVL22:'+rbsx] += rpd2*(dFdu[0]*X[2]**2+dFdu[2]*X[0]**2-dFdu[4]*X[0]*X[2]) |
---|
213 | dFdvDict[pfx+'RBVL33:'+rbsx] += rpd2*(dFdu[0]*X[1]**2+dFdu[1]*X[0]**2-dFdu[3]*X[0]*X[1]) |
---|
214 | dFdvDict[pfx+'RBVL12:'+rbsx] += rpd2*(-dFdu[3]*X[2]**2-2.*dFdu[2]*X[0]*X[1]+ |
---|
215 | dFdu[4]*X[1]*X[2]+dFdu[5]*X[0]*X[2]) |
---|
216 | dFdvDict[pfx+'RBVL13:'+rbsx] += rpd2*(-dFdu[4]*X[1]**2-2.*dFdu[1]*X[0]*X[2]+ |
---|
217 | dFdu[3]*X[1]*X[2]+dFdu[5]*X[0]*X[1]) |
---|
218 | dFdvDict[pfx+'RBVL23:'+rbsx] += rpd2*(-dFdu[5]*X[0]**2-2.*dFdu[0]*X[1]*X[2]+ |
---|
219 | dFdu[3]*X[0]*X[2]+dFdu[4]*X[0]*X[1]) |
---|
220 | if 'S' in RBObj['ThermalMotion'][0]: |
---|
221 | dFdvDict[pfx+'RBVS12:'+rbsx] += rpd*(dFdu[5]*X[1]-2.*dFdu[1]*X[2]) |
---|
222 | dFdvDict[pfx+'RBVS13:'+rbsx] += rpd*(-dFdu[5]*X[2]+2.*dFdu[2]*X[1]) |
---|
223 | dFdvDict[pfx+'RBVS21:'+rbsx] += rpd*(-dFdu[4]*X[0]+2.*dFdu[0]*X[2]) |
---|
224 | dFdvDict[pfx+'RBVS23:'+rbsx] += rpd*(dFdu[4]*X[2]-2.*dFdu[2]*X[0]) |
---|
225 | dFdvDict[pfx+'RBVS31:'+rbsx] += rpd*(dFdu[3]*X[0]-2.*dFdu[0]*X[1]) |
---|
226 | dFdvDict[pfx+'RBVS32:'+rbsx] += rpd*(-dFdu[3]*X[1]+2.*dFdu[1]*X[0]) |
---|
227 | dFdvDict[pfx+'RBVSAA:'+rbsx] += rpd*(dFdu[4]*X[1]-dFdu[3]*X[2]) |
---|
228 | dFdvDict[pfx+'RBVSBB:'+rbsx] += rpd*(dFdu[5]*X[0]-dFdu[3]*X[2]) |
---|
229 | if 'U' in RBObj['ThermalMotion'][0]: |
---|
230 | dFdvDict[pfx+'RBVU:'+rbsx] += dFdvDict[pfx+'AUiso:'+str(AtLookup[atId])] |
---|
231 | |
---|
232 | |
---|
233 | for irb,RBObj in enumerate(RBModels.get('Residue',[])): |
---|
234 | Q = RBObj['Orient'][0] |
---|
235 | Pos = RBObj['Orig'][0] |
---|
236 | jrb = RRBIds.index(RBObj['RBId']) |
---|
237 | torData = RBData['Residue'][RBObj['RBId']]['rbSeq'] |
---|
238 | rbsx = str(irb)+':'+str(jrb) |
---|
239 | XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Residue') |
---|
240 | for itors,tors in enumerate(RBObj['Torsions']): #derivative error? |
---|
241 | tname = pfx+'RBRTr;'+str(itors)+':'+rbsx |
---|
242 | orId,pvId = torData[itors][:2] |
---|
243 | pivotVec = Cart[orId]-Cart[pvId] |
---|
244 | QA = G2mth.AVdeg2Q(-0.001,pivotVec) |
---|
245 | QB = G2mth.AVdeg2Q(0.001,pivotVec) |
---|
246 | for ir in torData[itors][3]: |
---|
247 | atNum = AtLookup[RBObj['Ids'][ir]] |
---|
248 | rVec = Cart[ir]-Cart[pvId] |
---|
249 | dR = G2mth.prodQVQ(QB,rVec)-G2mth.prodQVQ(QA,rVec) |
---|
250 | dRdT = np.inner(Bmat,G2mth.prodQVQ(Q,dR))/.002 |
---|
251 | for ix in [0,1,2]: |
---|
252 | dFdvDict[tname] += dRdT[ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)] |
---|
253 | for ia,atId in enumerate(RBObj['Ids']): |
---|
254 | atNum = AtLookup[atId] |
---|
255 | dx = 0.00001 |
---|
256 | for i,name in enumerate(['RBRPx:','RBRPy:','RBRPz:']): |
---|
257 | dFdvDict[pfx+name+rbsx] += dFdvDict[pfx+atxIds[i]+str(atNum)] |
---|
258 | for iv in range(4): |
---|
259 | Q[iv] -= dx |
---|
260 | XYZ1 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q)) |
---|
261 | Q[iv] += 2.*dx |
---|
262 | XYZ2 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q)) |
---|
263 | Q[iv] -= dx |
---|
264 | dXdO = (XYZ2[ia]-XYZ1[ia])/(2.*dx) |
---|
265 | for ix in [0,1,2]: |
---|
266 | dFdvDict[pfx+'RBR'+OIds[iv]+rbsx] += dXdO[ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)] |
---|
267 | X = G2mth.prodQVQ(Q,Cart[ia]) |
---|
268 | dFdu = np.array([dFdvDict[pfx+Uid+str(AtLookup[atId])] for Uid in atuIds]).T/gvec |
---|
269 | dFdu = G2lat.U6toUij(dFdu.T) |
---|
270 | dFdu = np.tensordot(Amat.T,np.tensordot(Amat,dFdu,([1,0])),([0,1])) |
---|
271 | dFdu = G2lat.UijtoU6(dFdu) |
---|
272 | atNum = AtLookup[atId] |
---|
273 | if 'T' in RBObj['ThermalMotion'][0]: |
---|
274 | for i,name in enumerate(['RBRT11:','RBRT22:','RBRT33:','RBRT12:','RBRT13:','RBRT23:']): |
---|
275 | dFdvDict[pfx+name+rbsx] += dFdu[i] |
---|
276 | if 'L' in RBObj['ThermalMotion'][0]: |
---|
277 | dFdvDict[pfx+'RBRL11:'+rbsx] += rpd2*(dFdu[1]*X[2]**2+dFdu[2]*X[1]**2-dFdu[5]*X[1]*X[2]) |
---|
278 | dFdvDict[pfx+'RBRL22:'+rbsx] += rpd2*(dFdu[0]*X[2]**2+dFdu[2]*X[0]**2-dFdu[4]*X[0]*X[2]) |
---|
279 | dFdvDict[pfx+'RBRL33:'+rbsx] += rpd2*(dFdu[0]*X[1]**2+dFdu[1]*X[0]**2-dFdu[3]*X[0]*X[1]) |
---|
280 | dFdvDict[pfx+'RBRL12:'+rbsx] += rpd2*(-dFdu[3]*X[2]**2-2.*dFdu[2]*X[0]*X[1]+ |
---|
281 | dFdu[4]*X[1]*X[2]+dFdu[5]*X[0]*X[2]) |
---|
282 | dFdvDict[pfx+'RBRL13:'+rbsx] += rpd2*(dFdu[4]*X[1]**2-2.*dFdu[1]*X[0]*X[2]+ |
---|
283 | dFdu[3]*X[1]*X[2]+dFdu[5]*X[0]*X[1]) |
---|
284 | dFdvDict[pfx+'RBRL23:'+rbsx] += rpd2*(dFdu[5]*X[0]**2-2.*dFdu[0]*X[1]*X[2]+ |
---|
285 | dFdu[3]*X[0]*X[2]+dFdu[4]*X[0]*X[1]) |
---|
286 | if 'S' in RBObj['ThermalMotion'][0]: |
---|
287 | dFdvDict[pfx+'RBRS12:'+rbsx] += rpd*(dFdu[5]*X[1]-2.*dFdu[1]*X[2]) |
---|
288 | dFdvDict[pfx+'RBRS13:'+rbsx] += rpd*(-dFdu[5]*X[2]+2.*dFdu[2]*X[1]) |
---|
289 | dFdvDict[pfx+'RBRS21:'+rbsx] += rpd*(-dFdu[4]*X[0]+2.*dFdu[0]*X[2]) |
---|
290 | dFdvDict[pfx+'RBRS23:'+rbsx] += rpd*(dFdu[4]*X[2]-2.*dFdu[2]*X[0]) |
---|
291 | dFdvDict[pfx+'RBRS31:'+rbsx] += rpd*(dFdu[3]*X[0]-2.*dFdu[0]*X[1]) |
---|
292 | dFdvDict[pfx+'RBRS32:'+rbsx] += rpd*(-dFdu[3]*X[1]+2.*dFdu[1]*X[0]) |
---|
293 | dFdvDict[pfx+'RBRSAA:'+rbsx] += rpd*(dFdu[4]*X[1]-dFdu[3]*X[2]) |
---|
294 | dFdvDict[pfx+'RBRSBB:'+rbsx] += rpd*(dFdu[5]*X[0]-dFdu[3]*X[2]) |
---|
295 | if 'U' in RBObj['ThermalMotion'][0]: |
---|
296 | dFdvDict[pfx+'RBRU:'+rbsx] += dFdvDict[pfx+'AUiso:'+str(AtLookup[atId])] |
---|
297 | |
---|
298 | ################################################################################ |
---|
299 | ##### Penalty & restraint functions |
---|
300 | ################################################################################ |
---|
301 | |
---|
302 | def penaltyFxn(HistoPhases,parmDict,varyList): |
---|
303 | 'Needs a doc string' |
---|
304 | Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases |
---|
305 | pNames = [] |
---|
306 | pVals = [] |
---|
307 | pWt = [] |
---|
308 | negWt = {} |
---|
309 | pWsum = {} |
---|
310 | for phase in Phases: |
---|
311 | pId = Phases[phase]['pId'] |
---|
312 | negWt[pId] = Phases[phase]['General']['Pawley neg wt'] |
---|
313 | General = Phases[phase]['General'] |
---|
314 | textureData = General['SH Texture'] |
---|
315 | SGData = General['SGData'] |
---|
316 | AtLookup = G2mth.FillAtomLookUp(Phases[phase]['Atoms']) |
---|
317 | cell = General['Cell'][1:7] |
---|
318 | Amat,Bmat = G2lat.cell2AB(cell) |
---|
319 | if phase not in restraintDict: |
---|
320 | continue |
---|
321 | phaseRest = restraintDict[phase] |
---|
322 | names = [['Bond','Bonds'],['Angle','Angles'],['Plane','Planes'], |
---|
323 | ['Chiral','Volumes'],['Torsion','Torsions'],['Rama','Ramas'], |
---|
324 | ['ChemComp','Sites'],['Texture','HKLs']] |
---|
325 | for name,rest in names: |
---|
326 | pWsum[name] = 0. |
---|
327 | itemRest = phaseRest[name] |
---|
328 | if itemRest[rest] and itemRest['Use']: |
---|
329 | wt = itemRest['wtFactor'] |
---|
330 | if name in ['Bond','Angle','Plane','Chiral']: |
---|
331 | for i,[indx,ops,obs,esd] in enumerate(itemRest[rest]): |
---|
332 | pNames.append(str(pId)+':'+name+':'+str(i)) |
---|
333 | XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx)) |
---|
334 | XYZ = G2mth.getSyXYZ(XYZ,ops,SGData) |
---|
335 | if name == 'Bond': |
---|
336 | calc = G2mth.getRestDist(XYZ,Amat) |
---|
337 | elif name == 'Angle': |
---|
338 | calc = G2mth.getRestAngle(XYZ,Amat) |
---|
339 | elif name == 'Plane': |
---|
340 | calc = G2mth.getRestPlane(XYZ,Amat) |
---|
341 | elif name == 'Chiral': |
---|
342 | calc = G2mth.getRestChiral(XYZ,Amat) |
---|
343 | pVals.append(obs-calc) |
---|
344 | pWt.append(wt/esd**2) |
---|
345 | pWsum[name] += wt*((obs-calc)/esd)**2 |
---|
346 | elif name in ['Torsion','Rama']: |
---|
347 | coeffDict = itemRest['Coeff'] |
---|
348 | for i,[indx,ops,cofName,esd] in enumerate(itemRest[rest]): |
---|
349 | pNames.append(str(pId)+':'+name+':'+str(i)) |
---|
350 | XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx)) |
---|
351 | XYZ = G2mth.getSyXYZ(XYZ,ops,SGData) |
---|
352 | if name == 'Torsion': |
---|
353 | tor = G2mth.getRestTorsion(XYZ,Amat) |
---|
354 | restr,calc = G2mth.calcTorsionEnergy(tor,coeffDict[cofName]) |
---|
355 | else: |
---|
356 | phi,psi = G2mth.getRestRama(XYZ,Amat) |
---|
357 | restr,calc = G2mth.calcRamaEnergy(phi,psi,coeffDict[cofName]) |
---|
358 | pVals.append(obs-calc) |
---|
359 | pWt.append(wt/esd**2) |
---|
360 | pWsum[name] += wt*((obs-calc)/esd)**2 |
---|
361 | elif name == 'ChemComp': |
---|
362 | for i,[indx,factors,obs,esd] in enumerate(itemRest[rest]): |
---|
363 | pNames.append(str(pId)+':'+name+':'+str(i)) |
---|
364 | mul = np.array(G2mth.GetAtomItemsById(Atoms,AtLookUp,indx,cs+1)) |
---|
365 | frac = np.array(G2mth.GetAtomItemsById(Atoms,AtLookUp,indx,cs-1)) |
---|
366 | calc = np.sum(mul*frac*factors) |
---|
367 | pVals.append(obs-calc) |
---|
368 | pWt.append(wt/esd**2) |
---|
369 | pWsum[name] += wt*((obs-calc)/esd)**2 |
---|
370 | elif name == 'Texture': |
---|
371 | SHkeys = textureData['SH Coeff'][1].keys() |
---|
372 | SHCoef = G2mth.GetSHCoeff(pId,parmDict,SHkeys) |
---|
373 | shModels = ['cylindrical','none','shear - 2/m','rolling - mmm'] |
---|
374 | SamSym = dict(zip(shModels,['0','-1','2/m','mmm'])) |
---|
375 | for i,[hkl,grid,esd1,ifesd2,esd2] in enumerate(itemRest[rest]): |
---|
376 | PH = np.array(hkl) |
---|
377 | phi,beta = G2lat.CrsAng(np.array(hkl),cell,SGData) |
---|
378 | ODFln = G2lat.Flnh(False,SHCoef,phi,beta,SGData) |
---|
379 | R,P,Z = G2mth.getRestPolefig(ODFln,SamSym[textureData['Model']],grid) |
---|
380 | Z1 = -ma.masked_greater(Z,0.0) |
---|
381 | IndZ1 = np.array(ma.nonzero(Z1)) |
---|
382 | for ind in IndZ1.T: |
---|
383 | pNames.append('%d:%s:%d:%.2f:%.2f'%(pId,name,i,R[ind[0],ind[1]],P[ind[0],ind[1]])) |
---|
384 | pVals.append(Z1[ind[0]][ind[1]]) |
---|
385 | pWt.append(wt/esd1**2) |
---|
386 | pWsum[name] += wt*((obs-calc)/esd)**2 |
---|
387 | if ifesd2: |
---|
388 | Z2 = 1.-Z |
---|
389 | for ind in np.ndindex(grid,grid): |
---|
390 | pNames.append('%d:%s:%d:%.2f:%.2f'%(pId,name+'-unit',i,R[ind[0],ind[1]],P[ind[0],ind[1]])) |
---|
391 | pVals.append(Z1[ind[0]][ind[1]]) |
---|
392 | pWt.append(wt/esd2**2) |
---|
393 | pWsum[name] += wt*((obs-calc)/esd)**2 |
---|
394 | |
---|
395 | for item in varyList: |
---|
396 | if 'PWLref' in item and parmDict[item] < 0.: |
---|
397 | pId = int(item.split(':')[0]) |
---|
398 | if negWt[pId]: |
---|
399 | pNames.append(item) |
---|
400 | pVals.append(-parmDict[item]) |
---|
401 | pWt.append(negWt[pId]) |
---|
402 | pWsum[name] += negWt[pId]*(-parmDict[item])**2 |
---|
403 | pVals = np.array(pVals) |
---|
404 | pWt = np.array(pWt) #should this be np.sqrt? |
---|
405 | return pNames,pVals,pWt,pWsum |
---|
406 | |
---|
407 | def penaltyDeriv(pNames,pVal,HistoPhases,parmDict,varyList): |
---|
408 | 'Needs a doc string' |
---|
409 | Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases |
---|
410 | pDerv = np.zeros((len(varyList),len(pVal))) |
---|
411 | for phase in Phases: |
---|
412 | # if phase not in restraintDict: |
---|
413 | # continue |
---|
414 | pId = Phases[phase]['pId'] |
---|
415 | General = Phases[phase]['General'] |
---|
416 | SGData = General['SGData'] |
---|
417 | AtLookup = G2mth.FillAtomLookUp(Phases[phase]['Atoms']) |
---|
418 | cell = General['Cell'][1:7] |
---|
419 | Amat,Bmat = G2lat.cell2AB(cell) |
---|
420 | textureData = General['SH Texture'] |
---|
421 | |
---|
422 | SHkeys = textureData['SH Coeff'][1].keys() |
---|
423 | SHCoef = G2mth.GetSHCoeff(pId,parmDict,SHkeys) |
---|
424 | shModels = ['cylindrical','none','shear - 2/m','rolling - mmm'] |
---|
425 | SamSym = dict(zip(shModels,['0','-1','2/m','mmm'])) |
---|
426 | sam = SamSym[textureData['Model']] |
---|
427 | phaseRest = restraintDict.get(phase,{}) |
---|
428 | names = {'Bond':'Bonds','Angle':'Angles','Plane':'Planes', |
---|
429 | 'Chiral':'Volumes','Torsion':'Torsions','Rama':'Ramas', |
---|
430 | 'ChemComp':'Sites','Texture':'HKLs'} |
---|
431 | lasthkl = np.array([0,0,0]) |
---|
432 | for ip,pName in enumerate(pNames): |
---|
433 | pnames = pName.split(':') |
---|
434 | if pId == int(pnames[0]): |
---|
435 | name = pnames[1] |
---|
436 | if 'PWL' in pName: |
---|
437 | pDerv[varyList.index(pName)][ip] += 1. |
---|
438 | continue |
---|
439 | id = int(pnames[2]) |
---|
440 | itemRest = phaseRest[name] |
---|
441 | if name in ['Bond','Angle','Plane','Chiral']: |
---|
442 | indx,ops,obs,esd = itemRest[names[name]][id] |
---|
443 | dNames = [] |
---|
444 | for ind in indx: |
---|
445 | dNames += [str(pId)+'::dA'+Xname+':'+str(AtLookup[ind]) for Xname in ['x','y','z']] |
---|
446 | XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx)) |
---|
447 | if name == 'Bond': |
---|
448 | deriv = G2mth.getRestDeriv(G2mth.getRestDist,XYZ,Amat,ops,SGData) |
---|
449 | elif name == 'Angle': |
---|
450 | deriv = G2mth.getRestDeriv(G2mth.getRestAngle,XYZ,Amat,ops,SGData) |
---|
451 | elif name == 'Plane': |
---|
452 | deriv = G2mth.getRestDeriv(G2mth.getRestPlane,XYZ,Amat,ops,SGData) |
---|
453 | elif name == 'Chiral': |
---|
454 | deriv = G2mth.getRestDeriv(G2mth.getRestChiral,XYZ,Amat,ops,SGData) |
---|
455 | elif name in ['Torsion','Rama']: |
---|
456 | coffDict = itemRest['Coeff'] |
---|
457 | indx,ops,cofName,esd = itemRest[names[name]][id] |
---|
458 | dNames = [] |
---|
459 | for ind in indx: |
---|
460 | dNames += [str(pId)+'::dA'+Xname+':'+str(AtLookup[ind]) for Xname in ['x','y','z']] |
---|
461 | XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx)) |
---|
462 | if name == 'Torsion': |
---|
463 | deriv = G2mth.getTorsionDeriv(XYZ,Amat,coffDict[cofName]) |
---|
464 | else: |
---|
465 | deriv = G2mth.getRamaDeriv(XYZ,Amat,coffDict[cofName]) |
---|
466 | elif name == 'ChemComp': |
---|
467 | indx,factors,obs,esd = itemRest[names[name]][id] |
---|
468 | dNames = [] |
---|
469 | for ind in indx: |
---|
470 | dNames += [str(pId)+'::Afrac:'+str(AtLookup[ind])] |
---|
471 | mul = np.array(G2mth.GetAtomItemsById(Atoms,AtLookUp,indx,cs+1)) |
---|
472 | deriv = mul*factors |
---|
473 | elif 'Texture' in name: |
---|
474 | deriv = [] |
---|
475 | dNames = [] |
---|
476 | hkl,grid,esd1,ifesd2,esd2 = itemRest[names[name]][id] |
---|
477 | hkl = np.array(hkl) |
---|
478 | if np.any(lasthkl-hkl): |
---|
479 | PH = np.array(hkl) |
---|
480 | phi,beta = G2lat.CrsAng(np.array(hkl),cell,SGData) |
---|
481 | ODFln = G2lat.Flnh(False,SHCoef,phi,beta,SGData) |
---|
482 | lasthkl = copy.copy(hkl) |
---|
483 | if 'unit' in name: |
---|
484 | pass |
---|
485 | else: |
---|
486 | gam = float(pnames[3]) |
---|
487 | psi = float(pnames[4]) |
---|
488 | for SHname in ODFln: |
---|
489 | l,m,n = eval(SHname[1:]) |
---|
490 | Ksl = G2lat.GetKsl(l,m,sam,psi,gam)[0] |
---|
491 | dNames += [str(pId)+'::'+SHname] |
---|
492 | deriv.append(-ODFln[SHname][0]*Ksl/SHCoef[SHname]) |
---|
493 | for dName,drv in zip(dNames,deriv): |
---|
494 | try: |
---|
495 | ind = varyList.index(dName) |
---|
496 | pDerv[ind][ip] += drv |
---|
497 | except ValueError: |
---|
498 | pass |
---|
499 | return pDerv |
---|
500 | |
---|
501 | ################################################################################ |
---|
502 | ##### Function & derivative calculations |
---|
503 | ################################################################################ |
---|
504 | |
---|
505 | def GetAtomFXU(pfx,calcControls,parmDict): |
---|
506 | 'Needs a doc string' |
---|
507 | Natoms = calcControls['Natoms'][pfx] |
---|
508 | Tdata = Natoms*[' ',] |
---|
509 | Mdata = np.zeros(Natoms) |
---|
510 | IAdata = Natoms*[' ',] |
---|
511 | Fdata = np.zeros(Natoms) |
---|
512 | FFdata = [] |
---|
513 | BLdata = [] |
---|
514 | Xdata = np.zeros((3,Natoms)) |
---|
515 | dXdata = np.zeros((3,Natoms)) |
---|
516 | Uisodata = np.zeros(Natoms) |
---|
517 | Uijdata = np.zeros((6,Natoms)) |
---|
518 | keys = {'Atype:':Tdata,'Amul:':Mdata,'Afrac:':Fdata,'AI/A:':IAdata, |
---|
519 | 'dAx:':dXdata[0],'dAy:':dXdata[1],'dAz:':dXdata[2], |
---|
520 | 'Ax:':Xdata[0],'Ay:':Xdata[1],'Az:':Xdata[2],'AUiso:':Uisodata, |
---|
521 | 'AU11:':Uijdata[0],'AU22:':Uijdata[1],'AU33:':Uijdata[2], |
---|
522 | 'AU12:':Uijdata[3],'AU13:':Uijdata[4],'AU23:':Uijdata[5]} |
---|
523 | for iatm in range(Natoms): |
---|
524 | for key in keys: |
---|
525 | parm = pfx+key+str(iatm) |
---|
526 | if parm in parmDict: |
---|
527 | keys[key][iatm] = parmDict[parm] |
---|
528 | return Tdata,Mdata,Fdata,Xdata,dXdata,IAdata,Uisodata,Uijdata |
---|
529 | |
---|
530 | def StructureFactor(refList,G,hfx,pfx,SGData,calcControls,parmDict): |
---|
531 | ''' Compute structure factors for all h,k,l for phase |
---|
532 | puts the result, F^2, in each ref[8] in refList |
---|
533 | input: |
---|
534 | |
---|
535 | :param list refList: [ref] where each ref = h,k,l,m,d,...,[equiv h,k,l],phase[equiv] |
---|
536 | :param np.array G: reciprocal metric tensor |
---|
537 | :param str pfx: phase id string |
---|
538 | :param dict SGData: space group info. dictionary output from SpcGroup |
---|
539 | :param dict calcControls: |
---|
540 | :param dict ParmDict: |
---|
541 | |
---|
542 | ''' |
---|
543 | twopi = 2.0*np.pi |
---|
544 | twopisq = 2.0*np.pi**2 |
---|
545 | phfx = pfx.split(':')[0]+hfx |
---|
546 | ast = np.sqrt(np.diag(G)) |
---|
547 | Mast = twopisq*np.multiply.outer(ast,ast) |
---|
548 | FFtables = calcControls['FFtables'] |
---|
549 | BLtables = calcControls['BLtables'] |
---|
550 | Tdata,Mdata,Fdata,Xdata,dXdata,IAdata,Uisodata,Uijdata = GetAtomFXU(pfx,calcControls,parmDict) |
---|
551 | FF = np.zeros(len(Tdata)) |
---|
552 | if 'N' in calcControls[hfx+'histType']: |
---|
553 | FP,FPP = G2el.BlenRes(Tdata,BLtables,parmDict[hfx+'Lam']) |
---|
554 | else: |
---|
555 | FP = np.array([FFtables[El][hfx+'FP'] for El in Tdata]) |
---|
556 | FPP = np.array([FFtables[El][hfx+'FPP'] for El in Tdata]) |
---|
557 | Uij = np.array(G2lat.U6toUij(Uijdata)) |
---|
558 | bij = Mast*Uij.T |
---|
559 | for refl in refList: |
---|
560 | fbs = np.array([0,0]) |
---|
561 | H = refl[:3] |
---|
562 | SQ = 1./(2.*refl[4])**2 |
---|
563 | SQfactor = 4.0*SQ*twopisq |
---|
564 | Bab = parmDict[phfx+'BabA']*np.exp(-parmDict[phfx+'BabU']*SQfactor) |
---|
565 | if not len(refl[-1]): #no form factors |
---|
566 | if 'N' in calcControls[hfx+'histType']: |
---|
567 | refl[-1] = G2el.getBLvalues(BLtables) |
---|
568 | else: #'X' |
---|
569 | refl[-1] = G2el.getFFvalues(FFtables,SQ) |
---|
570 | for i,El in enumerate(Tdata): |
---|
571 | FF[i] = refl[-1][El] |
---|
572 | Uniq = refl[11] |
---|
573 | phi = refl[12] |
---|
574 | phase = twopi*(np.inner(Uniq,(dXdata.T+Xdata.T))+phi[:,np.newaxis]) |
---|
575 | sinp = np.sin(phase) |
---|
576 | cosp = np.cos(phase) |
---|
577 | occ = Mdata*Fdata/len(Uniq) |
---|
578 | biso = -SQfactor*Uisodata |
---|
579 | Tiso = np.where(biso<1.,np.exp(biso),1.0) |
---|
580 | HbH = np.array([-np.inner(h,np.inner(bij,h)) for h in Uniq]) |
---|
581 | Tuij = np.where(HbH<1.,np.exp(HbH),1.0) |
---|
582 | Tcorr = Tiso*Tuij |
---|
583 | fa = np.array([(FF+FP-Bab)*occ*cosp*Tcorr,-FPP*occ*sinp*Tcorr]) |
---|
584 | fas = np.sum(np.sum(fa,axis=1),axis=1) #real |
---|
585 | if not SGData['SGInv']: |
---|
586 | fb = np.array([(FF+FP-Bab)*occ*sinp*Tcorr,FPP*occ*cosp*Tcorr]) |
---|
587 | fbs = np.sum(np.sum(fb,axis=1),axis=1) |
---|
588 | fasq = fas**2 |
---|
589 | fbsq = fbs**2 #imaginary |
---|
590 | refl[9] = np.sum(fasq)+np.sum(fbsq) |
---|
591 | refl[10] = atan2d(fbs[0],fas[0]) |
---|
592 | |
---|
593 | def StructureFactorDerv(refList,G,hfx,pfx,SGData,calcControls,parmDict): |
---|
594 | 'Needs a doc string' |
---|
595 | twopi = 2.0*np.pi |
---|
596 | twopisq = 2.0*np.pi**2 |
---|
597 | phfx = pfx.split(':')[0]+hfx |
---|
598 | ast = np.sqrt(np.diag(G)) |
---|
599 | Mast = twopisq*np.multiply.outer(ast,ast) |
---|
600 | FFtables = calcControls['FFtables'] |
---|
601 | BLtables = calcControls['BLtables'] |
---|
602 | Tdata,Mdata,Fdata,Xdata,dXdata,IAdata,Uisodata,Uijdata = GetAtomFXU(pfx,calcControls,parmDict) |
---|
603 | FF = np.zeros(len(Tdata)) |
---|
604 | if 'N' in calcControls[hfx+'histType']: |
---|
605 | FP = 0. |
---|
606 | FPP = 0. |
---|
607 | else: |
---|
608 | FP = np.array([FFtables[El][hfx+'FP'] for El in Tdata]) |
---|
609 | FPP = np.array([FFtables[El][hfx+'FPP'] for El in Tdata]) |
---|
610 | Uij = np.array(G2lat.U6toUij(Uijdata)) |
---|
611 | bij = Mast*Uij.T |
---|
612 | dFdvDict = {} |
---|
613 | dFdfr = np.zeros((len(refList),len(Mdata))) |
---|
614 | dFdx = np.zeros((len(refList),len(Mdata),3)) |
---|
615 | dFdui = np.zeros((len(refList),len(Mdata))) |
---|
616 | dFdua = np.zeros((len(refList),len(Mdata),6)) |
---|
617 | dFdbab = np.zeros((len(refList),2)) |
---|
618 | for iref,refl in enumerate(refList): |
---|
619 | H = np.array(refl[:3]) |
---|
620 | SQ = 1./(2.*refl[4])**2 # or (sin(theta)/lambda)**2 |
---|
621 | SQfactor = 8.0*SQ*np.pi**2 |
---|
622 | dBabdA = np.exp(-parmDict[phfx+'BabU']*SQfactor) |
---|
623 | Bab = parmDict[phfx+'BabA']*dBabdA |
---|
624 | for i,El in enumerate(Tdata): |
---|
625 | FF[i] = refl[-1][El] |
---|
626 | Uniq = refl[11] |
---|
627 | phi = refl[12] |
---|
628 | phase = twopi*(np.inner((dXdata.T+Xdata.T),Uniq)+phi[np.newaxis,:]) |
---|
629 | sinp = np.sin(phase) |
---|
630 | cosp = np.cos(phase) |
---|
631 | occ = Mdata*Fdata/len(Uniq) |
---|
632 | biso = -SQfactor*Uisodata |
---|
633 | Tiso = np.where(biso<1.,np.exp(biso),1.0) |
---|
634 | HbH = -np.inner(H,np.inner(bij,H)) |
---|
635 | Hij = np.array([Mast*np.multiply.outer(U,U) for U in Uniq]) |
---|
636 | Hij = np.array([G2lat.UijtoU6(Uij) for Uij in Hij]) |
---|
637 | Tuij = np.where(HbH<1.,np.exp(HbH),1.0) |
---|
638 | Tcorr = Tiso*Tuij |
---|
639 | fot = (FF+FP-Bab)*occ*Tcorr |
---|
640 | fotp = FPP*occ*Tcorr |
---|
641 | fa = np.array([fot[:,np.newaxis]*cosp,fotp[:,np.newaxis]*cosp]) #non positions |
---|
642 | fb = np.array([fot[:,np.newaxis]*sinp,-fotp[:,np.newaxis]*sinp]) |
---|
643 | |
---|
644 | fas = np.sum(np.sum(fa,axis=1),axis=1) |
---|
645 | fbs = np.sum(np.sum(fb,axis=1),axis=1) |
---|
646 | fax = np.array([-fot[:,np.newaxis]*sinp,-fotp[:,np.newaxis]*sinp]) #positions |
---|
647 | fbx = np.array([fot[:,np.newaxis]*cosp,-fot[:,np.newaxis]*cosp]) |
---|
648 | #sum below is over Uniq |
---|
649 | dfadfr = np.sum(fa/occ[:,np.newaxis],axis=2) |
---|
650 | dfadx = np.sum(twopi*Uniq*fax[:,:,:,np.newaxis],axis=2) |
---|
651 | dfadui = np.sum(-SQfactor*fa,axis=2) |
---|
652 | dfadua = np.sum(-Hij*fa[:,:,:,np.newaxis],axis=2) |
---|
653 | dfadba = np.sum(-cosp*(occ*Tcorr)[:,np.newaxis],axis=1) |
---|
654 | #NB: the above have been checked against PA(1:10,1:2) in strfctr.for |
---|
655 | dFdfr[iref] = 2.*(fas[0]*dfadfr[0]+fas[1]*dfadfr[1])*Mdata/len(Uniq) |
---|
656 | dFdx[iref] = 2.*(fas[0]*dfadx[0]+fas[1]*dfadx[1]) |
---|
657 | dFdui[iref] = 2.*(fas[0]*dfadui[0]+fas[1]*dfadui[1]) |
---|
658 | dFdua[iref] = 2.*(fas[0]*dfadua[0]+fas[1]*dfadua[1]) |
---|
659 | dFdbab[iref] = np.array([np.sum(dfadba*dBabdA),np.sum(-dfadba*parmDict[phfx+'BabA']*SQfactor*dBabdA)]).T |
---|
660 | if not SGData['SGInv']: |
---|
661 | dfbdfr = np.sum(fb/occ[:,np.newaxis],axis=2) #problem here if occ=0 for some atom |
---|
662 | dfbdx = np.sum(twopi*Uniq*fbx[:,:,:,np.newaxis],axis=2) |
---|
663 | dfbdui = np.sum(-SQfactor*fb,axis=2) |
---|
664 | dfbdua = np.sum(-Hij*fb[:,:,:,np.newaxis],axis=2) |
---|
665 | dfbdba = np.sum(-sinp*(occ*Tcorr)[:,np.newaxis],axis=1) |
---|
666 | dFdfr[iref] += 2.*(fbs[0]*dfbdfr[0]-fbs[1]*dfbdfr[1])*Mdata/len(Uniq) |
---|
667 | dFdx[iref] += 2.*(fbs[0]*dfbdx[0]+fbs[1]*dfbdx[1]) |
---|
668 | dFdui[iref] += 2.*(fbs[0]*dfbdui[0]-fbs[1]*dfbdui[1]) |
---|
669 | dFdua[iref] += 2.*(fbs[0]*dfbdua[0]+fbs[1]*dfbdua[1]) |
---|
670 | dFdbab[iref] += np.array([np.sum(dfbdba*dBabdA),np.sum(-dfbdba*parmDict[phfx+'BabA']*SQfactor*dBabdA)]).T |
---|
671 | #loop over atoms - each dict entry is list of derivatives for all the reflections |
---|
672 | for i in range(len(Mdata)): |
---|
673 | dFdvDict[pfx+'Afrac:'+str(i)] = dFdfr.T[i] |
---|
674 | dFdvDict[pfx+'dAx:'+str(i)] = dFdx.T[0][i] |
---|
675 | dFdvDict[pfx+'dAy:'+str(i)] = dFdx.T[1][i] |
---|
676 | dFdvDict[pfx+'dAz:'+str(i)] = dFdx.T[2][i] |
---|
677 | dFdvDict[pfx+'AUiso:'+str(i)] = dFdui.T[i] |
---|
678 | dFdvDict[pfx+'AU11:'+str(i)] = dFdua.T[0][i] |
---|
679 | dFdvDict[pfx+'AU22:'+str(i)] = dFdua.T[1][i] |
---|
680 | dFdvDict[pfx+'AU33:'+str(i)] = dFdua.T[2][i] |
---|
681 | dFdvDict[pfx+'AU12:'+str(i)] = 2.*dFdua.T[3][i] |
---|
682 | dFdvDict[pfx+'AU13:'+str(i)] = 2.*dFdua.T[4][i] |
---|
683 | dFdvDict[pfx+'AU23:'+str(i)] = 2.*dFdua.T[5][i] |
---|
684 | dFdvDict[pfx+'BabA'] = dFdbab.T[0] |
---|
685 | dFdvDict[pfx+'BabU'] = dFdbab.T[1] |
---|
686 | return dFdvDict |
---|
687 | |
---|
688 | def SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varyList): |
---|
689 | ''' Single crystal extinction function; puts correction in ref[13] and returns |
---|
690 | corrections needed for derivatives |
---|
691 | ''' |
---|
692 | ref[13] = 1.0 |
---|
693 | dervCor = 1.0 |
---|
694 | dervDict = {} |
---|
695 | if calcControls[phfx+'EType'] != 'None': |
---|
696 | cos2T = 1.0-0.5*(parmDict[hfx+'Lam']/ref[4])**2 #cos(2theta) |
---|
697 | if 'SXC' in parmDict[hfx+'Type']: |
---|
698 | AV = 7.9406e5/parmDict[pfx+'Vol']**2 |
---|
699 | PL = np.sqrt(1.0-cos2T**2)/parmDict[hfx+'Lam'] |
---|
700 | P12 = (calcControls[phfx+'Cos2TM']+cos2T**4)/(calcControls[phfx+'Cos2TM']+cos2T**2) |
---|
701 | elif 'SNT' in parmDict[hfx+'Type']: |
---|
702 | AV = 1.e7/parmDict[pfx+'Vol']**2 |
---|
703 | PL = 1./(4.*refl[4]**2) |
---|
704 | P12 = 1.0 |
---|
705 | elif 'SNC' in parmDict[hfx+'Type']: |
---|
706 | AV = 1.e7/parmDict[pfx+'Vol']**2 |
---|
707 | PL = np.sqrt(1.0-cos2T**2)/parmDict[hfx+'Lam'] |
---|
708 | P12 = 1.0 |
---|
709 | |
---|
710 | PLZ = AV*P12*parmDict[hfx+'Lam']**2*ref[7] |
---|
711 | if 'Primary' in calcControls[phfx+'EType']: |
---|
712 | PLZ *= 1.5 |
---|
713 | else: |
---|
714 | PLZ *= calcControls[phfx+'Tbar'] |
---|
715 | |
---|
716 | if 'Primary' in calcControls[phfx+'EType']: |
---|
717 | PSIG = parmDict[phfx+'Ep'] |
---|
718 | elif 'I & II' in calcControls[phfx+'EType']: |
---|
719 | PSIG = parmDict[phfx+'Eg']/np.sqrt(1.+(parmDict[phfx+'Es']*PL/parmDict[phfx+'Eg'])**2) |
---|
720 | elif 'Type II' in calcControls[phfx+'EType']: |
---|
721 | PSIG = parmDict[phfx+'Es'] |
---|
722 | else: # 'Secondary Type I' |
---|
723 | PSIG = parmDict[phfx+'Eg']/PL |
---|
724 | |
---|
725 | AG = 0.58+0.48*cos2T+0.24*cos2T**2 |
---|
726 | AL = 0.025+0.285*cos2T |
---|
727 | BG = 0.02-0.025*cos2T |
---|
728 | BL = 0.15-0.2*(0.75-cos2T)**2 |
---|
729 | if cos2T < 0.: |
---|
730 | BL = -0.45*cos2T |
---|
731 | CG = 2. |
---|
732 | CL = 2. |
---|
733 | PF = PLZ*PSIG |
---|
734 | |
---|
735 | if 'Gaussian' in calcControls[phfx+'EApprox']: |
---|
736 | PF4 = 1.+CG*PF+AG*PF**2/(1.+BG*PF) |
---|
737 | extCor = np.sqrt(PF4) |
---|
738 | PF3 = 0.5*(CG+2.*AG*PF/(1.+BG*PF)-AG*PF**2*BG/(1.+BG*PF)**2)/(PF4*extCor) |
---|
739 | else: |
---|
740 | PF4 = 1.+CL*PF+AL*PF**2/(1.+BL*PF) |
---|
741 | extCor = np.sqrt(PF4) |
---|
742 | PF3 = 0.5*(CL+2.*AL*PF/(1.+BL*PF)-AL*PF**2*BL/(1.+BL*PF)**2)/(PF4*extCor) |
---|
743 | |
---|
744 | dervCor = (1.+PF)*PF3 |
---|
745 | if 'Primary' in calcControls[phfx+'EType'] and phfx+'Ep' in varyList: |
---|
746 | dervDict[phfx+'Ep'] = -ref[7]*PLZ*PF3 |
---|
747 | if 'II' in calcControls[phfx+'EType'] and phfx+'Es' in varyList: |
---|
748 | dervDict[phfx+'Es'] = -ref[7]*PLZ*PF3*(PSIG/parmDict[phfx+'Es'])**3 |
---|
749 | if 'I' in calcControls[phfx+'EType'] and phfx+'Eg' in varyList: |
---|
750 | dervDict[phfx+'Eg'] = -ref[7]*PLZ*PF3*(PSIG/parmDict[phfx+'Eg'])**3*PL**2 |
---|
751 | |
---|
752 | ref[13] = 1./extCor |
---|
753 | return dervCor,dervDict |
---|
754 | |
---|
755 | |
---|
756 | def Dict2Values(parmdict, varylist): |
---|
757 | '''Use before call to leastsq to setup list of values for the parameters |
---|
758 | in parmdict, as selected by key in varylist''' |
---|
759 | return [parmdict[key] for key in varylist] |
---|
760 | |
---|
761 | def Values2Dict(parmdict, varylist, values): |
---|
762 | ''' Use after call to leastsq to update the parameter dictionary with |
---|
763 | values corresponding to keys in varylist''' |
---|
764 | parmdict.update(zip(varylist,values)) |
---|
765 | |
---|
766 | def GetNewCellParms(parmDict,varyList): |
---|
767 | 'Needs a doc string' |
---|
768 | newCellDict = {} |
---|
769 | Anames = ['A'+str(i) for i in range(6)] |
---|
770 | Ddict = dict(zip(['D11','D22','D33','D12','D13','D23'],Anames)) |
---|
771 | for item in varyList: |
---|
772 | keys = item.split(':') |
---|
773 | if keys[2] in Ddict: |
---|
774 | key = keys[0]+'::'+Ddict[keys[2]] #key is e.g. '0::A0' |
---|
775 | parm = keys[0]+'::'+keys[2] #parm is e.g. '0::D11' |
---|
776 | newCellDict[parm] = [key,parmDict[key]+parmDict[item]] |
---|
777 | return newCellDict # is e.g. {'0::D11':A0+D11} |
---|
778 | |
---|
779 | def ApplyXYZshifts(parmDict,varyList): |
---|
780 | ''' |
---|
781 | takes atom x,y,z shift and applies it to corresponding atom x,y,z value |
---|
782 | |
---|
783 | :param dict parmDict: parameter dictionary |
---|
784 | :param list varyList: list of variables |
---|
785 | :returns: newAtomDict - dictionary of new atomic coordinate names & values; key is parameter shift name |
---|
786 | |
---|
787 | ''' |
---|
788 | newAtomDict = {} |
---|
789 | for item in parmDict: |
---|
790 | if 'dA' in item: |
---|
791 | parm = ''.join(item.split('d')) |
---|
792 | parmDict[parm] += parmDict[item] |
---|
793 | newAtomDict[item] = [parm,parmDict[parm]] |
---|
794 | return newAtomDict |
---|
795 | |
---|
796 | def SHTXcal(refl,g,pfx,hfx,SGData,calcControls,parmDict): |
---|
797 | 'Spherical harmonics texture' |
---|
798 | IFCoup = 'Bragg' in calcControls[hfx+'instType'] |
---|
799 | odfCor = 1.0 |
---|
800 | H = refl[:3] |
---|
801 | cell = G2lat.Gmat2cell(g) |
---|
802 | Sangls = [parmDict[pfx+'SH omega'],parmDict[pfx+'SH chi'],parmDict[pfx+'SH phi']] |
---|
803 | Gangls = [parmDict[hfx+'Omega'],parmDict[hfx+'Chi'],parmDict[hfx+'Phi'],parmDict[hfx+'Azimuth']] |
---|
804 | phi,beta = G2lat.CrsAng(H,cell,SGData) |
---|
805 | psi,gam,x,x = G2lat.SamAng(refl[5]/2.,Gangls,Sangls,IFCoup) #ignore 2 sets of angle derivs. |
---|
806 | SHnames = G2lat.GenSHCoeff(SGData['SGLaue'],parmDict[pfx+'SHmodel'],parmDict[pfx+'SHorder']) |
---|
807 | for item in SHnames: |
---|
808 | L,M,N = eval(item.strip('C')) |
---|
809 | Kcl = G2lat.GetKcl(L,N,SGData['SGLaue'],phi,beta) |
---|
810 | Ksl,x,x = G2lat.GetKsl(L,M,parmDict[pfx+'SHmodel'],psi,gam) |
---|
811 | Lnorm = G2lat.Lnorm(L) |
---|
812 | odfCor += parmDict[pfx+item]*Lnorm*Kcl*Ksl |
---|
813 | return odfCor |
---|
814 | |
---|
815 | def SHTXcalDerv(refl,g,pfx,hfx,SGData,calcControls,parmDict): |
---|
816 | 'Spherical harmonics texture derivatives' |
---|
817 | FORPI = 4.0*np.pi |
---|
818 | IFCoup = 'Bragg' in calcControls[hfx+'instType'] |
---|
819 | odfCor = 1.0 |
---|
820 | dFdODF = {} |
---|
821 | dFdSA = [0,0,0] |
---|
822 | H = refl[:3] |
---|
823 | cell = G2lat.Gmat2cell(g) |
---|
824 | Sangls = [parmDict[pfx+'SH omega'],parmDict[pfx+'SH chi'],parmDict[pfx+'SH phi']] |
---|
825 | Gangls = [parmDict[hfx+'Omega'],parmDict[hfx+'Chi'],parmDict[hfx+'Phi'],parmDict[hfx+'Azimuth']] |
---|
826 | phi,beta = G2lat.CrsAng(H,cell,SGData) |
---|
827 | psi,gam,dPSdA,dGMdA = G2lat.SamAng(refl[5]/2.,Gangls,Sangls,IFCoup) |
---|
828 | SHnames = G2lat.GenSHCoeff(SGData['SGLaue'],parmDict[pfx+'SHmodel'],parmDict[pfx+'SHorder']) |
---|
829 | for item in SHnames: |
---|
830 | L,M,N = eval(item.strip('C')) |
---|
831 | Kcl = G2lat.GetKcl(L,N,SGData['SGLaue'],phi,beta) |
---|
832 | Ksl,dKsdp,dKsdg = G2lat.GetKsl(L,M,parmDict[pfx+'SHmodel'],psi,gam) |
---|
833 | Lnorm = G2lat.Lnorm(L) |
---|
834 | odfCor += parmDict[pfx+item]*Lnorm*Kcl*Ksl |
---|
835 | dFdODF[pfx+item] = Lnorm*Kcl*Ksl |
---|
836 | for i in range(3): |
---|
837 | dFdSA[i] += parmDict[pfx+item]*Lnorm*Kcl*(dKsdp*dPSdA[i]+dKsdg*dGMdA[i]) |
---|
838 | return odfCor,dFdODF,dFdSA |
---|
839 | |
---|
840 | def SHPOcal(refl,g,phfx,hfx,SGData,calcControls,parmDict): |
---|
841 | 'spherical harmonics preferred orientation (cylindrical symmetry only)' |
---|
842 | odfCor = 1.0 |
---|
843 | H = refl[:3] |
---|
844 | cell = G2lat.Gmat2cell(g) |
---|
845 | Sangl = [0.,0.,0.] |
---|
846 | if 'Bragg' in calcControls[hfx+'instType']: |
---|
847 | Gangls = [0.,90.,0.,parmDict[hfx+'Azimuth']] |
---|
848 | IFCoup = True |
---|
849 | else: |
---|
850 | Gangls = [0.,0.,0.,parmDict[hfx+'Azimuth']] |
---|
851 | IFCoup = False |
---|
852 | phi,beta = G2lat.CrsAng(H,cell,SGData) |
---|
853 | psi,gam,x,x = G2lat.SamAng(refl[5]/2.,Gangls,Sangl,IFCoup) #ignore 2 sets of angle derivs. |
---|
854 | SHnames = G2lat.GenSHCoeff(SGData['SGLaue'],'0',calcControls[phfx+'SHord'],False) |
---|
855 | for item in SHnames: |
---|
856 | L,N = eval(item.strip('C')) |
---|
857 | Kcsl,Lnorm = G2lat.GetKclKsl(L,N,SGData['SGLaue'],psi,phi,beta) |
---|
858 | odfCor += parmDict[phfx+item]*Lnorm*Kcsl |
---|
859 | return np.squeeze(odfCor) |
---|
860 | |
---|
861 | def SHPOcalDerv(refl,g,phfx,hfx,SGData,calcControls,parmDict): |
---|
862 | 'spherical harmonics preferred orientation derivatives (cylindrical symmetry only)' |
---|
863 | FORPI = 12.5663706143592 |
---|
864 | odfCor = 1.0 |
---|
865 | dFdODF = {} |
---|
866 | H = refl[:3] |
---|
867 | cell = G2lat.Gmat2cell(g) |
---|
868 | Sangl = [0.,0.,0.] |
---|
869 | if 'Bragg' in calcControls[hfx+'instType']: |
---|
870 | Gangls = [0.,90.,0.,parmDict[hfx+'Azimuth']] |
---|
871 | IFCoup = True |
---|
872 | else: |
---|
873 | Gangls = [0.,0.,0.,parmDict[hfx+'Azimuth']] |
---|
874 | IFCoup = False |
---|
875 | phi,beta = G2lat.CrsAng(H,cell,SGData) |
---|
876 | psi,gam,x,x = G2lat.SamAng(refl[5]/2.,Gangls,Sangl,IFCoup) #ignore 2 sets of angle derivs. |
---|
877 | SHnames = G2lat.GenSHCoeff(SGData['SGLaue'],'0',calcControls[phfx+'SHord'],False) |
---|
878 | for item in SHnames: |
---|
879 | L,N = eval(item.strip('C')) |
---|
880 | Kcsl,Lnorm = G2lat.GetKclKsl(L,N,SGData['SGLaue'],psi,phi,beta) |
---|
881 | odfCor += parmDict[phfx+item]*Lnorm*Kcsl |
---|
882 | dFdODF[phfx+item] = Kcsl*Lnorm |
---|
883 | return odfCor,dFdODF |
---|
884 | |
---|
885 | def GetPrefOri(refl,G,g,phfx,hfx,SGData,calcControls,parmDict): |
---|
886 | 'Needs a doc string' |
---|
887 | POcorr = 1.0 |
---|
888 | if calcControls[phfx+'poType'] == 'MD': |
---|
889 | MD = parmDict[phfx+'MD'] |
---|
890 | if MD != 1.0: |
---|
891 | MDAxis = calcControls[phfx+'MDAxis'] |
---|
892 | sumMD = 0 |
---|
893 | for H in refl[11]: |
---|
894 | cosP,sinP = G2lat.CosSinAngle(H,MDAxis,G) |
---|
895 | A = 1.0/np.sqrt((MD*cosP)**2+sinP**2/MD) |
---|
896 | sumMD += A**3 |
---|
897 | POcorr = sumMD/len(refl[11]) |
---|
898 | else: #spherical harmonics |
---|
899 | if calcControls[phfx+'SHord']: |
---|
900 | POcorr = SHPOcal(refl,g,phfx,hfx,SGData,calcControls,parmDict) |
---|
901 | return POcorr |
---|
902 | |
---|
903 | def GetPrefOriDerv(refl,G,g,phfx,hfx,SGData,calcControls,parmDict): |
---|
904 | 'Needs a doc string' |
---|
905 | POcorr = 1.0 |
---|
906 | POderv = {} |
---|
907 | if calcControls[phfx+'poType'] == 'MD': |
---|
908 | MD = parmDict[phfx+'MD'] |
---|
909 | MDAxis = calcControls[phfx+'MDAxis'] |
---|
910 | sumMD = 0 |
---|
911 | sumdMD = 0 |
---|
912 | for H in refl[11]: |
---|
913 | cosP,sinP = G2lat.CosSinAngle(H,MDAxis,G) |
---|
914 | A = 1.0/np.sqrt((MD*cosP)**2+sinP**2/MD) |
---|
915 | sumMD += A**3 |
---|
916 | sumdMD -= (1.5*A**5)*(2.0*MD*cosP**2-(sinP/MD)**2) |
---|
917 | POcorr = sumMD/len(refl[11]) |
---|
918 | POderv[phfx+'MD'] = sumdMD/len(refl[11]) |
---|
919 | else: #spherical harmonics |
---|
920 | if calcControls[phfx+'SHord']: |
---|
921 | POcorr,POderv = SHPOcalDerv(refl,g,phfx,hfx,SGData,calcControls,parmDict) |
---|
922 | return POcorr,POderv |
---|
923 | |
---|
924 | def GetAbsorb(refl,hfx,calcControls,parmDict): |
---|
925 | 'Needs a doc string' |
---|
926 | if 'Debye' in calcControls[hfx+'instType']: |
---|
927 | return G2pwd.Absorb('Cylinder',parmDict[hfx+'Absorption'],refl[5],0,0) |
---|
928 | else: |
---|
929 | return 1.0 |
---|
930 | |
---|
931 | def GetAbsorbDerv(refl,hfx,calcControls,parmDict): |
---|
932 | 'Needs a doc string' |
---|
933 | if 'Debye' in calcControls[hfx+'instType']: |
---|
934 | return G2pwd.AbsorbDerv('Cylinder',parmDict[hfx+'Absorption'],refl[5],0,0) |
---|
935 | else: |
---|
936 | return 0.0 |
---|
937 | |
---|
938 | def GetIntensityCorr(refl,G,g,pfx,phfx,hfx,SGData,calcControls,parmDict): |
---|
939 | 'Needs a doc string' |
---|
940 | Icorr = parmDict[phfx+'Scale']*parmDict[hfx+'Scale']*refl[3] #scale*multiplicity |
---|
941 | if 'X' in parmDict[hfx+'Type']: |
---|
942 | Icorr *= G2pwd.Polarization(parmDict[hfx+'Polariz.'],refl[5],parmDict[hfx+'Azimuth'])[0] |
---|
943 | Icorr *= GetPrefOri(refl,G,g,phfx,hfx,SGData,calcControls,parmDict) |
---|
944 | if pfx+'SHorder' in parmDict: |
---|
945 | Icorr *= SHTXcal(refl,g,pfx,hfx,SGData,calcControls,parmDict) |
---|
946 | Icorr *= GetAbsorb(refl,hfx,calcControls,parmDict) |
---|
947 | refl[13] = Icorr |
---|
948 | |
---|
949 | def GetIntensityDerv(refl,G,g,pfx,phfx,hfx,SGData,calcControls,parmDict): |
---|
950 | 'Needs a doc string' |
---|
951 | dIdsh = 1./parmDict[hfx+'Scale'] |
---|
952 | dIdsp = 1./parmDict[phfx+'Scale'] |
---|
953 | if 'X' in parmDict[hfx+'Type']: |
---|
954 | pola,dIdPola = G2pwd.Polarization(parmDict[hfx+'Polariz.'],refl[5],parmDict[hfx+'Azimuth']) |
---|
955 | dIdPola /= pola |
---|
956 | else: #'N' |
---|
957 | dIdPola = 0.0 |
---|
958 | POcorr,dIdPO = GetPrefOriDerv(refl,G,g,phfx,hfx,SGData,calcControls,parmDict) |
---|
959 | for iPO in dIdPO: |
---|
960 | dIdPO[iPO] /= POcorr |
---|
961 | dFdODF = {} |
---|
962 | dFdSA = [0,0,0] |
---|
963 | if pfx+'SHorder' in parmDict: |
---|
964 | odfCor,dFdODF,dFdSA = SHTXcalDerv(refl,g,pfx,hfx,SGData,calcControls,parmDict) |
---|
965 | for iSH in dFdODF: |
---|
966 | dFdODF[iSH] /= odfCor |
---|
967 | for i in range(3): |
---|
968 | dFdSA[i] /= odfCor |
---|
969 | dFdAb = GetAbsorbDerv(refl,hfx,calcControls,parmDict) |
---|
970 | return dIdsh,dIdsp,dIdPola,dIdPO,dFdODF,dFdSA,dFdAb |
---|
971 | |
---|
972 | def GetSampleSigGam(refl,wave,G,GB,phfx,calcControls,parmDict): |
---|
973 | 'Needs a doc string' |
---|
974 | costh = cosd(refl[5]/2.) |
---|
975 | #crystallite size |
---|
976 | if calcControls[phfx+'SizeType'] == 'isotropic': |
---|
977 | Sgam = 1.8*wave/(np.pi*parmDict[phfx+'Size;i']*costh) |
---|
978 | elif calcControls[phfx+'SizeType'] == 'uniaxial': |
---|
979 | H = np.array(refl[:3]) |
---|
980 | P = np.array(calcControls[phfx+'SizeAxis']) |
---|
981 | cosP,sinP = G2lat.CosSinAngle(H,P,G) |
---|
982 | Sgam = (1.8*wave/np.pi)/(parmDict[phfx+'Size;i']*parmDict[phfx+'Size;a']*costh) |
---|
983 | Sgam *= np.sqrt((sinP*parmDict[phfx+'Size;a'])**2+(cosP*parmDict[phfx+'Size;i'])**2) |
---|
984 | else: #ellipsoidal crystallites |
---|
985 | Sij =[parmDict[phfx+'Size:%d'%(i)] for i in range(6)] |
---|
986 | H = np.array(refl[:3]) |
---|
987 | lenR = G2pwd.ellipseSize(H,Sij,GB) |
---|
988 | Sgam = 1.8*wave/(np.pi*costh*lenR) |
---|
989 | #microstrain |
---|
990 | if calcControls[phfx+'MustrainType'] == 'isotropic': |
---|
991 | Mgam = 0.018*parmDict[phfx+'Mustrain;i']*tand(refl[5]/2.)/np.pi |
---|
992 | elif calcControls[phfx+'MustrainType'] == 'uniaxial': |
---|
993 | H = np.array(refl[:3]) |
---|
994 | P = np.array(calcControls[phfx+'MustrainAxis']) |
---|
995 | cosP,sinP = G2lat.CosSinAngle(H,P,G) |
---|
996 | Si = parmDict[phfx+'Mustrain;i'] |
---|
997 | Sa = parmDict[phfx+'Mustrain;a'] |
---|
998 | Mgam = 0.018*Si*Sa*tand(refl[5]/2.)/(np.pi*np.sqrt((Si*cosP)**2+(Sa*sinP)**2)) |
---|
999 | else: #generalized - P.W. Stephens model |
---|
1000 | pwrs = calcControls[phfx+'MuPwrs'] |
---|
1001 | sum = 0 |
---|
1002 | for i,pwr in enumerate(pwrs): |
---|
1003 | sum += parmDict[phfx+'Mustrain:'+str(i)]*refl[0]**pwr[0]*refl[1]**pwr[1]*refl[2]**pwr[2] |
---|
1004 | Mgam = 0.018*refl[4]**2*tand(refl[5]/2.)*sum |
---|
1005 | gam = Sgam*parmDict[phfx+'Size;mx']+Mgam*parmDict[phfx+'Mustrain;mx'] |
---|
1006 | sig = (Sgam*(1.-parmDict[phfx+'Size;mx']))**2+(Mgam*(1.-parmDict[phfx+'Mustrain;mx']))**2 |
---|
1007 | sig /= ateln2 |
---|
1008 | return sig,gam |
---|
1009 | |
---|
1010 | def GetSampleSigGamDerv(refl,wave,G,GB,phfx,calcControls,parmDict): |
---|
1011 | 'Needs a doc string' |
---|
1012 | gamDict = {} |
---|
1013 | sigDict = {} |
---|
1014 | costh = cosd(refl[5]/2.) |
---|
1015 | tanth = tand(refl[5]/2.) |
---|
1016 | #crystallite size derivatives |
---|
1017 | if calcControls[phfx+'SizeType'] == 'isotropic': |
---|
1018 | Sgam = 1.8*wave/(np.pi*parmDict[phfx+'Size;i']*costh) |
---|
1019 | gamDict[phfx+'Size;i'] = -1.8*wave*parmDict[phfx+'Size;mx']/(np.pi*costh) |
---|
1020 | sigDict[phfx+'Size;i'] = -3.6*Sgam*wave*(1.-parmDict[phfx+'Size;mx'])**2/(np.pi*costh*ateln2) |
---|
1021 | elif calcControls[phfx+'SizeType'] == 'uniaxial': |
---|
1022 | H = np.array(refl[:3]) |
---|
1023 | P = np.array(calcControls[phfx+'SizeAxis']) |
---|
1024 | cosP,sinP = G2lat.CosSinAngle(H,P,G) |
---|
1025 | Si = parmDict[phfx+'Size;i'] |
---|
1026 | Sa = parmDict[phfx+'Size;a'] |
---|
1027 | gami = (1.8*wave/np.pi)/(Si*Sa) |
---|
1028 | sqtrm = np.sqrt((sinP*Sa)**2+(cosP*Si)**2) |
---|
1029 | Sgam = gami*sqtrm |
---|
1030 | gam = Sgam/costh |
---|
1031 | dsi = (gami*Si*cosP**2/(sqtrm*costh)-gam/Si) |
---|
1032 | dsa = (gami*Sa*sinP**2/(sqtrm*costh)-gam/Sa) |
---|
1033 | gamDict[phfx+'Size;i'] = dsi*parmDict[phfx+'Size;mx'] |
---|
1034 | gamDict[phfx+'Size;a'] = dsa*parmDict[phfx+'Size;mx'] |
---|
1035 | sigDict[phfx+'Size;i'] = 2.*dsi*Sgam*(1.-parmDict[phfx+'Size;mx'])**2/ateln2 |
---|
1036 | sigDict[phfx+'Size;a'] = 2.*dsa*Sgam*(1.-parmDict[phfx+'Size;mx'])**2/ateln2 |
---|
1037 | else: #ellipsoidal crystallites |
---|
1038 | const = 1.8*wave/(np.pi*costh) |
---|
1039 | Sij =[parmDict[phfx+'Size:%d'%(i)] for i in range(6)] |
---|
1040 | H = np.array(refl[:3]) |
---|
1041 | lenR,dRdS = G2pwd.ellipseSizeDerv(H,Sij,GB) |
---|
1042 | Sgam = 1.8*wave/(np.pi*costh*lenR) |
---|
1043 | for i,item in enumerate([phfx+'Size:%d'%(j) for j in range(6)]): |
---|
1044 | gamDict[item] = -(const/lenR**2)*dRdS[i]*parmDict[phfx+'Size;mx'] |
---|
1045 | sigDict[item] = -2.*Sgam*(const/lenR**2)*dRdS[i]*(1.-parmDict[phfx+'Size;mx'])**2/ateln2 |
---|
1046 | gamDict[phfx+'Size;mx'] = Sgam |
---|
1047 | sigDict[phfx+'Size;mx'] = -2.*Sgam**2*(1.-parmDict[phfx+'Size;mx'])/ateln2 |
---|
1048 | |
---|
1049 | #microstrain derivatives |
---|
1050 | if calcControls[phfx+'MustrainType'] == 'isotropic': |
---|
1051 | Mgam = 0.018*parmDict[phfx+'Mustrain;i']*tand(refl[5]/2.)/np.pi |
---|
1052 | gamDict[phfx+'Mustrain;i'] = 0.018*tanth*parmDict[phfx+'Mustrain;mx']/np.pi |
---|
1053 | sigDict[phfx+'Mustrain;i'] = 0.036*Mgam*tanth*(1.-parmDict[phfx+'Mustrain;mx'])**2/(np.pi*ateln2) |
---|
1054 | elif calcControls[phfx+'MustrainType'] == 'uniaxial': |
---|
1055 | H = np.array(refl[:3]) |
---|
1056 | P = np.array(calcControls[phfx+'MustrainAxis']) |
---|
1057 | cosP,sinP = G2lat.CosSinAngle(H,P,G) |
---|
1058 | Si = parmDict[phfx+'Mustrain;i'] |
---|
1059 | Sa = parmDict[phfx+'Mustrain;a'] |
---|
1060 | gami = 0.018*Si*Sa*tanth/np.pi |
---|
1061 | sqtrm = np.sqrt((Si*cosP)**2+(Sa*sinP)**2) |
---|
1062 | Mgam = gami/sqtrm |
---|
1063 | dsi = -gami*Si*cosP**2/sqtrm**3 |
---|
1064 | dsa = -gami*Sa*sinP**2/sqtrm**3 |
---|
1065 | gamDict[phfx+'Mustrain;i'] = (Mgam/Si+dsi)*parmDict[phfx+'Mustrain;mx'] |
---|
1066 | gamDict[phfx+'Mustrain;a'] = (Mgam/Sa+dsa)*parmDict[phfx+'Mustrain;mx'] |
---|
1067 | sigDict[phfx+'Mustrain;i'] = 2*(Mgam/Si+dsi)*Mgam*(1.-parmDict[phfx+'Mustrain;mx'])**2/ateln2 |
---|
1068 | sigDict[phfx+'Mustrain;a'] = 2*(Mgam/Sa+dsa)*Mgam*(1.-parmDict[phfx+'Mustrain;mx'])**2/ateln2 |
---|
1069 | else: #generalized - P.W. Stephens model |
---|
1070 | pwrs = calcControls[phfx+'MuPwrs'] |
---|
1071 | const = 0.018*refl[4]**2*tanth |
---|
1072 | sum = 0 |
---|
1073 | for i,pwr in enumerate(pwrs): |
---|
1074 | term = refl[0]**pwr[0]*refl[1]**pwr[1]*refl[2]**pwr[2] |
---|
1075 | sum += parmDict[phfx+'Mustrain:'+str(i)]*term |
---|
1076 | gamDict[phfx+'Mustrain:'+str(i)] = const*term*parmDict[phfx+'Mustrain;mx'] |
---|
1077 | sigDict[phfx+'Mustrain:'+str(i)] = \ |
---|
1078 | 2.*const*term*(1.-parmDict[phfx+'Mustrain;mx'])**2/ateln2 |
---|
1079 | Mgam = 0.018*refl[4]**2*tand(refl[5]/2.)*sum |
---|
1080 | for i in range(len(pwrs)): |
---|
1081 | sigDict[phfx+'Mustrain:'+str(i)] *= Mgam |
---|
1082 | gamDict[phfx+'Mustrain;mx'] = Mgam |
---|
1083 | sigDict[phfx+'Mustrain;mx'] = -2.*Mgam**2*(1.-parmDict[phfx+'Mustrain;mx'])/ateln2 |
---|
1084 | return sigDict,gamDict |
---|
1085 | |
---|
1086 | def GetReflPos(refl,wave,G,hfx,calcControls,parmDict): |
---|
1087 | 'Needs a doc string' |
---|
1088 | h,k,l = refl[:3] |
---|
1089 | dsq = 1./G2lat.calc_rDsq2(np.array([h,k,l]),G) |
---|
1090 | d = np.sqrt(dsq) |
---|
1091 | |
---|
1092 | refl[4] = d |
---|
1093 | pos = 2.0*asind(wave/(2.0*d))+parmDict[hfx+'Zero'] |
---|
1094 | const = 9.e-2/(np.pi*parmDict[hfx+'Gonio. radius']) #shifts in microns |
---|
1095 | if 'Bragg' in calcControls[hfx+'instType']: |
---|
1096 | pos -= const*(4.*parmDict[hfx+'Shift']*cosd(pos/2.0)+ \ |
---|
1097 | parmDict[hfx+'Transparency']*sind(pos)*100.0) #trans(=1/mueff) in cm |
---|
1098 | else: #Debye-Scherrer - simple but maybe not right |
---|
1099 | pos -= const*(parmDict[hfx+'DisplaceX']*cosd(pos)+parmDict[hfx+'DisplaceY']*sind(pos)) |
---|
1100 | return pos |
---|
1101 | |
---|
1102 | def GetReflPosDerv(refl,wave,A,hfx,calcControls,parmDict): |
---|
1103 | 'Needs a doc string' |
---|
1104 | dpr = 180./np.pi |
---|
1105 | h,k,l = refl[:3] |
---|
1106 | dstsq = G2lat.calc_rDsq(np.array([h,k,l]),A) |
---|
1107 | dst = np.sqrt(dstsq) |
---|
1108 | pos = refl[5]-parmDict[hfx+'Zero'] |
---|
1109 | const = dpr/np.sqrt(1.0-wave**2*dstsq/4.0) |
---|
1110 | dpdw = const*dst |
---|
1111 | dpdA = np.array([h**2,k**2,l**2,h*k,h*l,k*l]) |
---|
1112 | dpdA *= const*wave/(2.0*dst) |
---|
1113 | dpdZ = 1.0 |
---|
1114 | const = 9.e-2/(np.pi*parmDict[hfx+'Gonio. radius']) #shifts in microns |
---|
1115 | if 'Bragg' in calcControls[hfx+'instType']: |
---|
1116 | dpdSh = -4.*const*cosd(pos/2.0) |
---|
1117 | dpdTr = -const*sind(pos)*100.0 |
---|
1118 | return dpdA,dpdw,dpdZ,dpdSh,dpdTr,0.,0. |
---|
1119 | else: #Debye-Scherrer - simple but maybe not right |
---|
1120 | dpdXd = -const*cosd(pos) |
---|
1121 | dpdYd = -const*sind(pos) |
---|
1122 | return dpdA,dpdw,dpdZ,0.,0.,dpdXd,dpdYd |
---|
1123 | |
---|
1124 | def GetHStrainShift(refl,SGData,phfx,parmDict): |
---|
1125 | 'Needs a doc string' |
---|
1126 | laue = SGData['SGLaue'] |
---|
1127 | uniq = SGData['SGUniq'] |
---|
1128 | h,k,l = refl[:3] |
---|
1129 | if laue in ['m3','m3m']: |
---|
1130 | Dij = parmDict[phfx+'D11']*(h**2+k**2+l**2)+ \ |
---|
1131 | refl[4]**2*parmDict[phfx+'eA']*((h*k)**2+(h*l)**2+(k*l)**2)/(h**2+k**2+l**2)**2 |
---|
1132 | elif laue in ['6/m','6/mmm','3m1','31m','3']: |
---|
1133 | Dij = parmDict[phfx+'D11']*(h**2+k**2+h*k)+parmDict[phfx+'D33']*l**2 |
---|
1134 | elif laue in ['3R','3mR']: |
---|
1135 | Dij = parmDict[phfx+'D11']*(h**2+k**2+l**2)+parmDict[phfx+'D12']*(h*k+h*l+k*l) |
---|
1136 | elif laue in ['4/m','4/mmm']: |
---|
1137 | Dij = parmDict[phfx+'D11']*(h**2+k**2)+parmDict[phfx+'D33']*l**2 |
---|
1138 | elif laue in ['mmm']: |
---|
1139 | Dij = parmDict[phfx+'D11']*h**2+parmDict[phfx+'D22']*k**2+parmDict[phfx+'D33']*l**2 |
---|
1140 | elif laue in ['2/m']: |
---|
1141 | Dij = parmDict[phfx+'D11']*h**2+parmDict[phfx+'D22']*k**2+parmDict[phfx+'D33']*l**2 |
---|
1142 | if uniq == 'a': |
---|
1143 | Dij += parmDict[phfx+'D23']*k*l |
---|
1144 | elif uniq == 'b': |
---|
1145 | Dij += parmDict[phfx+'D13']*h*l |
---|
1146 | elif uniq == 'c': |
---|
1147 | Dij += parmDict[phfx+'D12']*h*k |
---|
1148 | else: |
---|
1149 | Dij = parmDict[phfx+'D11']*h**2+parmDict[phfx+'D22']*k**2+parmDict[phfx+'D33']*l**2+ \ |
---|
1150 | parmDict[phfx+'D12']*h*k+parmDict[phfx+'D13']*h*l+parmDict[phfx+'D23']*k*l |
---|
1151 | return -Dij*refl[4]**2*tand(refl[5]/2.0) |
---|
1152 | |
---|
1153 | def GetHStrainShiftDerv(refl,SGData,phfx): |
---|
1154 | 'Needs a doc string' |
---|
1155 | laue = SGData['SGLaue'] |
---|
1156 | uniq = SGData['SGUniq'] |
---|
1157 | h,k,l = refl[:3] |
---|
1158 | if laue in ['m3','m3m']: |
---|
1159 | dDijDict = {phfx+'D11':h**2+k**2+l**2, |
---|
1160 | phfx+'eA':refl[4]**2*((h*k)**2+(h*l)**2+(k*l)**2)/(h**2+k**2+l**2)**2} |
---|
1161 | elif laue in ['6/m','6/mmm','3m1','31m','3']: |
---|
1162 | dDijDict = {phfx+'D11':h**2+k**2+h*k,phfx+'D33':l**2} |
---|
1163 | elif laue in ['3R','3mR']: |
---|
1164 | dDijDict = {phfx+'D11':h**2+k**2+l**2,phfx+'D12':h*k+h*l+k*l} |
---|
1165 | elif laue in ['4/m','4/mmm']: |
---|
1166 | dDijDict = {phfx+'D11':h**2+k**2,phfx+'D33':l**2} |
---|
1167 | elif laue in ['mmm']: |
---|
1168 | dDijDict = {phfx+'D11':h**2,phfx+'D22':k**2,phfx+'D33':l**2} |
---|
1169 | elif laue in ['2/m']: |
---|
1170 | dDijDict = {phfx+'D11':h**2,phfx+'D22':k**2,phfx+'D33':l**2} |
---|
1171 | if uniq == 'a': |
---|
1172 | dDijDict[phfx+'D23'] = k*l |
---|
1173 | elif uniq == 'b': |
---|
1174 | dDijDict[phfx+'D13'] = h*l |
---|
1175 | elif uniq == 'c': |
---|
1176 | dDijDict[phfx+'D12'] = h*k |
---|
1177 | names.append() |
---|
1178 | else: |
---|
1179 | dDijDict = {phfx+'D11':h**2,phfx+'D22':k**2,phfx+'D33':l**2, |
---|
1180 | phfx+'D12':h*k,phfx+'D13':h*l,phfx+'D23':k*l} |
---|
1181 | for item in dDijDict: |
---|
1182 | dDijDict[item] *= -refl[4]**2*tand(refl[5]/2.0) |
---|
1183 | return dDijDict |
---|
1184 | |
---|
1185 | def GetFobsSq(Histograms,Phases,parmDict,calcControls): |
---|
1186 | 'Needs a doc string' |
---|
1187 | histoList = Histograms.keys() |
---|
1188 | histoList.sort() |
---|
1189 | for histogram in histoList: |
---|
1190 | if 'PWDR' in histogram[:4]: |
---|
1191 | Histogram = Histograms[histogram] |
---|
1192 | hId = Histogram['hId'] |
---|
1193 | hfx = ':%d:'%(hId) |
---|
1194 | Limits = calcControls[hfx+'Limits'] |
---|
1195 | shl = max(parmDict[hfx+'SH/L'],0.0005) |
---|
1196 | Ka2 = False |
---|
1197 | kRatio = 0.0 |
---|
1198 | if hfx+'Lam1' in parmDict.keys(): |
---|
1199 | Ka2 = True |
---|
1200 | lamRatio = 360*(parmDict[hfx+'Lam2']-parmDict[hfx+'Lam1'])/(np.pi*parmDict[hfx+'Lam1']) |
---|
1201 | kRatio = parmDict[hfx+'I(L2)/I(L1)'] |
---|
1202 | x,y,w,yc,yb,yd = Histogram['Data'] |
---|
1203 | xB = np.searchsorted(x,Limits[0]) |
---|
1204 | xF = np.searchsorted(x,Limits[1]) |
---|
1205 | ymb = np.array(y-yb) |
---|
1206 | ymb = np.where(ymb,ymb,1.0) |
---|
1207 | ycmb = np.array(yc-yb) |
---|
1208 | ratio = 1./np.where(ycmb,ycmb/ymb,1.e10) |
---|
1209 | refLists = Histogram['Reflection Lists'] |
---|
1210 | for phase in refLists: |
---|
1211 | Phase = Phases[phase] |
---|
1212 | pId = Phase['pId'] |
---|
1213 | phfx = '%d:%d:'%(pId,hId) |
---|
1214 | refList = refLists[phase] |
---|
1215 | sumFo = 0.0 |
---|
1216 | sumdF = 0.0 |
---|
1217 | sumFosq = 0.0 |
---|
1218 | sumdFsq = 0.0 |
---|
1219 | for refl in refList: |
---|
1220 | if 'C' in calcControls[hfx+'histType']: |
---|
1221 | yp = np.zeros_like(yb) |
---|
1222 | Wd,fmin,fmax = G2pwd.getWidthsCW(refl[5],refl[6],refl[7],shl) |
---|
1223 | iBeg = max(xB,np.searchsorted(x,refl[5]-fmin)) |
---|
1224 | iFin = max(xB,min(np.searchsorted(x,refl[5]+fmax),xF)) |
---|
1225 | iFin2 = iFin |
---|
1226 | yp[iBeg:iFin] = refl[13]*refl[9]*G2pwd.getFCJVoigt3(refl[5],refl[6],refl[7],shl,x[iBeg:iFin]) #>90% of time spent here |
---|
1227 | if Ka2: |
---|
1228 | pos2 = refl[5]+lamRatio*tand(refl[5]/2.0) # + 360/pi * Dlam/lam * tan(th) |
---|
1229 | Wd,fmin,fmax = G2pwd.getWidthsCW(pos2,refl[6],refl[7],shl) |
---|
1230 | iBeg2 = max(xB,np.searchsorted(x,pos2-fmin)) |
---|
1231 | iFin2 = min(np.searchsorted(x,pos2+fmax),xF) |
---|
1232 | if not iBeg2+iFin2: #peak below low limit - skip peak |
---|
1233 | continue |
---|
1234 | elif not iBeg2-iFin2: #peak above high limit - done |
---|
1235 | break |
---|
1236 | yp[iBeg2:iFin2] += refl[13]*refl[9]*kRatio*G2pwd.getFCJVoigt3(pos2,refl[6],refl[7],shl,x[iBeg2:iFin2]) #and here |
---|
1237 | refl[8] = np.sum(np.where(ratio[iBeg:iFin2]>0.,yp[iBeg:iFin2]*ratio[iBeg:iFin2]/(refl[13]*(1.+kRatio)),0.0)) |
---|
1238 | elif 'T' in calcControls[hfx+'histType']: |
---|
1239 | print 'TOF Undefined at present' |
---|
1240 | raise Exception #no TOF yet |
---|
1241 | Fo = np.sqrt(np.abs(refl[8])) |
---|
1242 | Fc = np.sqrt(np.abs(refl[9])) |
---|
1243 | sumFo += Fo |
---|
1244 | sumFosq += refl[8]**2 |
---|
1245 | sumdF += np.abs(Fo-Fc) |
---|
1246 | sumdFsq += (refl[8]-refl[9])**2 |
---|
1247 | Histogram['Residuals'][phfx+'Rf'] = min(100.,(sumdF/sumFo)*100.) |
---|
1248 | Histogram['Residuals'][phfx+'Rf^2'] = min(100.,np.sqrt(sumdFsq/sumFosq)*100.) |
---|
1249 | Histogram['Residuals'][phfx+'Nref'] = len(refList) |
---|
1250 | Histogram['Residuals']['hId'] = hId |
---|
1251 | elif 'HKLF' in histogram[:4]: |
---|
1252 | Histogram = Histograms[histogram] |
---|
1253 | Histogram['Residuals']['hId'] = Histograms[histogram]['hId'] |
---|
1254 | |
---|
1255 | def getPowderProfile(parmDict,x,varylist,Histogram,Phases,calcControls,pawleyLookup): |
---|
1256 | 'Needs a doc string' |
---|
1257 | |
---|
1258 | def GetReflSigGam(refl,wave,G,GB,hfx,phfx,calcControls,parmDict): |
---|
1259 | U = parmDict[hfx+'U'] |
---|
1260 | V = parmDict[hfx+'V'] |
---|
1261 | W = parmDict[hfx+'W'] |
---|
1262 | X = parmDict[hfx+'X'] |
---|
1263 | Y = parmDict[hfx+'Y'] |
---|
1264 | tanPos = tand(refl[5]/2.0) |
---|
1265 | Ssig,Sgam = GetSampleSigGam(refl,wave,G,GB,phfx,calcControls,parmDict) |
---|
1266 | sig = U*tanPos**2+V*tanPos+W+Ssig #save peak sigma |
---|
1267 | sig = max(0.001,sig) |
---|
1268 | gam = X/cosd(refl[5]/2.0)+Y*tanPos+Sgam #save peak gamma |
---|
1269 | gam = max(0.001,gam) |
---|
1270 | return sig,gam |
---|
1271 | |
---|
1272 | hId = Histogram['hId'] |
---|
1273 | hfx = ':%d:'%(hId) |
---|
1274 | bakType = calcControls[hfx+'bakType'] |
---|
1275 | yb = G2pwd.getBackground(hfx,parmDict,bakType,x) |
---|
1276 | yc = np.zeros_like(yb) |
---|
1277 | |
---|
1278 | if 'C' in calcControls[hfx+'histType']: |
---|
1279 | shl = max(parmDict[hfx+'SH/L'],0.002) |
---|
1280 | Ka2 = False |
---|
1281 | if hfx+'Lam1' in parmDict.keys(): |
---|
1282 | wave = parmDict[hfx+'Lam1'] |
---|
1283 | Ka2 = True |
---|
1284 | lamRatio = 360*(parmDict[hfx+'Lam2']-parmDict[hfx+'Lam1'])/(np.pi*parmDict[hfx+'Lam1']) |
---|
1285 | kRatio = parmDict[hfx+'I(L2)/I(L1)'] |
---|
1286 | else: |
---|
1287 | wave = parmDict[hfx+'Lam'] |
---|
1288 | else: |
---|
1289 | print 'TOF Undefined at present' |
---|
1290 | raise ValueError |
---|
1291 | for phase in Histogram['Reflection Lists']: |
---|
1292 | refList = Histogram['Reflection Lists'][phase] |
---|
1293 | Phase = Phases[phase] |
---|
1294 | pId = Phase['pId'] |
---|
1295 | pfx = '%d::'%(pId) |
---|
1296 | phfx = '%d:%d:'%(pId,hId) |
---|
1297 | hfx = ':%d:'%(hId) |
---|
1298 | SGData = Phase['General']['SGData'] |
---|
1299 | A = [parmDict[pfx+'A%d'%(i)] for i in range(6)] |
---|
1300 | G,g = G2lat.A2Gmat(A) #recip & real metric tensors |
---|
1301 | GA,GB = G2lat.Gmat2AB(G) #Orthogonalization matricies |
---|
1302 | Vst = np.sqrt(nl.det(G)) #V* |
---|
1303 | if not Phase['General'].get('doPawley'): |
---|
1304 | time0 = time.time() |
---|
1305 | StructureFactor(refList,G,hfx,pfx,SGData,calcControls,parmDict) |
---|
1306 | # print 'sf calc time: %.3fs'%(time.time()-time0) |
---|
1307 | time0 = time.time() |
---|
1308 | for refl in refList: |
---|
1309 | if 'C' in calcControls[hfx+'histType']: |
---|
1310 | h,k,l = refl[:3] |
---|
1311 | refl[5] = GetReflPos(refl,wave,G,hfx,calcControls,parmDict) #corrected reflection position |
---|
1312 | Lorenz = 1./(2.*sind(refl[5]/2.)**2*cosd(refl[5]/2.)) #Lorentz correction |
---|
1313 | refl[5] += GetHStrainShift(refl,SGData,phfx,parmDict) #apply hydrostatic strain shift |
---|
1314 | refl[6:8] = GetReflSigGam(refl,wave,G,GB,hfx,phfx,calcControls,parmDict) #peak sig & gam |
---|
1315 | GetIntensityCorr(refl,G,g,pfx,phfx,hfx,SGData,calcControls,parmDict) #puts corrections in refl[13] |
---|
1316 | refl[13] *= Vst*Lorenz |
---|
1317 | if Phase['General'].get('doPawley'): |
---|
1318 | try: |
---|
1319 | pInd =pfx+'PWLref:%d'%(pawleyLookup[pfx+'%d,%d,%d'%(h,k,l)]) |
---|
1320 | refl[9] = parmDict[pInd] |
---|
1321 | except KeyError: |
---|
1322 | # print ' ***Error %d,%d,%d missing from Pawley reflection list ***'%(h,k,l) |
---|
1323 | continue |
---|
1324 | Wd,fmin,fmax = G2pwd.getWidthsCW(refl[5],refl[6],refl[7],shl) |
---|
1325 | iBeg = np.searchsorted(x,refl[5]-fmin) |
---|
1326 | iFin = np.searchsorted(x,refl[5]+fmax) |
---|
1327 | if not iBeg+iFin: #peak below low limit - skip peak |
---|
1328 | continue |
---|
1329 | elif not iBeg-iFin: #peak above high limit - done |
---|
1330 | break |
---|
1331 | yc[iBeg:iFin] += refl[13]*refl[9]*G2pwd.getFCJVoigt3(refl[5],refl[6],refl[7],shl,ma.getdata(x[iBeg:iFin])) #>90% of time spent here |
---|
1332 | if Ka2: |
---|
1333 | pos2 = refl[5]+lamRatio*tand(refl[5]/2.0) # + 360/pi * Dlam/lam * tan(th) |
---|
1334 | Wd,fmin,fmax = G2pwd.getWidthsCW(pos2,refl[6],refl[7],shl) |
---|
1335 | iBeg = np.searchsorted(x,pos2-fmin) |
---|
1336 | iFin = np.searchsorted(x,pos2+fmax) |
---|
1337 | if not iBeg+iFin: #peak below low limit - skip peak |
---|
1338 | continue |
---|
1339 | elif not iBeg-iFin: #peak above high limit - done |
---|
1340 | return yc,yb |
---|
1341 | yc[iBeg:iFin] += refl[13]*refl[9]*kRatio*G2pwd.getFCJVoigt3(pos2,refl[6],refl[7],shl,ma.getdata(x[iBeg:iFin])) #and here |
---|
1342 | elif 'T' in calcControls[hfx+'histType']: |
---|
1343 | print 'TOF Undefined at present' |
---|
1344 | raise Exception #no TOF yet |
---|
1345 | # print 'profile calc time: %.3fs'%(time.time()-time0) |
---|
1346 | return yc,yb |
---|
1347 | |
---|
1348 | def getPowderProfileDerv(parmDict,x,varylist,Histogram,Phases,rigidbodyDict,calcControls,pawleyLookup): |
---|
1349 | 'Needs a doc string' |
---|
1350 | |
---|
1351 | def cellVaryDerv(pfx,SGData,dpdA): |
---|
1352 | if SGData['SGLaue'] in ['-1',]: |
---|
1353 | return [[pfx+'A0',dpdA[0]],[pfx+'A1',dpdA[1]],[pfx+'A2',dpdA[2]], |
---|
1354 | [pfx+'A3',dpdA[3]],[pfx+'A4',dpdA[4]],[pfx+'A5',dpdA[5]]] |
---|
1355 | elif SGData['SGLaue'] in ['2/m',]: |
---|
1356 | if SGData['SGUniq'] == 'a': |
---|
1357 | return [[pfx+'A0',dpdA[0]],[pfx+'A1',dpdA[1]],[pfx+'A2',dpdA[2]],[pfx+'A3',dpdA[3]]] |
---|
1358 | elif SGData['SGUniq'] == 'b': |
---|
1359 | return [[pfx+'A0',dpdA[0]],[pfx+'A1',dpdA[1]],[pfx+'A2',dpdA[2]],[pfx+'A4',dpdA[4]]] |
---|
1360 | else: |
---|
1361 | return [[pfx+'A0',dpdA[0]],[pfx+'A1',dpdA[1]],[pfx+'A2',dpdA[2]],[pfx+'A5',dpdA[5]]] |
---|
1362 | elif SGData['SGLaue'] in ['mmm',]: |
---|
1363 | return [[pfx+'A0',dpdA[0]],[pfx+'A1',dpdA[1]],[pfx+'A2',dpdA[2]]] |
---|
1364 | elif SGData['SGLaue'] in ['4/m','4/mmm']: |
---|
1365 | return [[pfx+'A0',dpdA[0]],[pfx+'A2',dpdA[2]]] |
---|
1366 | elif SGData['SGLaue'] in ['6/m','6/mmm','3m1', '31m', '3']: |
---|
1367 | return [[pfx+'A0',dpdA[0]],[pfx+'A2',dpdA[2]]] |
---|
1368 | elif SGData['SGLaue'] in ['3R', '3mR']: |
---|
1369 | return [[pfx+'A0',dpdA[0]+dpdA[1]+dpdA[2]],[pfx+'A3',dpdA[3]+dpdA[4]+dpdA[5]]] |
---|
1370 | elif SGData['SGLaue'] in ['m3m','m3']: |
---|
1371 | return [[pfx+'A0',dpdA[0]]] |
---|
1372 | |
---|
1373 | # create a list of dependent variables and set up a dictionary to hold their derivatives |
---|
1374 | dependentVars = G2mv.GetDependentVars() |
---|
1375 | depDerivDict = {} |
---|
1376 | for j in dependentVars: |
---|
1377 | depDerivDict[j] = np.zeros(shape=(len(x))) |
---|
1378 | #print 'dependent vars',dependentVars |
---|
1379 | lenX = len(x) |
---|
1380 | hId = Histogram['hId'] |
---|
1381 | hfx = ':%d:'%(hId) |
---|
1382 | bakType = calcControls[hfx+'bakType'] |
---|
1383 | dMdv = np.zeros(shape=(len(varylist),len(x))) |
---|
1384 | dMdb,dMddb,dMdpk = G2pwd.getBackgroundDerv(hfx,parmDict,bakType,x) |
---|
1385 | if hfx+'Back:0' in varylist: # for now assume that Back:x vars to not appear in constraints |
---|
1386 | bBpos =varylist.index(hfx+'Back:0') |
---|
1387 | dMdv[bBpos:bBpos+len(dMdb)] = dMdb |
---|
1388 | names = [hfx+'DebyeA',hfx+'DebyeR',hfx+'DebyeU'] |
---|
1389 | for name in varylist: |
---|
1390 | if 'Debye' in name: |
---|
1391 | id = int(name.split(':')[-1]) |
---|
1392 | parm = name[:int(name.rindex(':'))] |
---|
1393 | ip = names.index(parm) |
---|
1394 | dMdv[varylist.index(name)] = dMddb[3*id+ip] |
---|
1395 | names = [hfx+'BkPkpos',hfx+'BkPkint',hfx+'BkPksig',hfx+'BkPkgam'] |
---|
1396 | for name in varylist: |
---|
1397 | if 'BkPk' in name: |
---|
1398 | id = int(name.split(':')[-1]) |
---|
1399 | parm = name[:int(name.rindex(':'))] |
---|
1400 | ip = names.index(parm) |
---|
1401 | dMdv[varylist.index(name)] = dMdpk[4*id+ip] |
---|
1402 | cw = np.diff(x) |
---|
1403 | cw = np.append(cw,cw[-1]) |
---|
1404 | if 'C' in calcControls[hfx+'histType']: |
---|
1405 | shl = max(parmDict[hfx+'SH/L'],0.002) |
---|
1406 | Ka2 = False |
---|
1407 | if hfx+'Lam1' in parmDict.keys(): |
---|
1408 | wave = parmDict[hfx+'Lam1'] |
---|
1409 | Ka2 = True |
---|
1410 | lamRatio = 360*(parmDict[hfx+'Lam2']-parmDict[hfx+'Lam1'])/(np.pi*parmDict[hfx+'Lam1']) |
---|
1411 | kRatio = parmDict[hfx+'I(L2)/I(L1)'] |
---|
1412 | else: |
---|
1413 | wave = parmDict[hfx+'Lam'] |
---|
1414 | else: |
---|
1415 | print 'TOF Undefined at present' |
---|
1416 | raise ValueError |
---|
1417 | for phase in Histogram['Reflection Lists']: |
---|
1418 | refList = Histogram['Reflection Lists'][phase] |
---|
1419 | Phase = Phases[phase] |
---|
1420 | SGData = Phase['General']['SGData'] |
---|
1421 | pId = Phase['pId'] |
---|
1422 | pfx = '%d::'%(pId) |
---|
1423 | phfx = '%d:%d:'%(pId,hId) |
---|
1424 | A = [parmDict[pfx+'A%d'%(i)] for i in range(6)] |
---|
1425 | G,g = G2lat.A2Gmat(A) #recip & real metric tensors |
---|
1426 | GA,GB = G2lat.Gmat2AB(G) #Orthogonalization matricies |
---|
1427 | if not Phase['General'].get('doPawley'): |
---|
1428 | time0 = time.time() |
---|
1429 | dFdvDict = StructureFactorDerv(refList,G,hfx,pfx,SGData,calcControls,parmDict) |
---|
1430 | # print 'sf-derv time %.3fs'%(time.time()-time0) |
---|
1431 | ApplyRBModelDervs(dFdvDict,parmDict,rigidbodyDict,Phase) |
---|
1432 | time0 = time.time() |
---|
1433 | for iref,refl in enumerate(refList): |
---|
1434 | if 'C' in calcControls[hfx+'histType']: #CW powder |
---|
1435 | h,k,l = refl[:3] |
---|
1436 | dIdsh,dIdsp,dIdpola,dIdPO,dFdODF,dFdSA,dFdAb = GetIntensityDerv(refl,G,g,pfx,phfx,hfx,SGData,calcControls,parmDict) |
---|
1437 | Wd,fmin,fmax = G2pwd.getWidthsCW(refl[5],refl[6],refl[7],shl) |
---|
1438 | iBeg = np.searchsorted(x,refl[5]-fmin) |
---|
1439 | iFin = np.searchsorted(x,refl[5]+fmax) |
---|
1440 | if not iBeg+iFin: #peak below low limit - skip peak |
---|
1441 | continue |
---|
1442 | elif not iBeg-iFin: #peak above high limit - done |
---|
1443 | break |
---|
1444 | pos = refl[5] |
---|
1445 | tanth = tand(pos/2.0) |
---|
1446 | costh = cosd(pos/2.0) |
---|
1447 | lenBF = iFin-iBeg |
---|
1448 | dMdpk = np.zeros(shape=(6,lenBF)) |
---|
1449 | dMdipk = G2pwd.getdFCJVoigt3(refl[5],refl[6],refl[7],shl,ma.getdata(x[iBeg:iFin])) |
---|
1450 | for i in range(5): |
---|
1451 | dMdpk[i] += 100.*cw[iBeg:iFin]*refl[13]*refl[9]*dMdipk[i] |
---|
1452 | dervDict = {'int':dMdpk[0],'pos':dMdpk[1],'sig':dMdpk[2],'gam':dMdpk[3],'shl':dMdpk[4],'L1/L2':np.zeros_like(dMdpk[0])} |
---|
1453 | if Ka2: |
---|
1454 | pos2 = refl[5]+lamRatio*tanth # + 360/pi * Dlam/lam * tan(th) |
---|
1455 | iBeg2 = np.searchsorted(x,pos2-fmin) |
---|
1456 | iFin2 = np.searchsorted(x,pos2+fmax) |
---|
1457 | if iBeg2-iFin2: |
---|
1458 | lenBF2 = iFin2-iBeg2 |
---|
1459 | dMdpk2 = np.zeros(shape=(6,lenBF2)) |
---|
1460 | dMdipk2 = G2pwd.getdFCJVoigt3(pos2,refl[6],refl[7],shl,ma.getdata(x[iBeg2:iFin2])) |
---|
1461 | for i in range(5): |
---|
1462 | dMdpk2[i] = 100.*cw[iBeg2:iFin2]*refl[13]*refl[9]*kRatio*dMdipk2[i] |
---|
1463 | dMdpk2[5] = 100.*cw[iBeg2:iFin2]*refl[13]*dMdipk2[0] |
---|
1464 | dervDict2 = {'int':dMdpk2[0],'pos':dMdpk2[1],'sig':dMdpk2[2],'gam':dMdpk2[3],'shl':dMdpk2[4],'L1/L2':dMdpk2[5]*refl[9]} |
---|
1465 | if Phase['General'].get('doPawley'): |
---|
1466 | dMdpw = np.zeros(len(x)) |
---|
1467 | try: |
---|
1468 | pIdx = pfx+'PWLref:'+str(pawleyLookup[pfx+'%d,%d,%d'%(h,k,l)]) |
---|
1469 | idx = varylist.index(pIdx) |
---|
1470 | dMdpw[iBeg:iFin] = dervDict['int']/refl[9] |
---|
1471 | if Ka2: |
---|
1472 | dMdpw[iBeg2:iFin2] += dervDict2['int']/refl[9] |
---|
1473 | dMdv[idx] = dMdpw |
---|
1474 | except: # ValueError: |
---|
1475 | pass |
---|
1476 | dpdA,dpdw,dpdZ,dpdSh,dpdTr,dpdX,dpdY = GetReflPosDerv(refl,wave,A,hfx,calcControls,parmDict) |
---|
1477 | names = {hfx+'Scale':[dIdsh,'int'],hfx+'Polariz.':[dIdpola,'int'],phfx+'Scale':[dIdsp,'int'], |
---|
1478 | hfx+'U':[tanth**2,'sig'],hfx+'V':[tanth,'sig'],hfx+'W':[1.0,'sig'], |
---|
1479 | hfx+'X':[1.0/costh,'gam'],hfx+'Y':[tanth,'gam'],hfx+'SH/L':[1.0,'shl'], |
---|
1480 | hfx+'I(L2)/I(L1)':[1.0,'L1/L2'],hfx+'Zero':[dpdZ,'pos'],hfx+'Lam':[dpdw,'pos'], |
---|
1481 | hfx+'Shift':[dpdSh,'pos'],hfx+'Transparency':[dpdTr,'pos'],hfx+'DisplaceX':[dpdX,'pos'], |
---|
1482 | hfx+'DisplaceY':[dpdY,'pos'],hfx+'Absorption':[dFdAb,'int'],} |
---|
1483 | for name in names: |
---|
1484 | item = names[name] |
---|
1485 | if name in varylist: |
---|
1486 | dMdv[varylist.index(name)][iBeg:iFin] += item[0]*dervDict[item[1]] |
---|
1487 | if Ka2: |
---|
1488 | dMdv[varylist.index(name)][iBeg2:iFin2] += item[0]*dervDict2[item[1]] |
---|
1489 | elif name in dependentVars: |
---|
1490 | depDerivDict[name][iBeg:iFin] += item[0]*dervDict[item[1]] |
---|
1491 | if Ka2: |
---|
1492 | depDerivDict[name][iBeg2:iFin2] += item[0]*dervDict2[item[1]] |
---|
1493 | for iPO in dIdPO: |
---|
1494 | if iPO in varylist: |
---|
1495 | dMdv[varylist.index(iPO)][iBeg:iFin] += dIdPO[iPO]*dervDict['int'] |
---|
1496 | if Ka2: |
---|
1497 | dMdv[varylist.index(iPO)][iBeg2:iFin2] += dIdPO[iPO]*dervDict2['int'] |
---|
1498 | elif iPO in dependentVars: |
---|
1499 | depDerivDict[iPO][iBeg:iFin] += dIdPO[iPO]*dervDict['int'] |
---|
1500 | if Ka2: |
---|
1501 | depDerivDict[iPO][iBeg2:iFin2] += dIdPO[iPO]*dervDict2['int'] |
---|
1502 | for i,name in enumerate(['omega','chi','phi']): |
---|
1503 | aname = pfx+'SH '+name |
---|
1504 | if aname in varylist: |
---|
1505 | dMdv[varylist.index(aname)][iBeg:iFin] += dFdSA[i]*dervDict['int'] |
---|
1506 | if Ka2: |
---|
1507 | dMdv[varylist.index(aname)][iBeg2:iFin2] += dFdSA[i]*dervDict2['int'] |
---|
1508 | elif aname in dependentVars: |
---|
1509 | depDerivDict[aname][iBeg:iFin] += dFdSA[i]*dervDict['int'] |
---|
1510 | if Ka2: |
---|
1511 | depDerivDict[aname][iBeg2:iFin2] += dFdSA[i]*dervDict2['int'] |
---|
1512 | for iSH in dFdODF: |
---|
1513 | if iSH in varylist: |
---|
1514 | dMdv[varylist.index(iSH)][iBeg:iFin] += dFdODF[iSH]*dervDict['int'] |
---|
1515 | if Ka2: |
---|
1516 | dMdv[varylist.index(iSH)][iBeg2:iFin2] += dFdODF[iSH]*dervDict2['int'] |
---|
1517 | elif iSH in dependentVars: |
---|
1518 | depDerivDict[iSH][iBeg:iFin] += dFdODF[iSH]*dervDict['int'] |
---|
1519 | if Ka2: |
---|
1520 | depDerivDict[iSH][iBeg2:iFin2] += dFdODF[iSH]*dervDict2['int'] |
---|
1521 | cellDervNames = cellVaryDerv(pfx,SGData,dpdA) |
---|
1522 | for name,dpdA in cellDervNames: |
---|
1523 | if name in varylist: |
---|
1524 | dMdv[varylist.index(name)][iBeg:iFin] += dpdA*dervDict['pos'] |
---|
1525 | if Ka2: |
---|
1526 | dMdv[varylist.index(name)][iBeg2:iFin2] += dpdA*dervDict2['pos'] |
---|
1527 | elif name in dependentVars: |
---|
1528 | depDerivDict[name][iBeg:iFin] += dpdA*dervDict['pos'] |
---|
1529 | if Ka2: |
---|
1530 | depDerivDict[name][iBeg2:iFin2] += dpdA*dervDict2['pos'] |
---|
1531 | dDijDict = GetHStrainShiftDerv(refl,SGData,phfx) |
---|
1532 | for name in dDijDict: |
---|
1533 | if name in varylist: |
---|
1534 | dMdv[varylist.index(name)][iBeg:iFin] += dDijDict[name]*dervDict['pos'] |
---|
1535 | if Ka2: |
---|
1536 | dMdv[varylist.index(name)][iBeg2:iFin2] += dDijDict[name]*dervDict2['pos'] |
---|
1537 | elif name in dependentVars: |
---|
1538 | depDerivDict[name][iBeg:iFin] += dDijDict[name]*dervDict['pos'] |
---|
1539 | if Ka2: |
---|
1540 | depDerivDict[name][iBeg2:iFin2] += dDijDict[name]*dervDict2['pos'] |
---|
1541 | sigDict,gamDict = GetSampleSigGamDerv(refl,wave,G,GB,phfx,calcControls,parmDict) |
---|
1542 | for name in gamDict: |
---|
1543 | if name in varylist: |
---|
1544 | dMdv[varylist.index(name)][iBeg:iFin] += gamDict[name]*dervDict['gam'] |
---|
1545 | if Ka2: |
---|
1546 | dMdv[varylist.index(name)][iBeg2:iFin2] += gamDict[name]*dervDict2['gam'] |
---|
1547 | elif name in dependentVars: |
---|
1548 | depDerivDict[name][iBeg:iFin] += gamDict[name]*dervDict['gam'] |
---|
1549 | if Ka2: |
---|
1550 | depDerivDict[name][iBeg2:iFin2] += gamDict[name]*dervDict2['gam'] |
---|
1551 | for name in sigDict: |
---|
1552 | if name in varylist: |
---|
1553 | dMdv[varylist.index(name)][iBeg:iFin] += sigDict[name]*dervDict['sig'] |
---|
1554 | if Ka2: |
---|
1555 | dMdv[varylist.index(name)][iBeg2:iFin2] += sigDict[name]*dervDict2['sig'] |
---|
1556 | elif name in dependentVars: |
---|
1557 | depDerivDict[name][iBeg:iFin] += sigDict[name]*dervDict['sig'] |
---|
1558 | if Ka2: |
---|
1559 | depDerivDict[name][iBeg2:iFin2] += sigDict[name]*dervDict2['sig'] |
---|
1560 | for name in ['BabA','BabU']: |
---|
1561 | if phfx+name in varylist: |
---|
1562 | dMdv[varylist.index(phfx+name)][iBeg:iFin] += dFdvDict[pfx+name][iref]*dervDict['int']*cw[iBeg:iFin] |
---|
1563 | if Ka2: |
---|
1564 | dMdv[varylist.index(phfx+name)][iBeg2:iFin2] += dFdvDict[pfx+name][iref]*dervDict2['int']*cw[iBeg2:iFin2] |
---|
1565 | elif phfx+name in dependentVars: |
---|
1566 | depDerivDict[phfx+name][iBeg:iFin] += dFdvDict[pfx+name][iref]*dervDict['int']*cw[iBeg:iFin] |
---|
1567 | if Ka2: |
---|
1568 | depDerivDict[phfx+name][iBeg2:iFin2] += dFdvDict[pfx+name][iref]*dervDict2['int']*cw[iBeg2:iFin2] |
---|
1569 | elif 'T' in calcControls[hfx+'histType']: |
---|
1570 | print 'TOF Undefined at present' |
---|
1571 | raise Exception #no TOF yet |
---|
1572 | if not Phase['General'].get('doPawley'): |
---|
1573 | #do atom derivatives - for RB,F,X & U so far |
---|
1574 | corr = dervDict['int']/refl[9] |
---|
1575 | if Ka2: |
---|
1576 | corr2 = dervDict2['int']/refl[9] |
---|
1577 | for name in varylist+dependentVars: |
---|
1578 | if '::RBV;' in name: |
---|
1579 | pass |
---|
1580 | else: |
---|
1581 | try: |
---|
1582 | aname = name.split(pfx)[1][:2] |
---|
1583 | if aname not in ['Af','dA','AU','RB']: continue # skip anything not an atom or rigid body param |
---|
1584 | except IndexError: |
---|
1585 | continue |
---|
1586 | if name in varylist: |
---|
1587 | dMdv[varylist.index(name)][iBeg:iFin] += dFdvDict[name][iref]*corr |
---|
1588 | if Ka2: |
---|
1589 | dMdv[varylist.index(name)][iBeg2:iFin2] += dFdvDict[name][iref]*corr2 |
---|
1590 | elif name in dependentVars: |
---|
1591 | depDerivDict[name][iBeg:iFin] += dFdvDict[name][iref]*corr |
---|
1592 | if Ka2: |
---|
1593 | depDerivDict[name][iBeg2:iFin2] += dFdvDict[name][iref]*corr2 |
---|
1594 | # print 'profile derv time: %.3fs'%(time.time()-time0) |
---|
1595 | # now process derivatives in constraints |
---|
1596 | G2mv.Dict2Deriv(varylist,depDerivDict,dMdv) |
---|
1597 | return dMdv |
---|
1598 | |
---|
1599 | def dervRefine(values,HistoPhases,parmDict,varylist,calcControls,pawleyLookup,dlg): |
---|
1600 | '''Loop over histograms and compute derivatives of the fitting |
---|
1601 | model (M) with respect to all parameters. Results are returned in |
---|
1602 | a Jacobian matrix (aka design matrix) of dimensions (n by m) where |
---|
1603 | n is the number of parameters and m is the number of data |
---|
1604 | points. This can exceed memory when m gets large. This routine is |
---|
1605 | used when refinement derivatives are selected as "analtytic |
---|
1606 | Jacobian" in Controls. |
---|
1607 | |
---|
1608 | :returns: Jacobian numpy.array dMdv for all histograms concatinated |
---|
1609 | ''' |
---|
1610 | parmDict.update(zip(varylist,values)) |
---|
1611 | G2mv.Dict2Map(parmDict,varylist) |
---|
1612 | Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases |
---|
1613 | nvar = len(varylist) |
---|
1614 | dMdv = np.empty(0) |
---|
1615 | histoList = Histograms.keys() |
---|
1616 | histoList.sort() |
---|
1617 | for histogram in histoList: |
---|
1618 | if 'PWDR' in histogram[:4]: |
---|
1619 | Histogram = Histograms[histogram] |
---|
1620 | hId = Histogram['hId'] |
---|
1621 | hfx = ':%d:'%(hId) |
---|
1622 | wtFactor = calcControls[hfx+'wtFactor'] |
---|
1623 | Limits = calcControls[hfx+'Limits'] |
---|
1624 | x,y,w,yc,yb,yd = Histogram['Data'] |
---|
1625 | W = wtFactor*w |
---|
1626 | xB = np.searchsorted(x,Limits[0]) |
---|
1627 | xF = np.searchsorted(x,Limits[1]) |
---|
1628 | dMdvh = np.sqrt(W[xB:xF])*getPowderProfileDerv(parmDict,x[xB:xF], |
---|
1629 | varylist,Histogram,Phases,rigidbodyDict,calcControls,pawleyLookup) |
---|
1630 | elif 'HKLF' in histogram[:4]: |
---|
1631 | Histogram = Histograms[histogram] |
---|
1632 | nobs = Histogram['Residuals']['Nobs'] |
---|
1633 | phase = Histogram['Reflection Lists'] |
---|
1634 | Phase = Phases[phase] |
---|
1635 | hId = Histogram['hId'] |
---|
1636 | hfx = ':%d:'%(hId) |
---|
1637 | wtFactor = calcControls[hfx+'wtFactor'] |
---|
1638 | pfx = '%d::'%(Phase['pId']) |
---|
1639 | phfx = '%d:%d:'%(Phase['pId'],hId) |
---|
1640 | SGData = Phase['General']['SGData'] |
---|
1641 | A = [parmDict[pfx+'A%d'%(i)] for i in range(6)] |
---|
1642 | G,g = G2lat.A2Gmat(A) #recip & real metric tensors |
---|
1643 | refList = Histogram['Data'] |
---|
1644 | dFdvDict = StructureFactorDerv(refList,G,hfx,pfx,SGData,calcControls,parmDict) |
---|
1645 | ApplyRBModelDervs(dFdvDict,parmDict,rigidbodyDict,Phase) |
---|
1646 | dMdvh = np.zeros((len(varylist),len(refList))) |
---|
1647 | dependentVars = G2mv.GetDependentVars() |
---|
1648 | depDerivDict = {} |
---|
1649 | for j in dependentVars: |
---|
1650 | depDerivDict[j] = np.zeros(shape=(len(refList))) |
---|
1651 | if calcControls['F**2']: |
---|
1652 | for iref,ref in enumerate(refList): |
---|
1653 | if ref[6] > 0: |
---|
1654 | dervCor,dervDict = SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varylist) #puts correction in refl[13] |
---|
1655 | w = 1.0/ref[6] |
---|
1656 | if w*ref[5] >= calcControls['minF/sig']: |
---|
1657 | for j,var in enumerate(varylist): |
---|
1658 | if var in dFdvDict: |
---|
1659 | dMdvh[j][iref] = w*dFdvDict[var][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1660 | for var in dependentVars: |
---|
1661 | if var in dFdvDict: |
---|
1662 | depDerivDict[var][iref] = w*dFdvDict[var][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1663 | if phfx+'Scale' in varylist: |
---|
1664 | dMdvh[varylist.index(phfx+'Scale')][iref] = w*ref[9]*dervCor |
---|
1665 | elif phfx+'Scale' in dependentVars: |
---|
1666 | depDerivDict[phfx+'Scale'][iref] = w*ref[9]*dervCor |
---|
1667 | for item in ['Ep','Es','Eg']: |
---|
1668 | if phfx+item in varylist: |
---|
1669 | dMdvh[varylist.index(phfx+item)][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1670 | elif phfx+item in dependentVars: |
---|
1671 | depDerivDict[phfx+item][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1672 | for item in ['BabA','BabU']: |
---|
1673 | if phfx+item in varylist: |
---|
1674 | dMdvh[varylist.index(phfx+item)][iref] = w*dervCor*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale'] |
---|
1675 | elif phfx+item in dependentVars: |
---|
1676 | depDerivDict[phfx+item][iref] = w*dervCor*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale'] |
---|
1677 | else: |
---|
1678 | for iref,ref in enumerate(refList): |
---|
1679 | if ref[5] > 0.: |
---|
1680 | dervCor,dervDict = SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varylist) #puts correction in refl[13] |
---|
1681 | Fo = np.sqrt(ref[5]) |
---|
1682 | Fc = np.sqrt(ref[7]) |
---|
1683 | w = 1.0/ref[6] |
---|
1684 | if 2.0*Fo*w*Fo >= calcControls['minF/sig']: |
---|
1685 | for j,var in enumerate(varylist): |
---|
1686 | if var in dFdvDict: |
---|
1687 | dMdvh[j][iref] = w*dFdvDict[var][iref]*dervCor*parmDict[phfx+'Scale'] |
---|
1688 | for var in dependentVars: |
---|
1689 | if var in dFdvDict: |
---|
1690 | depDerivDict[var][iref] = w*dFdvDict[var][iref]*dervCor*parmDict[phfx+'Scale'] |
---|
1691 | if phfx+'Scale' in varylist: |
---|
1692 | dMdvh[varylist.index(phfx+'Scale')][iref] = w*ref[9]*dervCor |
---|
1693 | elif phfx+'Scale' in dependentVars: |
---|
1694 | depDerivDict[phfx+'Scale'][iref] = w*ref[9]*dervCor |
---|
1695 | for item in ['Ep','Es','Eg']: |
---|
1696 | if phfx+item in varylist: |
---|
1697 | dMdvh[varylist.index(phfx+item)][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1698 | elif phfx+item in dependentVars: |
---|
1699 | depDerivDict[phfx+item][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1700 | for item in ['BabA','BabU']: |
---|
1701 | if phfx+item in varylist: |
---|
1702 | dMdvh[varylist.index(phfx+item)][iref] = w*dervCor*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale'] |
---|
1703 | elif phfx+item in dependentVars: |
---|
1704 | depDerivDict[phfx+item][iref] = w*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1705 | # now process derivatives in constraints |
---|
1706 | G2mv.Dict2Deriv(varylist,depDerivDict,dMdvh) |
---|
1707 | else: |
---|
1708 | continue #skip non-histogram entries |
---|
1709 | if len(dMdv): |
---|
1710 | dMdv = np.concatenate((dMdv.T,np.sqrt(wtFactor)*dMdvh.T)).T |
---|
1711 | else: |
---|
1712 | dMdv = np.sqrt(wtFactor)*dMdvh |
---|
1713 | |
---|
1714 | pNames,pVals,pWt,pWsum = penaltyFxn(HistoPhases,parmDict,varylist) |
---|
1715 | if np.any(pVals): |
---|
1716 | dpdv = penaltyDeriv(pNames,pVals,HistoPhases,parmDict,varylist) |
---|
1717 | dMdv = np.concatenate((dMdv.T,(np.sqrt(pWt)*dpdv).T)).T |
---|
1718 | |
---|
1719 | return dMdv |
---|
1720 | |
---|
1721 | def HessRefine(values,HistoPhases,parmDict,varylist,calcControls,pawleyLookup,dlg): |
---|
1722 | '''Loop over histograms and compute derivatives of the fitting |
---|
1723 | model (M) with respect to all parameters. For each histogram, the |
---|
1724 | Jacobian matrix, dMdv, with dimensions (n by m) where n is the |
---|
1725 | number of parameters and m is the number of data points *in the |
---|
1726 | histogram*. The (n by n) Hessian is computed from each Jacobian |
---|
1727 | and it is returned. This routine is used when refinement |
---|
1728 | derivatives are selected as "analtytic Hessian" in Controls. |
---|
1729 | |
---|
1730 | :returns: Vec,Hess where Vec is the least-squares vector and Hess is the Hessian |
---|
1731 | ''' |
---|
1732 | parmDict.update(zip(varylist,values)) |
---|
1733 | G2mv.Dict2Map(parmDict,varylist) |
---|
1734 | Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases |
---|
1735 | ApplyRBModels(parmDict,Phases,rigidbodyDict) #,Update=True?? |
---|
1736 | nvar = len(varylist) |
---|
1737 | Hess = np.empty(0) |
---|
1738 | histoList = Histograms.keys() |
---|
1739 | histoList.sort() |
---|
1740 | for histogram in histoList: |
---|
1741 | if 'PWDR' in histogram[:4]: |
---|
1742 | Histogram = Histograms[histogram] |
---|
1743 | hId = Histogram['hId'] |
---|
1744 | hfx = ':%d:'%(hId) |
---|
1745 | wtFactor = calcControls[hfx+'wtFactor'] |
---|
1746 | Limits = calcControls[hfx+'Limits'] |
---|
1747 | x,y,w,yc,yb,yd = Histogram['Data'] |
---|
1748 | W = wtFactor*w |
---|
1749 | dy = y-yc |
---|
1750 | xB = np.searchsorted(x,Limits[0]) |
---|
1751 | xF = np.searchsorted(x,Limits[1]) |
---|
1752 | dMdvh = getPowderProfileDerv(parmDict,x[xB:xF], |
---|
1753 | varylist,Histogram,Phases,rigidbodyDict,calcControls,pawleyLookup) |
---|
1754 | Wt = ma.sqrt(W[xB:xF])[np.newaxis,:] |
---|
1755 | Dy = dy[xB:xF][np.newaxis,:] |
---|
1756 | dMdvh *= Wt |
---|
1757 | if dlg: |
---|
1758 | dlg.Update(Histogram['Residuals']['wR'],newmsg='Hessian for histogram %d\nAll data Rw=%8.3f%s'%(hId,Histogram['Residuals']['wR'],'%'))[0] |
---|
1759 | if len(Hess): |
---|
1760 | Hess += np.inner(dMdvh,dMdvh) |
---|
1761 | dMdvh *= Wt*Dy |
---|
1762 | Vec += np.sum(dMdvh,axis=1) |
---|
1763 | else: |
---|
1764 | Hess = np.inner(dMdvh,dMdvh) |
---|
1765 | dMdvh *= Wt*Dy |
---|
1766 | Vec = np.sum(dMdvh,axis=1) |
---|
1767 | elif 'HKLF' in histogram[:4]: |
---|
1768 | Histogram = Histograms[histogram] |
---|
1769 | nobs = Histogram['Residuals']['Nobs'] |
---|
1770 | phase = Histogram['Reflection Lists'] |
---|
1771 | Phase = Phases[phase] |
---|
1772 | hId = Histogram['hId'] |
---|
1773 | hfx = ':%d:'%(hId) |
---|
1774 | wtFactor = calcControls[hfx+'wtFactor'] |
---|
1775 | pfx = '%d::'%(Phase['pId']) |
---|
1776 | phfx = '%d:%d:'%(Phase['pId'],hId) |
---|
1777 | SGData = Phase['General']['SGData'] |
---|
1778 | A = [parmDict[pfx+'A%d'%(i)] for i in range(6)] |
---|
1779 | G,g = G2lat.A2Gmat(A) #recip & real metric tensors |
---|
1780 | refList = Histogram['Data'] |
---|
1781 | time0 = time.time() |
---|
1782 | dFdvDict = StructureFactorDerv(refList,G,hfx,pfx,SGData,calcControls,parmDict) #accurate for powders! |
---|
1783 | ApplyRBModelDervs(dFdvDict,parmDict,rigidbodyDict,Phase) |
---|
1784 | dMdvh = np.zeros((len(varylist),len(refList))) |
---|
1785 | dependentVars = G2mv.GetDependentVars() |
---|
1786 | depDerivDict = {} |
---|
1787 | for j in dependentVars: |
---|
1788 | depDerivDict[j] = np.zeros(shape=(len(refList))) |
---|
1789 | wdf = np.zeros(len(refList)) |
---|
1790 | if calcControls['F**2']: |
---|
1791 | for iref,ref in enumerate(refList): |
---|
1792 | if ref[6] > 0: |
---|
1793 | dervCor,dervDict = SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varylist) #puts correction in refl[13] |
---|
1794 | w = 1.0/ref[6] |
---|
1795 | if w*ref[5] >= calcControls['minF/sig']: |
---|
1796 | wdf[iref] = w*(ref[5]-ref[7]) |
---|
1797 | for j,var in enumerate(varylist): |
---|
1798 | if var in dFdvDict: |
---|
1799 | dMdvh[j][iref] = w*dFdvDict[var][iref]*dervCor*parmDict[phfx+'Scale'] |
---|
1800 | for var in dependentVars: |
---|
1801 | if var in dFdvDict: |
---|
1802 | depDerivDict[var][iref] = w*dFdvDict[var][iref]*dervCor*parmDict[phfx+'Scale'] |
---|
1803 | if phfx+'Scale' in varylist: |
---|
1804 | dMdvh[varylist.index(phfx+'Scale')][iref] = w*ref[9]*dervCor |
---|
1805 | elif phfx+'Scale' in dependentVars: |
---|
1806 | depDerivDict[phfx+'Scale'][iref] = w*ref[9]*dervCor |
---|
1807 | for item in ['Ep','Es','Eg']: |
---|
1808 | if phfx+item in varylist: |
---|
1809 | dMdvh[varylist.index(phfx+item)][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1810 | elif phfx+item in dependentVars: |
---|
1811 | depDerivDict[phfx+item][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1812 | for item in ['BabA','BabU']: |
---|
1813 | if phfx+item in varylist: |
---|
1814 | dMdvh[varylist.index(phfx+item)][iref] = w*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1815 | elif phfx+item in dependentVars: |
---|
1816 | depDerivDict[phfx+item][iref] = w*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1817 | else: |
---|
1818 | for iref,ref in enumerate(refList): |
---|
1819 | if ref[5] > 0.: |
---|
1820 | dervCor,dervDict = SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varylist) #puts correction in refl[13] |
---|
1821 | Fo = np.sqrt(ref[5]) |
---|
1822 | Fc = np.sqrt(ref[7]) |
---|
1823 | w = 1.0/ref[6] |
---|
1824 | if 2.0*Fo*w*Fo >= calcControls['minF/sig']: |
---|
1825 | wdf[iref] = 2.0*Fo*w*(Fo-Fc) |
---|
1826 | for j,var in enumerate(varylist): |
---|
1827 | if var in dFdvDict: |
---|
1828 | dMdvh[j][iref] = w*dFdvDict[var][iref]*dervCor*parmDict[phfx+'Scale'] |
---|
1829 | for var in dependentVars: |
---|
1830 | if var in dFdvDict: |
---|
1831 | depDerivDict[var][iref] = w*dFdvDict[var][iref]*dervCor*parmDict[phfx+'Scale'] |
---|
1832 | if phfx+'Scale' in varylist: |
---|
1833 | dMdvh[varylist.index(phfx+'Scale')][iref] = w*ref[9]*dervCor |
---|
1834 | elif phfx+'Scale' in dependentVars: |
---|
1835 | depDerivDict[phfx+'Scale'][iref] = w*ref[9]*dervCor |
---|
1836 | for item in ['Ep','Es','Eg']: |
---|
1837 | if phfx+item in varylist: |
---|
1838 | dMdvh[varylist.index(phfx+item)][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1839 | elif phfx+item in dependentVars: |
---|
1840 | depDerivDict[phfx+item][iref] = w*dervDict[phfx+item]*parmDict[phfx+'Scale'] |
---|
1841 | for item in ['BabA','BabU']: |
---|
1842 | if phfx+item in varylist: |
---|
1843 | dMdvh[varylist.index(phfx+item)][iref] = w*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1844 | elif phfx+item in dependentVars: |
---|
1845 | depDerivDict[phfx+item][iref] = w*dFdvDict[pfx+item][iref]*parmDict[phfx+'Scale']*dervCor |
---|
1846 | # now process derivatives in constraints |
---|
1847 | G2mv.Dict2Deriv(varylist,depDerivDict,dMdvh) |
---|
1848 | |
---|
1849 | if dlg: |
---|
1850 | dlg.Update(Histogram['Residuals']['wR'],newmsg='Hessian for histogram %d Rw=%8.3f%s'%(hId,Histogram['Residuals']['wR'],'%'))[0] |
---|
1851 | if len(Hess): |
---|
1852 | Vec += wtFactor*np.sum(dMdvh*wdf,axis=1) |
---|
1853 | Hess += wtFactor*np.inner(dMdvh,dMdvh) |
---|
1854 | else: |
---|
1855 | Vec = wtFactor*np.sum(dMdvh*wdf,axis=1) |
---|
1856 | Hess = wtFactor*np.inner(dMdvh,dMdvh) |
---|
1857 | else: |
---|
1858 | continue #skip non-histogram entries |
---|
1859 | pNames,pVals,pWt,pWsum = penaltyFxn(HistoPhases,parmDict,varylist) |
---|
1860 | if np.any(pVals): |
---|
1861 | dpdv = penaltyDeriv(pNames,pVals,HistoPhases,parmDict,varylist) |
---|
1862 | Vec += np.sum(dpdv*pWt*pVals,axis=1) |
---|
1863 | Hess += np.inner(dpdv*pWt,dpdv) |
---|
1864 | return Vec,Hess |
---|
1865 | |
---|
1866 | def errRefine(values,HistoPhases,parmDict,varylist,calcControls,pawleyLookup,dlg): |
---|
1867 | 'Needs a doc string' |
---|
1868 | parmDict.update(zip(varylist,values)) |
---|
1869 | Values2Dict(parmDict, varylist, values) |
---|
1870 | G2mv.Dict2Map(parmDict,varylist) |
---|
1871 | Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases |
---|
1872 | M = np.empty(0) |
---|
1873 | SumwYo = 0 |
---|
1874 | Nobs = 0 |
---|
1875 | ApplyRBModels(parmDict,Phases,rigidbodyDict) |
---|
1876 | histoList = Histograms.keys() |
---|
1877 | histoList.sort() |
---|
1878 | for histogram in histoList: |
---|
1879 | if 'PWDR' in histogram[:4]: |
---|
1880 | Histogram = Histograms[histogram] |
---|
1881 | hId = Histogram['hId'] |
---|
1882 | hfx = ':%d:'%(hId) |
---|
1883 | wtFactor = calcControls[hfx+'wtFactor'] |
---|
1884 | Limits = calcControls[hfx+'Limits'] |
---|
1885 | x,y,w,yc,yb,yd = Histogram['Data'] |
---|
1886 | yc *= 0.0 #zero full calcd profiles |
---|
1887 | yb *= 0.0 |
---|
1888 | yd *= 0.0 |
---|
1889 | xB = np.searchsorted(x,Limits[0]) |
---|
1890 | xF = np.searchsorted(x,Limits[1]) |
---|
1891 | yc[xB:xF],yb[xB:xF] = getPowderProfile(parmDict,x[xB:xF], |
---|
1892 | varylist,Histogram,Phases,calcControls,pawleyLookup) |
---|
1893 | yc[xB:xF] += yb[xB:xF] |
---|
1894 | if not np.any(y): #fill dummy data |
---|
1895 | rv = st.poisson(yc[xB:xF]) |
---|
1896 | y[xB:xF] = rv.rvs() |
---|
1897 | Z = np.ones_like(yc[xB:xF]) |
---|
1898 | Z[1::2] *= -1 |
---|
1899 | y[xB:xF] = yc[xB:xF]+np.abs(y[xB:xF]-yc[xB:xF])*Z |
---|
1900 | w[xB:xF] = np.where(y[xB:xF]>0.,1./y[xB:xF],1.0) |
---|
1901 | yd[xB:xF] = y[xB:xF]-yc[xB:xF] |
---|
1902 | W = wtFactor*w |
---|
1903 | wdy = -ma.sqrt(W[xB:xF])*(yd[xB:xF]) |
---|
1904 | Histogram['Residuals']['Nobs'] = ma.count(x[xB:xF]) |
---|
1905 | Nobs += Histogram['Residuals']['Nobs'] |
---|
1906 | Histogram['Residuals']['sumwYo'] = ma.sum(W[xB:xF]*y[xB:xF]**2) |
---|
1907 | SumwYo += Histogram['Residuals']['sumwYo'] |
---|
1908 | Histogram['Residuals']['R'] = min(100.,ma.sum(ma.abs(yd[xB:xF]))/ma.sum(y[xB:xF])*100.) |
---|
1909 | Histogram['Residuals']['wR'] = min(100.,ma.sqrt(ma.sum(wdy**2)/Histogram['Residuals']['sumwYo'])*100.) |
---|
1910 | sumYmB = ma.sum(ma.where(yc[xB:xF]!=yb[xB:xF],ma.abs(y[xB:xF]-yb[xB:xF]),0.)) |
---|
1911 | sumwYmB2 = ma.sum(ma.where(yc[xB:xF]!=yb[xB:xF],W[xB:xF]*(y[xB:xF]-yb[xB:xF])**2,0.)) |
---|
1912 | sumYB = ma.sum(ma.where(yc[xB:xF]!=yb[xB:xF],ma.abs(y[xB:xF]-yc[xB:xF])*ma.abs(y[xB:xF]-yb[xB:xF])/y[xB:xF],0.)) |
---|
1913 | sumwYB2 = ma.sum(ma.where(yc[xB:xF]!=yb[xB:xF],W[xB:xF]*(ma.abs(y[xB:xF]-yc[xB:xF])*ma.abs(y[xB:xF]-yb[xB:xF])/y[xB:xF])**2,0.)) |
---|
1914 | Histogram['Residuals']['Rb'] = min(100.,100.*sumYB/sumYmB) |
---|
1915 | Histogram['Residuals']['wRb'] = min(100.,100.*ma.sqrt(sumwYB2/sumwYmB2)) |
---|
1916 | Histogram['Residuals']['wRmin'] = min(100.,100.*ma.sqrt(Histogram['Residuals']['Nobs']/Histogram['Residuals']['sumwYo'])) |
---|
1917 | if dlg: |
---|
1918 | dlg.Update(Histogram['Residuals']['wR'],newmsg='For histogram %d Rw=%8.3f%s'%(hId,Histogram['Residuals']['wR'],'%'))[0] |
---|
1919 | M = np.concatenate((M,wdy)) |
---|
1920 | #end of PWDR processing |
---|
1921 | elif 'HKLF' in histogram[:4]: |
---|
1922 | Histogram = Histograms[histogram] |
---|
1923 | Histogram['Residuals'] = {} |
---|
1924 | phase = Histogram['Reflection Lists'] |
---|
1925 | Phase = Phases[phase] |
---|
1926 | hId = Histogram['hId'] |
---|
1927 | hfx = ':%d:'%(hId) |
---|
1928 | wtFactor = calcControls[hfx+'wtFactor'] |
---|
1929 | pfx = '%d::'%(Phase['pId']) |
---|
1930 | phfx = '%d:%d:'%(Phase['pId'],hId) |
---|
1931 | SGData = Phase['General']['SGData'] |
---|
1932 | A = [parmDict[pfx+'A%d'%(i)] for i in range(6)] |
---|
1933 | G,g = G2lat.A2Gmat(A) #recip & real metric tensors |
---|
1934 | refList = Histogram['Data'] |
---|
1935 | time0 = time.time() |
---|
1936 | StructureFactor(refList,G,hfx,pfx,SGData,calcControls,parmDict) |
---|
1937 | # print 'sf-calc time: %.3f'%(time.time()-time0) |
---|
1938 | df = np.zeros(len(refList)) |
---|
1939 | sumwYo = 0 |
---|
1940 | sumFo = 0 |
---|
1941 | sumFo2 = 0 |
---|
1942 | sumdF = 0 |
---|
1943 | sumdF2 = 0 |
---|
1944 | nobs = 0 |
---|
1945 | if calcControls['F**2']: |
---|
1946 | for i,ref in enumerate(refList): |
---|
1947 | if ref[6] > 0: |
---|
1948 | SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varylist) #puts correction in refl[13] |
---|
1949 | w = 1.0/ref[6] |
---|
1950 | ref[7] = parmDict[phfx+'Scale']*ref[9] |
---|
1951 | ref[7] *= ref[13] #correct Fc^2 for extinction |
---|
1952 | ref[8] = ref[5]/parmDict[phfx+'Scale'] |
---|
1953 | if w*ref[5] >= calcControls['minF/sig']: |
---|
1954 | sumFo2 += ref[5] |
---|
1955 | Fo = np.sqrt(ref[5]) |
---|
1956 | sumFo += Fo |
---|
1957 | sumFo2 += ref[5] |
---|
1958 | sumdF += abs(Fo-np.sqrt(ref[7])) |
---|
1959 | sumdF2 += abs(ref[5]-ref[7]) |
---|
1960 | nobs += 1 |
---|
1961 | df[i] = -w*(ref[5]-ref[7]) |
---|
1962 | sumwYo += (w*ref[5])**2 |
---|
1963 | else: |
---|
1964 | for i,ref in enumerate(refList): |
---|
1965 | if ref[5] > 0.: |
---|
1966 | SCExtinction(ref,phfx,hfx,pfx,calcControls,parmDict,varylist) #puts correction in refl[13] |
---|
1967 | ref[7] = parmDict[phfx+'Scale']*ref[9] |
---|
1968 | ref[7] *= ref[13] #correct Fc^2 for extinction |
---|
1969 | Fo = np.sqrt(ref[5]) |
---|
1970 | Fc = np.sqrt(ref[7]) |
---|
1971 | w = 2.0*Fo/ref[6] |
---|
1972 | if w*Fo >= calcControls['minF/sig']: |
---|
1973 | sumFo += Fo |
---|
1974 | sumFo2 += ref[5] |
---|
1975 | sumdF += abs(Fo-Fc) |
---|
1976 | sumdF2 += abs(ref[5]-ref[7]) |
---|
1977 | nobs += 1 |
---|
1978 | df[i] = -w*(Fo-Fc) |
---|
1979 | sumwYo += (w*Fo)**2 |
---|
1980 | Histogram['Residuals']['Nobs'] = nobs |
---|
1981 | Histogram['Residuals']['sumwYo'] = sumwYo |
---|
1982 | SumwYo += sumwYo |
---|
1983 | Histogram['Residuals']['wR'] = min(100.,np.sqrt(np.sum(df**2)/Histogram['Residuals']['sumwYo'])*100.) |
---|
1984 | Histogram['Residuals'][phfx+'Rf'] = 100.*sumdF/sumFo |
---|
1985 | Histogram['Residuals'][phfx+'Rf^2'] = 100.*sumdF2/sumFo2 |
---|
1986 | Histogram['Residuals'][phfx+'Nref'] = nobs |
---|
1987 | Nobs += nobs |
---|
1988 | if dlg: |
---|
1989 | dlg.Update(Histogram['Residuals']['wR'],newmsg='For histogram %d Rw=%8.3f%s'%(hId,Histogram['Residuals']['wR'],'%'))[0] |
---|
1990 | M = np.concatenate((M,wtFactor*df)) |
---|
1991 | # end of HKLF processing |
---|
1992 | Histograms['sumwYo'] = SumwYo |
---|
1993 | Histograms['Nobs'] = Nobs |
---|
1994 | Rw = min(100.,np.sqrt(np.sum(M**2)/SumwYo)*100.) |
---|
1995 | if dlg: |
---|
1996 | GoOn = dlg.Update(Rw,newmsg='%s%8.3f%s'%('All data Rw =',Rw,'%'))[0] |
---|
1997 | if not GoOn: |
---|
1998 | parmDict['saved values'] = values |
---|
1999 | dlg.Destroy() |
---|
2000 | raise Exception #Abort!! |
---|
2001 | pDict,pVals,pWt,pWsum = penaltyFxn(HistoPhases,parmDict,varylist) |
---|
2002 | if np.any(pVals): |
---|
2003 | pSum = np.sum(pWt*pVals**2) |
---|
2004 | for name in pWsum: |
---|
2005 | print ' Penalty function for %s = %.3f'%(name,pWsum[name]) |
---|
2006 | print 'Total penalty function: %.3f on %d terms'%(pSum,len(pVals)) |
---|
2007 | Nobs += len(pVals) |
---|
2008 | M = np.concatenate((M,np.sqrt(pWt)*pVals)) |
---|
2009 | return M |
---|
2010 | |
---|