source: trunk/GSASIIpwd.py @ 2355

Last change on this file since 2355 was 2355, checked in by vondreele, 7 years ago

fix plot bug
add RDF calc - in progress
add background correction to pdf

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 87.2 KB
Line 
1#/usr/bin/env python
2# -*- coding: utf-8 -*-
3'''
4*GSASII powder calculation module*
5==================================
6
7'''
8########### SVN repository information ###################
9# $Date: 2016-06-30 13:36:29 +0000 (Thu, 30 Jun 2016) $
10# $Author: vondreele $
11# $Revision: 2355 $
12# $URL: trunk/GSASIIpwd.py $
13# $Id: GSASIIpwd.py 2355 2016-06-30 13:36:29Z vondreele $
14########### SVN repository information ###################
15import sys
16import math
17import time
18import os
19import subprocess as subp
20
21import numpy as np
22import scipy as sp
23import numpy.linalg as nl
24import numpy.ma as ma
25import random as rand
26from numpy.fft import ifft, fft, fftshift
27import scipy.interpolate as si
28import scipy.stats as st
29import scipy.optimize as so
30
31import GSASIIpath
32GSASIIpath.SetVersionNumber("$Revision: 2355 $")
33import GSASIIlattice as G2lat
34import GSASIIspc as G2spc
35import GSASIIElem as G2elem
36import GSASIIgrid as G2gd
37import GSASIIIO as G2IO
38import GSASIImath as G2mth
39import pypowder as pyd
40try:
41    import pydiffax as pyx
42except ImportError:
43    print 'pydiffax is not available for this platform - under develpment'
44
45   
46# trig functions in degrees
47sind = lambda x: math.sin(x*math.pi/180.)
48asind = lambda x: 180.*math.asin(x)/math.pi
49tand = lambda x: math.tan(x*math.pi/180.)
50atand = lambda x: 180.*math.atan(x)/math.pi
51atan2d = lambda y,x: 180.*math.atan2(y,x)/math.pi
52cosd = lambda x: math.cos(x*math.pi/180.)
53acosd = lambda x: 180.*math.acos(x)/math.pi
54rdsq2d = lambda x,p: round(1.0/math.sqrt(x),p)
55#numpy versions
56npsind = lambda x: np.sin(x*np.pi/180.)
57npasind = lambda x: 180.*np.arcsin(x)/math.pi
58npcosd = lambda x: np.cos(x*math.pi/180.)
59npacosd = lambda x: 180.*np.arccos(x)/math.pi
60nptand = lambda x: np.tan(x*math.pi/180.)
61npatand = lambda x: 180.*np.arctan(x)/np.pi
62npatan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
63npT2stl = lambda tth, wave: 2.0*npsind(tth/2.0)/wave
64npT2q = lambda tth,wave: 2.0*np.pi*npT2stl(tth,wave)
65ateln2 = 8.0*math.log(2.0)
66   
67#GSASII pdf calculation routines
68       
69def Transmission(Geometry,Abs,Diam):
70    '''
71    Calculate sample transmission
72
73    :param str Geometry: one of 'Cylinder','Bragg-Brentano','Tilting flat plate in transmission','Fixed flat plate'
74    :param float Abs: absorption coeff in cm-1
75    :param float Diam: sample thickness/diameter in mm
76    '''
77    if 'Cylinder' in Geometry:      #Lobanov & Alte da Veiga for 2-theta = 0; beam fully illuminates sample
78        MuR = Abs*Diam/20.0
79        if MuR <= 3.0:
80            T0 = 16/(3.*math.pi)
81            T1 = -0.045780
82            T2 = -0.02489
83            T3 = 0.003045
84            T = -T0*MuR-T1*MuR**2-T2*MuR**3-T3*MuR**4
85            if T < -20.:
86                return 2.06e-9
87            else:
88                return math.exp(T)
89        else:
90            T1 = 1.433902
91            T2 = 0.013869+0.337894
92            T3 = 1.933433+1.163198
93            T4 = 0.044365-0.04259
94            T = (T1-T4)/(1.0+T2*(MuR-3.0))**T3+T4
95            return T/100.
96    elif 'plate' in Geometry:
97        MuR = Abs*Diam/10.
98        return math.exp(-MuR)
99    elif 'Bragg' in Geometry:
100        return 0.0
101       
102def SurfaceRough(SRA,SRB,Tth):
103    ''' Suortti (J. Appl. Cryst, 5,325-331, 1972) surface roughness correction
104    :param float SRA: Suortti surface roughness parameter
105    :param float SRB: Suortti surface roughness parameter
106    :param float Tth: 2-theta(deg) - can be numpy array
107   
108    '''
109    sth = npsind(Tth/2.)
110    T1 = np.exp(-SRB/sth)
111    T2 = SRA+(1.-SRA)*np.exp(-SRB)
112    return (SRA+(1.-SRA)*T1)/T2
113   
114def SurfaceRoughDerv(SRA,SRB,Tth):
115    ''' Suortti surface roughness correction derivatives
116    :param float SRA: Suortti surface roughness parameter (dimensionless)
117    :param float SRB: Suortti surface roughness parameter (dimensionless)
118    :param float Tth: 2-theta(deg) - can be numpy array
119    :return list: [dydSRA,dydSRB] derivatives to be used for intensity derivative
120    '''
121    sth = npsind(Tth/2.)
122    T1 = np.exp(-SRB/sth)
123    T2 = SRA+(1.-SRA)*np.exp(-SRB)
124    Trans = (SRA+(1.-SRA)*T1)/T2
125    dydSRA = ((1.-T1)*T2-(1.-np.exp(-SRB))*Trans)/T2**2
126    dydSRB = ((SRA-1.)*T1*T2/sth-Trans*(SRA-T2))/T2**2
127    return [dydSRA,dydSRB]
128
129def Absorb(Geometry,MuR,Tth,Phi=0,Psi=0):
130    '''Calculate sample absorption
131    :param str Geometry: one of 'Cylinder','Bragg-Brentano','Tilting Flat Plate in transmission','Fixed flat plate'
132    :param float MuR: absorption coeff * sample thickness/2 or radius
133    :param Tth: 2-theta scattering angle - can be numpy array
134    :param float Phi: flat plate tilt angle - future
135    :param float Psi: flat plate tilt axis - future
136    '''
137   
138    def muRunder3(MuR,Sth2):
139        T0 = 16.0/(3.*np.pi)
140        T1 = (25.99978-0.01911*Sth2**0.25)*np.exp(-0.024551*Sth2)+ \
141            0.109561*np.sqrt(Sth2)-26.04556
142        T2 = -0.02489-0.39499*Sth2+1.219077*Sth2**1.5- \
143            1.31268*Sth2**2+0.871081*Sth2**2.5-0.2327*Sth2**3
144        T3 = 0.003045+0.018167*Sth2-0.03305*Sth2**2
145        Trns = -T0*MuR-T1*MuR**2-T2*MuR**3-T3*MuR**4
146        return np.exp(Trns)
147   
148    def muRover3(MuR,Sth2):
149        T1 = 1.433902+11.07504*Sth2-8.77629*Sth2*Sth2+ \
150            10.02088*Sth2**3-3.36778*Sth2**4
151        T2 = (0.013869-0.01249*Sth2)*np.exp(3.27094*Sth2)+ \
152            (0.337894+13.77317*Sth2)/(1.0+11.53544*Sth2)**1.555039
153        T3 = 1.933433/(1.0+23.12967*Sth2)**1.686715- \
154            0.13576*np.sqrt(Sth2)+1.163198
155        T4 = 0.044365-0.04259/(1.0+0.41051*Sth2)**148.4202
156        Trns = (T1-T4)/(1.0+T2*(MuR-3.0))**T3+T4
157        return Trns/100.
158       
159    Sth2 = npsind(Tth/2.0)**2
160    Cth2 = 1.-Sth2
161    if 'Cylinder' in Geometry:      #Lobanov & Alte da Veiga for 2-theta = 0; beam fully illuminates sample
162        if 'array' in str(type(MuR)):
163            MuRSTh2 = np.concatenate((MuR,Sth2))
164            AbsCr = np.where(MuRSTh2[0]<=3.0,muRunder3(MuRSTh2[0],MuRSTh2[1]),muRover3(MuRSTh2[0],MuRSTh2[1]))
165            return AbsCr
166        else:
167            if MuR <= 3.0:
168                return muRunder3(MuR,Sth2)
169            else:
170                return muRover3(MuR,Sth2)
171    elif 'Bragg' in Geometry:
172        return 1.0
173    elif 'Fixed' in Geometry: #assumes sample plane is perpendicular to incident beam
174        # and only defined for 2theta < 90
175        MuT = 2.*MuR
176        T1 = np.exp(-MuT)
177        T2 = np.exp(-MuT/npcosd(Tth))
178        Tb = MuT-MuT/npcosd(Tth)
179        return (T2-T1)/Tb
180    elif 'Tilting' in Geometry: #assumes symmetric tilt so sample plane is parallel to diffraction vector
181        MuT = 2.*MuR
182        cth = npcosd(Tth/2.0)
183        return np.exp(-MuT/cth)/cth
184       
185def AbsorbDerv(Geometry,MuR,Tth,Phi=0,Psi=0):
186    'needs a doc string'
187    dA = 0.001
188    AbsP = Absorb(Geometry,MuR+dA,Tth,Phi,Psi)
189    if MuR:
190        AbsM = Absorb(Geometry,MuR-dA,Tth,Phi,Psi)
191        return (AbsP-AbsM)/(2.0*dA)
192    else:
193        return (AbsP-1.)/dA
194       
195def Polarization(Pola,Tth,Azm=0.0):
196    """   Calculate angle dependent x-ray polarization correction (not scaled correctly!)
197
198    :param Pola: polarization coefficient e.g 1.0 fully polarized, 0.5 unpolarized
199    :param Azm: azimuthal angle e.g. 0.0 in plane of polarization
200    :param Tth: 2-theta scattering angle - can be numpy array
201      which (if either) of these is "right"?
202    :return: (pola, dpdPola)
203      * pola = ((1-Pola)*npcosd(Azm)**2+Pola*npsind(Azm)**2)*npcosd(Tth)**2+ \
204        (1-Pola)*npsind(Azm)**2+Pola*npcosd(Azm)**2
205      * dpdPola: derivative needed for least squares
206
207    """
208    pola = ((1.0-Pola)*npcosd(Azm)**2+Pola*npsind(Azm)**2)*npcosd(Tth)**2+   \
209        (1.0-Pola)*npsind(Azm)**2+Pola*npcosd(Azm)**2
210    dpdPola = -npsind(Tth)**2*(npsind(Azm)**2-npcosd(Azm)**2)
211    return pola,dpdPola
212   
213def Oblique(ObCoeff,Tth):
214    'currently assumes detector is normal to beam'
215    if ObCoeff:
216        return (1.-ObCoeff)/(1.0-np.exp(np.log(ObCoeff)/npcosd(Tth)))
217    else:
218        return 1.0
219               
220def Ruland(RulCoff,wave,Q,Compton):
221    'needs a doc string'
222    C = 2.9978e8
223    D = 1.5e-3
224    hmc = 0.024262734687
225    sinth2 = (Q*wave/(4.0*np.pi))**2
226    dlam = (wave**2)*Compton*Q/C
227    dlam_c = 2.0*hmc*sinth2-D*wave**2
228    return 1.0/((1.0+dlam/RulCoff)*(1.0+(np.pi*dlam_c/(dlam+RulCoff))**2))
229   
230def LorchWeight(Q):
231    'needs a doc string'
232    return np.sin(np.pi*(Q[-1]-Q)/(2.0*Q[-1]))
233           
234def GetAsfMean(ElList,Sthl2):
235    '''Calculate various scattering factor terms for PDF calcs
236
237    :param dict ElList: element dictionary contains scattering factor coefficients, etc.
238    :param np.array Sthl2: numpy array of sin theta/lambda squared values
239    :returns: mean(f^2), mean(f)^2, mean(compton)
240    '''
241    sumNoAtoms = 0.0
242    FF = np.zeros_like(Sthl2)
243    FF2 = np.zeros_like(Sthl2)
244    CF = np.zeros_like(Sthl2)
245    for El in ElList:
246        sumNoAtoms += ElList[El]['FormulaNo']
247    for El in ElList:
248        el = ElList[El]
249        ff2 = (G2elem.ScatFac(el,Sthl2)+el['fp'])**2+el['fpp']**2
250        cf = G2elem.ComptonFac(el,Sthl2)
251        FF += np.sqrt(ff2)*el['FormulaNo']/sumNoAtoms
252        FF2 += ff2*el['FormulaNo']/sumNoAtoms
253        CF += cf*el['FormulaNo']/sumNoAtoms
254    return FF2,FF**2,CF
255   
256def GetNumDensity(ElList,Vol):
257    'needs a doc string'
258    sumNoAtoms = 0.0
259    for El in ElList:
260        sumNoAtoms += ElList[El]['FormulaNo']
261    return sumNoAtoms/Vol
262           
263def CalcPDF(data,inst,xydata):
264    'needs a doc string'
265    auxPlot = []
266    import copy
267    import scipy.fftpack as ft
268    #subtract backgrounds - if any
269    xydata['IofQ'] = copy.deepcopy(xydata['Sample'])
270    if data['Sample Bkg.']['Name']:
271        xydata['IofQ'][1][1] += (xydata['Sample Bkg.'][1][1]+
272            data['Sample Bkg.']['Add'])*data['Sample Bkg.']['Mult']
273    if data['Container']['Name']:
274        xycontainer = (xydata['Container'][1][1]+data['Container']['Add'])*data['Container']['Mult']
275        if data['Container Bkg.']['Name']:
276            xycontainer += (xydata['Container Bkg.'][1][1]+
277                data['Container Bkg.']['Add'])*data['Container Bkg.']['Mult']
278        xydata['IofQ'][1][1] += xycontainer
279    #get element data & absorption coeff.
280    ElList = data['ElList']
281    Abs = G2lat.CellAbsorption(ElList,data['Form Vol'])
282    #Apply angle dependent corrections
283    Tth = xydata['Sample'][1][0]
284    dt = (Tth[1]-Tth[0])
285    MuR = Abs*data['Diam']/20.0
286    xydata['IofQ'][1][1] /= Absorb(data['Geometry'],MuR,Tth)
287    if 'X' in inst['Type'][0]:
288        xydata['IofQ'][1][1] /= Polarization(inst['Polariz.'][1],Tth,Azm=inst['Azimuth'][1])[0]
289    if data['DetType'] == 'Image plate':
290        xydata['IofQ'][1][1] *= Oblique(data['ObliqCoeff'],Tth)
291    XY = xydata['IofQ'][1]   
292    #convert to Q
293    hc = 12.397639
294    wave = G2mth.getWave(inst)
295    keV = hc/wave
296    minQ = npT2q(Tth[0],wave)
297    maxQ = npT2q(Tth[-1],wave)   
298    Qpoints = np.linspace(0.,maxQ,len(XY[0]),endpoint=True)
299    dq = Qpoints[1]-Qpoints[0]
300    XY[0] = npT2q(XY[0],wave)   
301#    Qdata = np.nan_to_num(si.griddata(XY[0],XY[1],Qpoints,method='linear')) #only OK for scipy 0.9!
302    T = si.interp1d(XY[0],XY[1],bounds_error=False,fill_value=XY[1][0])      #OK for scipy 0.8
303    Qdata = T(Qpoints)
304    Qdata -= np.min(Qdata)*data['BackRatio']
305   
306    qLimits = data['QScaleLim']
307    minQ = np.searchsorted(Qpoints,qLimits[0])
308    maxQ = np.searchsorted(Qpoints,qLimits[1])
309    newdata = []
310    xydata['IofQ'][1][0] = Qpoints
311    xydata['IofQ'][1][1] = Qdata
312    for item in xydata['IofQ'][1]:
313        newdata.append(item[:maxQ])
314    xydata['IofQ'][1] = newdata
315   
316    xydata['SofQ'] = copy.deepcopy(xydata['IofQ'])
317    FFSq,SqFF,CF = GetAsfMean(ElList,(xydata['SofQ'][1][0]/(4.0*np.pi))**2)  #these are <f^2>,<f>^2,Cf
318    Q = xydata['SofQ'][1][0]
319    ruland = Ruland(data['Ruland'],wave,Q,CF)
320#    auxPlot.append([Q,ruland,'Ruland'])     
321    CF *= ruland
322#    auxPlot.append([Q,CF,'CF-Corr'])
323    scale = np.sum((FFSq+CF)[minQ:maxQ])/np.sum(xydata['SofQ'][1][1][minQ:maxQ])
324    xydata['SofQ'][1][1] *= scale
325    xydata['SofQ'][1][1] -= CF
326    xydata['SofQ'][1][1] = xydata['SofQ'][1][1]/SqFF
327    scale = len(xydata['SofQ'][1][1][minQ:maxQ])/np.sum(xydata['SofQ'][1][1][minQ:maxQ])
328    xydata['SofQ'][1][1] *= scale
329#    if data.get('sinDamp',False):   #not the right thing but leave as place holder
330#        sinDamp = np.sin(np.pi*xydata['SofQ'][1][0]/qLimits[1])
331#        xydata['SofQ'][1][1] *= sinDamp
332    xydata['FofQ'] = copy.deepcopy(xydata['SofQ'])
333    xydata['FofQ'][1][1] = xydata['FofQ'][1][0]*(xydata['SofQ'][1][1]-1.0)
334    if data['Lorch']:
335        xydata['FofQ'][1][1] *= LorchWeight(Q)   
336    xydata['GofR'] = copy.deepcopy(xydata['FofQ'])
337    nR = len(xydata['GofR'][1][1])
338    xydata['GofR'][1][1] = -dq*np.imag(ft.fft(xydata['FofQ'][1][1],4*nR)[:nR])
339    xydata['GofR'][1][0] = 0.5*np.pi*np.linspace(0,nR,nR)/qLimits[1]
340    return auxPlot
341   
342def MakeRDF(RDFcontrols,background,inst,pwddata,xydata):
343#    GSASIIpath.IPyBreak()
344    auxPlot = []
345    if 'C' in inst['Type'][0]:
346        Tth = pwddata[0]
347        hc = 12.397639
348        wave = G2mth.getWave(inst)
349        keV = hc/wave
350        minQ = npT2q(Tth[0],wave)
351        maxQ = npT2q(Tth[-1],wave)   
352       
353
354    elif 'T' in inst['Type'][0]:
355        TOF = pwddata[0]
356        difC = inst['difC'][1]
357        minQ = 2.*np.pi*difC/TOF[-1]
358        maxQ = 2.*np.pi*difC/TOF[0]
359       
360    Qpoints = np.linspace(0.,maxQ,len(pwddata[0]),endpoint=True)
361    dq = Qpoints[1]-Qpoints[0]
362    T = si.interp1d(pwddata[0],pwddata[1],bounds_error=False,fill_value=0.0,kind=RDFcontrols['Smooth'])      #OK for scipy 0.8
363    Qdata = T(Qpoints)
364    auxPlot.append([Qpoints,Qdata,'interp1d:'+RDFcontrols['Smooth']])
365    return auxPlot
366    print 'make RDF'
367
368################################################################################       
369#GSASII peak fitting routines: Finger, Cox & Jephcoat model       
370################################################################################
371
372def factorize(num):
373    ''' Provide prime number factors for integer num
374    :returns: dictionary of prime factors (keys) & power for each (data)
375    '''
376    factors = {}
377    orig = num
378
379    # we take advantage of the fact that (i +1)**2 = i**2 + 2*i +1
380    i, sqi = 2, 4
381    while sqi <= num:
382        while not num%i:
383            num /= i
384            factors[i] = factors.get(i, 0) + 1
385
386        sqi += 2*i + 1
387        i += 1
388
389    if num != 1 and num != orig:
390        factors[num] = factors.get(num, 0) + 1
391
392    if factors:
393        return factors
394    else:
395        return {num:1}          #a prime number!
396           
397def makeFFTsizeList(nmin=1,nmax=1023,thresh=15):
398    ''' Provide list of optimal data sizes for FFT calculations
399
400    :param int nmin: minimum data size >= 1
401    :param int nmax: maximum data size > nmin
402    :param int thresh: maximum prime factor allowed
403    :Returns: list of data sizes where the maximum prime factor is < thresh
404    ''' 
405    plist = []
406    nmin = max(1,nmin)
407    nmax = max(nmin+1,nmax)
408    for p in range(nmin,nmax):
409        if max(factorize(p).keys()) < thresh:
410            plist.append(p)
411    return plist
412
413np.seterr(divide='ignore')
414
415# Normal distribution
416
417# loc = mu, scale = std
418_norm_pdf_C = 1./math.sqrt(2*math.pi)
419class norm_gen(st.rv_continuous):
420    'needs a doc string'
421     
422    def pdf(self,x,*args,**kwds):
423        loc,scale=kwds['loc'],kwds['scale']
424        x = (x-loc)/scale
425        return np.exp(-x**2/2.0) * _norm_pdf_C / scale
426       
427norm = norm_gen(name='norm',longname='A normal',extradoc="""
428
429Normal distribution
430
431The location (loc) keyword specifies the mean.
432The scale (scale) keyword specifies the standard deviation.
433
434normal.pdf(x) = exp(-x**2/2)/sqrt(2*pi)
435""")
436
437## Cauchy
438
439# median = loc
440
441class cauchy_gen(st.rv_continuous):
442    'needs a doc string'
443
444    def pdf(self,x,*args,**kwds):
445        loc,scale=kwds['loc'],kwds['scale']
446        x = (x-loc)/scale
447        return 1.0/np.pi/(1.0+x*x) / scale
448       
449cauchy = cauchy_gen(name='cauchy',longname='Cauchy',extradoc="""
450
451Cauchy distribution
452
453cauchy.pdf(x) = 1/(pi*(1+x**2))
454
455This is the t distribution with one degree of freedom.
456""")
457   
458   
459#GSASII peak fitting routine: Finger, Cox & Jephcoat model       
460
461
462class fcjde_gen(st.rv_continuous):
463    """
464    Finger-Cox-Jephcoat D(2phi,2th) function for S/L = H/L
465    Ref: J. Appl. Cryst. (1994) 27, 892-900.
466
467    :param x: array -1 to 1
468    :param t: 2-theta position of peak
469    :param s: sum(S/L,H/L); S: sample height, H: detector opening,
470      L: sample to detector opening distance
471    :param dx: 2-theta step size in deg
472
473    :returns: for fcj.pdf
474
475     * T = x*dx+t
476     * s = S/L+H/L
477     * if x < 0::
478
479        fcj.pdf = [1/sqrt({cos(T)**2/cos(t)**2}-1) - 1/s]/|cos(T)|
480
481     * if x >= 0: fcj.pdf = 0   
482    """
483    def _pdf(self,x,t,s,dx):
484        T = dx*x+t
485        ax2 = abs(npcosd(T))
486        ax = ax2**2
487        bx = npcosd(t)**2
488        bx = np.where(ax>bx,bx,ax)
489        fx = np.where(ax>bx,(np.sqrt(bx/(ax-bx))-1./s)/ax2,0.0)
490        fx = np.where(fx > 0.,fx,0.0)
491        return fx
492             
493    def pdf(self,x,*args,**kwds):
494        loc=kwds['loc']
495        return self._pdf(x-loc,*args)
496       
497fcjde = fcjde_gen(name='fcjde',shapes='t,s,dx')
498               
499def getWidthsCW(pos,sig,gam,shl):
500    '''Compute the peak widths used for computing the range of a peak
501    for constant wavelength data. On low-angle side, 10 FWHM are used,
502    on high-angle side 15 are used, low angle side extended by axial divergence
503    (for peaks above 90 deg, these are reversed.)
504    '''
505    widths = [np.sqrt(sig)/100.,gam/100.]
506    fwhm = 2.355*widths[0]+widths[1]
507    fmin = 10.*(fwhm+shl*abs(npcosd(pos)))
508    fmax = 15.0*fwhm
509    if pos > 90:
510        fmin,fmax = [fmax,fmin]         
511    return widths,fmin,fmax
512   
513def getWidthsTOF(pos,alp,bet,sig,gam):
514    '''Compute the peak widths used for computing the range of a peak
515    for constant wavelength data. 10 FWHM are used on both sides each
516    extended by exponential coeff.
517    '''
518    lnf = 3.3      # =log(0.001/2)
519    widths = [np.sqrt(sig),gam]
520    fwhm = 2.355*widths[0]+2.*widths[1]
521    fmin = 10.*fwhm*(1.+1./alp)   
522    fmax = 10.*fwhm*(1.+1./bet)
523    return widths,fmin,fmax
524   
525def getFWHM(pos,Inst):
526    '''Compute total FWHM from Thompson, Cox & Hastings (1987), J. Appl. Cryst. 20, 79-83
527    via getgamFW(g,s).
528   
529    :param pos: float peak position in deg 2-theta or tof in musec
530    :param Inst: dict instrument parameters
531   
532    :returns float: total FWHM of pseudoVoigt in deg or musec
533    ''' 
534   
535    sig = lambda Th,U,V,W: np.sqrt(max(0.001,U*tand(Th)**2+V*tand(Th)+W))
536    sigTOF = lambda dsp,S0,S1,S2,Sq:  S0+S1*dsp**2+S2*dsp**4+Sq/dsp**2
537    gam = lambda Th,X,Y: (X/cosd(Th)+Y*tand(Th))
538    gamTOF = lambda dsp,X,Y: X*dsp+Y*dsp**2
539    if 'C' in Inst['Type'][0]:
540        s = sig(pos/2.,Inst['U'][1],Inst['V'][1],Inst['W'][1])
541        g = gam(pos/2.,Inst['X'][1],Inst['Y'][1])
542        return getgamFW(g,s)/100.  #returns FWHM in deg
543    else:
544        dsp = pos/Inst['difC'][0]
545        s = sigTOF(dsp,Inst['sig-0'][1],Inst['sig-1'][1],Inst['sig-2'][1],Inst['sig-q'][1])
546        g = gamTOF(dsp,Inst['X'][1],Inst['Y'][1])
547        return getgamFW(g,s)
548   
549def getgamFW(g,s):
550    '''Compute total FWHM from Thompson, Cox & Hastings (1987), J. Appl. Cryst. 20, 79-83
551    lambda fxn needs FWHM for both Gaussian & Lorentzian components
552   
553    :param g: float Lorentzian gamma = FWHM(L)
554    :param s: float Gaussian sig
555   
556    :returns float: total FWHM of pseudoVoigt
557    ''' 
558    gamFW = lambda s,g: np.exp(np.log(s**5+2.69269*s**4*g+2.42843*s**3*g**2+4.47163*s**2*g**3+0.07842*s*g**4+g**5)/5.)
559    return gamFW(2.35482*s,g)   #sqrt(8ln2)*sig = FWHM(G)
560               
561def getFCJVoigt(pos,intens,sig,gam,shl,xdata):   
562    'needs a doc string'
563    DX = xdata[1]-xdata[0]
564    widths,fmin,fmax = getWidthsCW(pos,sig,gam,shl)
565    x = np.linspace(pos-fmin,pos+fmin,256)
566    dx = x[1]-x[0]
567    Norm = norm.pdf(x,loc=pos,scale=widths[0])
568    Cauchy = cauchy.pdf(x,loc=pos,scale=widths[1])
569    arg = [pos,shl/57.2958,dx,]
570    FCJ = fcjde.pdf(x,*arg,loc=pos)
571    if len(np.nonzero(FCJ)[0])>5:
572        z = np.column_stack([Norm,Cauchy,FCJ]).T
573        Z = fft(z)
574        Df = ifft(Z.prod(axis=0)).real
575    else:
576        z = np.column_stack([Norm,Cauchy]).T
577        Z = fft(z)
578        Df = fftshift(ifft(Z.prod(axis=0))).real
579    Df /= np.sum(Df)
580    Df = si.interp1d(x,Df,bounds_error=False,fill_value=0.0)
581    return intens*Df(xdata)*DX/dx
582
583def getBackground(pfx,parmDict,bakType,dataType,xdata):
584    'needs a doc string'
585    if 'T' in dataType:
586        q = 2.*np.pi*parmDict[pfx+'difC']/xdata
587    elif 'C' in dataType:
588        wave = parmDict.get(pfx+'Lam',parmDict.get(pfx+'Lam1',1.0))
589        q = 2.*np.pi*npsind(xdata/2.)/wave
590    yb = np.zeros_like(xdata)
591    nBak = 0
592    cw = np.diff(xdata)
593    cw = np.append(cw,cw[-1])
594    sumBk = [0.,0.,0]
595    while True:
596        key = pfx+'Back;'+str(nBak)
597        if key in parmDict:
598            nBak += 1
599        else:
600            break
601#empirical functions
602    if bakType in ['chebyschev','cosine']:
603        dt = xdata[-1]-xdata[0]   
604        for iBak in range(nBak):
605            key = pfx+'Back;'+str(iBak)
606            if bakType == 'chebyschev':
607                ybi = parmDict[key]*(2.*(xdata-xdata[0])/dt-1.)**iBak
608            elif bakType == 'cosine':
609                ybi = parmDict[key]*npcosd(xdata*iBak)
610            yb += ybi
611        sumBk[0] = np.sum(yb)
612    elif bakType in ['Q^2 power series','Q^-2 power series']:
613        QT = 1.
614        yb += np.ones_like(yb)*parmDict[pfx+'Back;0']
615        for iBak in range(nBak-1):
616            key = pfx+'Back;'+str(iBak+1)
617            if '-2' in bakType:
618                QT *= (iBak+1)*q**-2
619            else:
620                QT *= q**2/(iBak+1)
621            yb += QT*parmDict[key]
622        sumBk[0] = np.sum(yb)
623    elif bakType in ['lin interpolate','inv interpolate','log interpolate',]:
624        if nBak == 1:
625            yb = np.ones_like(xdata)*parmDict[pfx+'Back;0']
626        elif nBak == 2:
627            dX = xdata[-1]-xdata[0]
628            T2 = (xdata-xdata[0])/dX
629            T1 = 1.0-T2
630            yb = parmDict[pfx+'Back;0']*T1+parmDict[pfx+'Back;1']*T2
631        else:
632            if bakType == 'lin interpolate':
633                bakPos = np.linspace(xdata[0],xdata[-1],nBak,True)
634            elif bakType == 'inv interpolate':
635                bakPos = 1./np.linspace(1./xdata[-1],1./xdata[0],nBak,True)
636            elif bakType == 'log interpolate':
637                bakPos = np.exp(np.linspace(np.log(xdata[0]),np.log(xdata[-1]),nBak,True))
638            bakPos[0] = xdata[0]
639            bakPos[-1] = xdata[-1]
640            bakVals = np.zeros(nBak)
641            for i in range(nBak):
642                bakVals[i] = parmDict[pfx+'Back;'+str(i)]
643            bakInt = si.interp1d(bakPos,bakVals,'linear')
644            yb = bakInt(ma.getdata(xdata))
645        sumBk[0] = np.sum(yb)
646#Debye function       
647    if pfx+'difC' in parmDict:
648        ff = 1.
649    else:       
650        try:
651            wave = parmDict[pfx+'Lam']
652        except KeyError:
653            wave = parmDict[pfx+'Lam1']
654        SQ = (q/(4.*np.pi))**2
655        FF = G2elem.GetFormFactorCoeff('Si')[0]
656        ff = np.array(G2elem.ScatFac(FF,SQ)[0])**2
657    iD = 0       
658    while True:
659        try:
660            dbA = parmDict[pfx+'DebyeA;'+str(iD)]
661            dbR = parmDict[pfx+'DebyeR;'+str(iD)]
662            dbU = parmDict[pfx+'DebyeU;'+str(iD)]
663            ybi = ff*dbA*np.sin(q*dbR)*np.exp(-dbU*q**2)/(q*dbR)
664            yb += ybi
665            sumBk[1] += np.sum(ybi)
666            iD += 1       
667        except KeyError:
668            break
669#peaks
670    iD = 0
671    while True:
672        try:
673            pkP = parmDict[pfx+'BkPkpos;'+str(iD)]
674            pkI = parmDict[pfx+'BkPkint;'+str(iD)]
675            pkS = parmDict[pfx+'BkPksig;'+str(iD)]
676            pkG = parmDict[pfx+'BkPkgam;'+str(iD)]
677            if 'C' in dataType:
678                Wd,fmin,fmax = getWidthsCW(pkP,pkS,pkG,.002)
679            else: #'T'OF
680                Wd,fmin,fmax = getWidthsTOF(pkP,1.,1.,pkS,pkG)
681            iBeg = np.searchsorted(xdata,pkP-fmin)
682            iFin = np.searchsorted(xdata,pkP+fmax)
683            lenX = len(xdata)
684            if not iBeg:
685                iFin = np.searchsorted(xdata,pkP+fmax)
686            elif iBeg == lenX:
687                iFin = iBeg
688            else:
689                iFin = np.searchsorted(xdata,pkP+fmax)
690            if 'C' in dataType:
691                ybi = pkI*getFCJVoigt3(pkP,pkS,pkG,0.002,xdata[iBeg:iFin])
692                yb[iBeg:iFin] += ybi
693            else:   #'T'OF
694                ybi = pkI*getEpsVoigt(pkP,1.,1.,pkS,pkG,xdata[iBeg:iFin])
695                yb[iBeg:iFin] += ybi
696            sumBk[2] += np.sum(ybi)
697            iD += 1       
698        except KeyError:
699            break
700        except ValueError:
701            print '**** WARNING - backround peak '+str(iD)+' sigma is negative; fix & try again ****'
702            break       
703    return yb,sumBk
704   
705def getBackgroundDerv(hfx,parmDict,bakType,dataType,xdata):
706    'needs a doc string'
707    if 'T' in dataType:
708        q = 2.*np.pi*parmDict[hfx+'difC']/xdata
709    elif 'C' in dataType:
710        wave = parmDict.get(hfx+'Lam',parmDict.get(hfx+'Lam1',1.0))
711        q = 2.*np.pi*npsind(xdata/2.)/wave
712    nBak = 0
713    while True:
714        key = hfx+'Back;'+str(nBak)
715        if key in parmDict:
716            nBak += 1
717        else:
718            break
719    dydb = np.zeros(shape=(nBak,len(xdata)))
720    dyddb = np.zeros(shape=(3*parmDict[hfx+'nDebye'],len(xdata)))
721    dydpk = np.zeros(shape=(4*parmDict[hfx+'nPeaks'],len(xdata)))
722    cw = np.diff(xdata)
723    cw = np.append(cw,cw[-1])
724
725    if bakType in ['chebyschev','cosine']:
726        dt = xdata[-1]-xdata[0]   
727        for iBak in range(nBak):   
728            if bakType == 'chebyschev':
729                dydb[iBak] = (2.*(xdata-xdata[0])/dt-1.)**iBak
730            elif bakType == 'cosine':
731                dydb[iBak] = npcosd(xdata*iBak)
732    elif bakType in ['Q^2 power series','Q^-2 power series']:
733        QT = 1.
734        dydb[0] = np.ones_like(xdata)
735        for iBak in range(nBak-1):
736            if '-2' in bakType:
737                QT *= (iBak+1)*q**-2
738            else:
739                QT *= q**2/(iBak+1)
740            dydb[iBak+1] = QT
741    elif bakType in ['lin interpolate','inv interpolate','log interpolate',]:
742        if nBak == 1:
743            dydb[0] = np.ones_like(xdata)
744        elif nBak == 2:
745            dX = xdata[-1]-xdata[0]
746            T2 = (xdata-xdata[0])/dX
747            T1 = 1.0-T2
748            dydb = [T1,T2]
749        else:
750            if bakType == 'lin interpolate':
751                bakPos = np.linspace(xdata[0],xdata[-1],nBak,True)
752            elif bakType == 'inv interpolate':
753                bakPos = 1./np.linspace(1./xdata[-1],1./xdata[0],nBak,True)
754            elif bakType == 'log interpolate':
755                bakPos = np.exp(np.linspace(np.log(xdata[0]),np.log(xdata[-1]),nBak,True))
756            bakPos[0] = xdata[0]
757            bakPos[-1] = xdata[-1]
758            for i,pos in enumerate(bakPos):
759                if i == 0:
760                    dydb[0] = np.where(xdata<bakPos[1],(bakPos[1]-xdata)/(bakPos[1]-bakPos[0]),0.)
761                elif i == len(bakPos)-1:
762                    dydb[i] = np.where(xdata>bakPos[-2],(bakPos[-1]-xdata)/(bakPos[-1]-bakPos[-2]),0.)
763                else:
764                    dydb[i] = np.where(xdata>bakPos[i],
765                        np.where(xdata<bakPos[i+1],(bakPos[i+1]-xdata)/(bakPos[i+1]-bakPos[i]),0.),
766                        np.where(xdata>bakPos[i-1],(xdata-bakPos[i-1])/(bakPos[i]-bakPos[i-1]),0.))
767    if hfx+'difC' in parmDict:
768        ff = 1.
769    else:
770        try:
771            wave = parmDict[hfx+'Lam']
772        except KeyError:
773            wave = parmDict[hfx+'Lam1']
774        q = 4.0*np.pi*npsind(xdata/2.0)/wave
775        SQ = (q/(4*np.pi))**2
776        FF = G2elem.GetFormFactorCoeff('Si')[0]
777        ff = np.array(G2elem.ScatFac(FF,SQ)[0])
778    iD = 0       
779    while True:
780        try:
781            if hfx+'difC' in parmDict:
782                q = 2*np.pi*parmDict[hfx+'difC']/xdata
783            dbA = parmDict[hfx+'DebyeA;'+str(iD)]
784            dbR = parmDict[hfx+'DebyeR;'+str(iD)]
785            dbU = parmDict[hfx+'DebyeU;'+str(iD)]
786            sqr = np.sin(q*dbR)/(q*dbR)
787            cqr = np.cos(q*dbR)
788            temp = np.exp(-dbU*q**2)
789            dyddb[3*iD] = ff*sqr*temp/(np.pi*cw)
790            dyddb[3*iD+1] = ff*dbA*temp*(cqr-sqr)/(np.pi*dbR*cw)
791            dyddb[3*iD+2] = -ff*dbA*sqr*temp*q**2/(np.pi*cw)
792            iD += 1       #ff*dbA*np.sin(q*dbR)*np.exp(-dbU*q**2)/(q*dbR)
793        except KeyError:
794            break
795    iD = 0
796    while True:
797        try:
798            pkP = parmDict[hfx+'BkPkpos;'+str(iD)]
799            pkI = parmDict[hfx+'BkPkint;'+str(iD)]
800            pkS = parmDict[hfx+'BkPksig;'+str(iD)]
801            pkG = parmDict[hfx+'BkPkgam;'+str(iD)]
802            if 'C' in dataType:
803                Wd,fmin,fmax = getWidthsCW(pkP,pkS,pkG,.002)
804            else: #'T'OF
805                Wd,fmin,fmax = getWidthsTOF(pkP,1.,1.,pkS,pkG)
806            iBeg = np.searchsorted(xdata,pkP-fmin)
807            iFin = np.searchsorted(xdata,pkP+fmax)
808            lenX = len(xdata)
809            if not iBeg:
810                iFin = np.searchsorted(xdata,pkP+fmax)
811            elif iBeg == lenX:
812                iFin = iBeg
813            else:
814                iFin = np.searchsorted(xdata,pkP+fmax)
815            if 'C' in dataType:
816                Df,dFdp,dFds,dFdg,x = getdFCJVoigt3(pkP,pkS,pkG,.002,xdata[iBeg:iFin])
817                dydpk[4*iD][iBeg:iFin] += 100.*cw[iBeg:iFin]*pkI*dFdp
818                dydpk[4*iD+1][iBeg:iFin] += 100.*cw[iBeg:iFin]*Df
819                dydpk[4*iD+2][iBeg:iFin] += 100.*cw[iBeg:iFin]*pkI*dFds
820                dydpk[4*iD+3][iBeg:iFin] += 100.*cw[iBeg:iFin]*pkI*dFdg
821            else:   #'T'OF
822                Df,dFdp,x,x,dFds,dFdg = getdEpsVoigt(pkP,1.,1.,pkS,pkG,xdata[iBeg:iFin])
823                dydpk[4*iD][iBeg:iFin] += pkI*dFdp
824                dydpk[4*iD+1][iBeg:iFin] += Df
825                dydpk[4*iD+2][iBeg:iFin] += pkI*dFds
826                dydpk[4*iD+3][iBeg:iFin] += pkI*dFdg
827            iD += 1       
828        except KeyError:
829            break
830        except ValueError:
831            print '**** WARNING - backround peak '+str(iD)+' sigma is negative; fix & try again ****'
832            break       
833    return dydb,dyddb,dydpk
834
835#use old fortran routine
836def getFCJVoigt3(pos,sig,gam,shl,xdata):
837    'needs a doc string'
838   
839    Df = pyd.pypsvfcj(len(xdata),xdata-pos,pos,sig,gam,shl)
840#    Df = pyd.pypsvfcjo(len(xdata),xdata-pos,pos,sig,gam,shl)
841    Df /= np.sum(Df)
842    return Df
843
844def getdFCJVoigt3(pos,sig,gam,shl,xdata):
845    'needs a doc string'
846   
847    Df,dFdp,dFds,dFdg,dFdsh = pyd.pydpsvfcj(len(xdata),xdata-pos,pos,sig,gam,shl)
848#    Df,dFdp,dFds,dFdg,dFdsh = pyd.pydpsvfcjo(len(xdata),xdata-pos,pos,sig,gam,shl)
849    sumDf = np.sum(Df)
850    return Df,dFdp,dFds,dFdg,dFdsh
851
852def getPsVoigt(pos,sig,gam,xdata):
853    'needs a doc string'
854   
855    Df = pyd.pypsvoigt(len(xdata),xdata-pos,sig,gam)
856    Df /= np.sum(Df)
857    return Df
858
859def getdPsVoigt(pos,sig,gam,xdata):
860    'needs a doc string'
861   
862    Df,dFdp,dFds,dFdg = pyd.pydpsvoigt(len(xdata),xdata-pos,sig,gam)
863    sumDf = np.sum(Df)
864    return Df,dFdp,dFds,dFdg
865
866def getEpsVoigt(pos,alp,bet,sig,gam,xdata):
867    'needs a doc string'
868    Df = pyd.pyepsvoigt(len(xdata),xdata-pos,alp,bet,sig,gam)
869    Df /= np.sum(Df)
870    return Df 
871   
872def getdEpsVoigt(pos,alp,bet,sig,gam,xdata):
873    'needs a doc string'
874    Df,dFdp,dFda,dFdb,dFds,dFdg = pyd.pydepsvoigt(len(xdata),xdata-pos,alp,bet,sig,gam)
875    sumDf = np.sum(Df)
876    return Df,dFdp,dFda,dFdb,dFds,dFdg   
877
878def ellipseSize(H,Sij,GB):
879    'needs a doc string'
880    HX = np.inner(H.T,GB)
881    lenHX = np.sqrt(np.sum(HX**2))
882    Esize,Rsize = nl.eigh(G2lat.U6toUij(Sij))           
883    R = np.inner(HX/lenHX,Rsize)*Esize         #want column length for hkl in crystal
884    lenR = np.sqrt(np.sum(R**2))
885    return lenR
886
887def ellipseSizeDerv(H,Sij,GB):
888    'needs a doc string'
889    lenR = ellipseSize(H,Sij,GB)
890    delt = 0.001
891    dRdS = np.zeros(6)
892    for i in range(6):
893        Sij[i] -= delt
894        lenM = ellipseSize(H,Sij,GB)
895        Sij[i] += 2.*delt
896        lenP = ellipseSize(H,Sij,GB)
897        Sij[i] -= delt
898        dRdS[i] = (lenP-lenM)/(2.*delt)
899    return lenR,dRdS
900
901def getHKLpeak(dmin,SGData,A,Inst=None):
902    'needs a doc string'
903    HKL = G2lat.GenHLaue(dmin,SGData,A)       
904    HKLs = []
905    for h,k,l,d in HKL:
906        ext = G2spc.GenHKLf([h,k,l],SGData)[0]
907        if not ext:
908            if Inst == None:
909                HKLs.append([h,k,l,d,0,-1])
910            else:
911                HKLs.append([h,k,l,d,G2lat.Dsp2pos(Inst,d),-1])
912    return HKLs
913
914def getHKLMpeak(dmin,Inst,SGData,SSGData,Vec,maxH,A):
915    'needs a doc string'
916    HKLs = []
917    vec = np.array(Vec)
918    vstar = np.sqrt(G2lat.calc_rDsq(vec,A))     #find extra needed for -n SS reflections
919    dvec = 1./(maxH*vstar+1./dmin)
920    HKL = G2lat.GenHLaue(dvec,SGData,A)       
921    SSdH = [vec*h for h in range(-maxH,maxH+1)]
922    SSdH = dict(zip(range(-maxH,maxH+1),SSdH))
923    for h,k,l,d in HKL:
924        ext = G2spc.GenHKLf([h,k,l],SGData)[0]
925        if not ext and d >= dmin:
926            HKLs.append([h,k,l,0,d,G2lat.Dsp2pos(Inst,d),-1])
927        for dH in SSdH:
928            if dH:
929                DH = SSdH[dH]
930                H = [h+DH[0],k+DH[1],l+DH[2]]
931                d = float(1/np.sqrt(G2lat.calc_rDsq(H,A)))
932                if d >= dmin:
933                    HKLM = np.array([h,k,l,dH])
934                    if G2spc.checkSSextc(HKLM,SSGData):
935                        HKLs.append([h,k,l,dH,d,G2lat.Dsp2pos(Inst,d),-1])   
936#    GSASIIpath.IPyBreak()
937    return G2lat.sortHKLd(HKLs,True,True,True)
938
939def getPeakProfile(dataType,parmDict,xdata,varyList,bakType):
940    'needs a doc string'
941   
942    yb = getBackground('',parmDict,bakType,dataType,xdata)[0]
943    yc = np.zeros_like(yb)
944    cw = np.diff(xdata)
945    cw = np.append(cw,cw[-1])
946    if 'C' in dataType:
947        shl = max(parmDict['SH/L'],0.002)
948        Ka2 = False
949        if 'Lam1' in parmDict.keys():
950            Ka2 = True
951            lamRatio = 360*(parmDict['Lam2']-parmDict['Lam1'])/(np.pi*parmDict['Lam1'])
952            kRatio = parmDict['I(L2)/I(L1)']
953        iPeak = 0
954        while True:
955            try:
956                pos = parmDict['pos'+str(iPeak)]
957                tth = (pos-parmDict['Zero'])
958                intens = parmDict['int'+str(iPeak)]
959                sigName = 'sig'+str(iPeak)
960                if sigName in varyList:
961                    sig = parmDict[sigName]
962                else:
963                    sig = G2mth.getCWsig(parmDict,tth)
964                sig = max(sig,0.001)          #avoid neg sigma^2
965                gamName = 'gam'+str(iPeak)
966                if gamName in varyList:
967                    gam = parmDict[gamName]
968                else:
969                    gam = G2mth.getCWgam(parmDict,tth)
970                gam = max(gam,0.001)             #avoid neg gamma
971                Wd,fmin,fmax = getWidthsCW(pos,sig,gam,shl)
972                iBeg = np.searchsorted(xdata,pos-fmin)
973                iFin = np.searchsorted(xdata,pos+fmin)
974                if not iBeg+iFin:       #peak below low limit
975                    iPeak += 1
976                    continue
977                elif not iBeg-iFin:     #peak above high limit
978                    return yb+yc
979                yc[iBeg:iFin] += intens*getFCJVoigt3(pos,sig,gam,shl,xdata[iBeg:iFin])
980                if Ka2:
981                    pos2 = pos+lamRatio*tand(pos/2.0)       # + 360/pi * Dlam/lam * tan(th)
982                    iBeg = np.searchsorted(xdata,pos2-fmin)
983                    iFin = np.searchsorted(xdata,pos2+fmin)
984                    if iBeg-iFin:
985                        yc[iBeg:iFin] += intens*kRatio*getFCJVoigt3(pos2,sig,gam,shl,xdata[iBeg:iFin])
986                iPeak += 1
987            except KeyError:        #no more peaks to process
988                return yb+yc
989    else:
990        Pdabc = parmDict['Pdabc']
991        difC = parmDict['difC']
992        iPeak = 0
993        while True:
994            try:
995                pos = parmDict['pos'+str(iPeak)]               
996                tof = pos-parmDict['Zero']
997                dsp = tof/difC
998                intens = parmDict['int'+str(iPeak)]
999                alpName = 'alp'+str(iPeak)
1000                if alpName in varyList:
1001                    alp = parmDict[alpName]
1002                else:
1003                    if len(Pdabc):
1004                        alp = np.interp(dsp,Pdabc[0],Pdabc[1])
1005                    else:
1006                        alp = G2mth.getTOFalpha(parmDict,dsp)
1007                alp = max(0.0001,alp)
1008                betName = 'bet'+str(iPeak)
1009                if betName in varyList:
1010                    bet = parmDict[betName]
1011                else:
1012                    if len(Pdabc):
1013                        bet = np.interp(dsp,Pdabc[0],Pdabc[2])
1014                    else:
1015                        bet = G2mth.getTOFbeta(parmDict,dsp)
1016                bet = max(0.0001,bet)
1017                sigName = 'sig'+str(iPeak)
1018                if sigName in varyList:
1019                    sig = parmDict[sigName]
1020                else:
1021                    sig = G2mth.getTOFsig(parmDict,dsp)
1022                gamName = 'gam'+str(iPeak)
1023                if gamName in varyList:
1024                    gam = parmDict[gamName]
1025                else:
1026                    gam = G2mth.getTOFgamma(parmDict,dsp)
1027                gam = max(gam,0.001)             #avoid neg gamma
1028                Wd,fmin,fmax = getWidthsTOF(pos,alp,bet,sig,gam)
1029                iBeg = np.searchsorted(xdata,pos-fmin)
1030                iFin = np.searchsorted(xdata,pos+fmax)
1031                lenX = len(xdata)
1032                if not iBeg:
1033                    iFin = np.searchsorted(xdata,pos+fmax)
1034                elif iBeg == lenX:
1035                    iFin = iBeg
1036                else:
1037                    iFin = np.searchsorted(xdata,pos+fmax)
1038                if not iBeg+iFin:       #peak below low limit
1039                    iPeak += 1
1040                    continue
1041                elif not iBeg-iFin:     #peak above high limit
1042                    return yb+yc
1043                yc[iBeg:iFin] += intens*getEpsVoigt(pos,alp,bet,sig,gam,xdata[iBeg:iFin])
1044                iPeak += 1
1045            except KeyError:        #no more peaks to process
1046                return yb+yc
1047           
1048def getPeakProfileDerv(dataType,parmDict,xdata,varyList,bakType):
1049    'needs a doc string'
1050# needs to return np.array([dMdx1,dMdx2,...]) in same order as varylist = backVary,insVary,peakVary order
1051    dMdv = np.zeros(shape=(len(varyList),len(xdata)))
1052    dMdb,dMddb,dMdpk = getBackgroundDerv('',parmDict,bakType,dataType,xdata)
1053    if 'Back;0' in varyList:            #background derivs are in front if present
1054        dMdv[0:len(dMdb)] = dMdb
1055    names = ['DebyeA','DebyeR','DebyeU']
1056    for name in varyList:
1057        if 'Debye' in name:
1058            parm,id = name.split(';')
1059            ip = names.index(parm)
1060            dMdv[varyList.index(name)] = dMddb[3*int(id)+ip]
1061    names = ['BkPkpos','BkPkint','BkPksig','BkPkgam']
1062    for name in varyList:
1063        if 'BkPk' in name:
1064            parm,id = name.split(';')
1065            ip = names.index(parm)
1066            dMdv[varyList.index(name)] = dMdpk[4*int(id)+ip]
1067    cw = np.diff(xdata)
1068    cw = np.append(cw,cw[-1])
1069    if 'C' in dataType:
1070        shl = max(parmDict['SH/L'],0.002)
1071        Ka2 = False
1072        if 'Lam1' in parmDict.keys():
1073            Ka2 = True
1074            lamRatio = 360*(parmDict['Lam2']-parmDict['Lam1'])/(np.pi*parmDict['Lam1'])
1075            kRatio = parmDict['I(L2)/I(L1)']
1076        iPeak = 0
1077        while True:
1078            try:
1079                pos = parmDict['pos'+str(iPeak)]
1080                tth = (pos-parmDict['Zero'])
1081                intens = parmDict['int'+str(iPeak)]
1082                sigName = 'sig'+str(iPeak)
1083                if sigName in varyList:
1084                    sig = parmDict[sigName]
1085                    dsdU = dsdV = dsdW = 0
1086                else:
1087                    sig = G2mth.getCWsig(parmDict,tth)
1088                    dsdU,dsdV,dsdW = G2mth.getCWsigDeriv(tth)
1089                sig = max(sig,0.001)          #avoid neg sigma
1090                gamName = 'gam'+str(iPeak)
1091                if gamName in varyList:
1092                    gam = parmDict[gamName]
1093                    dgdX = dgdY = 0
1094                else:
1095                    gam = G2mth.getCWgam(parmDict,tth)
1096                    dgdX,dgdY = G2mth.getCWgamDeriv(tth)
1097                gam = max(gam,0.001)             #avoid neg gamma
1098                Wd,fmin,fmax = getWidthsCW(pos,sig,gam,shl)
1099                iBeg = np.searchsorted(xdata,pos-fmin)
1100                iFin = np.searchsorted(xdata,pos+fmin)
1101                if not iBeg+iFin:       #peak below low limit
1102                    iPeak += 1
1103                    continue
1104                elif not iBeg-iFin:     #peak above high limit
1105                    break
1106                dMdpk = np.zeros(shape=(6,len(xdata)))
1107                dMdipk = getdFCJVoigt3(pos,sig,gam,shl,xdata[iBeg:iFin])
1108                for i in range(1,5):
1109                    dMdpk[i][iBeg:iFin] += 100.*cw[iBeg:iFin]*intens*dMdipk[i]
1110                dMdpk[0][iBeg:iFin] += 100.*cw[iBeg:iFin]*dMdipk[0]
1111                dervDict = {'int':dMdpk[0],'pos':dMdpk[1],'sig':dMdpk[2],'gam':dMdpk[3],'shl':dMdpk[4]}
1112                if Ka2:
1113                    pos2 = pos+lamRatio*tand(pos/2.0)       # + 360/pi * Dlam/lam * tan(th)
1114                    iBeg = np.searchsorted(xdata,pos2-fmin)
1115                    iFin = np.searchsorted(xdata,pos2+fmin)
1116                    if iBeg-iFin:
1117                        dMdipk2 = getdFCJVoigt3(pos2,sig,gam,shl,xdata[iBeg:iFin])
1118                        for i in range(1,5):
1119                            dMdpk[i][iBeg:iFin] += 100.*cw[iBeg:iFin]*intens*kRatio*dMdipk2[i]
1120                        dMdpk[0][iBeg:iFin] += 100.*cw[iBeg:iFin]*kRatio*dMdipk2[0]
1121                        dMdpk[5][iBeg:iFin] += 100.*cw[iBeg:iFin]*dMdipk2[0]
1122                        dervDict = {'int':dMdpk[0],'pos':dMdpk[1],'sig':dMdpk[2],'gam':dMdpk[3],'shl':dMdpk[4],'L1/L2':dMdpk[5]*intens}
1123                for parmName in ['pos','int','sig','gam']:
1124                    try:
1125                        idx = varyList.index(parmName+str(iPeak))
1126                        dMdv[idx] = dervDict[parmName]
1127                    except ValueError:
1128                        pass
1129                if 'U' in varyList:
1130                    dMdv[varyList.index('U')] += dsdU*dervDict['sig']
1131                if 'V' in varyList:
1132                    dMdv[varyList.index('V')] += dsdV*dervDict['sig']
1133                if 'W' in varyList:
1134                    dMdv[varyList.index('W')] += dsdW*dervDict['sig']
1135                if 'X' in varyList:
1136                    dMdv[varyList.index('X')] += dgdX*dervDict['gam']
1137                if 'Y' in varyList:
1138                    dMdv[varyList.index('Y')] += dgdY*dervDict['gam']
1139                if 'SH/L' in varyList:
1140                    dMdv[varyList.index('SH/L')] += dervDict['shl']         #problem here
1141                if 'I(L2)/I(L1)' in varyList:
1142                    dMdv[varyList.index('I(L2)/I(L1)')] += dervDict['L1/L2']
1143                iPeak += 1
1144            except KeyError:        #no more peaks to process
1145                break
1146    else:
1147        Pdabc = parmDict['Pdabc']
1148        difC = parmDict['difC']
1149        iPeak = 0
1150        while True:
1151            try:
1152                pos = parmDict['pos'+str(iPeak)]               
1153                tof = pos-parmDict['Zero']
1154                dsp = tof/difC
1155                intens = parmDict['int'+str(iPeak)]
1156                alpName = 'alp'+str(iPeak)
1157                if alpName in varyList:
1158                    alp = parmDict[alpName]
1159                else:
1160                    if len(Pdabc):
1161                        alp = np.interp(dsp,Pdabc[0],Pdabc[1])
1162                        dad0 = 0
1163                    else:
1164                        alp = G2mth.getTOFalpha(parmDict,dsp)
1165                        dada0 = G2mth.getTOFalphaDeriv(dsp)
1166                betName = 'bet'+str(iPeak)
1167                if betName in varyList:
1168                    bet = parmDict[betName]
1169                else:
1170                    if len(Pdabc):
1171                        bet = np.interp(dsp,Pdabc[0],Pdabc[2])
1172                        dbdb0 = dbdb1 = dbdb2 = 0
1173                    else:
1174                        bet = G2mth.getTOFbeta(parmDict,dsp)
1175                        dbdb0,dbdb1,dbdb2 = G2mth.getTOFbetaDeriv(dsp)
1176                sigName = 'sig'+str(iPeak)
1177                if sigName in varyList:
1178                    sig = parmDict[sigName]
1179                    dsds0 = dsds1 = dsds2 = dsds3 = 0
1180                else:
1181                    sig = G2mth.getTOFsig(parmDict,dsp)
1182                    dsds0,dsds1,dsds2,dsds3 = G2mth.getTOFsigDeriv(dsp)
1183                gamName = 'gam'+str(iPeak)
1184                if gamName in varyList:
1185                    gam = parmDict[gamName]
1186                    dsdX = dsdY = 0
1187                else:
1188                    gam = G2mth.getTOFgamma(parmDict,dsp)
1189                    dsdX,dsdY = G2mth.getTOFgammaDeriv(dsp)
1190                gam = max(gam,0.001)             #avoid neg gamma
1191                Wd,fmin,fmax = getWidthsTOF(pos,alp,bet,sig,gam)
1192                iBeg = np.searchsorted(xdata,pos-fmin)
1193                lenX = len(xdata)
1194                if not iBeg:
1195                    iFin = np.searchsorted(xdata,pos+fmax)
1196                elif iBeg == lenX:
1197                    iFin = iBeg
1198                else:
1199                    iFin = np.searchsorted(xdata,pos+fmax)
1200                if not iBeg+iFin:       #peak below low limit
1201                    iPeak += 1
1202                    continue
1203                elif not iBeg-iFin:     #peak above high limit
1204                    break
1205                dMdpk = np.zeros(shape=(7,len(xdata)))
1206                dMdipk = getdEpsVoigt(pos,alp,bet,sig,gam,xdata[iBeg:iFin])
1207                for i in range(1,6):
1208                    dMdpk[i][iBeg:iFin] += intens*cw[iBeg:iFin]*dMdipk[i]
1209                dMdpk[0][iBeg:iFin] += cw[iBeg:iFin]*dMdipk[0]
1210                dervDict = {'int':dMdpk[0],'pos':dMdpk[1],'alp':dMdpk[2],'bet':dMdpk[3],'sig':dMdpk[4],'gam':dMdpk[5]}
1211                for parmName in ['pos','int','alp','bet','sig','gam']:
1212                    try:
1213                        idx = varyList.index(parmName+str(iPeak))
1214                        dMdv[idx] = dervDict[parmName]
1215                    except ValueError:
1216                        pass
1217                if 'alpha' in varyList:
1218                    dMdv[varyList.index('alpha')] += dada0*dervDict['alp']
1219                if 'beta-0' in varyList:
1220                    dMdv[varyList.index('beta-0')] += dbdb0*dervDict['bet']
1221                if 'beta-1' in varyList:
1222                    dMdv[varyList.index('beta-1')] += dbdb1*dervDict['bet']
1223                if 'beta-q' in varyList:
1224                    dMdv[varyList.index('beta-q')] += dbdb2*dervDict['bet']
1225                if 'sig-0' in varyList:
1226                    dMdv[varyList.index('sig-0')] += dsds0*dervDict['sig']
1227                if 'sig-1' in varyList:
1228                    dMdv[varyList.index('sig-1')] += dsds1*dervDict['sig']
1229                if 'sig-2' in varyList:
1230                    dMdv[varyList.index('sig-2')] += dsds2*dervDict['sig']
1231                if 'sig-q' in varyList:
1232                    dMdv[varyList.index('sig-q')] += dsds3*dervDict['sig']
1233                if 'X' in varyList:
1234                    dMdv[varyList.index('X')] += dsdX*dervDict['gam']
1235                if 'Y' in varyList:
1236                    dMdv[varyList.index('Y')] += dsdY*dervDict['gam']         #problem here
1237                iPeak += 1
1238            except KeyError:        #no more peaks to process
1239                break
1240    return dMdv
1241       
1242def Dict2Values(parmdict, varylist):
1243    '''Use before call to leastsq to setup list of values for the parameters
1244    in parmdict, as selected by key in varylist'''
1245    return [parmdict[key] for key in varylist] 
1246   
1247def Values2Dict(parmdict, varylist, values):
1248    ''' Use after call to leastsq to update the parameter dictionary with
1249    values corresponding to keys in varylist'''
1250    parmdict.update(zip(varylist,values))
1251   
1252def SetBackgroundParms(Background):
1253    'needs a doc string'
1254    if len(Background) == 1:            # fix up old backgrounds
1255        BackGround.append({'nDebye':0,'debyeTerms':[]})
1256    bakType,bakFlag = Background[0][:2]
1257    backVals = Background[0][3:]
1258    backNames = ['Back;'+str(i) for i in range(len(backVals))]
1259    Debye = Background[1]           #also has background peaks stuff
1260    backDict = dict(zip(backNames,backVals))
1261    backVary = []
1262    if bakFlag:
1263        backVary = backNames
1264
1265    backDict['nDebye'] = Debye['nDebye']
1266    debyeDict = {}
1267    debyeList = []
1268    for i in range(Debye['nDebye']):
1269        debyeNames = ['DebyeA;'+str(i),'DebyeR;'+str(i),'DebyeU;'+str(i)]
1270        debyeDict.update(dict(zip(debyeNames,Debye['debyeTerms'][i][::2])))
1271        debyeList += zip(debyeNames,Debye['debyeTerms'][i][1::2])
1272    debyeVary = []
1273    for item in debyeList:
1274        if item[1]:
1275            debyeVary.append(item[0])
1276    backDict.update(debyeDict)
1277    backVary += debyeVary
1278
1279    backDict['nPeaks'] = Debye['nPeaks']
1280    peaksDict = {}
1281    peaksList = []
1282    for i in range(Debye['nPeaks']):
1283        peaksNames = ['BkPkpos;'+str(i),'BkPkint;'+str(i),'BkPksig;'+str(i),'BkPkgam;'+str(i)]
1284        peaksDict.update(dict(zip(peaksNames,Debye['peaksList'][i][::2])))
1285        peaksList += zip(peaksNames,Debye['peaksList'][i][1::2])
1286    peaksVary = []
1287    for item in peaksList:
1288        if item[1]:
1289            peaksVary.append(item[0])
1290    backDict.update(peaksDict)
1291    backVary += peaksVary   
1292    return bakType,backDict,backVary
1293   
1294def DoCalibInst(IndexPeaks,Inst):
1295   
1296    def SetInstParms():
1297        dataType = Inst['Type'][0]
1298        insVary = []
1299        insNames = []
1300        insVals = []
1301        for parm in Inst:
1302            insNames.append(parm)
1303            insVals.append(Inst[parm][1])
1304            if parm in ['Lam','difC','difA','difB','Zero',]:
1305                if Inst[parm][2]:
1306                    insVary.append(parm)
1307        instDict = dict(zip(insNames,insVals))
1308        return dataType,instDict,insVary
1309       
1310    def GetInstParms(parmDict,Inst,varyList):
1311        for name in Inst:
1312            Inst[name][1] = parmDict[name]
1313       
1314    def InstPrint(Inst,sigDict):
1315        print 'Instrument Parameters:'
1316        if 'C' in Inst['Type'][0]:
1317            ptfmt = "%12.6f"
1318        else:
1319            ptfmt = "%12.3f"
1320        ptlbls = 'names :'
1321        ptstr =  'values:'
1322        sigstr = 'esds  :'
1323        for parm in Inst:
1324            if parm in  ['Lam','difC','difA','difB','Zero',]:
1325                ptlbls += "%s" % (parm.center(12))
1326                ptstr += ptfmt % (Inst[parm][1])
1327                if parm in sigDict:
1328                    sigstr += ptfmt % (sigDict[parm])
1329                else:
1330                    sigstr += 12*' '
1331        print ptlbls
1332        print ptstr
1333        print sigstr
1334       
1335    def errPeakPos(values,peakDsp,peakPos,peakWt,dataType,parmDict,varyList):
1336        parmDict.update(zip(varyList,values))
1337        return np.sqrt(peakWt)*(G2lat.getPeakPos(dataType,parmDict,peakDsp)-peakPos)
1338
1339    peakPos = []
1340    peakDsp = []
1341    peakWt = []
1342    for peak,sig in zip(IndexPeaks[0],IndexPeaks[1]):
1343        if peak[2] and peak[3] and sig > 0.:
1344            peakPos.append(peak[0])
1345            peakDsp.append(peak[-1])    #d-calc
1346            peakWt.append(peak[-1]**2/sig**2)   #weight by d**2
1347    peakPos = np.array(peakPos)
1348    peakDsp = np.array(peakDsp)
1349    peakWt = np.array(peakWt)
1350    dataType,insDict,insVary = SetInstParms()
1351    parmDict = {}
1352    parmDict.update(insDict)
1353    varyList = insVary
1354    if not len(varyList):
1355        print '**** ERROR - nothing to refine! ****'
1356        return False
1357    while True:
1358        begin = time.time()
1359        values =  np.array(Dict2Values(parmDict, varyList))
1360        result = so.leastsq(errPeakPos,values,full_output=True,ftol=0.000001,
1361            args=(peakDsp,peakPos,peakWt,dataType,parmDict,varyList))
1362        ncyc = int(result[2]['nfev']/2)
1363        runtime = time.time()-begin   
1364        chisq = np.sum(result[2]['fvec']**2)
1365        Values2Dict(parmDict, varyList, result[0])
1366        GOF = chisq/(len(peakPos)-len(varyList))       #reduced chi^2
1367        print 'Number of function calls:',result[2]['nfev'],' Number of observations: ',len(peakPos),' Number of parameters: ',len(varyList)
1368        print 'calib time = %8.3fs, %8.3fs/cycle'%(runtime,runtime/ncyc)
1369        print 'chi**2 = %12.6g, reduced chi**2 = %6.2f'%(chisq,GOF)
1370        try:
1371            sig = np.sqrt(np.diag(result[1])*GOF)
1372            if np.any(np.isnan(sig)):
1373                print '*** Least squares aborted - some invalid esds possible ***'
1374            break                   #refinement succeeded - finish up!
1375        except ValueError:          #result[1] is None on singular matrix
1376            print '**** Refinement failed - singular matrix ****'
1377       
1378    sigDict = dict(zip(varyList,sig))
1379    GetInstParms(parmDict,Inst,varyList)
1380    InstPrint(Inst,sigDict)
1381    return True
1382           
1383def DoPeakFit(FitPgm,Peaks,Background,Limits,Inst,Inst2,data,prevVaryList=[],oneCycle=False,controls=None,dlg=None):
1384    '''Called to perform a peak fit, refining the selected items in the peak
1385    table as well as selected items in the background.
1386
1387    :param str FitPgm: type of fit to perform. At present this is ignored.
1388    :param list Peaks: a list of peaks. Each peak entry is a list with 8 values:
1389      four values followed by a refine flag where the values are: position, intensity,
1390      sigma (Gaussian width) and gamma (Lorentzian width). From the Histogram/"Peak List"
1391      tree entry, dict item "peaks"
1392    :param list Background: describes the background. List with two items.
1393      Item 0 specifies a background model and coefficients. Item 1 is a dict.
1394      From the Histogram/Background tree entry.
1395    :param list Limits: min and max x-value to use
1396    :param dict Inst: Instrument parameters
1397    :param dict Inst2: more Instrument parameters
1398    :param numpy.array data: a 5xn array. data[0] is the x-values,
1399      data[1] is the y-values, data[2] are weight values, data[3], [4] and [5] are
1400      calc, background and difference intensities, respectively.
1401    :param list prevVaryList: Used in sequential refinements to override the
1402      variable list. Defaults as an empty list.
1403    :param bool oneCycle: True if only one cycle of fitting should be performed
1404    :param dict controls: a dict specifying two values, Ftol = controls['min dM/M']
1405      and derivType = controls['deriv type']. If None default values are used.
1406    :param wx.Dialog dlg: A dialog box that is updated with progress from the fit.
1407      Defaults to None, which means no updates are done.
1408    '''
1409    def GetBackgroundParms(parmList,Background):
1410        iBak = 0
1411        while True:
1412            try:
1413                bakName = 'Back;'+str(iBak)
1414                Background[0][iBak+3] = parmList[bakName]
1415                iBak += 1
1416            except KeyError:
1417                break
1418        iDb = 0
1419        while True:
1420            names = ['DebyeA;','DebyeR;','DebyeU;']
1421            try:
1422                for i,name in enumerate(names):
1423                    val = parmList[name+str(iDb)]
1424                    Background[1]['debyeTerms'][iDb][2*i] = val
1425                iDb += 1
1426            except KeyError:
1427                break
1428        iDb = 0
1429        while True:
1430            names = ['BkPkpos;','BkPkint;','BkPksig;','BkPkgam;']
1431            try:
1432                for i,name in enumerate(names):
1433                    val = parmList[name+str(iDb)]
1434                    Background[1]['peaksList'][iDb][2*i] = val
1435                iDb += 1
1436            except KeyError:
1437                break
1438               
1439    def BackgroundPrint(Background,sigDict):
1440        print 'Background coefficients for',Background[0][0],'function'
1441        ptfmt = "%12.5f"
1442        ptstr =  'value: '
1443        sigstr = 'esd  : '
1444        for i,back in enumerate(Background[0][3:]):
1445            ptstr += ptfmt % (back)
1446            if Background[0][1]:
1447                prm = 'Back;'+str(i)
1448                if prm in sigDict:
1449                    sigstr += ptfmt % (sigDict[prm])
1450                else:
1451                    sigstr += " "*12
1452            if len(ptstr) > 75:
1453                print ptstr
1454                if Background[0][1]: print sigstr
1455                ptstr =  'value: '
1456                sigstr = 'esd  : '
1457        if len(ptstr) > 8:
1458            print ptstr
1459            if Background[0][1]: print sigstr
1460
1461        if Background[1]['nDebye']:
1462            parms = ['DebyeA;','DebyeR;','DebyeU;']
1463            print 'Debye diffuse scattering coefficients'
1464            ptfmt = "%12.5f"
1465            print ' term       DebyeA       esd        DebyeR       esd        DebyeU        esd'
1466            for term in range(Background[1]['nDebye']):
1467                line = ' term %d'%(term)
1468                for ip,name in enumerate(parms):
1469                    line += ptfmt%(Background[1]['debyeTerms'][term][2*ip])
1470                    if name+str(term) in sigDict:
1471                        line += ptfmt%(sigDict[name+str(term)])
1472                    else:
1473                        line += " "*12
1474                print line
1475        if Background[1]['nPeaks']:
1476            print 'Coefficients for Background Peaks'
1477            ptfmt = "%15.3f"
1478            for j,pl in enumerate(Background[1]['peaksList']):
1479                names =  'peak %3d:'%(j+1)
1480                ptstr =  'values  :'
1481                sigstr = 'esds    :'
1482                for i,lbl in enumerate(['BkPkpos','BkPkint','BkPksig','BkPkgam']):
1483                    val = pl[2*i]
1484                    prm = lbl+";"+str(j)
1485                    names += '%15s'%(prm)
1486                    ptstr += ptfmt%(val)
1487                    if prm in sigDict:
1488                        sigstr += ptfmt%(sigDict[prm])
1489                    else:
1490                        sigstr += " "*15
1491                print names
1492                print ptstr
1493                print sigstr
1494                           
1495    def SetInstParms(Inst):
1496        dataType = Inst['Type'][0]
1497        insVary = []
1498        insNames = []
1499        insVals = []
1500        for parm in Inst:
1501            insNames.append(parm)
1502            insVals.append(Inst[parm][1])
1503            if parm in ['U','V','W','X','Y','SH/L','I(L2)/I(L1)','alpha',
1504                'beta-0','beta-1','beta-q','sig-0','sig-1','sig-2','sig-q',] and Inst[parm][2]:
1505                    insVary.append(parm)
1506        instDict = dict(zip(insNames,insVals))
1507        instDict['X'] = max(instDict['X'],0.01)
1508        instDict['Y'] = max(instDict['Y'],0.01)
1509        if 'SH/L' in instDict:
1510            instDict['SH/L'] = max(instDict['SH/L'],0.002)
1511        return dataType,instDict,insVary
1512       
1513    def GetInstParms(parmDict,Inst,varyList,Peaks):
1514        for name in Inst:
1515            Inst[name][1] = parmDict[name]
1516        iPeak = 0
1517        while True:
1518            try:
1519                sigName = 'sig'+str(iPeak)
1520                pos = parmDict['pos'+str(iPeak)]
1521                if sigName not in varyList:
1522                    if 'C' in Inst['Type'][0]:
1523                        parmDict[sigName] = G2mth.getCWsig(parmDict,pos)
1524                    else:
1525                        dsp = G2lat.Pos2dsp(Inst,pos)
1526                        parmDict[sigName] = G2mth.getTOFsig(parmDict,dsp)
1527                gamName = 'gam'+str(iPeak)
1528                if gamName not in varyList:
1529                    if 'C' in Inst['Type'][0]:
1530                        parmDict[gamName] = G2mth.getCWgam(parmDict,pos)
1531                    else:
1532                        dsp = G2lat.Pos2dsp(Inst,pos)
1533                        parmDict[gamName] = G2mth.getTOFgamma(parmDict,dsp)
1534                iPeak += 1
1535            except KeyError:
1536                break
1537       
1538    def InstPrint(Inst,sigDict):
1539        print 'Instrument Parameters:'
1540        ptfmt = "%12.6f"
1541        ptlbls = 'names :'
1542        ptstr =  'values:'
1543        sigstr = 'esds  :'
1544        for parm in Inst:
1545            if parm in  ['U','V','W','X','Y','SH/L','I(L2)/I(L1)','alpha',
1546                'beta-0','beta-1','beta-q','sig-0','sig-1','sig-2','sig-q',]:
1547                ptlbls += "%s" % (parm.center(12))
1548                ptstr += ptfmt % (Inst[parm][1])
1549                if parm in sigDict:
1550                    sigstr += ptfmt % (sigDict[parm])
1551                else:
1552                    sigstr += 12*' '
1553        print ptlbls
1554        print ptstr
1555        print sigstr
1556
1557    def SetPeaksParms(dataType,Peaks):
1558        peakNames = []
1559        peakVary = []
1560        peakVals = []
1561        if 'C' in dataType:
1562            names = ['pos','int','sig','gam']
1563        else:
1564            names = ['pos','int','alp','bet','sig','gam']
1565        for i,peak in enumerate(Peaks):
1566            for j,name in enumerate(names):
1567                peakVals.append(peak[2*j])
1568                parName = name+str(i)
1569                peakNames.append(parName)
1570                if peak[2*j+1]:
1571                    peakVary.append(parName)
1572        return dict(zip(peakNames,peakVals)),peakVary
1573               
1574    def GetPeaksParms(Inst,parmDict,Peaks,varyList):
1575        if 'C' in Inst['Type'][0]:
1576            names = ['pos','int','sig','gam']
1577        else:   #'T'
1578            names = ['pos','int','alp','bet','sig','gam']
1579        for i,peak in enumerate(Peaks):
1580            pos = parmDict['pos'+str(i)]
1581            if 'difC' in Inst:
1582                dsp = pos/Inst['difC'][1]
1583            for j in range(len(names)):
1584                parName = names[j]+str(i)
1585                if parName in varyList:
1586                    peak[2*j] = parmDict[parName]
1587                elif 'alpha' in parName:
1588                    peak[2*j] = parmDict['alpha']/dsp
1589                elif 'beta' in parName:
1590                    peak[2*j] = G2mth.getTOFbeta(parmDict,dsp)
1591                elif 'sig' in parName:
1592                    if 'C' in Inst['Type'][0]:
1593                        peak[2*j] = G2mth.getCWsig(parmDict,pos)
1594                    else:
1595                        peak[2*j] = G2mth.getTOFsig(parmDict,dsp)
1596                elif 'gam' in parName:
1597                    if 'C' in Inst['Type'][0]:
1598                        peak[2*j] = G2mth.getCWgam(parmDict,pos)
1599                    else:
1600                        peak[2*j] = G2mth.getTOFgamma(parmDict,dsp)
1601                       
1602    def PeaksPrint(dataType,parmDict,sigDict,varyList):
1603        print 'Peak coefficients:'
1604        if 'C' in dataType:
1605            names = ['pos','int','sig','gam']
1606        else:   #'T'
1607            names = ['pos','int','alp','bet','sig','gam']           
1608        head = 13*' '
1609        for name in names:
1610            if name in ['alp','bet']:
1611                head += name.center(8)+'esd'.center(8)
1612            else:
1613                head += name.center(10)+'esd'.center(10)
1614        print head
1615        if 'C' in dataType:
1616            ptfmt = {'pos':"%10.5f",'int':"%10.1f",'sig':"%10.3f",'gam':"%10.3f"}
1617        else:
1618            ptfmt = {'pos':"%10.2f",'int':"%10.4f",'alp':"%8.3f",'bet':"%8.5f",'sig':"%10.3f",'gam':"%10.3f"}
1619        for i,peak in enumerate(Peaks):
1620            ptstr =  ':'
1621            for j in range(len(names)):
1622                name = names[j]
1623                parName = name+str(i)
1624                ptstr += ptfmt[name] % (parmDict[parName])
1625                if parName in varyList:
1626                    ptstr += ptfmt[name] % (sigDict[parName])
1627                else:
1628                    if name in ['alp','bet']:
1629                        ptstr += 8*' '
1630                    else:
1631                        ptstr += 10*' '
1632            print '%s'%(('Peak'+str(i+1)).center(8)),ptstr
1633               
1634    def devPeakProfile(values,xdata,ydata, weights,dataType,parmdict,varylist,bakType,dlg):
1635        parmdict.update(zip(varylist,values))
1636        return np.sqrt(weights)*getPeakProfileDerv(dataType,parmdict,xdata,varylist,bakType)
1637           
1638    def errPeakProfile(values,xdata,ydata,weights,dataType,parmdict,varylist,bakType,dlg):       
1639        parmdict.update(zip(varylist,values))
1640        M = np.sqrt(weights)*(getPeakProfile(dataType,parmdict,xdata,varylist,bakType)-ydata)
1641        Rwp = min(100.,np.sqrt(np.sum(M**2)/np.sum(weights*ydata**2))*100.)
1642        if dlg:
1643            GoOn = dlg.Update(Rwp,newmsg='%s%8.3f%s'%('Peak fit Rwp =',Rwp,'%'))[0]
1644            if not GoOn:
1645                return -M           #abort!!
1646        return M
1647       
1648    if controls:
1649        Ftol = controls['min dM/M']
1650        derivType = controls['deriv type']
1651    else:
1652        Ftol = 0.0001
1653        derivType = 'analytic'
1654    if oneCycle:
1655        Ftol = 1.0
1656    x,y,w,yc,yb,yd = data               #these are numpy arrays!
1657    yc *= 0.                            #set calcd ones to zero
1658    yb *= 0.
1659    yd *= 0.
1660    xBeg = np.searchsorted(x,Limits[0])
1661    xFin = np.searchsorted(x,Limits[1])
1662    bakType,bakDict,bakVary = SetBackgroundParms(Background)
1663    dataType,insDict,insVary = SetInstParms(Inst)
1664    peakDict,peakVary = SetPeaksParms(Inst['Type'][0],Peaks)
1665    parmDict = {}
1666    parmDict.update(bakDict)
1667    parmDict.update(insDict)
1668    parmDict.update(peakDict)
1669    parmDict['Pdabc'] = []      #dummy Pdabc
1670    parmDict.update(Inst2)      #put in real one if there
1671    if prevVaryList:
1672        varyList = prevVaryList[:]
1673    else:
1674        varyList = bakVary+insVary+peakVary
1675    fullvaryList = varyList[:]
1676    while True:
1677        begin = time.time()
1678        values =  np.array(Dict2Values(parmDict, varyList))
1679        Rvals = {}
1680        badVary = []
1681        result = so.leastsq(errPeakProfile,values,Dfun=devPeakProfile,full_output=True,ftol=Ftol,col_deriv=True,
1682               args=(x[xBeg:xFin],y[xBeg:xFin],w[xBeg:xFin],dataType,parmDict,varyList,bakType,dlg))
1683        ncyc = int(result[2]['nfev']/2)
1684        runtime = time.time()-begin   
1685        chisq = np.sum(result[2]['fvec']**2)
1686        Values2Dict(parmDict, varyList, result[0])
1687        Rvals['Rwp'] = np.sqrt(chisq/np.sum(w[xBeg:xFin]*y[xBeg:xFin]**2))*100.      #to %
1688        Rvals['GOF'] = chisq/(xFin-xBeg-len(varyList))       #reduced chi^2
1689        print 'Number of function calls:',result[2]['nfev'],' Number of observations: ',xFin-xBeg,' Number of parameters: ',len(varyList)
1690        if ncyc:
1691            print 'fitpeak time = %8.3fs, %8.3fs/cycle'%(runtime,runtime/ncyc)
1692        print 'Rwp = %7.2f%%, chi**2 = %12.6g, reduced chi**2 = %6.2f'%(Rvals['Rwp'],chisq,Rvals['GOF'])
1693        sig = [0]*len(varyList)
1694        if len(varyList) == 0: break  # if nothing was refined
1695        try:
1696            sig = np.sqrt(np.diag(result[1])*Rvals['GOF'])
1697            if np.any(np.isnan(sig)):
1698                print '*** Least squares aborted - some invalid esds possible ***'
1699            break                   #refinement succeeded - finish up!
1700        except ValueError:          #result[1] is None on singular matrix
1701            print '**** Refinement failed - singular matrix ****'
1702            Ipvt = result[2]['ipvt']
1703            for i,ipvt in enumerate(Ipvt):
1704                if not np.sum(result[2]['fjac'],axis=1)[i]:
1705                    print 'Removing parameter: ',varyList[ipvt-1]
1706                    badVary.append(varyList[ipvt-1])
1707                    del(varyList[ipvt-1])
1708                    break
1709            else: # nothing removed
1710                break
1711    if dlg: dlg.Destroy()
1712    sigDict = dict(zip(varyList,sig))
1713    yb[xBeg:xFin] = getBackground('',parmDict,bakType,dataType,x[xBeg:xFin])[0]
1714    yc[xBeg:xFin] = getPeakProfile(dataType,parmDict,x[xBeg:xFin],varyList,bakType)
1715    yd[xBeg:xFin] = y[xBeg:xFin]-yc[xBeg:xFin]
1716    GetBackgroundParms(parmDict,Background)
1717    if bakVary: BackgroundPrint(Background,sigDict)
1718    GetInstParms(parmDict,Inst,varyList,Peaks)
1719    if insVary: InstPrint(Inst,sigDict)
1720    GetPeaksParms(Inst,parmDict,Peaks,varyList)   
1721    if peakVary: PeaksPrint(dataType,parmDict,sigDict,varyList)
1722    return sigDict,result,sig,Rvals,varyList,parmDict,fullvaryList,badVary
1723
1724def calcIncident(Iparm,xdata):
1725    'needs a doc string'
1726
1727    def IfunAdv(Iparm,xdata):
1728        Itype = Iparm['Itype']
1729        Icoef = Iparm['Icoeff']
1730        DYI = np.ones((12,xdata.shape[0]))
1731        YI = np.ones_like(xdata)*Icoef[0]
1732       
1733        x = xdata/1000.                 #expressions are in ms
1734        if Itype == 'Exponential':
1735            for i in [1,3,5,7,9]:
1736                Eterm = np.exp(-Icoef[i+1]*x**((i+1)/2))
1737                YI += Icoef[i]*Eterm
1738                DYI[i] *= Eterm
1739                DYI[i+1] *= -Icoef[i]*Eterm*x**((i+1)/2)           
1740        elif 'Maxwell'in Itype:
1741            Eterm = np.exp(-Icoef[2]/x**2)
1742            DYI[1] = Eterm/x**5
1743            DYI[2] = -Icoef[1]*DYI[1]/x**2
1744            YI += (Icoef[1]*Eterm/x**5)
1745            if 'Exponential' in Itype:
1746                for i in range(3,12,2):
1747                    Eterm = np.exp(-Icoef[i+1]*x**((i+1)/2))
1748                    YI += Icoef[i]*Eterm
1749                    DYI[i] *= Eterm
1750                    DYI[i+1] *= -Icoef[i]*Eterm*x**((i+1)/2)
1751            else:   #Chebyschev
1752                T = (2./x)-1.
1753                Ccof = np.ones((12,xdata.shape[0]))
1754                Ccof[1] = T
1755                for i in range(2,12):
1756                    Ccof[i] = 2*T*Ccof[i-1]-Ccof[i-2]
1757                for i in range(1,10):
1758                    YI += Ccof[i]*Icoef[i+2]
1759                    DYI[i+2] =Ccof[i]
1760        return YI,DYI
1761       
1762    Iesd = np.array(Iparm['Iesd'])
1763    Icovar = Iparm['Icovar']
1764    YI,DYI = IfunAdv(Iparm,xdata)
1765    YI = np.where(YI>0,YI,1.)
1766    WYI = np.zeros_like(xdata)
1767    vcov = np.zeros((12,12))
1768    k = 0
1769    for i in range(12):
1770        for j in range(i,12):
1771            vcov[i][j] = Icovar[k]*Iesd[i]*Iesd[j]
1772            vcov[j][i] = Icovar[k]*Iesd[i]*Iesd[j]
1773            k += 1
1774    M = np.inner(vcov,DYI.T)
1775    WYI = np.sum(M*DYI,axis=0)
1776    WYI = np.where(WYI>0.,WYI,0.)
1777    return YI,WYI
1778   
1779################################################################################
1780# Stacking fault simulation codes
1781################################################################################
1782
1783def GetStackParms(Layers):
1784   
1785    Parms = []
1786#cell parms
1787    cell = Layers['Cell'] 
1788    if Layers['Laue'] in ['-3','-3m','4/m','4/mmm','6/m','6/mmm']:
1789        Parms.append('cellA')
1790        Parms.append('cellC')
1791    else:
1792        Parms.append('cellA')
1793        Parms.append('cellB')
1794        Parms.append('cellC')
1795        if Layers['Laue'] != 'mmm':
1796            Parms.append('cellG')
1797#Transition parms
1798    Trans = Layers['Transitions']
1799    for iY in range(len(Layers['Layers'])):
1800        for iX in range(len(Layers['Layers'])):
1801            Parms.append('TransP;%d;%d'%(iY,iX))
1802            Parms.append('TransX;%d;%d'%(iY,iX))
1803            Parms.append('TransY;%d;%d'%(iY,iX))
1804            Parms.append('TransZ;%d;%d'%(iY,iX))
1805    return Parms
1806
1807def StackSim(Layers,ctrls,scale=0.,background={},limits=[],inst={},profile=[]):
1808    '''Simulate powder or selected area diffraction pattern from stacking faults using DIFFaX
1809   
1810    param: Layers dict: 'Laue':'-1','Cell':[False,1.,1.,1.,90.,90.,90,1.],
1811                        'Width':[[10.,10.],[False,False]],'Toler':0.01,'AtInfo':{},
1812                        'Layers':[],'Stacking':[],'Transitions':[]}
1813    param: ctrls string: controls string to be written on DIFFaX controls.dif file
1814    param: scale float: scale factor
1815    param: background dict: background parameters
1816    param: limits list: min/max 2-theta to be calculated
1817    param: inst dict: instrument parameters dictionary
1818    param: profile list: powder pattern data
1819   
1820    all updated in place   
1821    '''
1822    import atmdata
1823    path = sys.path
1824    for name in path:
1825        if 'bin' in name:
1826            DIFFaX = name+'/DIFFaX.exe'
1827            print ' Execute ',DIFFaX
1828            break
1829    # make form factor file that DIFFaX wants - atom types are GSASII style
1830    sf = open('data.sfc','w')
1831    sf.write('GSASII special form factor file for DIFFaX\n\n')
1832    atTypes = Layers['AtInfo'].keys()
1833    if 'H' not in atTypes:
1834        atTypes.insert(0,'H')
1835    for atType in atTypes:
1836        if atType == 'H': 
1837            blen = -.3741
1838        else:
1839            blen = Layers['AtInfo'][atType]['Isotopes']['Nat. Abund.']['SL'][0]
1840        Adat = atmdata.XrayFF[atType]
1841        text = '%4s'%(atType.ljust(4))
1842        for i in range(4):
1843            text += '%11.6f%11.6f'%(Adat['fa'][i],Adat['fb'][i])
1844        text += '%11.6f%11.6f'%(Adat['fc'],blen)
1845        text += '%3d\n'%(Adat['Z'])
1846        sf.write(text)
1847    sf.close()
1848    #make DIFFaX control.dif file - future use GUI to set some of these flags
1849    cf = open('control.dif','w')
1850    if ctrls == '0\n0\n3\n' or ctrls == '0\n1\n3\n': 
1851        x0 = profile[0]
1852        iBeg = np.searchsorted(x0,limits[0])
1853        iFin = np.searchsorted(x0,limits[1])
1854        if iFin-iBeg > 20000:
1855            iFin = iBeg+20000
1856        Dx = (x0[iFin]-x0[iBeg])/(iFin-iBeg)
1857        cf.write('GSASII-DIFFaX.dat\n'+ctrls)
1858        cf.write('%.6f %.6f %.6f\n1\n1\nend\n'%(x0[iBeg],x0[iFin],Dx))
1859    else:
1860        cf.write('GSASII-DIFFaX.dat\n'+ctrls)
1861        inst = {'Type':['XSC','XSC',]}
1862    cf.close()
1863    #make DIFFaX data file
1864    df = open('GSASII-DIFFaX.dat','w')
1865    df.write('INSTRUMENTAL\n')
1866    if 'X' in inst['Type'][0]:
1867        df.write('X-RAY\n')
1868    elif 'N' in inst['Type'][0]:
1869        df.write('NEUTRON\n')
1870    if ctrls == '0\n0\n3\n' or ctrls == '0\n1\n3\n': 
1871        df.write('%.4f\n'%(G2mth.getMeanWave(inst)))
1872        U = ateln2*inst['U'][1]/10000.
1873        V = ateln2*inst['V'][1]/10000.
1874        W = ateln2*inst['W'][1]/10000.
1875        HWHM = U*nptand(x0[iBeg:iFin]/2.)**2+V*nptand(x0[iBeg:iFin]/2.)+W
1876        HW = np.sqrt(np.mean(HWHM))
1877    #    df.write('PSEUDO-VOIGT 0.015 -0.0036 0.009 0.605 TRIM\n')
1878        if 'Mean' in Layers['selInst']:
1879            df.write('GAUSSIAN %.6f TRIM\n'%(HW))     #fast option - might not really matter
1880        elif 'Gaussian' in Layers['selInst']:
1881            df.write('GAUSSIAN %.6f %.6f %.6f TRIM\n'%(U,V,W))    #slow - make a GUI option?
1882        else:
1883            df.write('None\n')
1884    else:
1885        df.write('0.10\nNone\n')
1886    df.write('STRUCTURAL\n')
1887    a,b,c = Layers['Cell'][1:4]
1888    gam = Layers['Cell'][6]
1889    df.write('%.4f %.4f %.4f %.3f\n'%(a,b,c,gam))
1890    laue = Layers['Laue']
1891    if laue == '2/m(ab)':
1892        laue = '2/m(1)'
1893    elif laue == '2/m(c)':
1894        laue = '2/m(2)'
1895    if 'unknown' in Layers['Laue']:
1896        df.write('%s %.3f\n'%(laue,Layers['Toler']))
1897    else:   
1898        df.write('%s\n'%(laue))
1899    df.write('%d\n'%(len(Layers['Layers'])))
1900    if Layers['Width'][0][0] < 1. or Layers['Width'][0][1] < 1.:
1901        df.write('%.1f %.1f\n'%(Layers['Width'][0][0]*10000.,Layers['Width'][0][0]*10000.))    #mum to A
1902    layerNames = []
1903    for layer in Layers['Layers']:
1904        layerNames.append(layer['Name'])
1905    for il,layer in enumerate(Layers['Layers']):
1906        if layer['SameAs']:
1907            df.write('LAYER %d = %d\n'%(il+1,layerNames.index(layer['SameAs'])+1))
1908            continue
1909        df.write('LAYER %d\n'%(il+1))
1910        if '-1' in layer['Symm']:
1911            df.write('CENTROSYMMETRIC\n')
1912        else:
1913            df.write('NONE\n')
1914        for ia,atom in enumerate(layer['Atoms']):
1915            [name,atype,x,y,z,frac,Uiso] = atom
1916            if '-1' in layer['Symm'] and [x,y,z] == [0.,0.,0.]:
1917                frac /= 2.
1918            df.write('%4s %3d %.5f %.5f %.5f %.4f %.2f\n'%(atype.ljust(6),ia,x,y,z,78.9568*Uiso,frac))
1919    df.write('STACKING\n')
1920    df.write('%s\n'%(Layers['Stacking'][0]))
1921    if 'recursive' in Layers['Stacking'][0]:
1922        df.write('%s\n'%Layers['Stacking'][1])
1923    else:
1924        if 'list' in Layers['Stacking'][1]:
1925            Slen = len(Layers['Stacking'][2])
1926            iB = 0
1927            iF = 0
1928            while True:
1929                iF += 68
1930                if iF >= Slen:
1931                    break
1932                iF = min(iF,Slen)
1933                df.write('%s\n'%(Layers['Stacking'][2][iB:iF]))
1934                iB = iF
1935        else:
1936            df.write('%s\n'%Layers['Stacking'][1])   
1937    df.write('TRANSITIONS\n')
1938    for iY in range(len(Layers['Layers'])):
1939        sumPx = 0.
1940        for iX in range(len(Layers['Layers'])):
1941            p,dx,dy,dz = Layers['Transitions'][iY][iX][:4]
1942            p = round(p,3)
1943            df.write('%.3f %.5f %.5f %.5f\n'%(p,dx,dy,dz))
1944            sumPx += p
1945        if sumPx != 1.0:    #this has to be picky since DIFFaX is.
1946            print 'ERROR - Layer probabilities sum to ',sumPx,' DIFFaX will insist it = 1.0'
1947            df.close()
1948            os.remove('data.sfc')
1949            os.remove('control.dif')
1950            os.remove('GSASII-DIFFaX.dat')
1951            return       
1952    df.close()   
1953    time0 = time.time()
1954    try:
1955        subp.call(DIFFaX)
1956    except OSError:
1957        print ' DIFFax.exe is not available for this platform - under development'
1958    print ' DIFFaX time = %.2fs'%(time.time()-time0)
1959    if os.path.exists('GSASII-DIFFaX.spc'):
1960        Xpat = np.loadtxt('GSASII-DIFFaX.spc').T
1961        iFin = iBeg+Xpat.shape[1]
1962        bakType,backDict,backVary = SetBackgroundParms(background)
1963        backDict['Lam1'] = G2mth.getWave(inst)
1964        profile[4][iBeg:iFin] = getBackground('',backDict,bakType,inst['Type'][0],profile[0][iBeg:iFin])[0]   
1965        profile[3][iBeg:iFin] = Xpat[-1]*scale+profile[4][iBeg:iFin]
1966        if not np.any(profile[1]):                   #fill dummy data x,y,w,yc,yb,yd
1967            rv = st.poisson(profile[3][iBeg:iFin])
1968            profile[1][iBeg:iFin] = rv.rvs()
1969            Z = np.ones_like(profile[3][iBeg:iFin])
1970            Z[1::2] *= -1
1971            profile[1][iBeg:iFin] = profile[3][iBeg:iFin]+np.abs(profile[1][iBeg:iFin]-profile[3][iBeg:iFin])*Z
1972            profile[2][iBeg:iFin] = np.where(profile[1][iBeg:iFin]>0.,1./profile[1][iBeg:iFin],1.0)
1973        profile[5][iBeg:iFin] = profile[1][iBeg:iFin]-profile[3][iBeg:iFin]
1974    #cleanup files..
1975        os.remove('GSASII-DIFFaX.spc')
1976    elif os.path.exists('GSASII-DIFFaX.sadp'):
1977        Sadp = np.fromfile('GSASII-DIFFaX.sadp','>u2')
1978        Sadp = np.reshape(Sadp,(256,-1))
1979        Layers['Sadp']['Img'] = Sadp
1980        os.remove('GSASII-DIFFaX.sadp')
1981    os.remove('data.sfc')
1982    os.remove('control.dif')
1983    os.remove('GSASII-DIFFaX.dat')
1984   
1985def SetPWDRscan(inst,limits,profile):
1986   
1987    wave = G2mth.getMeanWave(inst)
1988    x0 = profile[0]
1989    iBeg = np.searchsorted(x0,limits[0])
1990    iFin = np.searchsorted(x0,limits[1])
1991    if iFin-iBeg > 20000:
1992        iFin = iBeg+20000
1993    Dx = (x0[iFin]-x0[iBeg])/(iFin-iBeg)
1994    pyx.pygetinst(wave,x0[iBeg],x0[iFin],Dx)
1995    return iFin-iBeg
1996       
1997def SetStackingSF(Layers,debug):
1998# Load scattering factors into DIFFaX arrays
1999    import atmdata
2000    atTypes = Layers['AtInfo'].keys()
2001    aTypes = []
2002    for atype in atTypes:
2003        aTypes.append('%4s'%(atype.ljust(4)))
2004    SFdat = []
2005    for atType in atTypes:
2006        if atType == 'H': 
2007            blen = -.3741
2008        else:
2009            blen = Layers['AtInfo'][atType]['Isotopes']['Nat. Abund.']['SL'][0]
2010        Adat = atmdata.XrayFF[atType]
2011        SF = np.zeros(9)
2012        SF[:8:2] = Adat['fa']
2013        SF[1:8:2] = Adat['fb']
2014        SF[8] = Adat['fc']
2015        SFdat.append(SF)
2016    SFdat = np.array(SFdat)
2017    pyx.pyloadscf(len(atTypes),aTypes,SFdat.T,debug)
2018   
2019def SetStackingClay(Layers,Type):
2020# Controls
2021    rand.seed()
2022    ranSeed = rand.randint(1,2**16-1)
2023    try:   
2024        laueId = ['-1','2/m(ab)','2/m(c)','mmm','-3','-3m','4/m','4/mmm',
2025            '6/m','6/mmm'].index(Layers['Laue'])+1
2026    except ValueError:  #for 'unknown'
2027        laueId = -1
2028    if 'SADP' in Type:
2029        planeId = ['h0l','0kl','hhl','h-hl'].index(Layers['Sadp']['Plane'])+1
2030        lmax = int(Layers['Sadp']['Lmax'])
2031    else:
2032        planeId = 0
2033        lmax = 0
2034# Sequences
2035    StkType = ['recursive','explicit'].index(Layers['Stacking'][0])
2036    try:
2037        StkParm = ['infinite','random','list'].index(Layers['Stacking'][1])
2038    except ValueError:
2039        StkParm = -1
2040    if StkParm == 2:    #list
2041        StkSeq = [int(val) for val in Layers['Stacking'][2].split()]
2042        Nstk = len(StkSeq)
2043    else:
2044        Nstk = 1
2045        StkSeq = [0,]
2046    if StkParm == -1:
2047        StkParm = int(Layers['Stacking'][1])
2048    Wdth = Layers['Width'][0]
2049    mult = 1
2050    controls = [laueId,planeId,lmax,mult,StkType,StkParm,ranSeed]
2051    LaueSym = Layers['Laue'].ljust(12)
2052    pyx.pygetclay(controls,LaueSym,Wdth,Nstk,StkSeq)
2053    return laueId,controls
2054   
2055def SetCellAtoms(Layers):
2056    Cell = Layers['Cell'][1:4]+Layers['Cell'][6:7]
2057# atoms in layers
2058    atTypes = Layers['AtInfo'].keys()
2059    AtomXOU = []
2060    AtomTp = []
2061    LayerSymm = []
2062    LayerNum = []
2063    layerNames = []
2064    Natm = 0
2065    Nuniq = 0
2066    for layer in Layers['Layers']:
2067        layerNames.append(layer['Name'])
2068    for il,layer in enumerate(Layers['Layers']):
2069        if layer['SameAs']:
2070            LayerNum.append(layerNames.index(layer['SameAs'])+1)
2071            continue
2072        else:
2073            LayerNum.append(il+1)
2074            Nuniq += 1
2075        if '-1' in layer['Symm']:
2076            LayerSymm.append(1)
2077        else:
2078            LayerSymm.append(0)
2079        for ia,atom in enumerate(layer['Atoms']):
2080            [name,atype,x,y,z,frac,Uiso] = atom
2081            Natm += 1
2082            AtomTp.append('%4s'%(atype.ljust(4)))
2083            Ta = atTypes.index(atype)+1
2084            AtomXOU.append([float(Nuniq),float(ia+1),float(Ta),x,y,z,frac,Uiso*78.9568])
2085    AtomXOU = np.array(AtomXOU)
2086    Nlayers = len(layerNames)
2087    pyx.pycellayer(Cell,Natm,AtomTp,AtomXOU.T,Nuniq,LayerSymm,Nlayers,LayerNum)
2088    return Nlayers
2089   
2090def SetStackingTrans(Layers,Nlayers):
2091# Transitions
2092    TransX = []
2093    TransP = []
2094    for Ytrans in Layers['Transitions']:
2095        TransP.append([trans[0] for trans in Ytrans])   #get just the numbers
2096        TransX.append([trans[1:4] for trans in Ytrans])   #get just the numbers
2097    TransP = np.array(TransP,dtype='float').T
2098    TransX = np.array(TransX,dtype='float')
2099#    GSASIIpath.IPyBreak()
2100    pyx.pygettrans(Nlayers,TransP,TransX)
2101   
2102def CalcStackingPWDR(Layers,scale,background,limits,inst,profile,debug):
2103# Scattering factors
2104    SetStackingSF(Layers,debug)
2105# Controls & sequences
2106    laueId,controls = SetStackingClay(Layers,'PWDR')
2107# cell & atoms
2108    Nlayers = SetCellAtoms(Layers)
2109    Volume = Layers['Cell'][7]   
2110# Transitions
2111    SetStackingTrans(Layers,Nlayers)
2112# PWDR scan
2113    Nsteps = SetPWDRscan(inst,limits,profile)
2114# result as Spec
2115    x0 = profile[0]
2116    profile[3] = np.zeros(len(profile[0]))
2117    profile[4] = np.zeros(len(profile[0]))
2118    profile[5] = np.zeros(len(profile[0]))
2119    iBeg = np.searchsorted(x0,limits[0])
2120    iFin = np.searchsorted(x0,limits[1])
2121    if iFin-iBeg > 20000:
2122        iFin = iBeg+20000
2123    Nspec = 20001       
2124    spec = np.zeros(Nspec,dtype='double')   
2125    time0 = time.time()
2126    pyx.pygetspc(controls,Nspec,spec)
2127    print ' GETSPC time = %.2fs'%(time.time()-time0)
2128    time0 = time.time()
2129    U = ateln2*inst['U'][1]/10000.
2130    V = ateln2*inst['V'][1]/10000.
2131    W = ateln2*inst['W'][1]/10000.
2132    HWHM = U*nptand(x0[iBeg:iFin]/2.)**2+V*nptand(x0[iBeg:iFin]/2.)+W
2133    HW = np.sqrt(np.mean(HWHM))
2134    BrdSpec = np.zeros(Nsteps)
2135    if 'Mean' in Layers['selInst']:
2136        pyx.pyprofile(U,V,W,HW,1,Nsteps,BrdSpec)
2137    elif 'Gaussian' in Layers['selInst']:
2138        pyx.pyprofile(U,V,W,HW,4,Nsteps,BrdSpec)
2139    else:
2140        BrdSpec = spec[:Nsteps]
2141    BrdSpec /= Volume
2142    iFin = iBeg+Nsteps
2143    bakType,backDict,backVary = SetBackgroundParms(background)
2144    backDict['Lam1'] = G2mth.getWave(inst)
2145    profile[4][iBeg:iFin] = getBackground('',backDict,bakType,inst['Type'][0],profile[0][iBeg:iFin])[0]   
2146    profile[3][iBeg:iFin] = BrdSpec*scale+profile[4][iBeg:iFin]
2147    if not np.any(profile[1]):                   #fill dummy data x,y,w,yc,yb,yd
2148        rv = st.poisson(profile[3][iBeg:iFin])
2149        profile[1][iBeg:iFin] = rv.rvs()
2150        Z = np.ones_like(profile[3][iBeg:iFin])
2151        Z[1::2] *= -1
2152        profile[1][iBeg:iFin] = profile[3][iBeg:iFin]+np.abs(profile[1][iBeg:iFin]-profile[3][iBeg:iFin])*Z
2153        profile[2][iBeg:iFin] = np.where(profile[1][iBeg:iFin]>0.,1./profile[1][iBeg:iFin],1.0)
2154    profile[5][iBeg:iFin] = profile[1][iBeg:iFin]-profile[3][iBeg:iFin]
2155    print ' Broadening time = %.2fs'%(time.time()-time0)
2156   
2157def CalcStackingSADP(Layers,debug):
2158   
2159# Scattering factors
2160    SetStackingSF(Layers,debug)
2161# Controls & sequences
2162    laueId,controls = SetStackingClay(Layers,'SADP')
2163# cell & atoms
2164    Nlayers = SetCellAtoms(Layers)   
2165# Transitions
2166    SetStackingTrans(Layers,Nlayers)
2167# result as Sadp
2168    mirror = laueId in [-1,2,3,7,8,9,10]
2169    Nspec = 20001       
2170    spec = np.zeros(Nspec,dtype='double')   
2171    time0 = time.time()
2172    hkLim,Incr,Nblk = pyx.pygetsadp(controls,Nspec,spec)
2173    Sapd = np.zeros((256,256))
2174    iB = 0
2175    for i in range(hkLim):
2176        iF = iB+Nblk
2177        p1 = 127+int(i*Incr)
2178        p2 = 128-int(i*Incr)
2179        if Nblk == 128:
2180            if i:
2181                Sapd[128:,p1] = spec[iB:iF]
2182                Sapd[:128,p1] = spec[iF:iB:-1]
2183            Sapd[128:,p2] = spec[iB:iF]
2184            Sapd[:128,p2] = spec[iF:iB:-1]
2185        else:
2186            if i:
2187                Sapd[:,p1] = spec[iB:iF]
2188            Sapd[:,p2] = spec[iB:iF]
2189        iB += Nblk
2190    Layers['Sadp']['Img'] = Sapd
2191    print ' GETSAD time = %.2fs'%(time.time()-time0)
2192#    GSASIIpath.IPyBreak()
2193   
2194#testing data
2195NeedTestData = True
2196def TestData():
2197    'needs a doc string'
2198#    global NeedTestData
2199    NeedTestData = False
2200    global bakType
2201    bakType = 'chebyschev'
2202    global xdata
2203    xdata = np.linspace(4.0,40.0,36000)
2204    global parmDict0
2205    parmDict0 = {
2206        'pos0':5.6964,'int0':8835.8,'sig0':1.0,'gam0':1.0,
2207        'pos1':11.4074,'int1':3922.3,'sig1':1.0,'gam1':1.0,
2208        'pos2':20.6426,'int2':1573.7,'sig2':1.0,'gam2':1.0,
2209        'pos3':26.9568,'int3':925.1,'sig3':1.0,'gam3':1.0,
2210        'U':1.163,'V':-0.605,'W':0.093,'X':0.0,'Y':2.183,'SH/L':0.002,
2211        'Back0':5.384,'Back1':-0.015,'Back2':.004,
2212        }
2213    global parmDict1
2214    parmDict1 = {
2215        'pos0':13.4924,'int0':48697.6,'sig0':1.0,'gam0':1.0,
2216        'pos1':23.4360,'int1':43685.5,'sig1':1.0,'gam1':1.0,
2217        'pos2':27.1152,'int2':123712.6,'sig2':1.0,'gam2':1.0,
2218        'pos3':33.7196,'int3':65349.4,'sig3':1.0,'gam3':1.0,
2219        'pos4':36.1119,'int4':115829.8,'sig4':1.0,'gam4':1.0,
2220        'pos5':39.0122,'int5':6916.9,'sig5':1.0,'gam5':1.0,
2221        'U':22.75,'V':-17.596,'W':10.594,'X':1.577,'Y':5.778,'SH/L':0.002,
2222        'Back0':36.897,'Back1':-0.508,'Back2':.006,
2223        'Lam1':1.540500,'Lam2':1.544300,'I(L2)/I(L1)':0.5,
2224        }
2225    global parmDict2
2226    parmDict2 = {
2227        'pos0':5.7,'int0':1000.0,'sig0':0.5,'gam0':0.5,
2228        'U':2.,'V':-2.,'W':5.,'X':0.5,'Y':0.5,'SH/L':0.02,
2229        'Back0':5.,'Back1':-0.02,'Back2':.004,
2230#        'Lam1':1.540500,'Lam2':1.544300,'I(L2)/I(L1)':0.5,
2231        }
2232    global varyList
2233    varyList = []
2234
2235def test0():
2236    if NeedTestData: TestData()
2237    msg = 'test '
2238    gplot = plotter.add('FCJ-Voigt, 11BM').gca()
2239    gplot.plot(xdata,getBackground('',parmDict0,bakType,'PXC',xdata)[0])   
2240    gplot.plot(xdata,getPeakProfile(parmDict0,xdata,varyList,bakType))
2241    fplot = plotter.add('FCJ-Voigt, Ka1+2').gca()
2242    fplot.plot(xdata,getBackground('',parmDict1,bakType,'PXC',xdata)[0])   
2243    fplot.plot(xdata,getPeakProfile(parmDict1,xdata,varyList,bakType))
2244   
2245def test1():
2246    if NeedTestData: TestData()
2247    time0 = time.time()
2248    for i in range(100):
2249        y = getPeakProfile(parmDict1,xdata,varyList,bakType)
2250    print '100+6*Ka1-2 peaks=1200 peaks',time.time()-time0
2251   
2252def test2(name,delt):
2253    if NeedTestData: TestData()
2254    varyList = [name,]
2255    xdata = np.linspace(5.6,5.8,400)
2256    hplot = plotter.add('derivatives test for '+name).gca()
2257    hplot.plot(xdata,getPeakProfileDerv(parmDict2,xdata,varyList,bakType)[0])
2258    y0 = getPeakProfile(parmDict2,xdata,varyList,bakType)
2259    parmDict2[name] += delt
2260    y1 = getPeakProfile(parmDict2,xdata,varyList,bakType)
2261    hplot.plot(xdata,(y1-y0)/delt,'r+')
2262   
2263def test3(name,delt):
2264    if NeedTestData: TestData()
2265    names = ['pos','sig','gam','shl']
2266    idx = names.index(name)
2267    myDict = {'pos':parmDict2['pos0'],'sig':parmDict2['sig0'],'gam':parmDict2['gam0'],'shl':parmDict2['SH/L']}
2268    xdata = np.linspace(5.6,5.8,800)
2269    dx = xdata[1]-xdata[0]
2270    hplot = plotter.add('derivatives test for '+name).gca()
2271    hplot.plot(xdata,100.*dx*getdFCJVoigt3(myDict['pos'],myDict['sig'],myDict['gam'],myDict['shl'],xdata)[idx+1])
2272    y0 = getFCJVoigt3(myDict['pos'],myDict['sig'],myDict['gam'],myDict['shl'],xdata)
2273    myDict[name] += delt
2274    y1 = getFCJVoigt3(myDict['pos'],myDict['sig'],myDict['gam'],myDict['shl'],xdata)
2275    hplot.plot(xdata,(y1-y0)/delt,'r+')
2276
2277if __name__ == '__main__':
2278    import GSASIItestplot as plot
2279    global plotter
2280    plotter = plot.PlotNotebook()
2281#    test0()
2282#    for name in ['int0','pos0','sig0','gam0','U','V','W','X','Y','SH/L','I(L2)/I(L1)']:
2283    for name,shft in [['int0',0.1],['pos0',0.0001],['sig0',0.01],['gam0',0.00001],
2284        ['U',0.1],['V',0.01],['W',0.01],['X',0.0001],['Y',0.0001],['SH/L',0.00005]]:
2285        test2(name,shft)
2286    for name,shft in [['pos',0.0001],['sig',0.01],['gam',0.0001],['shl',0.00005]]:
2287        test3(name,shft)
2288    print "OK"
2289    plotter.StartEventLoop()
Note: See TracBrowser for help on using the repository browser.