source: trunk/GSASIIpwd.py @ 2334

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

fix mask array issue with si.interp1d in scipy 17.0+; use ma.getdata to use nonmasked data
fix problem with update of anomalous scattering density on copy substances for small angle stuff
fix this in AddSubstance? also.

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