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
/
GSASIIsasd.py
1695 lines (1525 loc) · 65.8 KB
/
GSASIIsasd.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
#/usr/bin/env python
# -*- coding: utf-8 -*-
########### SVN repository information ###################
# $Date: 2023-11-06 09:07:21 -0600 (Mon, 06 Nov 2023) $
# $Author: vondreele $
# $Revision: 5698 $
# $URL: https://subversion.xray.aps.anl.gov/pyGSAS/trunk/GSASIIsasd.py $
# $Id: GSASIIsasd.py 5698 2023-11-06 15:07:21Z vondreele $
########### SVN repository information ###################
'''
Classes and routines defined in :mod:`GSASIIsasd` follow.
'''
from __future__ import division, print_function
import os
import math
import numpy as np
import scipy.special as scsp
import scipy.optimize as so
#import pdb
import GSASIIpath
GSASIIpath.SetVersionNumber("$Revision: 5698 $")
import GSASIIpwd as G2pwd
# trig functions in degrees
sind = lambda x: math.sin(x*math.pi/180.)
asind = lambda x: 180.*math.asin(x)/math.pi
tand = lambda x: math.tan(x*math.pi/180.)
atand = lambda x: 180.*math.atan(x)/math.pi
atan2d = lambda y,x: 180.*math.atan2(y,x)/math.pi
cosd = lambda x: math.cos(x*math.pi/180.)
acosd = lambda x: 180.*math.acos(x)/math.pi
rdsq2d = lambda x,p: round(1.0/math.sqrt(x),p)
#numpy versions
npsind = lambda x: np.sin(x*np.pi/180.)
npasind = lambda x: 180.*np.arcsin(x)/math.pi
npcosd = lambda x: np.cos(x*math.pi/180.)
npacosd = lambda x: 180.*np.arccos(x)/math.pi
nptand = lambda x: np.tan(x*math.pi/180.)
npatand = lambda x: 180.*np.arctan(x)/np.pi
npatan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
# npT2stl = lambda tth, wave: 2.0*npsind(tth/2.0)/wave
# npT2q = lambda tth,wave: 2.0*np.pi*npT2stl(tth,wave)
###############################################################################
#### Particle form factors
###############################################################################
def SphereFF(Q,R,args=()):
''' Compute hard sphere form factor - can use numpy arrays
:param float Q: Q value array (usually in A-1)
:param float R: sphere radius (Usually in A - must match Q-1 units)
:param array args: ignored
:returns: form factors as array as needed (float)
'''
QR = Q[:,np.newaxis]*R
return (3./(QR**3))*(np.sin(QR)-(QR*np.cos(QR)))
def SphericalShellFF(Q,R,args=()):
''' Compute spherical shell form factor - can use numpy arrays
:param float Q: Q value array (usually in A-1)
:param float R: sphere radius (Usually in A - must match Q-1 units)
:param array args: [float r]: controls the shell thickness: R_inner = min(r*R,R), R_outer = max(r*R,R)
:returns float: form factors as array as needed
Contributed by: L.A. Avakyan, Southern Federal University, Russia
'''
r = args[0]
if r < 0: # truncate to positive value
r = 0
if r < 1: # r controls inner sphere radius
return SphereFF(Q,R) - SphereFF(Q,R*r)
else: # r controls outer sphere radius
return SphereFF(Q,R*r) - SphereFF(Q,R)
def SpheroidFF(Q,R,args):
''' Compute form factor of cylindrically symmetric ellipsoid (spheroid)
- can use numpy arrays for R & AR; will return corresponding numpy array
param float Q : Q value array (usually in A-1)
param float R: radius along 2 axes of spheroid
param array args: [float AR]: aspect ratio so 3rd axis = R*AR
returns float: form factors as array as needed
'''
NP = 50
AR = args[0]
if 0.99 < AR < 1.01:
return SphereFF(Q,R,0)
else:
cth = np.linspace(0,1.,NP)
try:
Rct = R[:,np.newaxis]*np.sqrt(1.+(AR**2-1.)*cth**2)
except:
Rct = R*np.sqrt(1.+(AR**2-1.)*cth**2)
return np.sqrt(np.sum(SphereFF(Q[:,np.newaxis],Rct,0)**2,axis=2)/NP)
def CylinderFF(Q,R,args):
''' Compute form factor for cylinders - can use numpy arrays
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: [float L]: cylinder length (A)
returns float: form factor
'''
L = args[0]
NP = 200
alp = np.linspace(0,np.pi/2.,NP)
QL = Q[:,np.newaxis]*L
LBessArg = 0.5*QL[:,:,np.newaxis]**np.cos(alp)
LBess = np.where(LBessArg<1.e-6,1.,np.sin(LBessArg)/LBessArg)
QR = Q[:,np.newaxis]*R
SBessArg = QR[:,:,np.newaxis]*np.sin(alp)
SBess = np.where(SBessArg<1.e-6,0.5,scsp.jv(1,SBessArg)/SBessArg)
LSBess = LBess*SBess
return np.sqrt(2.*np.pi*np.sum(np.sin(alp)*LSBess**2,axis=2)/NP)
def CylinderDFF(Q,L,args):
''' Compute form factor for cylinders - can use numpy arrays
param float Q: Q value array (A-1)
param float L: cylinder half length (A)
param array args: [float R]: cylinder radius (A)
returns float: form factor
'''
R = args[0]
return CylinderFF(Q,R,args=[2.*L])
def CylinderARFF(Q,R,args):
''' Compute form factor for cylinders - can use numpy arrays
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: [float AR]: cylinder aspect ratio = L/D = L/2R
returns float: form factor
'''
AR = args[0]
return CylinderFF(Q,R,args=[2.*R*AR])
def UniSphereFF(Q,R,args=0):
''' Compute form factor for unified sphere - can use numpy arrays
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: ignored
returns float: form factor
'''
Rg = np.sqrt(3./5.)*R
B = np.pi*1.62/(Rg**4) #are we missing *np.pi? 1.62 = 6*(3/5)**2/(4/3) sense?
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*Rg/np.sqrt(6)))**3
return np.sqrt(np.exp((-Q[:,np.newaxis]**2*Rg**2)/3.)+(B/QstV**4))
def UniRodFF(Q,R,args):
''' Compute form factor for unified rod - can use numpy arrays
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: [float R]: cylinder radius (A)
returns float: form factor
'''
L = args[0]
Rg2 = np.sqrt(R**2/2+L**2/12)
B2 = np.pi/L
Rg1 = np.sqrt(3.)*R/2.
G1 = (2./3.)*R/L
B1 = 4.*(L+R)/(R**3*L**2)
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*Rg2/np.sqrt(6)))**3
FF = np.exp(-Q[:,np.newaxis]**2*Rg2**2/3.)+(B2/QstV)*np.exp(-Rg1**2*Q[:,np.newaxis]**2/3.)
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*Rg1/np.sqrt(6)))**3
FF += G1*np.exp(-Q[:,np.newaxis]**2*Rg1**2/3.)+(B1/QstV**4)
return np.sqrt(FF)
def UniRodARFF(Q,R,args):
''' Compute form factor for unified rod of fixed aspect ratio - can use numpy arrays
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: [float AR]: cylinder aspect ratio = L/D = L/2R
returns float: form factor
'''
AR = args[0]
return UniRodFF(Q,R,args=[2.*AR*R,])
def UniDiskFF(Q,R,args):
''' Compute form factor for unified disk - can use numpy arrays
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: [float T]: disk thickness (A)
returns float: form factor
'''
T = args[0]
Rg2 = np.sqrt(R**2/2.+T**2/12.)
B2 = 2./R**2
Rg1 = np.sqrt(3.)*T/2.
RgC2 = 1.1*Rg1
G1 = (2./3.)*(T/R)**2
B1 = 4.*(T+R)/(R**3*T**2)
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*Rg2/np.sqrt(6)))**3
FF = np.exp(-Q[:,np.newaxis]**2*Rg2**2/3.)+(B2/QstV**2)*np.exp(-RgC2**2*Q[:,np.newaxis]**2/3.)
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*Rg1/np.sqrt(6)))**3
FF += G1*np.exp(-Q[:,np.newaxis]**2*Rg1**2/3.)+(B1/QstV**4)
return np.sqrt(FF)
def UniTubeFF(Q,R,args):
''' Compute form factor for unified tube - can use numpy arrays
assumes that core of tube is same as the matrix/solvent so contrast
is from tube wall vs matrix
param float Q: Q value array (A-1)
param float R: cylinder radius (A)
param array args: [float L,T]: tube length & wall thickness(A)
returns float: form factor
'''
L,T = args[:2]
Ri = R-T
DR2 = R**2-Ri**2
Vt = np.pi*DR2*L
Rg3 = np.sqrt(DR2/2.+L**2/12.)
B1 = 4.*np.pi**2*(DR2+L*(R+Ri))/Vt**2
B2 = np.pi**2*T/Vt
B3 = np.pi/L
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*Rg3/np.sqrt(6)))**3
FF = np.exp(-Q[:,np.newaxis]**2*Rg3**2/3.)+(B3/QstV)*np.exp(-Q[:,np.newaxis]**2*R**2/3.)
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*R/np.sqrt(6)))**3
FF += (B2/QstV**2)*np.exp(-Q[:,np.newaxis]**2*T**2/3.)
QstV = Q[:,np.newaxis]/(scsp.erf(Q[:,np.newaxis]*T/np.sqrt(6)))**3
FF += B1/QstV**4
return np.sqrt(FF)
###############################################################################
#### Particle volumes
###############################################################################
def SphereVol(R,args=()):
''' Compute volume of sphere
- numpy array friendly
param float R: sphere radius
param array args: ignored
returns float: volume
'''
return (4./3.)*np.pi*R**3
def SphericalShellVol(R,args):
''' Compute volume of spherical shell
- numpy array friendly
param float R: sphere radius
param array args: [float r]: controls shell thickness, see SphericalShellFF description
returns float: volume
'''
r = args[0]
if r < 0:
r = 0
if r < 1:
return SphereVol(R) - SphereVol(R*r)
else:
return SphereVol(R*r) - SphereVol(R)
def SpheroidVol(R,args):
''' Compute volume of cylindrically symmetric ellipsoid (spheroid)
- numpy array friendly
param float R: radius along 2 axes of spheroid
param array args: [float AR]: aspect ratio so radius of 3rd axis = R*AR
returns float: volume
'''
AR = args[0]
return AR*SphereVol(R)
def CylinderVol(R,args):
''' Compute cylinder volume for radius & length
- numpy array friendly
param float R: diameter (A)
param array args: [float L]: length (A)
returns float:volume (A^3)
'''
L = args[0]
return np.pi*L*R**2
def CylinderDVol(L,args):
''' Compute cylinder volume for length & diameter
- numpy array friendly
param float: L half length (A)
param array args: [float D]: diameter (A)
returns float:volume (A^3)
'''
D = args[0]
return CylinderVol(D/2.,[2.*L,])
def CylinderARVol(R,args):
''' Compute cylinder volume for radius & aspect ratio = L/D
- numpy array friendly
param float: R radius (A)
param array args: [float AR]: =L/D=L/2R aspect ratio
returns float:volume
'''
AR = args[0]
return CylinderVol(R,[2.*R*AR,])
def UniSphereVol(R,args=()):
''' Compute volume of sphere
- numpy array friendly
param float R: sphere radius
param array args: ignored
returns float: volume
'''
return SphereVol(R)
def UniRodVol(R,args):
''' Compute cylinder volume for radius & length
- numpy array friendly
param float R: diameter (A)
param array args: [float L]: length (A)
returns float:volume (A^3)
'''
L = args[0]
return CylinderVol(R,[L,])
def UniRodARVol(R,args):
''' Compute rod volume for radius & aspect ratio
- numpy array friendly
param float R: diameter (A)
param array args: [float AR]: =L/D=L/2R aspect ratio
returns float:volume (A^3)
'''
AR = args[0]
return CylinderARVol(R,[AR,])
def UniDiskVol(R,args):
''' Compute disk volume for radius & thickness
- numpy array friendly
param float R: diameter (A)
param array args: [float T]: thickness
returns float:volume (A^3)
'''
T = args[0]
return CylinderVol(R,[T,])
def UniTubeVol(R,args):
''' Compute tube volume for radius, length & wall thickness
- numpy array friendly
param float R: diameter (A)
param array args: [float L,T]: tube length & wall thickness(A)
returns float: volume (A^3) of tube wall
'''
L,T = args[:2]
return CylinderVol(R,[L,])-CylinderVol(R-T,[L,])
################################################################################
#### Distribution functions & their cumulative fxns
################################################################################
def LogNormalDist(x,pos,args):
''' Standard LogNormal distribution - numpy friendly on x axis
ref: http://www.itl.nist.gov/div898/handbook/index.htm 1.3.6.6.9
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: width of distribution (m)
param float shape: shape - (sigma of log(LogNormal))
returns float: LogNormal distribution
'''
scale,shape = args
return np.exp(-np.log((x-pos)/scale)**2/(2.*shape**2))/(np.sqrt(2.*np.pi)*(x-pos)*shape)
def GaussDist(x,pos,args):
''' Standard Normal distribution - numpy friendly on x axis
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: width of distribution (sigma)
param float shape: not used
returns float: Normal distribution
'''
scale = args[0]
return (1./(scale*np.sqrt(2.*np.pi)))*np.exp(-(x-pos)**2/(2.*scale**2))
def LSWDist(x,pos,args=[]):
''' Lifshitz-Slyozov-Wagner Ostwald ripening distribution - numpy friendly on x axis
ref:
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: not used
param float shape: not used
returns float: LSW distribution
'''
redX = x/pos
result = (81.*2**(-5/3.))*(redX**2*np.exp(-redX/(1.5-redX)))/((1.5-redX)**(11/3.)*(3.-redX)**(7/3.))
return np.nan_to_num(result/pos)
def SchulzZimmDist(x,pos,args):
''' Schulz-Zimm macromolecule distribution - numpy friendly on x axis
ref: http://goldbook.iupac.org/S05502.html
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: width of distribution (sigma)
param float shape: not used
returns float: Schulz-Zimm distribution
'''
scale = args[0]
b = (2.*pos/scale)**2
a = b/pos
if b<70: #why bother?
return (a**(b+1.))*x**b*np.exp(-a*x)/scsp.gamma(b+1.)
else:
return np.exp((b+1.)*np.log(a)-scsp.gammaln(b+1.)+b*np.log(x)-(a*x))
def LogNormalCume(x,pos,args):
''' Standard LogNormal cumulative distribution - numpy friendly on x axis
ref: http://www.itl.nist.gov/div898/handbook/index.htm 1.3.6.6.9
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: width of distribution (sigma)
param float shape: shape parameter
returns float: LogNormal cumulative distribution
'''
scale,shape = args
lredX = np.log((x-pos)/scale)
return (scsp.erf((lredX/shape)/np.sqrt(2.))+1.)/2.
def GaussCume(x,pos,args):
''' Standard Normal cumulative distribution - numpy friendly on x axis
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: width of distribution (sigma)
param float shape: not used
returns float: Normal cumulative distribution
'''
scale = args[0]
redX = (x-pos)/scale
return (scsp.erf(redX/np.sqrt(2.))+1.)/2.
def LSWCume(x,pos,args=[]):
''' Lifshitz-Slyozov-Wagner Ostwald ripening cumulative distribution - numpy friendly on x axis
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: not used
param float shape: not used
returns float: LSW cumulative distribution
'''
nP = 500
xMin,xMax = [0.,2.*pos]
X = np.linspace(xMin,xMax,nP,True)
fxn = LSWDist(X,pos)
mat = np.outer(np.ones(nP),fxn)
cume = np.sum(np.tril(mat),axis=1)/np.sum(fxn)
return np.interp(x,X,cume,0,1)
def SchulzZimmCume(x,pos,args):
''' Schulz-Zimm cumulative distribution - numpy friendly on x axis
param float x: independent axis (can be numpy array)
param float pos: location of distribution
param float scale: width of distribution (sigma)
param float shape: not used
returns float: Normal distribution
'''
scale = args[0]
nP = 500
xMin = np.fmax(0.,pos-4.*scale)
xMax = np.fmin(pos+4.*scale,1.e5)
X = np.linspace(xMin,xMax,nP,True)
fxn = LSWDist(X,pos)
mat = np.outer(np.ones(nP),fxn)
cume = np.sum(np.tril(mat),axis=1)/np.sum(fxn)
return np.interp(x,X,cume,0,1)
return []
################################################################################
#### Structure factors for condensed systems
################################################################################
def DiluteSF(Q,args=[]):
''' Default: no structure factor correction for dilute system
'''
return np.ones_like(Q) #or return 1.0
def HardSpheresSF(Q,args):
''' Computes structure factor for not dilute monodisperse hard spheres
Refs.: PERCUS,YEVICK PHYS. REV. 110 1 (1958),THIELE J. CHEM PHYS. 39 474 (1968),
WERTHEIM PHYS. REV. LETT. 47 1462 (1981)
param float Q: Q value array (A-1)
param array args: [float R, float VolFrac]: interparticle distance & volume fraction
returns numpy array S(Q)
'''
R,VolFr = args
denom = (1.-VolFr)**4
num = (1.+2.*VolFr)**2
alp = num/denom
bet = -6.*VolFr*(1.+VolFr/2.)**2/denom
gamm = 0.5*VolFr*num/denom
A = 2.*Q*R
A2 = A**2
A3 = A**3
A4 = A**4
Rca = np.cos(A)
Rsa = np.sin(A)
calp = alp*(Rsa/A2-Rca/A)
cbet = bet*(2.*Rsa/A2-(A2-2.)*Rca/A3-2./A3)
cgam = gamm*(-Rca/A+(4./A)*((3.*A2-6.)*Rca/A4+(A2-6.)*Rsa/A3+6./A4))
pfac = -24.*VolFr/A
C = pfac*(calp+cbet+cgam)
return 1./(1.-C)
def SquareWellSF(Q,args):
'''Computes structure factor for not dilute monodisperse hard sphere with a
square well potential interaction.
Refs.: SHARMA,SHARMA, PHYSICA 89A,(1977),213-
:param float Q: Q value array (A-1)
:param array args: [float R, float VolFrac, float depth, float width]:
interparticle distance, volume fraction (<0.08), well depth (e/kT<1.5kT),
well width
:returns: numpy array S(Q)
well depth > 0 attractive & values outside above limits nonphysical cf.
Monte Carlo simulations
'''
R,VolFr,Depth,Width = args
eta = VolFr
eta2 = eta*eta
eta3 = eta*eta2
eta4 = eta*eta3
etam1 = 1. - eta
etam14 = etam1**4
alp = ( ( (1. + 2.*eta)**2 ) + eta3*( eta-4.0 ) )/etam14
bet = -(eta/3.0) * ( 18. + 20.*eta - 12.*eta2 + eta4 )/etam14
gam = 0.5*eta*( (1. + 2.*eta)**2 + eta3*(eta-4.) )/etam14
SK = 2.*Q*R
SK2 = SK*SK
SK3 = SK*SK2
SK4 = SK3*SK
T1 = alp * SK3 * ( np.sin(SK) - SK * np.cos(SK) )
T2 = bet * SK2 * ( 2.*SK*np.sin(SK) - (SK2-2.)*np.cos(SK) - 2.0 )
T3 = ( 4.0*SK3 - 24.*SK ) * np.sin(SK)
T3 = T3 - ( SK4 - 12.0*SK2 + 24.0 )*np.cos(SK) + 24.0
T3 = gam*T3
T4 = -Depth*SK3*(np.sin(Width*SK) - Width*SK*np.cos(Width*SK)+ SK*np.cos(SK) - np.sin(SK) )
CK = -24.0*eta*( T1 + T2 + T3 + T4 )/SK3/SK3
return 1./(1.-CK)
def StickyHardSpheresSF(Q,args):
''' Computes structure factor for not dilute monodisperse hard spheres
Refs.: PERCUS,YEVICK PHYS. REV. 110 1 (1958),THIELE J. CHEM PHYS. 39 474 (1968),
WERTHEIM PHYS. REV. LETT. 47 1462 (1981)
param float Q: Q value array (A-1)
param array args: [float R, float VolFrac]: sphere radius & volume fraction
returns numpy array S(Q)
'''
R,VolFr,epis,sticky = args
eta = VolFr/(1.0-epis)/(1.0-epis)/(1.0-epis)
sig = 2.0 * R
aa = sig/(1.0 - epis)
etam1 = 1.0 - eta
# SOLVE QUADRATIC FOR LAMBDA
qa = eta/12.0
qb = -1.0*(sticky + eta/(etam1))
qc = (1.0 + eta/2.0)/(etam1*etam1)
radic = qb*qb - 4.0*qa*qc
# KEEP THE SMALLER ROOT, THE LARGER ONE IS UNPHYSICAL
lam1 = (-1.0*qb-np.sqrt(radic))/(2.0*qa)
lam2 = (-1.0*qb+np.sqrt(radic))/(2.0*qa)
lam = min(lam1,lam2)
mu = lam*eta*etam1
alp = (1.0 + 2.0*eta - mu)/(etam1*etam1)
bet = (mu - 3.0*eta)/(2.0*etam1*etam1)
# CALCULATE THE STRUCTURE FACTOR<P></P>
kk = Q*aa
k2 = kk*kk
k3 = kk*k2
ds = np.sin(kk)
dc = np.cos(kk)
aq1 = ((ds - kk*dc)*alp)/k3
aq2 = (bet*(1.0-dc))/k2
aq3 = (lam*ds)/(12.0*kk)
aq = 1.0 + 12.0*eta*(aq1+aq2-aq3)
bq1 = alp*(0.5/kk - ds/k2 + (1.0 - dc)/k3)
bq2 = bet*(1.0/kk - ds/k2)
bq3 = (lam/12.0)*((1.0 - dc)/kk)
bq = 12.0*eta*(bq1+bq2-bq3)
sq = 1.0/(aq*aq +bq*bq)
return sq
#def HayterPenfoldSF(Q,args): #big & ugly function - do later (if ever)
# pass
def InterPrecipitateSF(Q,args):
''' Computes structure factor for precipitates in a matrix
Refs.: E-Wen Huang, Peter K. Liaw, Lionel Porcar, Yun Liu, Yee-Lang Liu,
Ji-Jung Kai, and Wei-Ren Chen,APPLIED PHYSICS LETTERS 93, 161904 (2008)
R. Giordano, A. Grasso, and J. Teixeira, Phys. Rev. A 43, 6894 1991
param float Q: Q value array (A-1)
param array args: [float R, float VolFr]: "radius" & volume fraction
returns numpy array S(Q)
'''
R,VolFr = args
QV2 = Q**2*VolFr**2
top = 1 - np.exp(-QV2/4)*np.cos(2.*Q*R)
bot = 1-2*np.exp(-QV2/4)*np.cos(2.*Q*R) + np.exp(-QV2/2)
return 2*(top/bot) - 1
################################################################################
##### SB-MaxEnt
################################################################################
def G_matrix(q,r,contrast,FFfxn,Volfxn,args=()):
'''Calculates the response matrix :math:`G(Q,r)`
:param float q: :math:`Q`
:param float r: :math:`r`
:param float contrast: :math:`|\\Delta\\rho|^2`, the scattering contrast
:param function FFfxn: form factor function FF(q,r,args)
:param function Volfxn: volume function Vol(r,args)
:returns float: G(Q,r)
'''
FF = FFfxn(q,r,args)
Vol = Volfxn(r,args)
return 1.e-4*(contrast*Vol*FF**2).T #10^-20 vs 10^-24
'''
sbmaxent
Entropy maximization routine as described in the article
J Skilling and RK Bryan; MNRAS 211 (1984) 111 - 124.
("MNRAS": "Monthly Notices of the Royal Astronomical Society")
:license: Copyright (c) 2013, UChicago Argonne, LLC
:license: This file is distributed subject to a Software License Agreement found
in the file LICENSE that is included with this distribution.
References:
1. J Skilling and RK Bryan; MON NOT R ASTR SOC 211 (1984) 111 - 124.
2. JA Potton, GJ Daniell, and BD Rainford; Proc. Workshop
Neutron Scattering Data Analysis, Rutherford
Appleton Laboratory, UK, 1986; ed. MW Johnson,
IOP Conference Series 81 (1986) 81 - 86, Institute
of Physics, Bristol, UK.
3. ID Culverwell and GP Clarke; Ibid. 87 - 96.
4. JA Potton, GK Daniell, & BD Rainford,
J APPL CRYST 21 (1988) 663 - 668.
5. JA Potton, GJ Daniell, & BD Rainford,
J APPL CRYST 21 (1988) 891 - 897.
'''
class MaxEntException(Exception):
'''Any exception from this module'''
pass
def MaxEnt_SB(datum, sigma, G, base, IterMax, image_to_data=None, data_to_image=None, report=False):
'''
do the complete Maximum Entropy algorithm of Skilling and Bryan
:param float datum[]:
:param float sigma[]:
:param float[][] G: transformation matrix
:param float base[]:
:param int IterMax:
:param obj image_to_data: opus function (defaults to opus)
:param obj data_to_image: tropus function (defaults to tropus)
:returns float[]: :math:`f(r) dr`
'''
TEST_LIMIT = 0.05 # for convergence
CHI_SQR_LIMIT = 0.01 # maximum difference in ChiSqr for a solution
SEARCH_DIRECTIONS = 3 # <10. This code requires value = 3
RESET_STRAYS = 1 # was 0.001, correction of stray negative values
DISTANCE_LIMIT_FACTOR = 0.1 # limitation on df to constrain runaways
MAX_MOVE_LOOPS = 5000 # for no solution in routine: move,
MOVE_PASSES = 0.001 # convergence test in routine: move
def tropus (data, G):
'''
tropus: transform data-space -> solution-space: [G] * data
default definition, caller can use this definition or provide an alternative
:param float[M] data: observations, ndarray of shape (M)
:param float[M][N] G: transformation matrix, ndarray of shape (M,N)
:returns float[N]: calculated image, ndarray of shape (N)
'''
return G.dot(data)
def opus (image, G):
'''
opus: transform solution-space -> data-space: [G]^tr * image
default definition, caller can use this definition or provide an alternative
:param float[N] image: solution, ndarray of shape (N)
:param float[M][N] G: transformation matrix, ndarray of shape (M,N)
:returns float[M]: calculated data, ndarray of shape (M)
'''
return np.dot(G.T,image) #G.transpose().dot(image)
def Dist(s2, beta):
'''measure the distance of this possible solution'''
w = 0
n = beta.shape[0]
for k in range(n):
z = -sum(s2[k] * beta)
w += beta[k] * z
return w
def ChiNow(ax, c1, c2, s1, s2):
'''
ChiNow
:returns tuple: (ChiNow computation of ``w``, beta)
'''
bx = 1 - ax
a = bx * c2 - ax * s2
b = -(bx * c1 - ax * s1)
beta = ChoSol(a, b)
w = 1.0
for k in range(SEARCH_DIRECTIONS):
w += beta[k] * (c1[k] + 0.5*sum(c2[k] * beta))
return w, beta
def ChoSol(a, b):
'''
ChoSol: ? chop the solution vectors ?
:returns: new vector beta
'''
n = b.shape[0]
fl = np.zeros((n,n))
bl = np.zeros_like(b)
#print_arr("ChoSol: a", a)
#print_vec("ChoSol: b", b)
if (a[0][0] <= 0):
msg = "ChoSol: a[0][0] = "
msg += str(a[0][0])
msg += ' Value must be positive'
raise MaxEntException(msg)
# first, compute fl from a
# note fl is a lower triangular matrix
fl[0][0] = math.sqrt (a[0][0])
for i in (1, 2):
fl[i][0] = a[i][0] / fl[0][0]
for j in range(1, i+1):
z = 0.0
for k in range(j):
z += fl[i][k] * fl[j][k]
#print "ChoSol: %d %d %d z = %lg" % ( i, j, k, z)
z = a[i][j] - z
if j == i:
y = math.sqrt(max(0.,z))
else:
y = z / fl[j][j]
fl[i][j] = y
#print_arr("ChoSol: fl", fl)
# next, compute bl from fl and b
bl[0] = b[0] / fl[0][0]
for i in (1, 2):
z = 0.0
for k in range(i):
z += fl[i][k] * bl[k]
#print "\t", i, k, z
bl[i] = (b[i] - z) / fl[i][i]
#print_vec("ChoSol: bl", bl)
# last, compute beta from bl and fl
beta = np.empty((n))
beta[-1] = bl[-1] / fl[-1][-1]
for i in (1, 0):
z = 0.0
for k in range(i+1, n):
z += fl[k][i] * beta[k]
#print "\t\t", i, k, 'z=', z
beta[i] = (bl[i] - z) / fl[i][i]
#print_vec("ChoSol: beta", beta)
return beta
def MaxEntMove(fSum, blank, chisq, chizer, c1, c2, s1, s2):
'''
move beta one step closer towards the solution
'''
a_lower, a_upper = 0., 1. # bracket "a"
cmin, beta = ChiNow (a_lower, c1, c2, s1, s2)
#print "MaxEntMove: cmin = %g" % cmin
if cmin*chisq > chizer:
ctarg = (1.0 + cmin)/2
else:
ctarg = chizer/chisq
f_lower = cmin - ctarg
c_upper, beta = ChiNow (a_upper, c1, c2, s1, s2)
f_upper = c_upper - ctarg
fx = 2*MOVE_PASSES # just to start off
loop = 1
while abs(fx) >= MOVE_PASSES and loop <= MAX_MOVE_LOOPS:
a_new = (a_lower + a_upper) * 0.5 # search by bisection
c_new, beta = ChiNow (a_new, c1, c2, s1, s2)
fx = c_new - ctarg
# tighten the search range for the next pass
if f_lower*fx > 0:
a_lower, f_lower = a_new, fx
if f_upper*fx > 0:
a_upper, f_upper = a_new, fx
loop += 1
if abs(fx) >= MOVE_PASSES or loop > MAX_MOVE_LOOPS:
msg = "MaxEntMove: Loop counter = "
msg += str(MAX_MOVE_LOOPS)
msg += ' No convergence in alpha chop'
raise MaxEntException(msg)
w = Dist (s2, beta);
m = SEARCH_DIRECTIONS
if (w > DISTANCE_LIMIT_FACTOR*fSum/blank): # invoke the distance penalty, SB eq. 17
for k in range(m):
beta[k] *= math.sqrt (fSum/(blank*w))
chtarg = ctarg * chisq
return w, chtarg, loop, a_new, fx, beta
#MaxEnt_SB starts here
if image_to_data == None:
image_to_data = opus
if data_to_image == None:
data_to_image = tropus
n = len(base)
npt = len(datum)
# Note that the order of subscripts for
# "xi" and "eta" has been reversed from
# the convention used in the FORTRAN version
# to enable parts of them to be passed as
# as vectors to "image_to_data" and "data_to_image".
xi = np.zeros((SEARCH_DIRECTIONS, n))
eta = np.zeros((SEARCH_DIRECTIONS, npt))
beta = np.zeros((SEARCH_DIRECTIONS))
# s1 = np.zeros((SEARCH_DIRECTIONS))
# c1 = np.zeros((SEARCH_DIRECTIONS))
s2 = np.zeros((SEARCH_DIRECTIONS, SEARCH_DIRECTIONS))
c2 = np.zeros((SEARCH_DIRECTIONS, SEARCH_DIRECTIONS))
# TODO: replace blank (scalar) with base (vector)
blank = sum(base) / len(base) # use the average value of base
chizer, chtarg = npt*1.0, npt*1.0
f = base * 1.0 # starting distribution is base
fSum = sum(f) # find the sum of the f-vector
z = (datum - image_to_data (f, G)) / sigma # standardized residuals, SB eq. 3
chisq = sum(z*z) # Chi^2, SB eq. 4
for iter in range(IterMax):
ox = -2 * z / sigma # gradient of Chi^2
cgrad = data_to_image (ox, G) # cgrad[i] = del(C)/del(f[i]), SB eq. 8
sgrad = -np.log(f/base) / (blank*math.exp (1.0)) # sgrad[i] = del(S)/del(f[i])
snorm = math.sqrt(sum(f * sgrad*sgrad)) # entropy term, SB eq. 22
cnorm = math.sqrt(sum(f * cgrad*cgrad)) # ChiSqr term, SB eq. 22
tnorm = sum(f * sgrad * cgrad) # norm for gradient term TEST
a = 1.0
b = 1.0 / cnorm
if iter == 0:
test = 0.0 # mismatch between entropy and ChiSquared gradients
else:
test = math.sqrt( ( 1.0 - tnorm/(snorm*cnorm) )/2 ) # SB eq. 37?
a = 0.5 / (snorm * test)
b *= 0.5 / test
xi[0] = f * cgrad / cnorm
xi[1] = f * (a * sgrad - b * cgrad)
eta[0] = image_to_data (xi[0], G); # image --> data
eta[1] = image_to_data (xi[1], G); # image --> data
ox = eta[1] / (sigma * sigma)
xi[2] = data_to_image (ox, G); # data --> image
a = 1.0 / math.sqrt(sum(f * xi[2]*xi[2]))
xi[2] = f * xi[2] * a
eta[2] = image_to_data (xi[2], G) # image --> data
# print_arr("MaxEnt: eta.transpose()", eta.transpose())
# print_arr("MaxEnt: xi.transpose()", xi.transpose())
# prepare the search directions for the conjugate gradient technique
c1 = xi.dot(cgrad) / chisq # C_mu, SB eq. 24
s1 = xi.dot(sgrad) # S_mu, SB eq. 24
# print_vec("MaxEnt: c1", c1)
# print_vec("MaxEnt: s1", s1)
for k in range(SEARCH_DIRECTIONS):
for l in range(k+1):
c2[k][l] = 2 * sum(eta[k] * eta[l] / sigma/sigma) / chisq
s2[k][l] = -sum(xi[k] * xi[l] / f) / blank
# print_arr("MaxEnt: c2", c2)
# print_arr("MaxEnt: s2", s2)
# reflect across the body diagonal
for k, l in ((0,1), (0,2), (1,2)):
c2[k][l] = c2[l][k] # M_(mu,nu)
s2[k][l] = s2[l][k] # g_(mu,nu)
beta[0] = -0.5 * c1[0] / c2[0][0]
beta[1] = 0.0
beta[2] = 0.0
if (iter > 0):
w, chtarg, loop, a_new, fx, beta = MaxEntMove(fSum, blank, chisq, chizer, c1, c2, s1, s2)
f_old = f.copy() # preserve the last image
f += xi.transpose().dot(beta) # move the image towards the solution, SB eq. 25
# As mentioned at the top of p.119,
# need to protect against stray negative values.
# In this case, set them to RESET_STRAYS * base[i]
#f = f.clip(RESET_STRAYS * blank, f.max())
for i in range(n):
if f[i] <= 0.0:
f[i] = RESET_STRAYS * base[i]
df = f - f_old
fSum = sum(f)
fChange = sum(df)
# calculate the normalized entropy
S = sum((f/fSum) * np.log(f/fSum)) # normalized entropy, S&B eq. 1
z = (datum - image_to_data (f, G)) / sigma # standardized residuals
chisq = sum(z*z) # report this ChiSq
if report:
print (" MaxEnt trial/max: %3d/%3d" % ((iter+1), IterMax))
print (" Residual: %5.2lf%% Entropy: %8lg" % (100*test, S))
print (" Function sum: %.6lg Change from last: %.2lf%%\n" % (fSum,100*fChange/fSum))
# See if we have finished our task.
# do the hardest test first
if (abs(chisq/chizer-1.0) < CHI_SQR_LIMIT) and (test < TEST_LIMIT):
print (' Convergence achieved.')
return chisq,f,image_to_data(f, G) # solution FOUND returns here
print (' No convergence! Try increasing Error multiplier.')
return chisq,f,image_to_data(f, G) # no solution after IterMax iterations
###############################################################################
#### IPG/TNNLS Routines
###############################################################################
def IPG(datum,sigma,G,Bins,Dbins,IterMax,Qvec=[],approach=0.8,Power=-1,report=False):
''' An implementation of the Interior-Point Gradient method of
Michael Merritt & Yin Zhang, Technical Report TR04-08, Dept. of Comp. and
Appl. Math., Rice Univ., Houston, Texas 77005, U.S.A. found on the web at
http://www.caam.rice.edu/caam/trs/2004/TR04-08.pdf
Problem addressed: Total Non-Negative Least Squares (TNNLS)
:param float datum[]:
:param float sigma[]:
:param float[][] G: transformation matrix
:param int IterMax:
:param float Qvec: data positions for Power = 0-4
:param float approach: 0.8 default fitting parameter
:param int Power: 0-4 for Q^Power weighting, -1 to use input sigma
'''
if Power < 0:
GmatE = G/sigma[:np.newaxis]
IntE = datum/sigma
pwr = 0
QvecP = np.ones_like(datum)
else:
GmatE = G[:]
IntE = datum[:]
pwr = Power
QvecP = Qvec**pwr
Amat = GmatE*QvecP[:np.newaxis]
AAmat = np.inner(Amat,Amat)
Bvec = IntE*QvecP
Xw = np.ones_like(Bins)*1.e-32
calc = np.dot(G.T,Xw)
nIter = 0
err = 10.
while (nIter<IterMax) and (err > 1.):
#Step 1 in M&Z paper:
Qk = np.inner(AAmat,Xw)-np.inner(Amat,Bvec)
Dk = Xw/np.inner(AAmat,Xw)
Pk = -Dk*Qk
#Step 2 in M&Z paper:
alpSt = -np.inner(Pk,Qk)/np.inner(Pk,np.inner(AAmat,Pk))
alpWv = -Xw/Pk
AkSt = [np.where(Pk[i]<0,np.min((approach*alpWv[i],alpSt)),alpSt) for i in range(Pk.shape[0])]
#Step 3 in M&Z paper:
shift = AkSt*Pk
Xw += shift
#done IPG; now results
nIter += 1
calc = np.dot(G.T,Xw)
chisq = np.sum(((datum-calc)/sigma)**2)
err = chisq/len(datum)
if report:
print (' Iteration: %d, chisq: %.3g, sum(shift^2): %.3g'%(nIter,chisq,np.sum(shift**2)))
return chisq,Xw,calc
###############################################################################
#### SASD Utilities
###############################################################################