This repository has been archived by the owner on Sep 16, 2024. It is now read-only.
-
Notifications
You must be signed in to change notification settings - Fork 5
/
Copy pathGSASIIstrMath.py
5052 lines (4895 loc) · 265 KB
/
GSASIIstrMath.py
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
# -*- coding: utf-8 -*-
########### SVN repository information ###################
# $Date: 2023-12-12 12:48:03 -0600 (Tue, 12 Dec 2023) $
# $Author: vondreele $
# $Revision: 5706 $
# $URL: https://subversion.xray.aps.anl.gov/pyGSAS/trunk/GSASIIstrMath.py $
# $Id: GSASIIstrMath.py 5706 2023-12-12 18:48:03Z vondreele $
########### SVN repository information ###################
'''
:mod:`GSASIIstrMath` routines, used for refinement computations
are found below.
'''
from __future__ import division, print_function
import time
import copy
import numpy as np
import numpy.ma as ma
import numpy.linalg as nl
import scipy.stats as st
import scipy.special as sp
import multiprocessing as mp
import pickle
import GSASIIpath
GSASIIpath.SetVersionNumber("$Revision: 5706 $")
import GSASIIElem as G2el
import GSASIIlattice as G2lat
import GSASIIspc as G2spc
import GSASIIpwd as G2pwd
import GSASIImapvars as G2mv
import GSASIImath as G2mth
import GSASIIobj as G2obj
import GSASIImpsubs as G2mp
#G2mp.InitMP(False) # This disables multiprocessing
sind = lambda x: np.sin(x*np.pi/180.)
cosd = lambda x: np.cos(x*np.pi/180.)
tand = lambda x: np.tan(x*np.pi/180.)
asind = lambda x: 180.*np.arcsin(x)/np.pi
acosd = lambda x: 180.*np.arccos(x)/np.pi
atan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
try: # fails on doc build
ateln2 = 8.0*np.log(2.0)
twopi = 2.0*np.pi
twopisq = 2.0*np.pi**2
except TypeError:
pass
nxs = np.newaxis
keV = 12.397639
################################################################################
##### Rigid Body Models
################################################################################
def ApplyRBModels(parmDict,Phases,rigidbodyDict,Update=False):
''' Takes RB info from RBModels in Phase and RB data in rigidbodyDict along with
current RB values in parmDict & modifies atom contents (fxyz & Uij) of parmDict
'''
atxIds = ['Ax:','Ay:','Az:']
atuIds = ['AU11:','AU22:','AU33:','AU12:','AU13:','AU23:']
RBIds = rigidbodyDict.get('RBIds',{'Vector':[],'Residue':[],'Spin':[]}) #these are lists of rbIds
RBIds['Spin'] = RBIds.get('Spin',[]) #patch
if not RBIds['Vector'] and not RBIds['Residue'] and not RBIds['Spin']:
return
VRBIds = RBIds['Vector']
RRBIds = RBIds['Residue']
SRBIds = RBIds['Spin']
if Update:
RBData = rigidbodyDict
else:
RBData = copy.deepcopy(rigidbodyDict) # don't mess with original!
if RBIds['Vector']: # first update the vector magnitudes
VRBData = RBData['Vector']
for i,rbId in enumerate(VRBIds):
if VRBData[rbId]['useCount']:
for j in range(len(VRBData[rbId]['VectMag'])):
name = '::RBV;'+str(j)+':'+str(i)
VRBData[rbId]['VectMag'][j] = parmDict[name]
for phase in Phases:
Phase = Phases[phase]
General = Phase['General']
cx,ct,cs,cia = General['AtomPtrs']
cell = General['Cell'][1:7]
Amat,Bmat = G2lat.cell2AB(cell)
AtLookup = G2mth.FillAtomLookUp(Phase['Atoms'],cia+8)
pfx = str(Phase['pId'])+'::'
if Update:
RBModels = Phase['RBModels']
else:
RBModels = copy.deepcopy(Phase['RBModels']) # again don't mess with original!
for irb,RBObj in enumerate(RBModels.get('Vector',[])):
jrb = VRBIds.index(RBObj['RBId'])
rbsx = str(irb)+':'+str(jrb)
for i,px in enumerate(['RBVPx:','RBVPy:','RBVPz:']):
RBObj['Orig'][0][i] = parmDict[pfx+px+rbsx]
for i,po in enumerate(['RBVOa:','RBVOi:','RBVOj:','RBVOk:']):
RBObj['Orient'][0][i] = parmDict[pfx+po+rbsx]
RBObj['Orient'][0] = G2mth.normQ(RBObj['Orient'][0])
RBObj['AtomFrac'][0] = parmDict[pfx+'RBVf:'+rbsx]
TLS = RBObj['ThermalMotion']
if 'T' in TLS[0]:
for i,pt in enumerate(['RBVT11:','RBVT22:','RBVT33:','RBVT12:','RBVT13:','RBVT23:']):
TLS[1][i] = parmDict[pfx+pt+rbsx]
if 'L' in TLS[0]:
for i,pt in enumerate(['RBVL11:','RBVL22:','RBVL33:','RBVL12:','RBVL13:','RBVL23:']):
TLS[1][i+6] = parmDict[pfx+pt+rbsx]
if 'S' in TLS[0]:
for i,pt in enumerate(['RBVS12:','RBVS13:','RBVS21:','RBVS23:','RBVS31:','RBVS32:','RBVSAA:','RBVSBB:']):
TLS[1][i+12] = parmDict[pfx+pt+rbsx]
if 'U' in TLS[0]:
TLS[1][0] = parmDict[pfx+'RBVU:'+rbsx]
XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Vector')
UIJ = G2mth.UpdateRBUIJ(Bmat,Cart,RBObj)
for i,x in enumerate(XYZ):
atId = RBObj['Ids'][i]
if parmDict[pfx+'Afrac:'+str(AtLookup[atId])]:
parmDict[pfx+'Afrac:'+str(AtLookup[atId])] = RBObj['AtomFrac'][0]
for j in [0,1,2]:
parmDict[pfx+atxIds[j]+str(AtLookup[atId])] = x[j]
if UIJ[i][0] == 'A':
for j in range(6):
parmDict[pfx+atuIds[j]+str(AtLookup[atId])] = UIJ[i][j+2]
elif UIJ[i][0] == 'I':
parmDict[pfx+'AUiso:'+str(AtLookup[atId])] = UIJ[i][1]
for irb,RBObj in enumerate(RBModels.get('Residue',[])):
jrb = RRBIds.index(RBObj['RBId'])
rbsx = str(irb)+':'+str(jrb)
for i,px in enumerate(['RBRPx:','RBRPy:','RBRPz:']):
RBObj['Orig'][0][i] = parmDict[pfx+px+rbsx]
for i,po in enumerate(['RBROa:','RBROi:','RBROj:','RBROk:']):
RBObj['Orient'][0][i] = parmDict[pfx+po+rbsx]
RBObj['Orient'][0] = G2mth.normQ(RBObj['Orient'][0])
RBObj['AtomFrac'][0] = parmDict[pfx+'RBRf:'+rbsx]
TLS = RBObj['ThermalMotion']
if 'T' in TLS[0]:
for i,pt in enumerate(['RBRT11:','RBRT22:','RBRT33:','RBRT12:','RBRT13:','RBRT23:']):
RBObj['ThermalMotion'][1][i] = parmDict[pfx+pt+rbsx]
if 'L' in TLS[0]:
for i,pt in enumerate(['RBRL11:','RBRL22:','RBRL33:','RBRL12:','RBRL13:','RBRL23:']):
RBObj['ThermalMotion'][1][i+6] = parmDict[pfx+pt+rbsx]
if 'S' in TLS[0]:
for i,pt in enumerate(['RBRS12:','RBRS13:','RBRS21:','RBRS23:','RBRS31:','RBRS32:','RBRSAA:','RBRSBB:']):
RBObj['ThermalMotion'][1][i+12] = parmDict[pfx+pt+rbsx]
if 'U' in TLS[0]:
RBObj['ThermalMotion'][1][0] = parmDict[pfx+'RBRU:'+rbsx]
for itors,tors in enumerate(RBObj['Torsions']):
tors[0] = parmDict[pfx+'RBRTr;'+str(itors)+':'+rbsx]
XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Residue')
UIJ = G2mth.UpdateRBUIJ(Bmat,Cart,RBObj)
for i,x in enumerate(XYZ):
atId = RBObj['Ids'][i]
if parmDict[pfx+'Afrac:'+str(AtLookup[atId])]:
parmDict[pfx+'Afrac:'+str(AtLookup[atId])] = RBObj['AtomFrac'][0]
for j in [0,1,2]:
parmDict[pfx+atxIds[j]+str(AtLookup[atId])] = x[j]
if UIJ[i][0] == 'A':
for j in range(6):
parmDict[pfx+atuIds[j]+str(AtLookup[atId])] = UIJ[i][j+2]
elif UIJ[i][0] == 'I':
parmDict[pfx+'AUiso:'+str(AtLookup[atId])] = UIJ[i][1]
for irb,RBObj in enumerate(RBModels.get('Spin',[])):
iAt = AtLookup[RBObj['Ids'][0]]
jrb = SRBIds.index(RBObj['RBId'][0])
name = pfx+'RBSOa:%d:%d'%(iAt,jrb)
for i,po in enumerate(['RBSOa:','RBSOi:','RBSOj:','RBSOk:']):
name = pfx+'%s%d:%d'%(po,iAt,jrb)
RBObj['Orient'][0][i] = parmDict[name]
for ish in range(len(RBObj['RBId'])):
jrb = SRBIds.index(RBObj['RBId'][ish])
if 'Q' not in RBObj['atType']:
name = pfx+'RBSSh;%d;Radius:%d:%d'%(ish,iAt,jrb)
RBObj['Radius'][ish][0] = parmDict[name]
for item in RBObj['SHC'][ish]:
name = pfx+'RBSSh;%d;%s:%d:%d'%(ish,item,iAt,jrb)
RBObj['SHC'][ish][item][0] = parmDict[name]
def ApplyRBModelDervs(dFdvDict,parmDict,rigidbodyDict,Phase):
'Computes rigid body derivatives; there are none for Spin RBs'
atxIds = ['dAx:','dAy:','dAz:']
atuIds = ['AU11:','AU22:','AU33:','AU12:','AU13:','AU23:']
OIds = ['Oa:','Oi:','Oj:','Ok:']
RBIds = rigidbodyDict.get('RBIds',{'Vector':[],'Residue':[]}) #these are lists of rbIds
if not RBIds['Vector'] and not RBIds['Residue']:
return
VRBIds = RBIds['Vector']
RRBIds = RBIds['Residue']
RBData = rigidbodyDict
for item in parmDict:
if 'RBV' in item or 'RBR' in item:
dFdvDict[item] = 0. #NB: this is a vector which is no. refl. long & must be filled!
General = Phase['General']
cx,ct,cs,cia = General['AtomPtrs']
cell = General['Cell'][1:7]
Amat,Bmat = G2lat.cell2AB(cell)
rpd = np.pi/180.
rpd2 = rpd**2
g = nl.inv(np.inner(Bmat,Bmat))
gvec = np.sqrt(np.array([g[0][0]**2,g[1][1]**2,g[2][2]**2,
g[0][0]*g[1][1],g[0][0]*g[2][2],g[1][1]*g[2][2]]))
AtLookup = G2mth.FillAtomLookUp(Phase['Atoms'],cia+8)
pfx = str(Phase['pId'])+'::'
RBModels = Phase['RBModels']
for irb,RBObj in enumerate(RBModels.get('Vector',[])):
symAxis = RBObj.get('symAxis')
VModel = RBData['Vector'][RBObj['RBId']]
Q = RBObj['Orient'][0]
jrb = VRBIds.index(RBObj['RBId'])
rbsx = str(irb)+':'+str(jrb)
dXdv = []
for iv in range(len(VModel['VectMag'])):
dCdv = []
for vec in VModel['rbVect'][iv]:
dCdv.append(G2mth.prodQVQ(Q,vec))
dXdv.append(np.inner(Bmat,np.array(dCdv)).T)
XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Vector')
for ia,atId in enumerate(RBObj['Ids']):
atNum = AtLookup[atId]
if parmDict[pfx+'Afrac:'+str(AtLookup[atId])]:
dFdvDict[pfx+'RBVf:'+rbsx] += dFdvDict[pfx+'Afrac:'+str(atNum)]
dx = 0.00001
for iv in range(len(VModel['VectMag'])):
for ix in [0,1,2]:
dFdvDict['::RBV;'+str(iv)+':'+str(jrb)] += dXdv[iv][ia][ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)]
for i,name in enumerate(['RBVPx:','RBVPy:','RBVPz:']):
dFdvDict[pfx+name+rbsx] += dFdvDict[pfx+atxIds[i]+str(atNum)]
for iv in range(4):
Q[iv] -= dx
XYZ1 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q),symAxis)
Q[iv] += 2.*dx
XYZ2 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q),symAxis)
Q[iv] -= dx
dXdO = (XYZ2[ia]-XYZ1[ia])/(2.*dx)
for ix in [0,1,2]:
dFdvDict[pfx+'RBV'+OIds[iv]+rbsx] += dXdO[ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)]
X = G2mth.prodQVQ(Q,Cart[ia])
dFdu = np.array([dFdvDict[pfx+Uid+str(AtLookup[atId])] for Uid in atuIds]).T/gvec
dFdu = G2lat.U6toUij(dFdu.T)
dFdu = np.tensordot(Amat,np.tensordot(Amat,dFdu,([1,0])),([0,1]))
dFdu = G2lat.UijtoU6(dFdu)
atNum = AtLookup[atId]
if 'T' in RBObj['ThermalMotion'][0]:
for i,name in enumerate(['RBVT11:','RBVT22:','RBVT33:','RBVT12:','RBVT13:','RBVT23:']):
dFdvDict[pfx+name+rbsx] += dFdu[i]
if 'L' in RBObj['ThermalMotion'][0]:
dFdvDict[pfx+'RBVL11:'+rbsx] += rpd2*(dFdu[1]*X[2]**2+dFdu[2]*X[1]**2-dFdu[5]*X[1]*X[2])
dFdvDict[pfx+'RBVL22:'+rbsx] += rpd2*(dFdu[0]*X[2]**2+dFdu[2]*X[0]**2-dFdu[4]*X[0]*X[2])
dFdvDict[pfx+'RBVL33:'+rbsx] += rpd2*(dFdu[0]*X[1]**2+dFdu[1]*X[0]**2-dFdu[3]*X[0]*X[1])
dFdvDict[pfx+'RBVL12:'+rbsx] += rpd2*(-dFdu[3]*X[2]**2-2.*dFdu[2]*X[0]*X[1]+
dFdu[4]*X[1]*X[2]+dFdu[5]*X[0]*X[2])
dFdvDict[pfx+'RBVL13:'+rbsx] += rpd2*(-dFdu[4]*X[1]**2-2.*dFdu[1]*X[0]*X[2]+
dFdu[3]*X[1]*X[2]+dFdu[5]*X[0]*X[1])
dFdvDict[pfx+'RBVL23:'+rbsx] += rpd2*(-dFdu[5]*X[0]**2-2.*dFdu[0]*X[1]*X[2]+
dFdu[3]*X[0]*X[2]+dFdu[4]*X[0]*X[1])
if 'S' in RBObj['ThermalMotion'][0]:
dFdvDict[pfx+'RBVS12:'+rbsx] += rpd*(dFdu[5]*X[1]-2.*dFdu[1]*X[2])
dFdvDict[pfx+'RBVS13:'+rbsx] += rpd*(-dFdu[5]*X[2]+2.*dFdu[2]*X[1])
dFdvDict[pfx+'RBVS21:'+rbsx] += rpd*(-dFdu[4]*X[0]+2.*dFdu[0]*X[2])
dFdvDict[pfx+'RBVS23:'+rbsx] += rpd*(dFdu[4]*X[2]-2.*dFdu[2]*X[0])
dFdvDict[pfx+'RBVS31:'+rbsx] += rpd*(dFdu[3]*X[0]-2.*dFdu[0]*X[1])
dFdvDict[pfx+'RBVS32:'+rbsx] += rpd*(-dFdu[3]*X[1]+2.*dFdu[1]*X[0])
dFdvDict[pfx+'RBVSAA:'+rbsx] += rpd*(dFdu[4]*X[1]-dFdu[3]*X[2])
dFdvDict[pfx+'RBVSBB:'+rbsx] += rpd*(dFdu[5]*X[0]-dFdu[3]*X[2])
if 'U' in RBObj['ThermalMotion'][0]:
dFdvDict[pfx+'RBVU:'+rbsx] += dFdvDict[pfx+'AUiso:'+str(AtLookup[atId])]
for irb,RBObj in enumerate(RBModels.get('Residue',[])):
symAxis = RBObj.get('symAxis')
Q = RBObj['Orient'][0]
jrb = RRBIds.index(RBObj['RBId'])
torData = RBData['Residue'][RBObj['RBId']]['rbSeq']
rbsx = str(irb)+':'+str(jrb)
XYZ,Cart = G2mth.UpdateRBXYZ(Bmat,RBObj,RBData,'Residue')
for itors,tors in enumerate(RBObj['Torsions']): #derivative error?
tname = pfx+'RBRTr;'+str(itors)+':'+rbsx
orId,pvId = torData[itors][:2]
pivotVec = Cart[orId]-Cart[pvId]
QA = G2mth.AVdeg2Q(-0.001,pivotVec)
QB = G2mth.AVdeg2Q(0.001,pivotVec)
for ir in torData[itors][3]:
atNum = AtLookup[RBObj['Ids'][ir]]
rVec = Cart[ir]-Cart[pvId]
dR = G2mth.prodQVQ(QB,rVec)-G2mth.prodQVQ(QA,rVec)
dRdT = np.inner(Bmat,G2mth.prodQVQ(Q,dR))/.002
for ix in [0,1,2]:
dFdvDict[tname] += dRdT[ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)]
for ia,atId in enumerate(RBObj['Ids']):
atNum = AtLookup[atId]
if parmDict[pfx+'Afrac:'+str(AtLookup[atId])]:
dFdvDict[pfx+'RBRf:'+rbsx] += dFdvDict[pfx+'Afrac:'+str(atNum)]
dx = 0.00001
for i,name in enumerate(['RBRPx:','RBRPy:','RBRPz:']):
dFdvDict[pfx+name+rbsx] += dFdvDict[pfx+atxIds[i]+str(atNum)]
for iv in range(4):
Q[iv] -= dx
XYZ1 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q),symAxis)
Q[iv] += 2.*dx
XYZ2 = G2mth.RotateRBXYZ(Bmat,Cart,G2mth.normQ(Q),symAxis)
Q[iv] -= dx
dXdO = (XYZ2[ia]-XYZ1[ia])/(2.*dx)
for ix in [0,1,2]:
dFdvDict[pfx+'RBR'+OIds[iv]+rbsx] += dXdO[ix]*dFdvDict[pfx+atxIds[ix]+str(atNum)]
X = G2mth.prodQVQ(Q,Cart[ia])
dFdu = np.array([dFdvDict[pfx+Uid+str(AtLookup[atId])] for Uid in atuIds]).T/gvec
dFdu = G2lat.U6toUij(dFdu.T)
dFdu = np.tensordot(Amat.T,np.tensordot(Amat,dFdu,([1,0])),([0,1]))
dFdu = G2lat.UijtoU6(dFdu)
atNum = AtLookup[atId]
if 'T' in RBObj['ThermalMotion'][0]:
for i,name in enumerate(['RBRT11:','RBRT22:','RBRT33:','RBRT12:','RBRT13:','RBRT23:']):
dFdvDict[pfx+name+rbsx] += dFdu[i]
if 'L' in RBObj['ThermalMotion'][0]:
dFdvDict[pfx+'RBRL11:'+rbsx] += rpd2*(dFdu[1]*X[2]**2+dFdu[2]*X[1]**2-dFdu[5]*X[1]*X[2])
dFdvDict[pfx+'RBRL22:'+rbsx] += rpd2*(dFdu[0]*X[2]**2+dFdu[2]*X[0]**2-dFdu[4]*X[0]*X[2])
dFdvDict[pfx+'RBRL33:'+rbsx] += rpd2*(dFdu[0]*X[1]**2+dFdu[1]*X[0]**2-dFdu[3]*X[0]*X[1])
dFdvDict[pfx+'RBRL12:'+rbsx] += rpd2*(-dFdu[3]*X[2]**2-2.*dFdu[2]*X[0]*X[1]+
dFdu[4]*X[1]*X[2]+dFdu[5]*X[0]*X[2])
dFdvDict[pfx+'RBRL13:'+rbsx] += rpd2*(dFdu[4]*X[1]**2-2.*dFdu[1]*X[0]*X[2]+
dFdu[3]*X[1]*X[2]+dFdu[5]*X[0]*X[1])
dFdvDict[pfx+'RBRL23:'+rbsx] += rpd2*(dFdu[5]*X[0]**2-2.*dFdu[0]*X[1]*X[2]+
dFdu[3]*X[0]*X[2]+dFdu[4]*X[0]*X[1])
if 'S' in RBObj['ThermalMotion'][0]:
dFdvDict[pfx+'RBRS12:'+rbsx] += rpd*(dFdu[5]*X[1]-2.*dFdu[1]*X[2])
dFdvDict[pfx+'RBRS13:'+rbsx] += rpd*(-dFdu[5]*X[2]+2.*dFdu[2]*X[1])
dFdvDict[pfx+'RBRS21:'+rbsx] += rpd*(-dFdu[4]*X[0]+2.*dFdu[0]*X[2])
dFdvDict[pfx+'RBRS23:'+rbsx] += rpd*(dFdu[4]*X[2]-2.*dFdu[2]*X[0])
dFdvDict[pfx+'RBRS31:'+rbsx] += rpd*(dFdu[3]*X[0]-2.*dFdu[0]*X[1])
dFdvDict[pfx+'RBRS32:'+rbsx] += rpd*(-dFdu[3]*X[1]+2.*dFdu[1]*X[0])
dFdvDict[pfx+'RBRSAA:'+rbsx] += rpd*(dFdu[4]*X[1]-dFdu[3]*X[2])
dFdvDict[pfx+'RBRSBB:'+rbsx] += rpd*(dFdu[5]*X[0]-dFdu[3]*X[2])
if 'U' in RBObj['ThermalMotion'][0]:
dFdvDict[pfx+'RBRU:'+rbsx] += dFdvDict[pfx+'AUiso:'+str(AtLookup[atId])]
def MakeSpHarmFF(HKL,Bmat,SHCdict,Tdata,hType,FFtables,ORBtables,BLtables,FF,SQ,ifDeriv=False):
''' Computes hkl dependent form factors & derivatives from spinning rigid bodies
:param array HKL: reflection hkl set to be considered
:param array Bmat: inv crystal to Cartesian transfomation matrix
:param dict SHCdict: RB spin/deformation data
:param array Tdata: atom type info
:param str hType: histogram type
:param dict FFtables: x-ray form factor tables
:param dict ORBtables: x-ray orbital form factor tables
:param dict BLtables: neutron scattering lenghts
:param array FF: form factors - will be modified by adding the spin/deformation RB spherical harmonics terms
:param array SQ: 1/4d^2 for the HKL set
:param bool ifDeriv: True if dFF/dcoff to be returned
:returns: dict dFFdS of derivatives if ifDeriv = True
'''
def MakePolar(Orient,QB):
QA = G2mth.invQ(Orient) #rotates about chosen axis
Q = G2mth.prodQQ(QB,QA) #might be switched? QB,QA is order for plotting
return G2lat.H2ThPh(np.reshape(HKL,(-1,3)),Bmat,Q)
dFFdS = {}
atFlg = []
Th,Ph = G2lat.H2ThPh(np.reshape(HKL,(-1,3)),Bmat,[1.,0.,0.,1.])
SQR = np.repeat(SQ,HKL.shape[1])
for iAt,Atype in enumerate(Tdata):
if 'Q' in Atype:
Th,Ph = G2lat.H2ThPh(np.reshape(HKL,(-1,3)),Bmat,[1.,0.,0.,1.])
atFlg.append(1.0)
SHdat = SHCdict[iAt]
symAxis = np.array(SHdat['symAxis'])
QB = G2mth.make2Quat(np.array([0,0,1.]),symAxis)[0] #position obj polar axis
Th,Ph = MakePolar([SHdat['Oa'],SHdat['Oi'],SHdat['Oj'],SHdat['Ok']],QB)
ThP,PhP = MakePolar([SHdat['Oa']+.0001,SHdat['Oi'],SHdat['Oj'],SHdat['Ok']],QB)
dp = 0.00001
ThPi,PhPi = MakePolar([SHdat['Oa'],SHdat['Oi']+dp,SHdat['Oj'],SHdat['Ok']],QB)
ThPj,PhPj = MakePolar([SHdat['Oa'],SHdat['Oi'],SHdat['Oj']+dp,SHdat['Ok']],QB)
ThPk,PhPk = MakePolar([SHdat['Oa'],SHdat['Oi'],SHdat['Oj'],SHdat['Ok']+dp],QB)
ThM,PhM = MakePolar([SHdat['Oa']-.0001,SHdat['Oi'],SHdat['Oj'],SHdat['Ok']],QB)
ThMi,PhMi = MakePolar([SHdat['Oa'],SHdat['Oi']-dp,SHdat['Oj'],SHdat['Ok']],QB)
ThMj,PhMj = MakePolar([SHdat['Oa'],SHdat['Oi'],SHdat['Oj']-dp,SHdat['Ok']],QB)
ThMk,PhMk = MakePolar([SHdat['Oa'],SHdat['Oi'],SHdat['Oj'],SHdat['Ok']-dp],QB)
QR = np.repeat(twopi*np.sqrt(4.*SQ),HKL.shape[1]) #refl Q for Bessel fxn
FF[:,iAt] = 0.
ishl = 0
dSHdO = np.zeros(HKL.shape[0]*HKL.shape[1])
dSHdOi = np.zeros(HKL.shape[0]*HKL.shape[1])
dSHdOj = np.zeros(HKL.shape[0]*HKL.shape[1])
dSHdOk = np.zeros(HKL.shape[0]*HKL.shape[1])
if '0' not in SHdat: #no spin RB for atom Q??
break
Shell = SHdat['0']
Irb = Shell['ShR']
Oname = 'Oa:%d:%s'%(iAt,Irb)
Oiname = 'Oi:%d:%s'%(iAt,Irb)
Ojname = 'Oj:%d:%s'%(iAt,Irb)
Okname = 'Ok:%d:%s'%(iAt,Irb)
while True:
shl = '%d'%ishl
if shl not in SHdat:
break
Shell = SHdat[shl]
Atm = Shell['AtType']
Nat = Shell['Natoms']
Irb = Shell['ShR']
if 'X' in hType:
if 'Q' in Atm:
SFF = 0.0
else:
SFF = G2el.ScatFac(FFtables[Atm],SQR)
elif 'N' in hType:
SFF = G2el.getBLvalues(BLtables)[Atm]
Rname = 'Sh;%s;Radius:%d:%s'%(shl,iAt,Irb)
if 'Q' in Atm:
dBSdR= 0.0
FF[:,iAt] = 0.0
else:
R = Shell['Radius']
R0 = sp.spherical_jn(0,QR*R)/(4.*np.pi)
R0P = sp.spherical_jn(0,QR*(R+0.01))/(4.*np.pi)
R0M = sp.spherical_jn(0,QR*(R-0.01))/(4.*np.pi)
dBSdR = Nat*SFF*(R0P-R0M)/0.02
FF[:,iAt] += Nat*SFF*R0 #Bessel function; L=0 term
for item in Shell:
if 'C(' in item:
l,m = eval(item.strip('C').strip('c'))
SH = G2lat.KslCalc(item,Th,Ph)
SHP = G2lat.KslCalc(item,ThP,PhP)
SHPi = G2lat.KslCalc(item,ThPi,PhPi)
SHPj = G2lat.KslCalc(item,ThPj,PhPj)
SHPk = G2lat.KslCalc(item,ThPk,PhPk)
SHM = G2lat.KslCalc(item,ThM,PhM)
SHMi = G2lat.KslCalc(item,ThMi,PhMi)
SHMj = G2lat.KslCalc(item,ThMj,PhMj)
SHMk = G2lat.KslCalc(item,ThMk,PhMk)
BS = 1.0
if 'Q' in Atm:
BS = sp.spherical_jn(l,1.0) #Slater term here?
else:
BS = sp.spherical_jn(l,QR*R)/(4.*np.pi) #Bessel function
BSP = sp.spherical_jn(l,QR*(R+0.01))/(4.*np.pi)
BSM = sp.spherical_jn(l,QR*(R-0.01))/(4.*np.pi)
dBSdR += Nat*SFF*SH*Shell[item]*(BSP-BSM)/0.02
dSHdO += Nat*SFF*BS*Shell[item]*(SHP-SHM)/0.0002
dSHdOi += Nat*SFF*BS*Shell[item]*(SHPi-SHMi)/(2.*dp)
dSHdOj += Nat*SFF*BS*Shell[item]*(SHPj-SHMj)/(2.*dp)
dSHdOk += Nat*SFF*BS*Shell[item]*(SHPk-SHMk)/(2.*dp)
FF[:,iAt] += Nat*SFF*BS*SH*Shell[item]
name = 'Sh;%s;%s:%d:%s'%(shl,item,iAt,Irb)
dFFdS[name] = Nat*SFF*BS*SH
if 'Q' not in Atm:
dFFdS[Rname] = dBSdR
ishl += 1
dFFdS[Oname] = dSHdO
dFFdS[Oiname] = dSHdOi
dFFdS[Ojname] = dSHdOj
dFFdS[Okname] = dSHdOk
elif iAt in SHCdict and 'X' in hType:
orKeys = [item for item in ORBtables[Atype] if item not in ['ZSlater','NSlater','SZE','popCore','popVal']]
orbs = SHCdict[iAt]
UVmat = np.inner(nl.inv(SHCdict[-iAt]['UVmat']),Bmat)
Th,Ph = G2lat.H2ThPh(np.reshape(HKL,(-1,3)),UVmat,[1.,0.,0.,1.])
atFlg.append(1.0)
orbTable = ORBtables[Atype][orKeys[0]]
ffOrb = {item:orbTable[item] for item in orbTable if item not in ['ZSlater','NSlater','SZE','popCore','popVal']}
FFcore = G2el.ScatFac(ffOrb,SQR) #core
FFtot = np.zeros_like(FFcore)
for orb in orbs:
if 'UVmat' in orb:
continue
Ne = orbs[orb].get('Ne',1.0) # not there for non <j0> orbs
if 'kappa' in orbs[orb]:
kappa = orbs[orb]['kappa']
SQk = SQR/kappa**2
korb = orb
orbTable = ORBtables[Atype][orKeys[int(orb)+1]]
ffOrb = {item:orbTable[item] for item in orbTable if item not in ['ZSlater','NSlater','SZE','popCore','popVal']}
ff = Ne*G2el.ScatFac(ffOrb,SQk)
dffdk = G2el.ScatFacDer(ffOrb,SQk)
dSH = 0.0
if '<j0>' in orKeys[int(orb)+1]:
dSH = 1.0
for term in orbs[orb]:
if 'D(' in term:
item = term.replace('D','C')
SH = G2lat.KslCalc(item,Th,Ph)
FFtot += SH*orbs[orb][term]*ff
name = 'A%s%s:%d'%(term,orb,iAt)
dFFdS[name] = SH*ff
dSH += SH*orbs[orb][term]
elif 'Ne' in term:
name = 'ANe%s:%d'%(orb,iAt)
dFFdS[name] = ff/Ne
if 'j0' in orKeys[int(orb)+1]:
FFtot += ff
name = 'Akappa%s:%d'%(korb,iAt)
if name in dFFdS:
dFFdS[name] += -2.0*Ne*SQk*dSH*dffdk/kappa
else:
dFFdS[name] = -2.0*Ne*SQk*dSH*dffdk/kappa
FF[:,iAt] = FFcore+FFtot
else:
atFlg.append(0.)
if ifDeriv:
return dFFdS,atFlg
def GetSHC(pfx,parmDict):
SHCdict = {}
for parm in parmDict:
if pfx+'RBS' in parm and 'RBS;' not in parm: #skips radii parms
items = parm.split(':')
atid = int(items[-2])
name = items[2][3:] #strip 'RBS'
if atid not in SHCdict:
SHCdict[atid] = {}
if ';' not in name: # will get Oa, Oi ,Oj, Ok
if name not in ['AtNo','Px','Py','Pz','SytSym']:
SHCdict[atid][name] = parmDict[parm]
continue
bits = name.split(';')
shno = bits[1]
if shno not in SHCdict[atid]:
SHCdict[atid][shno] = {}
if 'AtType' in bits[0] or 'Natoms' in bits[0] or 'ShR' in bits[0]:
SHCdict[atid][shno][bits[0]] = parmDict[parm]
elif 'Sh' in name:
cof = bits[2]
SHCdict[atid][shno][cof] = parmDict[parm]
if pfx+'AD(' in parm or pfx+'Akappa' in parm or pfx+'ANe' in parm: #atom deformation parms
items = parm.split(':')
atid = int(items[-1])
name = items[2][1:] #strip 'A'
if atid not in SHCdict:
SHCdict[atid] = {}
orb = name[-1]
if orb not in SHCdict[atid]:
SHCdict[atid][orb] = {}
SHCdict[atid][orb][name[:-1]] = parmDict[parm] #[atom id][orb no.][deform. coef]
if pfx+'UVmat' in parm:
items = parm.split(':')
atid = int(items[-1])
if -atid not in SHCdict:
SHCdict[-atid] = {}
SHCdict[-atid]['UVmat'] = parmDict[parm]
if len(SHCdict):
return {pfx:SHCdict,}
else: return {}
################################################################################
##### Penalty & restraint functions
################################################################################
def penaltyFxn(HistoPhases,calcControls,parmDict,varyList):
'Compute user-supplied and built-in restraint functions'
Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases
pNames = []
pVals = []
pWt = []
negWt = {}
pWsum = {}
pWnum = {}
for phase in Phases:
pId = Phases[phase]['pId']
negWt[pId] = Phases[phase]['General']['Pawley neg wt']
General = Phases[phase]['General']
cx,ct,cs,cia = General['AtomPtrs']
textureData = General['SH Texture']
SGData = General['SGData']
Atoms = Phases[phase]['Atoms']
AtLookup = G2mth.FillAtomLookUp(Phases[phase]['Atoms'],cia+8)
cell = General['Cell'][1:7]
Amat,Bmat = G2lat.cell2AB(cell)
if phase not in restraintDict:
continue
phaseRest = restraintDict[phase]
names = [['Bond','Bonds'],['Angle','Angles'],['Plane','Planes'],
['Chiral','Volumes'],['Torsion','Torsions'],['Rama','Ramas'],
['ChemComp','Sites'],['Texture','HKLs'],['General','General'],
['Moments','Moments']]
for name,rest in names:
pWsum[name] = 0.
pWnum[name] = 0
if name not in phaseRest:
continue
itemRest = phaseRest[name]
if itemRest[rest] and itemRest['Use']:
wt = itemRest.get('wtFactor',1.)
if name in ['Bond','Angle','Plane','Chiral']:
for i,[indx,ops,obs,esd] in enumerate(itemRest[rest]):
pNames.append(str(pId)+':'+name+':'+str(i))
XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx))
XYZ = G2mth.getSyXYZ(XYZ,ops,SGData)
if name == 'Bond':
calc = G2mth.getRestDist(XYZ,Amat)
elif name == 'Angle':
calc = G2mth.getRestAngle(XYZ,Amat)
elif name == 'Plane':
calc = G2mth.getRestPlane(XYZ,Amat)
elif name == 'Chiral':
calc = G2mth.getRestChiral(XYZ,Amat)
pVals.append(obs-calc)
pWt.append(wt/esd**2)
pWsum[name] += wt*((obs-calc)/esd)**2
pWnum[name] += 1
elif name in ['Torsion','Rama']:
coeffDict = itemRest['Coeff']
for i,[indx,ops,cofName,esd] in enumerate(itemRest[rest]):
pNames.append(str(pId)+':'+name+':'+str(i))
XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx))
XYZ = G2mth.getSyXYZ(XYZ,ops,SGData)
if name == 'Torsion':
tor = G2mth.getRestTorsion(XYZ,Amat)
restr,calc = G2mth.calcTorsionEnergy(tor,coeffDict[cofName])
else:
phi,psi = G2mth.getRestRama(XYZ,Amat)
restr,calc = G2mth.calcRamaEnergy(phi,psi,coeffDict[cofName])
pVals.append(restr)
pWt.append(wt/esd**2)
pWsum[name] += wt*(restr/esd)**2
pWnum[name] += 1
elif name == 'ChemComp':
for i,[indx,factors,obs,esd] in enumerate(itemRest[rest]):
pNames.append(str(pId)+':'+name+':'+str(i))
mul = np.array(G2mth.GetAtomItemsById(Atoms,AtLookup,indx,cs+1))
frac = np.array(G2mth.GetAtomFracByID(pId,parmDict,AtLookup,indx))
calc = np.sum(mul*frac*factors)
pVals.append(obs-calc)
pWt.append(wt/esd**2)
pWsum[name] += wt*((obs-calc)/esd)**2
pWnum[name] += 1
elif name == 'Moments':
for i,[indx,obs,esd] in enumerate(itemRest[rest]):
pNames.append(str(pId)+':'+name+':'+str(i))
moms = G2mth.GetAtomMomsByID(pId,parmDict,AtLookup,indx)
obs = 0.
calcs = []
for i,mom in enumerate(moms):
calcs.append(G2mth.GetMag(mom,cell))
obs += calcs[-1]
obs /= len(indx)
for calc in calcs:
pVals.append(obs-calc)
pWt.append(wt/esd**2)
pWsum[name] += wt*((obs-calc)/esd)**2
pWnum[name] += 1
elif name == 'Texture':
SHkeys = list(textureData['SH Coeff'][1].keys())
SHCoef = G2mth.GetSHCoeff(pId,parmDict,SHkeys)
shModels = ['cylindrical','none','shear - 2/m','rolling - mmm']
SamSym = dict(zip(shModels,['0','-1','2/m','mmm']))
for i,[hkl,grid,esd1,ifesd2,esd2] in enumerate(itemRest[rest]):
PH = np.array(hkl)
phi,beta = G2lat.CrsAng(np.array(hkl),cell,SGData)
ODFln = G2lat.Flnh(SHCoef,phi,beta,SGData)
R,P,Z = G2mth.getRestPolefig(ODFln,SamSym[textureData['Model']],grid)
Z1 = ma.masked_greater(Z,0.0) #is this + or -?
IndZ1 = np.array(ma.nonzero(Z1))
for ind in IndZ1.T:
pNames.append('%d:%s:%d:%.2f:%.2f'%(pId,name,i,R[ind[0],ind[1]],P[ind[0],ind[1]]))
pVals.append(Z1[ind[0]][ind[1]])
pWt.append(wt/esd1**2)
pWsum[name] += wt*(-Z1[ind[0]][ind[1]]/esd1)**2
pWnum[name] += 1
if ifesd2:
Z2 = 1.-Z
for ind in np.ndindex(grid,grid):
pNames.append('%d:%s:%d:%.2f:%.2f'%(pId,name+'-unit',i,R[ind[0],ind[1]],P[ind[0],ind[1]]))
pVals.append(Z2[ind[0]][ind[1]])
pWt.append(wt/esd2**2)
pWsum[name] += wt*(Z2/esd2)**2
pWnum[name] += 1
elif name == 'General':
for i,(eq,obs,esd) in enumerate(itemRest[rest]):
calcobj = G2obj.ExpressionCalcObj(eq)
calcobj.SetupCalc(parmDict)
calc = calcobj.EvalExpression()
try:
pVals.append(obs-calc)
pWt.append(wt/esd**2)
pWsum[name] += wt*((obs-calc)/esd)**2
pWnum[name] += 1
pNames.append(str(pId)+':'+name+':'+str(i))
except:
print('Error computing General restraint #{}'.format(i+1))
for phase in Phases:
name = 'SH-Pref.Ori.'
pId = Phases[phase]['pId']
General = Phases[phase]['General']
SGData = General['SGData']
cell = General['Cell'][1:7]
pWsum[name] = 0.0
pWnum[name] = 0
for hist in Phases[phase]['Histograms']:
if not Phases[phase]['Histograms'][hist]['Use']:
continue
if hist in Histograms and 'PWDR' in hist:
hId = Histograms[hist]['hId']
phfx = '%d:%d:'%(pId,hId)
if calcControls.get(phfx+'poType','') == 'SH':
toler = calcControls[phfx+'SHtoler']
wt = 1./toler**2
HKLs = np.array(calcControls[phfx+'SHhkl'])
SHnames = calcControls[phfx+'SHnames']
SHcof = dict(zip(SHnames,[parmDict[phfx+cof] for cof in SHnames]))
for i,PH in enumerate(HKLs):
phi,beta = G2lat.CrsAng(PH,cell,SGData)
SH3Coef = {}
for item in SHcof:
L,N = eval(item.strip('C'))
SH3Coef['C%d,0,%d'%(L,N)] = SHcof[item]
ODFln = G2lat.Flnh(SH3Coef,phi,beta,SGData)
X = np.linspace(0,90.0,26)
Y = ma.masked_greater(G2lat.polfcal(ODFln,'0',X,0.0),0.0) #+ or -?
IndY = ma.nonzero(Y)
for ind in IndY[0]:
pNames.append('%d:%d:%s:%d:%.2f'%(pId,hId,name,i,X[ind]))
pVals.append(Y[ind])
pWt.append(wt)
pWsum[name] += wt*(Y[ind])**2
pWnum[name] += 1
pWsum['PWLref'] = 0.
pWnum['PWLref'] = 0
for item in varyList:
if 'PWLref' in item and parmDict[item] < 0.:
pId = int(item.split(':')[0])
if negWt[pId]:
pNames.append(item)
pVals.append(parmDict[item])
pWt.append(negWt[pId])
pWsum['PWLref'] += negWt[pId]*(parmDict[item])**2
pWnum['PWLref'] += 1
pVals = np.array(pVals)
pWt = np.array(pWt) #should this be np.sqrt?
return pNames,pVals,pWt,pWsum,pWnum
def penaltyDeriv(pNames,pVal,HistoPhases,calcControls,parmDict,varyList):
'''Compute derivatives on user-supplied and built-in restraint
(penalty) functions
where pNames is list of restraint labels
:returns: array pDerv: partial derivatives by variable# in varList and
restraint# in pNames (pDerv[variable#][restraint#])
'''
Histograms,Phases,restraintDict,rigidbodyDict = HistoPhases
pDerv = np.zeros((len(varyList),len(pVal)))
for pName in pNames: # loop over restraints
if 'General' == pName.split(':')[1]:
# initialize for General restraint(s) here
parmDict0 = parmDict.copy()
# setup steps for each parameter
stepDict = {}
for parm in varyList:
stepDict[parm] = G2obj.getVarStep(parm,parmDict)
break
for phase in Phases:
# if phase not in restraintDict:
# continue
pId = Phases[phase]['pId']
General = Phases[phase]['General']
cx,ct,cs,cia = General['AtomPtrs']
SGData = General['SGData']
Atoms = Phases[phase]['Atoms']
AtLookup = G2mth.FillAtomLookUp(Phases[phase]['Atoms'],cia+8)
cell = General['Cell'][1:7]
Amat,Bmat = G2lat.cell2AB(cell)
textureData = General['SH Texture']
SHkeys = list(textureData['SH Coeff'][1].keys())
SHCoef = G2mth.GetSHCoeff(pId,parmDict,SHkeys)
shModels = ['cylindrical','none','shear - 2/m','rolling - mmm']
SamSym = dict(zip(shModels,['0','-1','2/m','mmm']))
sam = SamSym[textureData['Model']]
phaseRest = restraintDict.get(phase,{})
names = {'Bond':'Bonds','Angle':'Angles','Plane':'Planes',
'Chiral':'Volumes','Torsion':'Torsions','Rama':'Ramas',
'ChemComp':'Sites','Texture':'HKLs','Moments':'Moments'}
lasthkl = np.array([0,0,0])
for ip,pName in enumerate(pNames): # loop over restraints
pnames = pName.split(':')
if pId == int(pnames[0]):
name = pnames[1]
if 'PWL' in pName:
pDerv[varyList.index(pName)][ip] += 1.
continue
elif 'SH-' in pName:
continue
Id = int(pnames[2])
itemRest = phaseRest[name]
if name in ['Bond','Angle','Plane','Chiral']:
indx,ops,obs,esd = itemRest[names[name]][Id]
dNames = []
for ind in indx:
dNames += [str(pId)+'::dA'+Xname+':'+str(AtLookup[ind]) for Xname in ['x','y','z']]
XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx))
if name == 'Bond':
deriv = G2mth.getRestDeriv(G2mth.getRestDist,XYZ,Amat,ops,SGData)
elif name == 'Angle':
deriv = G2mth.getRestDeriv(G2mth.getRestAngle,XYZ,Amat,ops,SGData)
elif name == 'Plane':
deriv = G2mth.getRestDeriv(G2mth.getRestPlane,XYZ,Amat,ops,SGData)
elif name == 'Chiral':
deriv = G2mth.getRestDeriv(G2mth.getRestChiral,XYZ,Amat,ops,SGData)
elif name in ['Torsion','Rama']:
coffDict = itemRest['Coeff']
indx,ops,cofName,esd = itemRest[names[name]][Id]
dNames = []
for ind in indx:
dNames += [str(pId)+'::dA'+Xname+':'+str(AtLookup[ind]) for Xname in ['x','y','z']]
XYZ = np.array(G2mth.GetAtomCoordsByID(pId,parmDict,AtLookup,indx))
if name == 'Torsion':
deriv = G2mth.getTorsionDeriv(XYZ,Amat,coffDict[cofName])
else:
deriv = G2mth.getRamaDeriv(XYZ,Amat,coffDict[cofName])
elif name == 'ChemComp':
indx,factors,obs,esd = itemRest[names[name]][Id]
dNames = []
for ind in indx:
dNames += [str(pId)+'::Afrac:'+str(AtLookup[ind])]
mul = np.array(G2mth.GetAtomItemsById(Atoms,AtLookup,indx,cs+1))
deriv = mul*factors
elif name == 'Moments':
indx,obs,esd = itemRest[names[name]][Id]
dNames = []
deriv = []
moms = G2mth.GetAtomMomsByID(pId,parmDict,AtLookup,indx)
for i,ind in enumerate(indx):
calc = G2mth.GetMag(moms[i],cell)
dNames += [str(pId)+'::'+Xname+':'+str(AtLookup[ind]) for Xname in ['AMx','AMy','AMz']]
deriv += list(G2mth.GetMagDerv(moms[i],cell)*np.sign((obs-calc)))
elif 'Texture' in name:
deriv = []
dNames = []
hkl,grid,esd1,ifesd2,esd2 = itemRest[names[name]][Id]
hkl = np.array(hkl)
if np.any(lasthkl-hkl):
phi,beta = G2lat.CrsAng(np.array(hkl),cell,SGData)
ODFln = G2lat.Flnh(SHCoef,phi,beta,SGData)
lasthkl = copy.copy(hkl)
if 'unit' in name:
pass
else:
gam = float(pnames[3])
psi = float(pnames[4])
for SHname in ODFln:
l,m,n = eval(SHname[1:])
Ksl = G2lat.GetKsl(l,m,sam,psi,gam)[0]
dNames += [str(pId)+'::'+SHname]
deriv.append(-ODFln[SHname][0]*Ksl/SHCoef[SHname])
elif name == 'General':
deriv = []
dNames = []
eq,obs,esd = itemRest[name][Id]
calcobj = G2obj.ExpressionCalcObj(eq)
parmlist = list(eq.assgnVars.values()) # parameters used in this expression
for parm in parmlist: # expand list if any parms are determined by constraints
if parm in G2mv.GetDependentVars():
parmlist += G2mv.GetIndependentVars()
break
for ind,var in enumerate(varyList):
drv = 0
if var in parmlist:
step = stepDict.get(var,1e-5)
calc = []
# apply step to parameter
oneparm = True
for s in -step,2*step:
parmDict[var] += s
# extend shift if needed to other parameters
if var in G2mv.independentVars:
G2mv.Dict2Map(parmDict)
oneparm = False
elif var in G2mv.dependentVars:
G2mv.Map2Dict(parmDict,[])
oneparm = False
if 'RB' in var:
ApplyRBModels(parmDict,Phases,rigidbodyDict)
# test
oneparm = False
calcobj.SetupCalc(parmDict)
calc.append(calcobj.EvalExpression())
drv = (calc[1]-calc[0])*.5/step
# restore the dict
if oneparm:
parmDict[var] = parmDict0[var]
else:
parmDict = parmDict0.copy()
else:
drv = 0
pDerv[ind][ip] = drv
# Add derivatives into matrix, if needed
for dName,drv in zip(dNames,deriv):
try:
ind = varyList.index(dName)
pDerv[ind][ip] += drv
except ValueError:
pass
lasthkl = np.array([0,0,0])
for ip,pName in enumerate(pNames):
deriv = []
dNames = []
pnames = pName.split(':')
if 'SH-' in pName and pId == int(pnames[0]):
hId = int(pnames[1])
phfx = '%d:%d:'%(pId,hId)
psi = float(pnames[4])
HKLs = calcControls[phfx+'SHhkl']
SHnames = calcControls[phfx+'SHnames']
SHcof = dict(zip(SHnames,[parmDict[phfx+cof] for cof in SHnames]))
hkl = np.array(HKLs[int(pnames[3])])
if np.any(lasthkl-hkl):
phi,beta = G2lat.CrsAng(np.array(hkl),cell,SGData)
SH3Coef = {}
for item in SHcof:
L,N = eval(item.strip('C'))
SH3Coef['C%d,0,%d'%(L,N)] = SHcof[item]
ODFln = G2lat.Flnh(SH3Coef,phi,beta,SGData)
lasthkl = copy.copy(hkl)
for SHname in SHnames:
l,n = eval(SHname[1:])
SH3name = 'C%d,0,%d'%(l,n)
Ksl = G2lat.GetKsl(l,0,'0',psi,0.0)[0]
dNames += [phfx+SHname]
deriv.append(ODFln[SH3name][0]*Ksl/SHcof[SHname])
for dName,drv in zip(dNames,deriv):
try:
ind = varyList.index(dName)
pDerv[ind][ip] += drv
except ValueError:
pass
return pDerv
################################################################################
##### Function & derivative calculations
################################################################################
def GetAtomFXU(pfx,calcControls,parmDict):
'Needs a doc string'
Natoms = calcControls['Natoms'][pfx]
Tdata = Natoms*[' ',]
Mdata = np.zeros(Natoms)
IAdata = Natoms*[' ',]
Fdata = np.zeros(Natoms)
Xdata = np.zeros((3,Natoms))
dXdata = np.zeros((3,Natoms))
Uisodata = np.zeros(Natoms)
Uijdata = np.zeros((6,Natoms))
Gdata = np.zeros((3,Natoms))
keys = {'Atype:':Tdata,'Amul:':Mdata,'Afrac:':Fdata,'AI/A:':IAdata,
'dAx:':dXdata[0],'dAy:':dXdata[1],'dAz:':dXdata[2],
'Ax:':Xdata[0],'Ay:':Xdata[1],'Az:':Xdata[2],'AUiso:':Uisodata,
'AU11:':Uijdata[0],'AU22:':Uijdata[1],'AU33:':Uijdata[2],
'AU12:':Uijdata[3],'AU13:':Uijdata[4],'AU23:':Uijdata[5],
'AMx:':Gdata[0],'AMy:':Gdata[1],'AMz:':Gdata[2],}
for iatm in range(Natoms):
for key in keys:
parm = pfx+key+str(iatm)
if parm in parmDict:
keys[key][iatm] = parmDict[parm]
Fdata = np.where(Fdata,Fdata,1.e-8) #avoid divide by zero in derivative calc.
Gdata = np.where(Gdata,Gdata,1.e-8) #avoid divide by zero in derivative calc.
return Tdata,Mdata,Fdata,Xdata,dXdata,IAdata,Uisodata,Uijdata,Gdata
def GetAtomSSFXU(pfx,calcControls,parmDict):
'Needs a doc string'
Natoms = calcControls['Natoms'][pfx]
maxSSwave = calcControls['maxSSwave'][pfx]
Nwave = {'F':maxSSwave['Sfrac'],'X':maxSSwave['Spos'],'Y':maxSSwave['Spos'],'Z':maxSSwave['Spos'],
'U':maxSSwave['Sadp'],'M':maxSSwave['Smag'],'T':maxSSwave['Spos']}
XSSdata = np.zeros((6,maxSSwave['Spos'],Natoms))
FSSdata = np.zeros((2,maxSSwave['Sfrac'],Natoms))
USSdata = np.zeros((12,maxSSwave['Sadp'],Natoms))
MSSdata = np.zeros((6,maxSSwave['Smag'],Natoms))
waveTypes = []
keys = {'Fsin:':FSSdata[0],'Fcos:':FSSdata[1],'Fzero:':FSSdata[0],'Fwid:':FSSdata[1],
'Tmin:':XSSdata[0],'Tmax:':XSSdata[1],'Xmax:':XSSdata[2],'Ymax:':XSSdata[3],'Zmax:':XSSdata[4],
'Xsin:':XSSdata[0],'Ysin:':XSSdata[1],'Zsin:':XSSdata[2],'Xcos:':XSSdata[3],'Ycos:':XSSdata[4],'Zcos:':XSSdata[5],
'U11sin:':USSdata[0],'U22sin:':USSdata[1],'U33sin:':USSdata[2],'U12sin:':USSdata[3],'U13sin:':USSdata[4],'U23sin:':USSdata[5],
'U11cos:':USSdata[6],'U22cos:':USSdata[7],'U33cos:':USSdata[8],'U12cos:':USSdata[9],'U13cos:':USSdata[10],'U23cos:':USSdata[11],
'MXsin:':MSSdata[0],'MYsin:':MSSdata[1],'MZsin:':MSSdata[2],'MXcos:':MSSdata[3],'MYcos:':MSSdata[4],'MZcos:':MSSdata[5]}
for iatm in range(Natoms):
wavetype = [parmDict.get(pfx+kind+'waveType:'+str(iatm),'') for kind in ['F','P','A','M']]
waveTypes.append(wavetype)
for key in keys:
for m in range(Nwave[key[0]]):
parm = pfx+key+str(iatm)+':%d'%(m)
if parm in parmDict:
keys[key][m][iatm] = parmDict[parm]
return waveTypes,FSSdata,XSSdata,USSdata,MSSdata
def StructureFactor2(refDict,G,hfx,pfx,SGData,calcControls,parmDict):
''' Compute structure factors for all h,k,l for phase
puts the result, F^2, in each ref[8] in refList
operates on blocks of 100 reflections for speed
input:
:param dict refDict: where
'RefList' list where each ref = h,k,l,it,d,...
'FF' dict of form factors - filed in below
:param np.array G: reciprocal metric tensor