source: trunk/GSASIImath.py @ 1244

Last change on this file since 1244 was 1244, checked in by vondreele, 9 years ago

fix geometric correction in integrate - too many 1/cos(2-theta)
plot of size distribution from SASD
MaxEnt? size distribution in operation (some tuning/errors)

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 105.0 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASIImath - major mathematics routines
3########### SVN repository information ###################
4# $Date: 2014-03-11 21:49:15 +0000 (Tue, 11 Mar 2014) $
5# $Author: vondreele $
6# $Revision: 1244 $
7# $URL: trunk/GSASIImath.py $
8# $Id: GSASIImath.py 1244 2014-03-11 21:49:15Z vondreele $
9########### SVN repository information ###################
10'''
11*GSASIImath: computation module*
12================================
13
14Routines for least-squares minimization and other stuff
15
16'''
17import sys
18import os
19import os.path as ospath
20import random as rn
21import numpy as np
22import numpy.linalg as nl
23import numpy.ma as ma
24import cPickle
25import time
26import math
27import copy
28import GSASIIpath
29GSASIIpath.SetVersionNumber("$Revision: 1244 $")
30import GSASIIElem as G2el
31import GSASIIlattice as G2lat
32import GSASIIspc as G2spc
33import GSASIIpwd as G2pwd
34import numpy.fft as fft
35import pypowder as pyd
36
37sind = lambda x: np.sin(x*np.pi/180.)
38cosd = lambda x: np.cos(x*np.pi/180.)
39tand = lambda x: np.tan(x*np.pi/180.)
40asind = lambda x: 180.*np.arcsin(x)/np.pi
41acosd = lambda x: 180.*np.arccos(x)/np.pi
42atand = lambda x: 180.*np.arctan(x)/np.pi
43atan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
44
45def HessianLSQ(func,x0,Hess,args=(),ftol=1.49012e-8,xtol=1.49012e-8, maxcyc=0,Print=False):
46   
47    """
48    Minimize the sum of squares of a function (:math:`f`) evaluated on a series of
49    values (y): :math:`\sum_{y=0}^{N_{obs}} f(y,{args})`
50   
51    ::
52   
53                    Nobs
54        x = arg min(sum(func(y)**2,axis=0))
55                    y=0
56
57    :param function func: callable method or function
58        should take at least one (possibly length N vector) argument and
59        returns M floating point numbers.
60    :param np.ndarray x0: The starting estimate for the minimization of length N
61    :param function Hess: callable method or function
62        A required function or method to compute the weighted vector and Hessian for func.
63        It must be a symmetric NxN array
64    :param tuple args: Any extra arguments to func are placed in this tuple.
65    :param float ftol: Relative error desired in the sum of squares.
66    :param float xtol: Relative error desired in the approximate solution.
67    :param int maxcyc: The maximum number of cycles of refinement to execute, if -1 refine
68        until other limits are met (ftol, xtol)
69    :param bool Print: True for printing results (residuals & times) by cycle
70
71    :returns: (x,cov_x,infodict) where
72
73      * x : ndarray
74        The solution (or the result of the last iteration for an unsuccessful
75        call).
76      * cov_x : ndarray
77        Uses the fjac and ipvt optional outputs to construct an
78        estimate of the jacobian around the solution.  ``None`` if a
79        singular matrix encountered (indicates very flat curvature in
80        some direction).  This matrix must be multiplied by the
81        residual standard deviation to get the covariance of the
82        parameter estimates -- see curve_fit.
83      * infodict : dict
84        a dictionary of optional outputs with the keys:
85
86         * 'fvec' : the function evaluated at the output
87         * 'num cyc':
88         * 'nfev':
89         * 'lamMax':
90         * 'psing':
91           
92    """
93               
94    x0 = np.array(x0, ndmin=1)      #might be redundant?
95    n = len(x0)
96    if type(args) != type(()):
97        args = (args,)
98       
99    icycle = 0
100    One = np.ones((n,n))
101    lam = 0.001
102    lamMax = lam
103    nfev = 0
104    if Print:
105        print ' Hessian refinement on %d variables:'%(n)
106    Lam = np.zeros((n,n))
107    while icycle < maxcyc:
108        time0 = time.time()
109        M = func(x0,*args)
110        nfev += 1
111        chisq0 = np.sum(M**2)
112        Yvec,Amat = Hess(x0,*args)
113        Adiag = np.sqrt(np.diag(Amat))
114        psing = np.where(np.abs(Adiag) < 1.e-14,True,False)
115        if np.any(psing):                #hard singularity in matrix
116            return [x0,None,{'num cyc':icycle,'fvec':M,'nfev':nfev,'lamMax':lamMax,'psing':psing}]
117        Anorm = np.outer(Adiag,Adiag)
118        Yvec /= Adiag
119        Amat /= Anorm
120        while True:
121            Lam = np.eye(Amat.shape[0])*lam
122            Amatlam = Amat*(One+Lam)
123            try:
124                Xvec = nl.solve(Amatlam,Yvec)
125            except nl.LinAlgError:
126                print 'ouch #1'
127                psing = list(np.where(np.diag(nl.qr(Amatlam)[1]) < 1.e-14)[0])
128                return [x0,None,{'num cyc':icycle,'fvec':M,'nfev':nfev,'lamMax':lamMax,'psing':psing}]
129            Xvec /= Adiag
130            M2 = func(x0+Xvec,*args)
131            nfev += 1
132            chisq1 = np.sum(M2**2)
133            if chisq1 > chisq0:
134                lam *= 10.
135            else:
136                x0 += Xvec
137                lam /= 10.
138                break
139            if lam > 10.e5:
140                print 'ouch #3 chisq1 ',chisq1,' stuck > chisq0 ',chisq0
141                break
142        lamMax = max(lamMax,lam)
143        if Print:
144            print ' Cycle: %d, Time: %.2fs, Chi**2: %.3g, Lambda: %.3g'%(icycle,time.time()-time0,chisq1,lam)
145        if (chisq0-chisq1)/chisq0 < ftol:
146            break
147        icycle += 1
148    else:       #after last cycle or if zero cycles
149        M = func(x0,*args)
150        nfev += 1
151        Yvec,Amat = Hess(x0,*args)
152        Adiag = np.sqrt(np.diag(Amat))
153        Anorm = np.outer(Adiag,Adiag)
154        Amatlam = Amat*(One+Lam)/Anorm              #scale Amat to Marquardt array       
155    try:
156        Bmat = nl.inv(Amatlam)*(One+Lam)/Anorm      #rescale Bmat to Marquardt array
157        return [x0,Bmat,{'num cyc':icycle,'fvec':M,'nfev':nfev,'lamMax':lamMax,'psing':[]}]
158    except nl.LinAlgError:
159        print 'ouch #2 linear algebra error in LS'
160        psing = []
161        if maxcyc:
162            psing = list(np.where(np.diag(nl.qr(Amat)[1]) < 1.e-14)[0])
163        return [x0,None,{'num cyc':icycle,'fvec':M,'nfev':nfev,'lamMax':lamMax,'psing':psing}]
164   
165def getVCov(varyNames,varyList,covMatrix):
166    '''obtain variance-covariance terms for a set of variables. NB: the varyList
167    and covMatrix were saved by the last least squares refinement so they must match.
168   
169    :param list varyNames: variable names to find v-cov matric for
170    :param list varyList: full list of all variables in v-cov matrix
171    :param nparray covMatrix: full variance-covariance matrix from the last
172     least squares refinement
173   
174    :returns: nparray vcov: variance-covariance matrix for the variables given
175     in varyNames
176   
177    '''
178    vcov = np.zeros((len(varyNames),len(varyNames)))
179    for i1,name1 in enumerate(varyNames):
180        for i2,name2 in enumerate(varyNames):
181            try:
182                vcov[i1][i2] = covMatrix[varyList.index(name1)][varyList.index(name2)]
183            except ValueError:
184                vcov[i1][i2] = 0.0
185    return vcov
186
187def FindAtomIndexByIDs(atomData,IDs,Draw=True):
188    '''finds the set of atom array indices for a list of atom IDs. Will search
189    either the Atom table or the drawAtom table.
190   
191    :param list atomData: Atom or drawAtom table containting coordinates, etc.
192    :param list IDs: atom IDs to be found
193    :param bool Draw: True if drawAtom table to be searched; False if Atom table
194      is searched
195   
196    :returns: list indx: atom (or drawAtom) indices
197   
198    '''
199    indx = []
200    for i,atom in enumerate(atomData):
201        if Draw and atom[-3] in IDs:
202            indx.append(i)
203        elif atom[-1] in IDs:
204            indx.append(i)
205    return indx
206
207def FillAtomLookUp(atomData):
208    '''create a dictionary of atom indexes with atom IDs as keys
209   
210    :param list atomData: Atom table to be used
211   
212    :returns: dict atomLookUp: dictionary of atom indexes with atom IDs as keys
213   
214    '''
215    atomLookUp = {}
216    for iatm,atom in enumerate(atomData):
217        atomLookUp[atom[-1]] = iatm
218    return atomLookUp
219
220def GetAtomsById(atomData,atomLookUp,IdList):
221    '''gets a list of atoms from Atom table that match a set of atom IDs
222   
223    :param list atomData: Atom table to be used
224    :param dict atomLookUp: dictionary of atom indexes with atom IDs as keys
225    :param list IdList: atom IDs to be found
226   
227    :returns: list atoms: list of atoms found
228   
229    '''
230    atoms = []
231    for id in IdList:
232        atoms.append(atomData[atomLookUp[id]])
233    return atoms
234   
235def GetAtomItemsById(atomData,atomLookUp,IdList,itemLoc,numItems=1):
236    '''gets atom parameters for atoms using atom IDs
237   
238    :param list atomData: Atom table to be used
239    :param dict atomLookUp: dictionary of atom indexes with atom IDs as keys
240    :param list IdList: atom IDs to be found
241    :param int itemLoc: pointer to desired 1st item in an atom table entry
242    :param int numItems: number of items to be retrieved
243   
244    :returns: type name: description
245   
246    '''
247    Items = []
248    if not isinstance(IdList,list):
249        IdList = [IdList,]
250    for id in IdList:
251        if numItems == 1:
252            Items.append(atomData[atomLookUp[id]][itemLoc])
253        else:
254            Items.append(atomData[atomLookUp[id]][itemLoc:itemLoc+numItems])
255    return Items
256   
257def GetAtomCoordsByID(pId,parmDict,AtLookup,indx):
258    '''default doc string
259   
260    :param type name: description
261   
262    :returns: type name: description
263   
264    '''
265    pfx = [str(pId)+'::A'+i+':' for i in ['x','y','z']]
266    dpfx = [str(pId)+'::dA'+i+':' for i in ['x','y','z']]
267    XYZ = []
268    for ind in indx:
269        names = [pfx[i]+str(AtLookup[ind]) for i in range(3)]
270        dnames = [dpfx[i]+str(AtLookup[ind]) for i in range(3)]
271        XYZ.append([parmDict[name]+parmDict[dname] for name,dname in zip(names,dnames)])
272    return XYZ
273
274def AtomUij2TLS(atomData,atPtrs,Amat,Bmat,rbObj):   #unfinished & not used
275    '''default doc string
276   
277    :param type name: description
278   
279    :returns: type name: description
280   
281    '''
282    for atom in atomData:
283        XYZ = np.inner(Amat,atom[cx:cx+3])
284        if atom[cia] == 'A':
285            UIJ = atom[cia+2:cia+8]
286               
287def TLS2Uij(xyz,g,Amat,rbObj):    #not used anywhere, but could be?
288    '''default doc string
289   
290    :param type name: description
291   
292    :returns: type name: description
293   
294    '''
295    TLStype,TLS = rbObj['ThermalMotion'][:2]
296    Tmat = np.zeros((3,3))
297    Lmat = np.zeros((3,3))
298    Smat = np.zeros((3,3))
299    gvec = np.sqrt(np.array([g[0][0]**2,g[1][1]**2,g[2][2]**2,
300        g[0][0]*g[1][1],g[0][0]*g[2][2],g[1][1]*g[2][2]]))
301    if 'T' in TLStype:
302        Tmat = G2lat.U6toUij(TLS[:6])
303    if 'L' in TLStype:
304        Lmat = G2lat.U6toUij(TLS[6:12])
305    if 'S' in TLStype:
306        Smat = np.array([[TLS[18],TLS[12],TLS[13]],[TLS[14],TLS[19],TLS[15]],[TLS[16],TLS[17],0] ])
307    XYZ = np.inner(Amat,xyz)
308    Axyz = np.array([[ 0,XYZ[2],-XYZ[1]], [-XYZ[2],0,XYZ[0]], [XYZ[1],-XYZ[0],0]] )
309    Umat = Tmat+np.inner(Axyz,Smat)+np.inner(Smat.T,Axyz.T)+np.inner(np.inner(Axyz,Lmat),Axyz.T)
310    beta = np.inner(np.inner(g,Umat),g)
311    return G2lat.UijtoU6(beta)*gvec   
312       
313def AtomTLS2UIJ(atomData,atPtrs,Amat,rbObj):    #not used anywhere, but could be?
314    '''default doc string
315   
316    :param type name: description
317   
318    :returns: type name: description
319   
320    '''
321    cx,ct,cs,cia = atPtrs
322    TLStype,TLS = rbObj['ThermalMotion'][:2]
323    Tmat = np.zeros((3,3))
324    Lmat = np.zeros((3,3))
325    Smat = np.zeros((3,3))
326    G,g = G2lat.A2Gmat(Amat)
327    gvec = 1./np.sqrt(np.array([g[0][0],g[1][1],g[2][2],g[0][1],g[0][2],g[1][2]]))
328    if 'T' in TLStype:
329        Tmat = G2lat.U6toUij(TLS[:6])
330    if 'L' in TLStype:
331        Lmat = G2lat.U6toUij(TLS[6:12])
332    if 'S' in TLStype:
333        Smat = np.array([ [TLS[18],TLS[12],TLS[13]], [TLS[14],TLS[19],TLS[15]], [TLS[16],TLS[17],0] ])
334    for atom in atomData:
335        XYZ = np.inner(Amat,atom[cx:cx+3])
336        Axyz = np.array([ 0,XYZ[2],-XYZ[1], -XYZ[2],0,XYZ[0], XYZ[1],-XYZ[0],0],ndmin=2 )
337        if 'U' in TSLtype:
338            atom[cia+1] = TLS[0]
339            atom[cia] = 'I'
340        else:
341            atom[cia] = 'A'
342            Umat = Tmat+np.inner(Axyz,Smat)+np.inner(Smat.T,Axyz.T)+np.inner(np.inner(Axyz,Lmat),Axyz.T)
343            beta = np.inner(np.inner(g,Umat),g)
344            atom[cia+2:cia+8] = G2spc.U2Uij(beta/gvec)
345
346def GetXYZDist(xyz,XYZ,Amat):
347    '''gets distance from position xyz to all XYZ, xyz & XYZ are np.array
348        and are in crystal coordinates; Amat is crystal to Cart matrix
349   
350    :param type name: description
351   
352    :returns: type name: description
353   
354    '''
355    return np.sqrt(np.sum(np.inner(Amat,XYZ-xyz)**2,axis=0))
356
357def getAtomXYZ(atoms,cx):
358    '''default doc string
359   
360    :param type name: description
361   
362    :returns: type name: description
363   
364    '''
365    XYZ = []
366    for atom in atoms:
367        XYZ.append(atom[cx:cx+3])
368    return np.array(XYZ)
369
370def RotateRBXYZ(Bmat,Cart,oriQ):
371    '''rotate & transform cartesian coordinates to crystallographic ones
372    no translation applied. To be used for numerical derivatives
373   
374    :param type name: description
375   
376    :returns: type name: description
377   
378    '''
379    ''' returns crystal coordinates for atoms described by RBObj
380    '''
381    XYZ = np.zeros_like(Cart)
382    for i,xyz in enumerate(Cart):
383        XYZ[i] = np.inner(Bmat,prodQVQ(oriQ,xyz))
384    return XYZ
385
386def UpdateRBXYZ(Bmat,RBObj,RBData,RBType):
387    '''default doc string
388   
389    :param type name: description
390   
391    :returns: type name: description
392   
393    '''
394    ''' returns crystal coordinates for atoms described by RBObj
395    '''
396    RBRes = RBData[RBType][RBObj['RBId']]
397    if RBType == 'Vector':
398        vecs = RBRes['rbVect']
399        mags = RBRes['VectMag']
400        Cart = np.zeros_like(vecs[0])
401        for vec,mag in zip(vecs,mags):
402            Cart += vec*mag
403    elif RBType == 'Residue':
404        Cart = np.array(RBRes['rbXYZ'])
405        for tor,seq in zip(RBObj['Torsions'],RBRes['rbSeq']):
406            QuatA = AVdeg2Q(tor[0],Cart[seq[0]]-Cart[seq[1]])
407            Cart[seq[3]] = prodQVQ(QuatA,(Cart[seq[3]]-Cart[seq[1]]))+Cart[seq[1]]
408    XYZ = np.zeros_like(Cart)
409    for i,xyz in enumerate(Cart):
410        XYZ[i] = np.inner(Bmat,prodQVQ(RBObj['Orient'][0],xyz))+RBObj['Orig'][0]
411    return XYZ,Cart
412
413def UpdateMCSAxyz(Bmat,MCSA):
414    '''default doc string
415   
416    :param type name: description
417   
418    :returns: type name: description
419   
420    '''
421    xyz = []
422    atTypes = []
423    iatm = 0
424    for model in MCSA['Models'][1:]:        #skip the MD model
425        if model['Type'] == 'Atom':
426            xyz.append(model['Pos'][0])
427            atTypes.append(model['atType'])
428            iatm += 1
429        else:
430            RBRes = MCSA['rbData'][model['Type']][model['RBId']]
431            Pos = np.array(model['Pos'][0])
432            Ori = np.array(model['Ori'][0])
433            Qori = AVdeg2Q(Ori[0],Ori[1:])
434            if model['Type'] == 'Vector':
435                vecs = RBRes['rbVect']
436                mags = RBRes['VectMag']
437                Cart = np.zeros_like(vecs[0])
438                for vec,mag in zip(vecs,mags):
439                    Cart += vec*mag
440            elif model['Type'] == 'Residue':
441                Cart = np.array(RBRes['rbXYZ'])
442                for itor,seq in enumerate(RBRes['rbSeq']):
443                    QuatA = AVdeg2Q(model['Tor'][0][itor],Cart[seq[0]]-Cart[seq[1]])
444                    Cart[seq[3]] = prodQVQ(QuatA,(Cart[seq[3]]-Cart[seq[1]]))+Cart[seq[1]]
445            if model['MolCent'][1]:
446                Cart -= model['MolCent'][0]
447            for i,x in enumerate(Cart):
448                xyz.append(np.inner(Bmat,prodQVQ(Qori,x))+Pos)
449                atType = RBRes['rbTypes'][i]
450                atTypes.append(atType)
451                iatm += 1
452    return np.array(xyz),atTypes
453   
454def SetMolCent(model,RBData):
455    '''default doc string
456   
457    :param type name: description
458   
459    :returns: type name: description
460   
461    '''
462    rideList = []
463    RBRes = RBData[model['Type']][model['RBId']]
464    if model['Type'] == 'Vector':
465        vecs = RBRes['rbVect']
466        mags = RBRes['VectMag']
467        Cart = np.zeros_like(vecs[0])
468        for vec,mag in zip(vecs,mags):
469            Cart += vec*mag
470    elif model['Type'] == 'Residue':
471        Cart = np.array(RBRes['rbXYZ'])
472        for seq in RBRes['rbSeq']:
473            rideList += seq[3]
474    centList = set(range(len(Cart)))-set(rideList)
475    cent = np.zeros(3)
476    for i in centList:
477        cent += Cart[i]
478    model['MolCent'][0] = cent/len(centList)   
479   
480def UpdateRBUIJ(Bmat,Cart,RBObj):
481    '''default doc string
482   
483    :param type name: description
484   
485    :returns: type name: description
486   
487    '''
488    ''' returns atom I/A, Uiso or UIJ for atoms at XYZ as described by RBObj
489    '''
490    TLStype,TLS = RBObj['ThermalMotion'][:2]
491    T = np.zeros(6)
492    L = np.zeros(6)
493    S = np.zeros(8)
494    if 'T' in TLStype:
495        T = TLS[:6]
496    if 'L' in TLStype:
497        L = np.array(TLS[6:12])*(np.pi/180.)**2
498    if 'S' in TLStype:
499        S = np.array(TLS[12:])*(np.pi/180.)
500    g = nl.inv(np.inner(Bmat,Bmat))
501    gvec = np.sqrt(np.array([g[0][0]**2,g[1][1]**2,g[2][2]**2,
502        g[0][0]*g[1][1],g[0][0]*g[2][2],g[1][1]*g[2][2]]))
503    Uout = []
504    Q = RBObj['Orient'][0]
505    for X in Cart:
506        X = prodQVQ(Q,X)
507        if 'U' in TLStype:
508            Uout.append(['I',TLS[0],0,0,0,0,0,0])
509        elif not 'N' in TLStype:
510            U = [0,0,0,0,0,0]
511            U[0] = T[0]+L[1]*X[2]**2+L[2]*X[1]**2-2.0*L[5]*X[1]*X[2]+2.0*(S[2]*X[2]-S[4]*X[1])
512            U[1] = T[1]+L[0]*X[2]**2+L[2]*X[0]**2-2.0*L[4]*X[0]*X[2]+2.0*(S[5]*X[0]-S[0]*X[2])
513            U[2] = T[2]+L[1]*X[0]**2+L[0]*X[1]**2-2.0*L[3]*X[1]*X[0]+2.0*(S[1]*X[1]-S[3]*X[0])
514            U[3] = T[3]+L[4]*X[1]*X[2]+L[5]*X[0]*X[2]-L[3]*X[2]**2-L[2]*X[0]*X[1]+  \
515                S[4]*X[0]-S[5]*X[1]-(S[6]+S[7])*X[2]
516            U[4] = T[4]+L[3]*X[1]*X[2]+L[5]*X[0]*X[1]-L[4]*X[1]**2-L[1]*X[0]*X[2]+  \
517                S[3]*X[2]-S[2]*X[0]+S[6]*X[1]
518            U[5] = T[5]+L[3]*X[0]*X[2]+L[4]*X[0]*X[1]-L[5]*X[0]**2-L[0]*X[2]*X[1]+  \
519                S[0]*X[1]-S[1]*X[2]+S[7]*X[0]
520            Umat = G2lat.U6toUij(U)
521            beta = np.inner(np.inner(Bmat.T,Umat),Bmat)
522            Uout.append(['A',0.0,]+list(G2lat.UijtoU6(beta)*gvec))
523        else:
524            Uout.append(['N',])
525    return Uout
526   
527def GetSHCoeff(pId,parmDict,SHkeys):
528    '''default doc string
529   
530    :param type name: description
531   
532    :returns: type name: description
533   
534    '''
535    SHCoeff = {}
536    for shkey in SHkeys:
537        shname = str(pId)+'::'+shkey
538        SHCoeff[shkey] = parmDict[shname]
539    return SHCoeff
540       
541def getMass(generalData):
542    '''Computes mass of unit cell contents
543   
544    :param dict generalData: The General dictionary in Phase
545   
546    :returns: float mass: Crystal unit cell mass in AMU.
547   
548    '''
549    mass = 0.
550    for i,elem in enumerate(generalData['AtomTypes']):
551        mass += generalData['NoAtoms'][elem]*generalData['AtomMass'][i]
552    return mass   
553
554def getDensity(generalData):
555    '''calculate crystal structure density
556   
557    :param dict generalData: The General dictionary in Phase
558   
559    :returns: float density: crystal density in gm/cm^3
560   
561    '''
562    mass = getMass(generalData)
563    Volume = generalData['Cell'][7]
564    density = mass/(0.6022137*Volume)
565    return density,Volume/mass
566   
567def getWave(Parms):
568    '''returns wavelength from Instrument parameters dictionary
569   
570    :param dict Parms: Instrument parameters;
571        must contain:
572        Lam: single wavelength
573        or
574        Lam1: Ka1 radiation wavelength
575   
576    :returns: float wave: wavelength
577   
578    '''
579    try:
580        return Parms['Lam'][1]
581    except KeyError:
582        return Parms['Lam1'][1]
583       
584def El2Mass(Elements):
585    '''compute molecular weight from Elements
586   
587    :param dict Elements: elements in molecular formula;
588        each must contain
589        Num: number of atoms in formula
590        Mass: at. wt.
591   
592    :returns: float mass: molecular weight.
593   
594    '''
595    mass = 0
596    for El in Elements:
597        mass += Elements[El]['Num']*Elements[El]['Mass']
598    return mass
599       
600def Den2Vol(Elements,density):
601    '''converts density to molecular volume
602   
603    :param dict Elements: elements in molecular formula;
604        each must contain
605        Num: number of atoms in formula
606        Mass: at. wt.
607    :param float density: material density in gm/cm^3
608   
609    :returns: float volume: molecular volume in A^3
610   
611    '''
612    return El2Mass(Elements)/(density*0.6022137)
613   
614def Vol2Den(Elements,volume):
615    '''converts volume to density
616   
617    :param dict Elements: elements in molecular formula;
618        each must contain
619        Num: number of atoms in formula
620        Mass: at. wt.
621    :param float volume: molecular volume in A^3
622   
623    :returns: float density: material density in gm/cm^3
624   
625    '''
626    return El2Mass(Elements)/(volume*0.6022137)
627   
628def El2EstVol(Elements):
629    '''Estimate volume from molecular formula; assumes atom volume = 10A^3
630   
631    :param dict Elements: elements in molecular formula;
632        each must contain
633        Num: number of atoms in formula
634   
635    :returns: float volume: estimate of molecular volume in A^3
636   
637    '''
638    vol = 0
639    for El in Elements:
640        vol += 10.*Elements[El]['Num']
641    return vol
642   
643def XScattDen(Elements,vol,wave=0.):
644    '''Estimate X-ray scattering density from molecular formula & volume;
645    ignores valence, but includes anomalous effects
646   
647    :param dict Elements: elements in molecular formula;
648        each element must contain
649        Num: number of atoms in formula
650        Z: atomic number
651    :param float vol: molecular volume in A^3
652    :param float wave: optional wavelength in A
653   
654    :returns: float rho: scattering density in 10^10cm^-2;
655        if wave > 0 the includes f' contribution
656    :returns: float mu: if wave>0 absorption coeff in cm^-1 ; otherwise 0
657   
658    '''
659    rho = 0
660    mu = 0
661    if wave:
662        Xanom = XAnomAbs(Elements,wave)
663    for El in Elements:
664        f0 = Elements[El]['Z']
665        if wave:
666            f0 += Xanom[El][0]
667            mu += Xanom[El][2]*Elements[El]['Num']
668        rho += Elements[El]['Num']*f0
669    return 28.179*rho/vol,0.1*mu/vol
670   
671def wavekE(wavekE):
672    '''Convert wavelength to energy & vise versa
673   
674    :param float waveKe:wavelength in A or energy in kE
675   
676    :returns float waveKe:the other one
677   
678    '''
679    return 12.397639/wavekE
680       
681def XAnomAbs(Elements,wave):
682    kE = wavekE(wave)
683    Xanom = {}
684    for El in Elements:
685        Orbs = G2el.GetXsectionCoeff(El)
686        Xanom[El] = G2el.FPcalc(Orbs, kE)
687    return Xanom
688   
689   
690################################################################################
691##### distance, angle, planes, torsion stuff stuff
692################################################################################
693
694def getSyXYZ(XYZ,ops,SGData):
695    '''default doc string
696   
697    :param type name: description
698   
699    :returns: type name: description
700   
701    '''
702    XYZout = np.zeros_like(XYZ)
703    for i,[xyz,op] in enumerate(zip(XYZ,ops)):
704        if op == '1':
705            XYZout[i] = xyz
706        else:
707            oprs = op.split('+')
708            unit = [0,0,0]
709            if oprs[1]:
710                unit = np.array(list(eval(oprs[1])))
711            syop =int(oprs[0])
712            inv = syop/abs(syop)
713            syop *= inv
714            cent = syop/100
715            syop %= 100
716            syop -= 1
717            M,T = SGData['SGOps'][syop]
718            XYZout[i] = (np.inner(M,xyz)+T)*inv+SGData['SGCen'][cent]+unit
719    return XYZout
720   
721def getRestDist(XYZ,Amat):
722    '''default doc string
723   
724    :param type name: description
725   
726    :returns: type name: description
727   
728    '''
729    return np.sqrt(np.sum(np.inner(Amat,(XYZ[1]-XYZ[0]))**2))
730   
731def getRestDeriv(Func,XYZ,Amat,ops,SGData):
732    '''default doc string
733   
734    :param type name: description
735   
736    :returns: type name: description
737   
738    '''
739    deriv = np.zeros((len(XYZ),3))
740    dx = 0.00001
741    for j,xyz in enumerate(XYZ):
742        for i,x in enumerate(np.array([[dx,0,0],[0,dx,0],[0,0,dx]])):
743            XYZ[j] -= x
744            d1 = Func(getSyXYZ(XYZ,ops,SGData),Amat)
745            XYZ[j] += 2*x
746            d2 = Func(getSyXYZ(XYZ,ops,SGData),Amat)
747            XYZ[j] -= x
748            deriv[j][i] = (d1-d2)/(2*dx)
749    return deriv.flatten()
750
751def getRestAngle(XYZ,Amat):
752    '''default doc string
753   
754    :param type name: description
755   
756    :returns: type name: description
757   
758    '''
759   
760    def calcVec(Ox,Tx,Amat):
761        return np.inner(Amat,(Tx-Ox))
762
763    VecA = calcVec(XYZ[1],XYZ[0],Amat)
764    VecA /= np.sqrt(np.sum(VecA**2))
765    VecB = calcVec(XYZ[1],XYZ[2],Amat)
766    VecB /= np.sqrt(np.sum(VecB**2))
767    edge = VecB-VecA
768    edge = np.sum(edge**2)
769    angle = (2.-edge)/2.
770    angle = max(angle,-1.)
771    return acosd(angle)
772   
773def getRestPlane(XYZ,Amat):
774    '''default doc string
775   
776    :param type name: description
777   
778    :returns: type name: description
779   
780    '''
781    sumXYZ = np.zeros(3)
782    for xyz in XYZ:
783        sumXYZ += xyz
784    sumXYZ /= len(XYZ)
785    XYZ = np.array(XYZ)-sumXYZ
786    XYZ = np.inner(Amat,XYZ).T
787    Zmat = np.zeros((3,3))
788    for i,xyz in enumerate(XYZ):
789        Zmat += np.outer(xyz.T,xyz)
790    Evec,Emat = nl.eig(Zmat)
791    Evec = np.sqrt(Evec)/(len(XYZ)-3)
792    Order = np.argsort(Evec)
793    return Evec[Order[0]]
794   
795def getRestChiral(XYZ,Amat):   
796    '''default doc string
797   
798    :param type name: description
799   
800    :returns: type name: description
801   
802    '''
803    VecA = np.empty((3,3))   
804    VecA[0] = np.inner(XYZ[1]-XYZ[0],Amat)
805    VecA[1] = np.inner(XYZ[2]-XYZ[0],Amat)
806    VecA[2] = np.inner(XYZ[3]-XYZ[0],Amat)
807    return nl.det(VecA)
808   
809def getRestTorsion(XYZ,Amat):
810    '''default doc string
811   
812    :param type name: description
813   
814    :returns: type name: description
815   
816    '''
817    VecA = np.empty((3,3))
818    VecA[0] = np.inner(XYZ[1]-XYZ[0],Amat)
819    VecA[1] = np.inner(XYZ[2]-XYZ[1],Amat)
820    VecA[2] = np.inner(XYZ[3]-XYZ[2],Amat)
821    D = nl.det(VecA)
822    Mag = np.sqrt(np.sum(VecA*VecA,axis=1))
823    P12 = np.sum(VecA[0]*VecA[1])/(Mag[0]*Mag[1])
824    P13 = np.sum(VecA[0]*VecA[2])/(Mag[0]*Mag[2])
825    P23 = np.sum(VecA[1]*VecA[2])/(Mag[1]*Mag[2])
826    Ang = 1.0
827    if abs(P12) < 1.0 and abs(P23) < 1.0:
828        Ang = (P12*P23-P13)/(np.sqrt(1.-P12**2)*np.sqrt(1.-P23**2))
829    TOR = (acosd(Ang)*D/abs(D)+720.)%360.
830    return TOR
831   
832def calcTorsionEnergy(TOR,Coeff=[]):
833    '''default doc string
834   
835    :param type name: description
836   
837    :returns: type name: description
838   
839    '''
840    sum = 0.
841    if len(Coeff):
842        cof = np.reshape(Coeff,(3,3)).T
843        delt = TOR-cof[1]
844        delt = np.where(delt<-180.,delt+360.,delt)
845        delt = np.where(delt>180.,delt-360.,delt)
846        term = -cof[2]*delt**2
847        val = cof[0]*np.exp(term/1000.0)
848        pMax = cof[0][np.argmin(val)]
849        Eval = np.sum(val)
850        sum = Eval-pMax
851    return sum,Eval
852
853def getTorsionDeriv(XYZ,Amat,Coeff):
854    '''default doc string
855   
856    :param type name: description
857   
858    :returns: type name: description
859   
860    '''
861    deriv = np.zeros((len(XYZ),3))
862    dx = 0.00001
863    for j,xyz in enumerate(XYZ):
864        for i,x in enumerate(np.array([[dx,0,0],[0,dx,0],[0,0,dx]])):
865            XYZ[j] -= x
866            tor = getRestTorsion(XYZ,Amat)
867            p1,d1 = calcTorsionEnergy(tor,Coeff)
868            XYZ[j] += 2*x
869            tor = getRestTorsion(XYZ,Amat)
870            p2,d2 = calcTorsionEnergy(tor,Coeff)           
871            XYZ[j] -= x
872            deriv[j][i] = (p2-p1)/(2*dx)
873    return deriv.flatten()
874
875def getRestRama(XYZ,Amat):
876    '''Computes a pair of torsion angles in a 5 atom string
877   
878    :param nparray XYZ: crystallographic coordinates of 5 atoms
879    :param nparray Amat: crystal to cartesian transformation matrix
880   
881    :returns: list (phi,psi) two torsion angles in degrees
882   
883    '''
884    phi = getRestTorsion(XYZ[:5],Amat)
885    psi = getRestTorsion(XYZ[1:],Amat)
886    return phi,psi
887   
888def calcRamaEnergy(phi,psi,Coeff=[]):
889    '''Computes pseudo potential energy from a pair of torsion angles and a
890    numerical description of the potential energy surface. Used to create
891    penalty function in LS refinement:     
892    :math:`Eval(\\phi,\\psi) = C[0]*exp(-V/1000)`
893
894    where :math:`V = -C[3] * (\\phi-C[1])^2 - C[4]*(\\psi-C[2])^2 - 2*(\\phi-C[1])*(\\psi-C[2])`
895   
896    :param float phi: first torsion angle (:math:`\\phi`)
897    :param float psi: second torsion angle (:math:`\\psi`)
898    :param list Coeff: pseudo potential coefficients
899   
900    :returns: list (sum,Eval): pseudo-potential difference from minimum & value;
901      sum is used for penalty function.
902   
903    '''
904    sum = 0.
905    Eval = 0.
906    if len(Coeff):
907        cof = Coeff.T
908        dPhi = phi-cof[1]
909        dPhi = np.where(dPhi<-180.,dPhi+360.,dPhi)
910        dPhi = np.where(dPhi>180.,dPhi-360.,dPhi)
911        dPsi = psi-cof[2]
912        dPsi = np.where(dPsi<-180.,dPsi+360.,dPsi)
913        dPsi = np.where(dPsi>180.,dPsi-360.,dPsi)
914        val = -cof[3]*dPhi**2-cof[4]*dPsi**2-2.0*cof[5]*dPhi*dPsi
915        val = cof[0]*np.exp(val/1000.)
916        pMax = cof[0][np.argmin(val)]
917        Eval = np.sum(val)
918        sum = Eval-pMax
919    return sum,Eval
920
921def getRamaDeriv(XYZ,Amat,Coeff):
922    '''Computes numerical derivatives of torsion angle pair pseudo potential
923    with respect of crystallographic atom coordinates of the 5 atom sequence
924   
925    :param nparray XYZ: crystallographic coordinates of 5 atoms
926    :param nparray Amat: crystal to cartesian transformation matrix
927    :param list Coeff: pseudo potential coefficients
928   
929    :returns: list (deriv) derivatives of pseudopotential with respect to 5 atom
930     crystallographic xyz coordinates.
931   
932    '''
933    deriv = np.zeros((len(XYZ),3))
934    dx = 0.00001
935    for j,xyz in enumerate(XYZ):
936        for i,x in enumerate(np.array([[dx,0,0],[0,dx,0],[0,0,dx]])):
937            XYZ[j] -= x
938            phi,psi = getRestRama(XYZ,Amat)
939            p1,d1 = calcRamaEnergy(phi,psi,Coeff)
940            XYZ[j] += 2*x
941            phi,psi = getRestRama(XYZ,Amat)
942            p2,d2 = calcRamaEnergy(phi,psi,Coeff)
943            XYZ[j] -= x
944            deriv[j][i] = (p2-p1)/(2*dx)
945    return deriv.flatten()
946
947def getRestPolefig(ODFln,SamSym,Grid):
948    '''default doc string
949   
950    :param type name: description
951   
952    :returns: type name: description
953   
954    '''
955    X,Y = np.meshgrid(np.linspace(1.,-1.,Grid),np.linspace(-1.,1.,Grid))
956    R,P = np.sqrt(X**2+Y**2).flatten(),atan2d(Y,X).flatten()
957    R = np.where(R <= 1.,2.*atand(R),0.0)
958    Z = np.zeros_like(R)
959    Z = G2lat.polfcal(ODFln,SamSym,R,P)
960    Z = np.reshape(Z,(Grid,Grid))
961    return np.reshape(R,(Grid,Grid)),np.reshape(P,(Grid,Grid)),Z
962
963def getRestPolefigDerv(HKL,Grid,SHCoeff):
964    '''default doc string
965   
966    :param type name: description
967   
968    :returns: type name: description
969   
970    '''
971    pass
972       
973def getDistDerv(Oxyz,Txyz,Amat,Tunit,Top,SGData):
974    '''default doc string
975   
976    :param type name: description
977   
978    :returns: type name: description
979   
980    '''
981    def calcDist(Ox,Tx,U,inv,C,M,T,Amat):
982        TxT = inv*(np.inner(M,Tx)+T)+C+U
983        return np.sqrt(np.sum(np.inner(Amat,(TxT-Ox))**2))
984       
985    inv = Top/abs(Top)
986    cent = abs(Top)/100
987    op = abs(Top)%100-1
988    M,T = SGData['SGOps'][op]
989    C = SGData['SGCen'][cent]
990    dx = .00001
991    deriv = np.zeros(6)
992    for i in [0,1,2]:
993        Oxyz[i] -= dx
994        d0 = calcDist(Oxyz,Txyz,Tunit,inv,C,M,T,Amat)
995        Oxyz[i] += 2*dx
996        deriv[i] = (calcDist(Oxyz,Txyz,Tunit,inv,C,M,T,Amat)-d0)/(2.*dx)
997        Oxyz[i] -= dx
998        Txyz[i] -= dx
999        d0 = calcDist(Oxyz,Txyz,Tunit,inv,C,M,T,Amat)
1000        Txyz[i] += 2*dx
1001        deriv[i+3] = (calcDist(Oxyz,Txyz,Tunit,inv,C,M,T,Amat)-d0)/(2.*dx)
1002        Txyz[i] -= dx
1003    return deriv
1004   
1005def getAngSig(VA,VB,Amat,SGData,covData={}):
1006    '''default doc string
1007   
1008    :param type name: description
1009   
1010    :returns: type name: description
1011   
1012    '''
1013    def calcVec(Ox,Tx,U,inv,C,M,T,Amat):
1014        TxT = inv*(np.inner(M,Tx)+T)+C
1015        TxT = G2spc.MoveToUnitCell(TxT)+U
1016        return np.inner(Amat,(TxT-Ox))
1017       
1018    def calcAngle(Ox,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat):
1019        VecA = calcVec(Ox,TxA,unitA,invA,CA,MA,TA,Amat)
1020        VecA /= np.sqrt(np.sum(VecA**2))
1021        VecB = calcVec(Ox,TxB,unitB,invB,CB,MB,TB,Amat)
1022        VecB /= np.sqrt(np.sum(VecB**2))
1023        edge = VecB-VecA
1024        edge = np.sum(edge**2)
1025        angle = (2.-edge)/2.
1026        angle = max(angle,-1.)
1027        return acosd(angle)
1028       
1029    OxAN,OxA,TxAN,TxA,unitA,TopA = VA
1030    OxBN,OxB,TxBN,TxB,unitB,TopB = VB
1031    invA = invB = 1
1032    invA = TopA/abs(TopA)
1033    invB = TopB/abs(TopB)
1034    centA = abs(TopA)/100
1035    centB = abs(TopB)/100
1036    opA = abs(TopA)%100-1
1037    opB = abs(TopB)%100-1
1038    MA,TA = SGData['SGOps'][opA]
1039    MB,TB = SGData['SGOps'][opB]
1040    CA = SGData['SGCen'][centA]
1041    CB = SGData['SGCen'][centB]
1042    if 'covMatrix' in covData:
1043        covMatrix = covData['covMatrix']
1044        varyList = covData['varyList']
1045        AngVcov = getVCov(OxAN+TxAN+TxBN,varyList,covMatrix)
1046        dx = .00001
1047        dadx = np.zeros(9)
1048        Ang = calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)
1049        for i in [0,1,2]:
1050            OxA[i] -= dx
1051            a0 = calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)
1052            OxA[i] += 2*dx
1053            dadx[i] = (calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)-a0)/(2*dx)
1054            OxA[i] -= dx
1055           
1056            TxA[i] -= dx
1057            a0 = calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)
1058            TxA[i] += 2*dx
1059            dadx[i+3] = (calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)-a0)/(2*dx)
1060            TxA[i] -= dx
1061           
1062            TxB[i] -= dx
1063            a0 = calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)
1064            TxB[i] += 2*dx
1065            dadx[i+6] = (calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat)-a0)/(2*dx)
1066            TxB[i] -= dx
1067           
1068        sigAng = np.sqrt(np.inner(dadx,np.inner(AngVcov,dadx)))
1069        if sigAng < 0.01:
1070            sigAng = 0.0
1071        return Ang,sigAng
1072    else:
1073        return calcAngle(OxA,TxA,TxB,unitA,unitB,invA,CA,MA,TA,invB,CB,MB,TB,Amat),0.0
1074
1075def GetDistSig(Oatoms,Atoms,Amat,SGData,covData={}):
1076    '''default doc string
1077   
1078    :param type name: description
1079   
1080    :returns: type name: description
1081   
1082    '''
1083    def calcDist(Atoms,SyOps,Amat):
1084        XYZ = []
1085        for i,atom in enumerate(Atoms):
1086            Inv,M,T,C,U = SyOps[i]
1087            XYZ.append(np.array(atom[1:4]))
1088            XYZ[-1] = Inv*(np.inner(M,np.array(XYZ[-1]))+T)+C+U
1089            XYZ[-1] = np.inner(Amat,XYZ[-1]).T
1090        V1 = XYZ[1]-XYZ[0]
1091        return np.sqrt(np.sum(V1**2))
1092       
1093    Inv = []
1094    SyOps = []
1095    names = []
1096    for i,atom in enumerate(Oatoms):
1097        names += atom[-1]
1098        Op,unit = Atoms[i][-1]
1099        inv = Op/abs(Op)
1100        m,t = SGData['SGOps'][abs(Op)%100-1]
1101        c = SGData['SGCen'][abs(Op)/100]
1102        SyOps.append([inv,m,t,c,unit])
1103    Dist = calcDist(Oatoms,SyOps,Amat)
1104   
1105    sig = -0.001
1106    if 'covMatrix' in covData:
1107        parmNames = []
1108        dx = .00001
1109        dadx = np.zeros(6)
1110        for i in range(6):
1111            ia = i/3
1112            ix = i%3
1113            Oatoms[ia][ix+1] += dx
1114            a0 = calcDist(Oatoms,SyOps,Amat)
1115            Oatoms[ia][ix+1] -= 2*dx
1116            dadx[i] = (calcDist(Oatoms,SyOps,Amat)-a0)/(2.*dx)
1117        covMatrix = covData['covMatrix']
1118        varyList = covData['varyList']
1119        DistVcov = getVCov(names,varyList,covMatrix)
1120        sig = np.sqrt(np.inner(dadx,np.inner(DistVcov,dadx)))
1121        if sig < 0.001:
1122            sig = -0.001
1123   
1124    return Dist,sig
1125
1126def GetAngleSig(Oatoms,Atoms,Amat,SGData,covData={}):
1127    '''default doc string
1128   
1129    :param type name: description
1130   
1131    :returns: type name: description
1132   
1133    '''
1134
1135    def calcAngle(Atoms,SyOps,Amat):
1136        XYZ = []
1137        for i,atom in enumerate(Atoms):
1138            Inv,M,T,C,U = SyOps[i]
1139            XYZ.append(np.array(atom[1:4]))
1140            XYZ[-1] = Inv*(np.inner(M,np.array(XYZ[-1]))+T)+C+U
1141            XYZ[-1] = np.inner(Amat,XYZ[-1]).T
1142        V1 = XYZ[1]-XYZ[0]
1143        V1 /= np.sqrt(np.sum(V1**2))
1144        V2 = XYZ[1]-XYZ[2]
1145        V2 /= np.sqrt(np.sum(V2**2))
1146        V3 = V2-V1
1147        cang = min(1.,max((2.-np.sum(V3**2))/2.,-1.))
1148        return acosd(cang)
1149
1150    Inv = []
1151    SyOps = []
1152    names = []
1153    for i,atom in enumerate(Oatoms):
1154        names += atom[-1]
1155        Op,unit = Atoms[i][-1]
1156        inv = Op/abs(Op)
1157        m,t = SGData['SGOps'][abs(Op)%100-1]
1158        c = SGData['SGCen'][abs(Op)/100]
1159        SyOps.append([inv,m,t,c,unit])
1160    Angle = calcAngle(Oatoms,SyOps,Amat)
1161   
1162    sig = -0.01
1163    if 'covMatrix' in covData:
1164        parmNames = []
1165        dx = .00001
1166        dadx = np.zeros(9)
1167        for i in range(9):
1168            ia = i/3
1169            ix = i%3
1170            Oatoms[ia][ix+1] += dx
1171            a0 = calcAngle(Oatoms,SyOps,Amat)
1172            Oatoms[ia][ix+1] -= 2*dx
1173            dadx[i] = (calcAngle(Oatoms,SyOps,Amat)-a0)/(2.*dx)
1174        covMatrix = covData['covMatrix']
1175        varyList = covData['varyList']
1176        AngVcov = getVCov(names,varyList,covMatrix)
1177        sig = np.sqrt(np.inner(dadx,np.inner(AngVcov,dadx)))
1178        if sig < 0.01:
1179            sig = -0.01
1180   
1181    return Angle,sig
1182
1183def GetTorsionSig(Oatoms,Atoms,Amat,SGData,covData={}):
1184    '''default doc string
1185   
1186    :param type name: description
1187   
1188    :returns: type name: description
1189   
1190    '''
1191
1192    def calcTorsion(Atoms,SyOps,Amat):
1193       
1194        XYZ = []
1195        for i,atom in enumerate(Atoms):
1196            Inv,M,T,C,U = SyOps[i]
1197            XYZ.append(np.array(atom[1:4]))
1198            XYZ[-1] = Inv*(np.inner(M,np.array(XYZ[-1]))+T)+C+U
1199            XYZ[-1] = np.inner(Amat,XYZ[-1]).T
1200        V1 = XYZ[1]-XYZ[0]
1201        V2 = XYZ[2]-XYZ[1]
1202        V3 = XYZ[3]-XYZ[2]
1203        V1 /= np.sqrt(np.sum(V1**2))
1204        V2 /= np.sqrt(np.sum(V2**2))
1205        V3 /= np.sqrt(np.sum(V3**2))
1206        M = np.array([V1,V2,V3])
1207        D = nl.det(M)
1208        Ang = 1.0
1209        P12 = np.dot(V1,V2)
1210        P13 = np.dot(V1,V3)
1211        P23 = np.dot(V2,V3)
1212        Tors = acosd((P12*P23-P13)/(np.sqrt(1.-P12**2)*np.sqrt(1.-P23**2)))*D/abs(D)
1213        return Tors
1214           
1215    Inv = []
1216    SyOps = []
1217    names = []
1218    for i,atom in enumerate(Oatoms):
1219        names += atom[-1]
1220        Op,unit = Atoms[i][-1]
1221        inv = Op/abs(Op)
1222        m,t = SGData['SGOps'][abs(Op)%100-1]
1223        c = SGData['SGCen'][abs(Op)/100]
1224        SyOps.append([inv,m,t,c,unit])
1225    Tors = calcTorsion(Oatoms,SyOps,Amat)
1226   
1227    sig = -0.01
1228    if 'covMatrix' in covData:
1229        parmNames = []
1230        dx = .00001
1231        dadx = np.zeros(12)
1232        for i in range(12):
1233            ia = i/3
1234            ix = i%3
1235            Oatoms[ia][ix+1] -= dx
1236            a0 = calcTorsion(Oatoms,SyOps,Amat)
1237            Oatoms[ia][ix+1] += 2*dx
1238            dadx[i] = (calcTorsion(Oatoms,SyOps,Amat)-a0)/(2.*dx)
1239            Oatoms[ia][ix+1] -= dx           
1240        covMatrix = covData['covMatrix']
1241        varyList = covData['varyList']
1242        TorVcov = getVCov(names,varyList,covMatrix)
1243        sig = np.sqrt(np.inner(dadx,np.inner(TorVcov,dadx)))
1244        if sig < 0.01:
1245            sig = -0.01
1246   
1247    return Tors,sig
1248       
1249def GetDATSig(Oatoms,Atoms,Amat,SGData,covData={}):
1250    '''default doc string
1251   
1252    :param type name: description
1253   
1254    :returns: type name: description
1255   
1256    '''
1257
1258    def calcDist(Atoms,SyOps,Amat):
1259        XYZ = []
1260        for i,atom in enumerate(Atoms):
1261            Inv,M,T,C,U = SyOps[i]
1262            XYZ.append(np.array(atom[1:4]))
1263            XYZ[-1] = Inv*(np.inner(M,np.array(XYZ[-1]))+T)+C+U
1264            XYZ[-1] = np.inner(Amat,XYZ[-1]).T
1265        V1 = XYZ[1]-XYZ[0]
1266        return np.sqrt(np.sum(V1**2))
1267       
1268    def calcAngle(Atoms,SyOps,Amat):
1269        XYZ = []
1270        for i,atom in enumerate(Atoms):
1271            Inv,M,T,C,U = SyOps[i]
1272            XYZ.append(np.array(atom[1:4]))
1273            XYZ[-1] = Inv*(np.inner(M,np.array(XYZ[-1]))+T)+C+U
1274            XYZ[-1] = np.inner(Amat,XYZ[-1]).T
1275        V1 = XYZ[1]-XYZ[0]
1276        V1 /= np.sqrt(np.sum(V1**2))
1277        V2 = XYZ[1]-XYZ[2]
1278        V2 /= np.sqrt(np.sum(V2**2))
1279        V3 = V2-V1
1280        cang = min(1.,max((2.-np.sum(V3**2))/2.,-1.))
1281        return acosd(cang)
1282
1283    def calcTorsion(Atoms,SyOps,Amat):
1284       
1285        XYZ = []
1286        for i,atom in enumerate(Atoms):
1287            Inv,M,T,C,U = SyOps[i]
1288            XYZ.append(np.array(atom[1:4]))
1289            XYZ[-1] = Inv*(np.inner(M,np.array(XYZ[-1]))+T)+C+U
1290            XYZ[-1] = np.inner(Amat,XYZ[-1]).T
1291        V1 = XYZ[1]-XYZ[0]
1292        V2 = XYZ[2]-XYZ[1]
1293        V3 = XYZ[3]-XYZ[2]
1294        V1 /= np.sqrt(np.sum(V1**2))
1295        V2 /= np.sqrt(np.sum(V2**2))
1296        V3 /= np.sqrt(np.sum(V3**2))
1297        M = np.array([V1,V2,V3])
1298        D = nl.det(M)
1299        Ang = 1.0
1300        P12 = np.dot(V1,V2)
1301        P13 = np.dot(V1,V3)
1302        P23 = np.dot(V2,V3)
1303        Tors = acosd((P12*P23-P13)/(np.sqrt(1.-P12**2)*np.sqrt(1.-P23**2)))*D/abs(D)
1304        return Tors
1305           
1306    Inv = []
1307    SyOps = []
1308    names = []
1309    for i,atom in enumerate(Oatoms):
1310        names += atom[-1]
1311        Op,unit = Atoms[i][-1]
1312        inv = Op/abs(Op)
1313        m,t = SGData['SGOps'][abs(Op)%100-1]
1314        c = SGData['SGCen'][abs(Op)/100]
1315        SyOps.append([inv,m,t,c,unit])
1316    M = len(Oatoms)
1317    if M == 2:
1318        Val = calcDist(Oatoms,SyOps,Amat)
1319    elif M == 3:
1320        Val = calcAngle(Oatoms,SyOps,Amat)
1321    else:
1322        Val = calcTorsion(Oatoms,SyOps,Amat)
1323   
1324    sigVals = [-0.001,-0.01,-0.01]
1325    sig = sigVals[M-3]
1326    if 'covMatrix' in covData:
1327        parmNames = []
1328        dx = .00001
1329        N = M*3
1330        dadx = np.zeros(N)
1331        for i in range(N):
1332            ia = i/3
1333            ix = i%3
1334            Oatoms[ia][ix+1] += dx
1335            if M == 2:
1336                a0 = calcDist(Oatoms,SyOps,Amat)
1337            elif M == 3:
1338                a0 = calcAngle(Oatoms,SyOps,Amat)
1339            else:
1340                a0 = calcTorsion(Oatoms,SyOps,Amat)
1341            Oatoms[ia][ix+1] -= 2*dx
1342            if M == 2:
1343                dadx[i] = (calcDist(Oatoms,SyOps,Amat)-a0)/(2.*dx)               
1344            elif M == 3:
1345                dadx[i] = (calcAngle(Oatoms,SyOps,Amat)-a0)/(2.*dx)               
1346            else:
1347                dadx[i] = (calcTorsion(Oatoms,SyOps,Amat)-a0)/(2.*dx)
1348        covMatrix = covData['covMatrix']
1349        varyList = covData['varyList']
1350        Vcov = getVCov(names,varyList,covMatrix)
1351        sig = np.sqrt(np.inner(dadx,np.inner(Vcov,dadx)))
1352        if sig < sigVals[M-3]:
1353            sig = sigVals[M-3]
1354   
1355    return Val,sig
1356       
1357def ValEsd(value,esd=0,nTZ=False):
1358    '''Format a floating point number with a given level of precision or
1359    with in crystallographic format with a "esd", as value(esd). If esd is
1360    negative the number is formatted with the level of significant figures
1361    appropriate if abs(esd) were the esd, but the esd is not included.
1362    if the esd is zero, approximately 6 significant figures are printed.
1363    nTZ=True causes "extra" zeros to be removed after the decimal place.
1364    for example:
1365
1366      * "1.235(3)" for value=1.2346 & esd=0.003
1367      * "1.235(3)e4" for value=12346. & esd=30
1368      * "1.235(3)e6" for value=0.12346e7 & esd=3000
1369      * "1.235" for value=1.2346 & esd=-0.003
1370      * "1.240" for value=1.2395 & esd=-0.003
1371      * "1.24" for value=1.2395 & esd=-0.003 with nTZ=True
1372      * "1.23460" for value=1.2346 & esd=0.0
1373
1374    :param float value: number to be formatted
1375    :param float esd: uncertainty or if esd < 0, specifies level of
1376      precision to be shown e.g. esd=-0.01 gives 2 places beyond decimal
1377    :param bool nTZ: True to remove trailing zeros (default is False)
1378    :returns: value(esd) or value as a string
1379
1380    '''
1381    # Note: this routine is Python 3 compatible -- I think
1382    if math.isnan(value): # invalid value, bail out
1383        return '?'
1384    if math.isnan(esd): # invalid esd, treat as zero
1385        esd = 0
1386        esdoff = 5
1387    elif esd != 0:
1388        # transform the esd to a one or two digit integer
1389        l = math.log10(abs(esd)) % 1
1390        # cut off of 19 1.9==>(19) but 1.95==>(2) N.B. log10(1.95) = 0.2900...
1391        if l < 0.290034611362518: l+= 1.       
1392        intesd = int(round(10**l)) # esd as integer
1393        # determine the number of digits offset for the esd
1394        esdoff = int(round(math.log10(intesd*1./abs(esd))))
1395    else:
1396        esdoff = 5
1397    valoff = 0
1398    if abs(value) < abs(esdoff): # value is effectively zero
1399        pass
1400    elif esdoff < 0 or abs(value) > 1.0e6 or abs(value) < 1.0e-4: # use scientific notation
1401        # where the digit offset is to the left of the decimal place or where too many
1402        # digits are needed
1403        if abs(value) > 1:
1404            valoff = int(math.log10(abs(value)))
1405        elif abs(value) > 0:
1406            valoff = int(math.log10(abs(value))-0.9999999)
1407        else:
1408            valoff = 0
1409    if esd != 0:
1410        if valoff+esdoff < 0:
1411            valoff = esdoff = 0
1412        out = ("{:."+str(valoff+esdoff)+"f}").format(value/10**valoff) # format the value
1413    elif valoff != 0: # esd = 0; exponential notation ==> esdoff decimal places
1414        out = ("{:."+str(esdoff)+"f}").format(value/10**valoff) # format the value
1415    else: # esd = 0; non-exponential notation ==> esdoff+1 significant digits
1416        if abs(value) > 0:           
1417            extra = -math.log10(abs(value))
1418        else:
1419            extra = 0
1420        if extra > 0: extra += 1
1421        out = ("{:."+str(max(0,esdoff+int(extra)))+"f}").format(value) # format the value
1422    if esd > 0:
1423        out += ("({:d})").format(intesd)  # add the esd
1424    elif nTZ and '.' in out:
1425        out = out.rstrip('0')  # strip zeros to right of decimal
1426        out = out.rstrip('.')  # and decimal place when not needed
1427    if valoff != 0:
1428        out += ("e{:d}").format(valoff) # add an exponent, when needed
1429    return out
1430
1431################################################################################
1432##### Fourier & charge flip stuff
1433################################################################################
1434
1435def adjHKLmax(SGData,Hmax):
1436    '''default doc string
1437   
1438    :param type name: description
1439   
1440    :returns: type name: description
1441   
1442    '''
1443    if SGData['SGLaue'] in ['3','3m1','31m','6/m','6/mmm']:
1444        Hmax[0] = ((Hmax[0]+3)/6)*6
1445        Hmax[1] = ((Hmax[1]+3)/6)*6
1446        Hmax[2] = ((Hmax[2]+1)/4)*4
1447    else:
1448        Hmax[0] = ((Hmax[0]+2)/4)*4
1449        Hmax[1] = ((Hmax[1]+2)/4)*4
1450        Hmax[2] = ((Hmax[2]+1)/4)*4
1451
1452def OmitMap(data,reflDict,pgbar=None):
1453    '''default doc string
1454   
1455    :param type name: description
1456   
1457    :returns: type name: description
1458   
1459    '''
1460    generalData = data['General']
1461    if not generalData['Map']['MapType']:
1462        print '**** ERROR - Fourier map not defined'
1463        return
1464    mapData = generalData['Map']
1465    dmin = mapData['Resolution']
1466    SGData = generalData['SGData']
1467    SGMT = np.array([ops[0].T for ops in SGData['SGOps']])
1468    SGT = np.array([ops[1] for ops in SGData['SGOps']])
1469    cell = generalData['Cell'][1:8]       
1470    A = G2lat.cell2A(cell[:6])
1471    Hmax = np.asarray(G2lat.getHKLmax(dmin,SGData,A),dtype='i')+1
1472    adjHKLmax(SGData,Hmax)
1473    Fhkl = np.zeros(shape=2*Hmax,dtype='c16')
1474    time0 = time.time()
1475    for iref,ref in enumerate(reflDict['RefList']):
1476        if ref[4] >= dmin:
1477            Fosq,Fcsq,ph = ref[8:11]
1478            Uniq = np.inner(ref[:3],SGMT)
1479            Phi = np.inner(ref[:3],SGT)
1480            for i,hkl in enumerate(Uniq):        #uses uniq
1481                hkl = np.asarray(hkl,dtype='i')
1482                dp = 360.*Phi[i]                #and phi
1483                a = cosd(ph+dp)
1484                b = sind(ph+dp)
1485                phasep = complex(a,b)
1486                phasem = complex(a,-b)
1487                Fo = np.sqrt(Fosq)
1488                if '2Fo-Fc' in mapData['MapType']:
1489                    F = 2.*np.sqrt(Fosq)-np.sqrt(Fcsq)
1490                else:
1491                    F = np.sqrt(Fosq)
1492                h,k,l = hkl+Hmax
1493                Fhkl[h,k,l] = F*phasep
1494                h,k,l = -hkl+Hmax
1495                Fhkl[h,k,l] = F*phasem
1496    rho0 = fft.fftn(fft.fftshift(Fhkl))/cell[6]
1497    M = np.mgrid[0:4,0:4,0:4]
1498    blkIds = np.array(zip(M[0].flatten(),M[1].flatten(),M[2].flatten()))
1499    iBeg = blkIds*rho0.shape/4
1500    iFin = (blkIds+1)*rho0.shape/4
1501    rho_omit = np.zeros_like(rho0)
1502    nBlk = 0
1503    for iB,iF in zip(iBeg,iFin):
1504        rho1 = np.copy(rho0)
1505        rho1[iB[0]:iF[0],iB[1]:iF[1],iB[2]:iF[2]] = 0.
1506        Fnew = fft.ifftshift(fft.ifftn(rho1))
1507        Fnew = np.where(Fnew,Fnew,1.0)           #avoid divide by zero
1508        phase = Fnew/np.absolute(Fnew)
1509        OFhkl = np.absolute(Fhkl)*phase
1510        rho1 = np.real(fft.fftn(fft.fftshift(OFhkl)))*(1.+0j)
1511        rho_omit[iB[0]:iF[0],iB[1]:iF[1],iB[2]:iF[2]] = np.copy(rho1[iB[0]:iF[0],iB[1]:iF[1],iB[2]:iF[2]])
1512        nBlk += 1
1513        pgbar.Update(nBlk)
1514    mapData['rho'] = np.real(rho_omit)/cell[6]
1515    mapData['rhoMax'] = max(np.max(mapData['rho']),-np.min(mapData['rho']))
1516    print 'Omit map time: %.4f'%(time.time()-time0),'no. elements: %d'%(Fhkl.size)
1517    return mapData
1518
1519def FourierMap(data,reflDict):
1520    '''default doc string
1521   
1522    :param type name: description
1523   
1524    :returns: type name: description
1525   
1526    '''
1527    generalData = data['General']
1528    if not generalData['Map']['MapType']:
1529        print '**** ERROR - Fourier map not defined'
1530        return
1531    mapData = generalData['Map']
1532    dmin = mapData['Resolution']
1533    SGData = generalData['SGData']
1534    SGMT = np.array([ops[0].T for ops in SGData['SGOps']])
1535    SGT = np.array([ops[1] for ops in SGData['SGOps']])
1536    cell = generalData['Cell'][1:8]       
1537    A = G2lat.cell2A(cell[:6])
1538    Hmax = np.asarray(G2lat.getHKLmax(dmin,SGData,A),dtype='i')+1
1539    adjHKLmax(SGData,Hmax)
1540    Fhkl = np.zeros(shape=2*Hmax,dtype='c16')
1541#    Fhkl[0,0,0] = generalData['F000X']
1542    time0 = time.time()
1543    for iref,ref in enumerate(reflDict['RefList']):
1544        if ref[4] >= dmin:
1545            Fosq,Fcsq,ph = ref[8:11]
1546            Uniq = np.inner(ref[:3],SGMT)
1547            Phi = np.inner(ref[:3],SGT)
1548            for i,hkl in enumerate(Uniq):        #uses uniq
1549                hkl = np.asarray(hkl,dtype='i')
1550                dp = 360.*Phi[i]                #and phi
1551                a = cosd(ph+dp)
1552                b = sind(ph+dp)
1553                phasep = complex(a,b)
1554                phasem = complex(a,-b)
1555                if 'Fobs' in mapData['MapType']:
1556                    F = np.where(Fosq>0.,np.sqrt(Fosq),0.)
1557                    h,k,l = hkl+Hmax
1558                    Fhkl[h,k,l] = F*phasep
1559                    h,k,l = -hkl+Hmax
1560                    Fhkl[h,k,l] = F*phasem
1561                elif 'Fcalc' in mapData['MapType']:
1562                    F = np.sqrt(Fcsq)
1563                    h,k,l = hkl+Hmax
1564                    Fhkl[h,k,l] = F*phasep
1565                    h,k,l = -hkl+Hmax
1566                    Fhkl[h,k,l] = F*phasem
1567                elif 'delt-F' in mapData['MapType']:
1568                    dF = np.where(Fosq>0.,np.sqrt(Fosq),0.)-np.sqrt(Fcsq)
1569                    h,k,l = hkl+Hmax
1570                    Fhkl[h,k,l] = dF*phasep
1571                    h,k,l = -hkl+Hmax
1572                    Fhkl[h,k,l] = dF*phasem
1573                elif '2*Fo-Fc' in mapData['MapType']:
1574                    F = 2.*np.where(Fosq>0.,np.sqrt(Fosq),0.)-np.sqrt(Fcsq)
1575                    h,k,l = hkl+Hmax
1576                    Fhkl[h,k,l] = F*phasep
1577                    h,k,l = -hkl+Hmax
1578                    Fhkl[h,k,l] = F*phasem
1579                elif 'Patterson' in mapData['MapType']:
1580                    h,k,l = hkl+Hmax
1581                    Fhkl[h,k,l] = complex(Fosq,0.)
1582                    h,k,l = -hkl+Hmax
1583                    Fhkl[h,k,l] = complex(Fosq,0.)
1584    rho = fft.fftn(fft.fftshift(Fhkl))/cell[6]
1585    print 'Fourier map time: %.4f'%(time.time()-time0),'no. elements: %d'%(Fhkl.size)
1586    mapData['rho'] = np.real(rho)
1587    mapData['rhoMax'] = max(np.max(mapData['rho']),-np.min(mapData['rho']))
1588    return mapData
1589   
1590# map printing for testing purposes
1591def printRho(SGLaue,rho,rhoMax):                         
1592    '''default doc string
1593   
1594    :param type name: description
1595   
1596    :returns: type name: description
1597   
1598    '''
1599    dim = len(rho.shape)
1600    if dim == 2:
1601        ix,jy = rho.shape
1602        for j in range(jy):
1603            line = ''
1604            if SGLaue in ['3','3m1','31m','6/m','6/mmm']:
1605                line += (jy-j)*'  '
1606            for i in range(ix):
1607                r = int(100*rho[i,j]/rhoMax)
1608                line += '%4d'%(r)
1609            print line+'\n'
1610    else:
1611        ix,jy,kz = rho.shape
1612        for k in range(kz):
1613            print 'k = ',k
1614            for j in range(jy):
1615                line = ''
1616                if SGLaue in ['3','3m1','31m','6/m','6/mmm']:
1617                    line += (jy-j)*'  '
1618                for i in range(ix):
1619                    r = int(100*rho[i,j,k]/rhoMax)
1620                    line += '%4d'%(r)
1621                print line+'\n'
1622## keep this
1623               
1624def findOffset(SGData,A,Fhkl):   
1625    '''default doc string
1626   
1627    :param type name: description
1628   
1629    :returns: type name: description
1630   
1631    '''
1632    if SGData['SpGrp'] == 'P 1':
1633        return [0,0,0]   
1634    hklShape = Fhkl.shape
1635    hklHalf = np.array(hklShape)/2
1636    sortHKL = np.argsort(Fhkl.flatten())
1637    Fdict = {}
1638    for hkl in sortHKL:
1639        HKL = np.unravel_index(hkl,hklShape)
1640        F = Fhkl[HKL[0]][HKL[1]][HKL[2]]
1641        if F == 0.:
1642            break
1643        Fdict['%.6f'%(np.absolute(F))] = hkl
1644    Flist = np.flipud(np.sort(Fdict.keys()))
1645    F = str(1.e6)
1646    i = 0
1647    DH = []
1648    Dphi = []
1649    Hmax = 2*np.asarray(G2lat.getHKLmax(3.5,SGData,A),dtype='i')
1650    for F in Flist:
1651        hkl = np.unravel_index(Fdict[F],hklShape)
1652        iabsnt,mulp,Uniq,Phi = G2spc.GenHKLf(list(hkl-hklHalf),SGData)
1653        Uniq = np.array(Uniq,dtype='i')
1654        Phi = np.array(Phi)
1655        Uniq = np.concatenate((Uniq,-Uniq))+hklHalf         # put in Friedel pairs & make as index to Farray
1656        Phi = np.concatenate((Phi,-Phi))                      # and their phase shifts
1657        Fh0 = Fhkl[hkl[0],hkl[1],hkl[2]]
1658        ang0 = np.angle(Fh0,deg=True)/360.
1659        for H,phi in zip(Uniq,Phi)[1:]:
1660            ang = (np.angle(Fhkl[H[0],H[1],H[2]],deg=True)/360.-phi)
1661            dH = H-hkl
1662            dang = ang-ang0
1663            if np.any(np.abs(dH)- Hmax > 0):    #keep low order DHs
1664                continue
1665            DH.append(dH)
1666            Dphi.append((dang+.5) % 1.0)
1667        if i > 20 or len(DH) > 30:
1668            break
1669        i += 1
1670    DH = np.array(DH)
1671    print ' map offset no.of terms: %d from %d reflections'%(len(DH),len(Flist))
1672    Dphi = np.array(Dphi)
1673    steps = np.array(hklShape)
1674    X,Y,Z = np.mgrid[0:1:1./steps[0],0:1:1./steps[1],0:1:1./steps[2]]
1675    XYZ = np.array(zip(X.flatten(),Y.flatten(),Z.flatten()))
1676    Dang = (np.dot(XYZ,DH.T)+.5)%1.-Dphi
1677    Mmap = np.reshape(np.sum((Dang)**2,axis=1),newshape=steps)/len(DH)
1678    hist,bins = np.histogram(Mmap,bins=1000)
1679#    for i,item in enumerate(hist[:10]):
1680#        print item,bins[i]
1681    chisq = np.min(Mmap)
1682    DX = -np.array(np.unravel_index(np.argmin(Mmap),Mmap.shape))
1683    print ' map offset chi**2: %.3f, map offset: %d %d %d'%(chisq,DX[0],DX[1],DX[2])
1684#    print (np.dot(DX,DH.T)+.5)%1.-Dphi
1685    return DX
1686   
1687def ChargeFlip(data,reflDict,pgbar):
1688    '''default doc string
1689   
1690    :param type name: description
1691   
1692    :returns: type name: description
1693   
1694    '''
1695    generalData = data['General']
1696    mapData = generalData['Map']
1697    flipData = generalData['Flip']
1698    FFtable = {}
1699    if 'None' not in flipData['Norm element']:
1700        normElem = flipData['Norm element'].upper()
1701        FFs = G2el.GetFormFactorCoeff(normElem.split('+')[0].split('-')[0])
1702        for ff in FFs:
1703            if ff['Symbol'] == normElem:
1704                FFtable.update(ff)
1705    dmin = flipData['Resolution']
1706    SGData = generalData['SGData']
1707    SGMT = np.array([ops[0].T for ops in SGData['SGOps']])
1708    SGT = np.array([ops[1] for ops in SGData['SGOps']])
1709    cell = generalData['Cell'][1:8]       
1710    A = G2lat.cell2A(cell[:6])
1711    Vol = cell[6]
1712    Hmax = np.asarray(G2lat.getHKLmax(dmin,SGData,A),dtype='i')+1
1713    adjHKLmax(SGData,Hmax)
1714    Ehkl = np.zeros(shape=2*Hmax,dtype='c16')       #2X64bits per complex no.
1715    time0 = time.time()
1716    for iref,ref in enumerate(reflDict['RefList']):
1717        dsp = ref[4]
1718        if dsp >= dmin:
1719            ff = 0.1*Vol    #est. no. atoms for ~10A**3/atom
1720            if FFtable:
1721                SQ = 0.25/dsp**2
1722                ff *= G2el.ScatFac(FFtable,SQ)[0]
1723            if ref[8] > 0.:         #use only +ve Fobs**2
1724                E = np.sqrt(ref[8])/ff
1725            else:
1726                E = 0.
1727            ph = ref[10]
1728            ph = rn.uniform(0.,360.)
1729            Uniq = np.inner(ref[:3],SGMT)
1730            Phi = np.inner(ref[:3],SGT)
1731            for i,hkl in enumerate(Uniq):        #uses uniq
1732                hkl = np.asarray(hkl,dtype='i')
1733                dp = 360.*Phi[i]                #and phi
1734                a = cosd(ph+dp)
1735                b = sind(ph+dp)
1736                phasep = complex(a,b)
1737                phasem = complex(a,-b)
1738                h,k,l = hkl+Hmax
1739                Ehkl[h,k,l] = E*phasep
1740                h,k,l = -hkl+Hmax       #Friedel pair refl.
1741                Ehkl[h,k,l] = E*phasem
1742#    Ehkl[Hmax] = 0.00001           #this to preserve F[0,0,0]
1743    CEhkl = copy.copy(Ehkl)
1744    MEhkl = ma.array(Ehkl,mask=(Ehkl==0.0))
1745    Emask = ma.getmask(MEhkl)
1746    sumE = np.sum(ma.array(np.absolute(CEhkl),mask=Emask))
1747    Ncyc = 0
1748    old = np.seterr(all='raise')
1749    while True:       
1750        CErho = np.real(fft.fftn(fft.fftshift(CEhkl)))*(1.+0j)
1751        CEsig = np.std(CErho)
1752        CFrho = np.where(np.real(CErho) >= flipData['k-factor']*CEsig,CErho,-CErho)
1753        CFrho = np.where(np.real(CErho) <= flipData['k-Max']*CEsig,CFrho,-CFrho)      #solves U atom problem!
1754        CFhkl = fft.ifftshift(fft.ifftn(CFrho))
1755        CFhkl = np.where(CFhkl,CFhkl,1.0)           #avoid divide by zero
1756        phase = CFhkl/np.absolute(CFhkl)
1757        CEhkl = np.absolute(Ehkl)*phase
1758        Ncyc += 1
1759        sumCF = np.sum(ma.array(np.absolute(CFhkl),mask=Emask))
1760        DEhkl = np.absolute(np.absolute(Ehkl)/sumE-np.absolute(CFhkl)/sumCF)
1761        Rcf = min(100.,np.sum(ma.array(DEhkl,mask=Emask)*100.))
1762        if Rcf < 5.:
1763            break
1764        GoOn = pgbar.Update(Rcf,newmsg='%s%8.3f%s\n%s %d'%('Residual Rcf =',Rcf,'%','No.cycles = ',Ncyc))[0]
1765        if not GoOn or Ncyc > 10000:
1766            break
1767    np.seterr(**old)
1768    print ' Charge flip time: %.4f'%(time.time()-time0),'no. elements: %d'%(Ehkl.size)
1769    CErho = np.real(fft.fftn(fft.fftshift(CEhkl)))
1770    print ' No.cycles = ',Ncyc,'Residual Rcf =%8.3f%s'%(Rcf,'%')+' Map size:',CErho.shape
1771    roll = findOffset(SGData,A,CEhkl)               #CEhkl needs to be just the observed set, not the full set!
1772       
1773    mapData['Rcf'] = Rcf
1774    mapData['rho'] = np.roll(np.roll(np.roll(CErho,roll[0],axis=0),roll[1],axis=1),roll[2],axis=2)
1775    mapData['rhoMax'] = max(np.max(mapData['rho']),-np.min(mapData['rho']))
1776    mapData['rollMap'] = [0,0,0]
1777    return mapData
1778   
1779def SearchMap(data):
1780    '''Does a search of a density map for peaks meeting the criterion of peak
1781    height is greater than mapData['cutOff']/100 of mapData['rhoMax'] where
1782    mapData is data['General']['mapData']; the map is also in mapData.
1783
1784    :param data: the phase data structure
1785
1786    :returns: (peaks,mags,dzeros) where
1787
1788        * peaks : ndarray
1789          x,y,z positions of the peaks found in the map
1790        * mags : ndarray
1791          the magnitudes of the peaks
1792        * dzeros : ndarray
1793          the distance of the peaks from  the unit cell origin
1794
1795    '''       
1796    rollMap = lambda rho,roll: np.roll(np.roll(np.roll(rho,roll[0],axis=0),roll[1],axis=1),roll[2],axis=2)
1797   
1798    norm = 1./(np.sqrt(3.)*np.sqrt(2.*np.pi)**3)
1799   
1800#    def noDuplicate(xyz,peaks,Amat):
1801#        XYZ = np.inner(Amat,xyz)
1802#        if True in [np.allclose(XYZ,np.inner(Amat,peak),atol=0.5) for peak in peaks]:
1803#            print ' Peak',xyz,' <0.5A from another peak'
1804#            return False
1805#        return True
1806#                           
1807    def fixSpecialPos(xyz,SGData,Amat):
1808        equivs = G2spc.GenAtom(xyz,SGData,Move=True)
1809        X = []
1810        xyzs = [equiv[0] for equiv in equivs]
1811        for x in xyzs:
1812            if np.sqrt(np.sum(np.inner(Amat,xyz-x)**2,axis=0))<0.5:
1813                X.append(x)
1814        if len(X) > 1:
1815            return np.average(X,axis=0)
1816        else:
1817            return xyz
1818       
1819    def rhoCalc(parms,rX,rY,rZ,res,SGLaue):
1820        Mag,x0,y0,z0,sig = parms
1821        z = -((x0-rX)**2+(y0-rY)**2+(z0-rZ)**2)/(2.*sig**2)
1822#        return norm*Mag*np.exp(z)/(sig*res**3)     #not slower but some faults in LS
1823        return norm*Mag*(1.+z+z**2/2.)/(sig*res**3)
1824       
1825    def peakFunc(parms,rX,rY,rZ,rho,res,SGLaue):
1826        Mag,x0,y0,z0,sig = parms
1827        M = rho-rhoCalc(parms,rX,rY,rZ,res,SGLaue)
1828        return M
1829       
1830    def peakHess(parms,rX,rY,rZ,rho,res,SGLaue):
1831        Mag,x0,y0,z0,sig = parms
1832        dMdv = np.zeros(([5,]+list(rX.shape)))
1833        delt = .01
1834        for i in range(5):
1835            parms[i] -= delt
1836            rhoCm = rhoCalc(parms,rX,rY,rZ,res,SGLaue)
1837            parms[i] += 2.*delt
1838            rhoCp = rhoCalc(parms,rX,rY,rZ,res,SGLaue)
1839            parms[i] -= delt
1840            dMdv[i] = (rhoCp-rhoCm)/(2.*delt)
1841        rhoC = rhoCalc(parms,rX,rY,rZ,res,SGLaue)
1842        Vec = np.sum(np.sum(np.sum(dMdv*(rho-rhoC),axis=3),axis=2),axis=1)
1843        dMdv = np.reshape(dMdv,(5,rX.size))
1844        Hess = np.inner(dMdv,dMdv)
1845       
1846        return Vec,Hess
1847       
1848    generalData = data['General']
1849    phaseName = generalData['Name']
1850    SGData = generalData['SGData']
1851    Amat,Bmat = G2lat.cell2AB(generalData['Cell'][1:7])
1852    drawingData = data['Drawing']
1853    peaks = []
1854    mags = []
1855    dzeros = []
1856    try:
1857        mapData = generalData['Map']
1858        contLevel = mapData['cutOff']*mapData['rhoMax']/100.
1859        rho = copy.copy(mapData['rho'])     #don't mess up original
1860        mapHalf = np.array(rho.shape)/2
1861        res = mapData['Resolution']
1862        incre = np.array(rho.shape,dtype=np.float)
1863        step = max(1.0,1./res)+1
1864        steps = np.array(3*[step,])
1865    except KeyError:
1866        print '**** ERROR - Fourier map not defined'
1867        return peaks,mags
1868    rhoMask = ma.array(rho,mask=(rho<contLevel))
1869    indices = (-1,0,1)
1870    rolls = np.array([[h,k,l] for h in indices for k in indices for l in indices])
1871    for roll in rolls:
1872        if np.any(roll):
1873            rhoMask = ma.array(rhoMask,mask=(rhoMask-rollMap(rho,roll)<=0.))
1874    indx = np.transpose(rhoMask.nonzero())
1875    peaks = indx/incre
1876    mags = rhoMask[rhoMask.nonzero()]
1877    for i,[ind,peak,mag] in enumerate(zip(indx,peaks,mags)):
1878        rho = rollMap(rho,ind)
1879        rMM = mapHalf-steps
1880        rMP = mapHalf+steps+1
1881        rhoPeak = rho[rMM[0]:rMP[0],rMM[1]:rMP[1],rMM[2]:rMP[2]]
1882        peakInt = np.sum(rhoPeak)*res**3
1883        rX,rY,rZ = np.mgrid[rMM[0]:rMP[0],rMM[1]:rMP[1],rMM[2]:rMP[2]]
1884        x0 = [peakInt,mapHalf[0],mapHalf[1],mapHalf[2],2.0]          #magnitude, position & width(sig)
1885        result = HessianLSQ(peakFunc,x0,Hess=peakHess,
1886            args=(rX,rY,rZ,rhoPeak,res,SGData['SGLaue']),ftol=.01,maxcyc=10)
1887        x1 = result[0]
1888        if not np.any(x1 < 0):
1889            mag = x1[0]
1890            peak = (np.array(x1[1:4])-ind)/incre
1891        peak = fixSpecialPos(peak,SGData,Amat)
1892        rho = rollMap(rho,-ind)       
1893    dzeros = np.sqrt(np.sum(np.inner(Amat,peaks)**2,axis=0))
1894    return np.array(peaks),np.array([mags,]).T,np.array([dzeros,]).T
1895   
1896def sortArray(data,pos,reverse=False):
1897    '''data is a list of items
1898    sort by pos in list; reverse if True
1899    '''
1900    T = []
1901    for i,M in enumerate(data):
1902        T.append((M[pos],i))
1903    D = dict(zip(T,data))
1904    T.sort()
1905    if reverse:
1906        T.reverse()
1907    X = []
1908    for key in T:
1909        X.append(D[key])
1910    return X
1911
1912def PeaksEquiv(data,Ind):
1913    '''Find the equivalent map peaks for those selected. Works on the
1914    contents of data['Map Peaks'].
1915
1916    :param data: the phase data structure
1917    :param list Ind: list of selected peak indices
1918    :returns: augmented list of peaks including those related by symmetry to the
1919      ones in Ind
1920
1921    '''       
1922    def Duplicate(xyz,peaks,Amat):
1923        if True in [np.allclose(np.inner(Amat,xyz),np.inner(Amat,peak),atol=0.5) for peak in peaks]:
1924            return True
1925        return False
1926                           
1927    generalData = data['General']
1928    cell = generalData['Cell'][1:7]
1929    Amat,Bmat = G2lat.cell2AB(generalData['Cell'][1:7])
1930    A = G2lat.cell2A(cell)
1931    SGData = generalData['SGData']
1932    mapPeaks = data['Map Peaks']
1933    XYZ = np.array([xyz[1:4] for xyz in mapPeaks])
1934    Indx = {}
1935    for ind in Ind:
1936        xyz = np.array(mapPeaks[ind][1:4])
1937        xyzs = np.array([equiv[0] for equiv in G2spc.GenAtom(xyz,SGData,Move=True)])
1938        for jnd,xyz in enumerate(XYZ):       
1939            Indx[jnd] = Duplicate(xyz,xyzs,Amat)
1940    Ind = []
1941    for ind in Indx:
1942        if Indx[ind]:
1943            Ind.append(ind)
1944    return Ind
1945               
1946def PeaksUnique(data,Ind):
1947    '''Finds the symmetry unique set of peaks from those selected. Works on the
1948    contents of data['Map Peaks'].
1949
1950    :param data: the phase data structure
1951    :param list Ind: list of selected peak indices
1952    :returns: the list of symmetry unique peaks from among those given in Ind
1953
1954    '''       
1955#    XYZE = np.array([[equiv[0] for equiv in G2spc.GenAtom(xyz[1:4],SGData,Move=True)] for xyz in mapPeaks]) #keep this!!
1956
1957    def noDuplicate(xyz,peaks,Amat):
1958        if True in [np.allclose(np.inner(Amat,xyz),np.inner(Amat,peak),atol=1.0) for peak in peaks]:
1959            return False
1960        return True
1961                           
1962    generalData = data['General']
1963    cell = generalData['Cell'][1:7]
1964    Amat,Bmat = G2lat.cell2AB(generalData['Cell'][1:7])
1965    A = G2lat.cell2A(cell)
1966    SGData = generalData['SGData']
1967    mapPeaks = data['Map Peaks']
1968    Indx = {}
1969    XYZ = {}
1970    for ind in Ind:
1971        XYZ[ind] = np.array(mapPeaks[ind][1:4])
1972        Indx[ind] = True
1973    for ind in Ind:
1974        if Indx[ind]:
1975            xyz = XYZ[ind]
1976            for jnd in Ind:
1977                if ind != jnd and Indx[jnd]:                       
1978                    Equiv = G2spc.GenAtom(XYZ[jnd],SGData,Move=True)
1979                    xyzs = np.array([equiv[0] for equiv in Equiv])
1980                    Indx[jnd] = noDuplicate(xyz,xyzs,Amat)
1981    Ind = []
1982    for ind in Indx:
1983        if Indx[ind]:
1984            Ind.append(ind)
1985    return Ind
1986   
1987################################################################################
1988##### single peak fitting profile fxn stuff
1989################################################################################
1990
1991def getCWsig(ins,pos):
1992    '''default doc string
1993   
1994    :param type name: description
1995   
1996    :returns: type name: description
1997   
1998    '''
1999    tp = tand(pos/2.0)
2000    return ins['U']*tp**2+ins['V']*tp+ins['W']
2001   
2002def getCWsigDeriv(pos):
2003    '''default doc string
2004   
2005    :param type name: description
2006   
2007    :returns: type name: description
2008   
2009    '''
2010    tp = tand(pos/2.0)
2011    return tp**2,tp,1.0
2012   
2013def getCWgam(ins,pos):
2014    '''default doc string
2015   
2016    :param type name: description
2017   
2018    :returns: type name: description
2019   
2020    '''
2021    return ins['X']/cosd(pos/2.0)+ins['Y']*tand(pos/2.0)
2022   
2023def getCWgamDeriv(pos):
2024    '''default doc string
2025   
2026    :param type name: description
2027   
2028    :returns: type name: description
2029   
2030    '''
2031    return 1./cosd(pos/2.0),tand(pos/2.0)
2032   
2033def getTOFsig(ins,dsp):
2034    '''default doc string
2035   
2036    :param type name: description
2037   
2038    :returns: type name: description
2039   
2040    '''
2041    return ins['sig-0']+ins['sig-1']*dsp**2+ins['sig-q']*dsp
2042   
2043def getTOFsigDeriv(dsp):
2044    '''default doc string
2045   
2046    :param type name: description
2047   
2048    :returns: type name: description
2049   
2050    '''
2051    return 1.0,dsp**2,dsp
2052   
2053def getTOFgamma(ins,dsp):
2054    '''default doc string
2055   
2056    :param type name: description
2057   
2058    :returns: type name: description
2059   
2060    '''
2061    return ins['X']*dsp+ins['Y']*dsp**2
2062   
2063def getTOFgammaDeriv(dsp):
2064    '''default doc string
2065   
2066    :param type name: description
2067   
2068    :returns: type name: description
2069   
2070    '''
2071    return dsp,dsp**2
2072   
2073def getTOFbeta(ins,dsp):
2074    '''default doc string
2075   
2076    :param type name: description
2077   
2078    :returns: type name: description
2079   
2080    '''
2081    return ins['beta-0']+ins['beta-1']/dsp**4+ins['beta-q']/dsp
2082   
2083def getTOFbetaDeriv(dsp):
2084    '''default doc string
2085   
2086    :param type name: description
2087   
2088    :returns: type name: description
2089   
2090    '''
2091    return 1.0,1./dsp**4,1./dsp
2092   
2093def getTOFalpha(ins,dsp):
2094    '''default doc string
2095   
2096    :param type name: description
2097   
2098    :returns: type name: description
2099   
2100    '''
2101    return ins['alpha']/dsp
2102   
2103def getTOFalphaDeriv(dsp):
2104    '''default doc string
2105   
2106    :param type name: description
2107   
2108    :returns: type name: description
2109   
2110    '''
2111    return 1./dsp
2112   
2113def setPeakparms(Parms,Parms2,pos,mag,ifQ=False,useFit=False):
2114    '''default doc string
2115   
2116    :param type name: description
2117   
2118    :returns: type name: description
2119   
2120    '''
2121    ind = 0
2122    if useFit:
2123        ind = 1
2124    ins = {}
2125    if 'C' in Parms['Type'][0]:                            #CW data - TOF later in an elif
2126        for x in ['U','V','W','X','Y']:
2127            ins[x] = Parms[x][ind]
2128        if ifQ:                              #qplot - convert back to 2-theta
2129            pos = 2.0*asind(pos*wave/(4*math.pi))
2130        sig = getCWsig(ins,pos)
2131        gam = getCWgam(ins,pos)           
2132        XY = [pos,0, mag,1, sig,0, gam,0]       #default refine intensity 1st
2133    else:
2134        if ifQ:
2135            dsp = 2.*np.pi/pos
2136            pos = Parms['difC']*dsp
2137        else:
2138            dsp = pos/Parms['difC'][1]
2139        if 'Pdabc' in Parms2:
2140            for x in ['sig-0','sig-1','sig-q','X','Y']:
2141                ins[x] = Parms[x][ind]
2142            Pdabc = Parms2['Pdabc'].T
2143            alp = np.interp(dsp,Pdabc[0],Pdabc[1])
2144            bet = np.interp(dsp,Pdabc[0],Pdabc[2])
2145        else:
2146            for x in ['alpha','beta-0','beta-1','beta-q','sig-0','sig-1','sig-q','X','Y']:
2147                ins[x] = Parms[x][ind]
2148            alp = getTOFalpha(ins,dsp)
2149            bet = getTOFbeta(ins,dsp)
2150        sig = getTOFsig(ins,dsp)
2151        gam = getTOFgamma(ins,dsp)
2152        XY = [pos,0,mag,1,alp,0,bet,0,sig,0,gam,0]
2153    return XY
2154   
2155################################################################################
2156##### MC/SA stuff
2157################################################################################
2158
2159#scipy/optimize/anneal.py code modified by R. Von Dreele 2013
2160# Original Author: Travis Oliphant 2002
2161# Bug-fixes in 2006 by Tim Leslie
2162
2163
2164import numpy
2165from numpy import asarray, tan, exp, ones, squeeze, sign, \
2166     all, log, sqrt, pi, shape, array, minimum, where
2167from numpy import random
2168
2169#__all__ = ['anneal']
2170
2171_double_min = numpy.finfo(float).min
2172_double_max = numpy.finfo(float).max
2173class base_schedule(object):
2174    def __init__(self):
2175        self.dwell = 20
2176        self.learn_rate = 0.5
2177        self.lower = -10
2178        self.upper = 10
2179        self.Ninit = 50
2180        self.accepted = 0
2181        self.tests = 0
2182        self.feval = 0
2183        self.k = 0
2184        self.T = None
2185
2186    def init(self, **options):
2187        self.__dict__.update(options)
2188        self.lower = asarray(self.lower)
2189        self.lower = where(self.lower == numpy.NINF, -_double_max, self.lower)
2190        self.upper = asarray(self.upper)
2191        self.upper = where(self.upper == numpy.PINF, _double_max, self.upper)
2192        self.k = 0
2193        self.accepted = 0
2194        self.feval = 0
2195        self.tests = 0
2196
2197    def getstart_temp(self, best_state):
2198        """ Find a matching starting temperature and starting parameters vector
2199        i.e. find x0 such that func(x0) = T0.
2200
2201        Parameters
2202        ----------
2203        best_state : _state
2204            A _state object to store the function value and x0 found.
2205
2206        returns
2207        -------
2208        x0 : array
2209            The starting parameters vector.
2210        """
2211
2212        assert(not self.dims is None)
2213        lrange = self.lower
2214        urange = self.upper
2215        fmax = _double_min
2216        fmin = _double_max
2217        for _ in range(self.Ninit):
2218            x0 = random.uniform(size=self.dims)*(urange-lrange) + lrange
2219            fval = self.func(x0, *self.args)
2220            self.feval += 1
2221            if fval > fmax:
2222                fmax = fval
2223            if fval < fmin:
2224                fmin = fval
2225                best_state.cost = fval
2226                best_state.x = array(x0)
2227
2228        self.T0 = (fmax-fmin)*1.5
2229        return best_state.x
2230       
2231    def set_range(self,x0,frac):
2232        delrange = frac*(self.upper-self.lower)
2233        self.upper = x0+delrange
2234        self.lower = x0-delrange
2235
2236    def accept_test(self, dE):
2237        T = self.T
2238        self.tests += 1
2239        if dE < 0:
2240            self.accepted += 1
2241            return 1
2242        p = exp(-dE*1.0/self.boltzmann/T)
2243        if (p > random.uniform(0.0, 1.0)):
2244            self.accepted += 1
2245            return 1
2246        return 0
2247
2248    def update_guess(self, x0):
2249        return np.squeeze(np.random.uniform(0.,1.,size=self.dims))*(self.upper-self.lower)+self.lower
2250
2251    def update_temp(self, x0):
2252        pass
2253
2254
2255#  A schedule due to Lester Ingber modified to use bounds - OK
2256class fast_sa(base_schedule):
2257    def init(self, **options):
2258        self.__dict__.update(options)
2259        if self.m is None:
2260            self.m = 1.0
2261        if self.n is None:
2262            self.n = 1.0
2263        self.c = self.m * exp(-self.n * self.quench)
2264
2265    def update_guess(self, x0):
2266        x0 = asarray(x0)
2267        u = squeeze(random.uniform(0.0, 1.0, size=self.dims))
2268        T = self.T
2269        xc = (sign(u-0.5)*T*((1+1.0/T)**abs(2*u-1)-1.0)+1.0)/2.0
2270        xnew = xc*(self.upper - self.lower)+self.lower
2271        return xnew
2272#        y = sign(u-0.5)*T*((1+1.0/T)**abs(2*u-1)-1.0)
2273#        xc = y*(self.upper - self.lower)
2274#        xnew = x0 + xc
2275#        return xnew
2276
2277    def update_temp(self):
2278        self.T = self.T0*exp(-self.c * self.k**(self.quench))
2279        self.k += 1
2280        return
2281
2282class cauchy_sa(base_schedule):     #modified to use bounds - not good
2283    def update_guess(self, x0):
2284        x0 = asarray(x0)
2285        numbers = squeeze(random.uniform(-pi/4, pi/4, size=self.dims))
2286        xc = (1.+(self.learn_rate * self.T * tan(numbers))%1.)
2287        xnew = xc*(self.upper - self.lower)+self.lower
2288        return xnew
2289#        numbers = squeeze(random.uniform(-pi/2, pi/2, size=self.dims))
2290#        xc = self.learn_rate * self.T * tan(numbers)
2291#        xnew = x0 + xc
2292#        return xnew
2293
2294    def update_temp(self):
2295        self.T = self.T0/(1+self.k)
2296        self.k += 1
2297        return
2298
2299class boltzmann_sa(base_schedule):
2300#    def update_guess(self, x0):
2301#        std = minimum(sqrt(self.T)*ones(self.dims), (self.upper-self.lower)/3.0/self.learn_rate)
2302#        x0 = asarray(x0)
2303#        xc = squeeze(random.normal(0, 1.0, size=self.dims))
2304#
2305#        xnew = x0 + xc*std*self.learn_rate
2306#        return xnew
2307
2308    def update_temp(self):
2309        self.k += 1
2310        self.T = self.T0 / log(self.k+1.0)
2311        return
2312
2313class log_sa(base_schedule):        #OK
2314
2315    def init(self,**options):
2316        self.__dict__.update(options)
2317       
2318    def update_guess(self,x0):     #same as default
2319        return np.squeeze(np.random.uniform(0.,1.,size=self.dims))*(self.upper-self.lower)+self.lower
2320       
2321    def update_temp(self):
2322        self.k += 1
2323        self.T = self.T0*self.slope**self.k
2324       
2325class _state(object):
2326    def __init__(self):
2327        self.x = None
2328        self.cost = None
2329
2330# TODO:
2331#     allow for general annealing temperature profile
2332#     in that case use update given by alpha and omega and
2333#     variation of all previous updates and temperature?
2334
2335# Simulated annealing
2336
2337def anneal(func, x0, args=(), schedule='fast', full_output=0,
2338           T0=None, Tf=1e-12, maxeval=None, maxaccept=None, maxiter=400,
2339           boltzmann=1.0, learn_rate=0.5, feps=1e-6, quench=1.0, m=1.0, n=1.0,
2340           lower=-100, upper=100, dwell=50, slope=0.9,ranStart=False,
2341           ranRange=0.10,autoRan=False,dlg=None):
2342    """Minimize a function using simulated annealing.
2343
2344    Schedule is a schedule class implementing the annealing schedule.
2345    Available ones are 'fast', 'cauchy', 'boltzmann'
2346
2347    :param callable func: f(x, \*args)
2348        Function to be optimized.
2349    :param ndarray x0:
2350        Initial guess.
2351    :param tuple args:
2352        Extra parameters to `func`.
2353    :param base_schedule schedule:
2354        Annealing schedule to use (a class).
2355    :param bool full_output:
2356        Whether to return optional outputs.
2357    :param float T0:
2358        Initial Temperature (estimated as 1.2 times the largest
2359        cost-function deviation over random points in the range).
2360    :param float Tf:
2361        Final goal temperature.
2362    :param int maxeval:
2363        Maximum function evaluations.
2364    :param int maxaccept:
2365        Maximum changes to accept.
2366    :param int maxiter:
2367        Maximum cooling iterations.
2368    :param float learn_rate:
2369        Scale constant for adjusting guesses.
2370    :param float boltzmann:
2371        Boltzmann constant in acceptance test
2372        (increase for less stringent test at each temperature).
2373    :param float feps:
2374        Stopping relative error tolerance for the function value in
2375        last four coolings.
2376    :param float quench,m,n:
2377        Parameters to alter fast_sa schedule.
2378    :param float/ndarray lower,upper:
2379        Lower and upper bounds on `x`.
2380    :param int dwell:
2381        The number of times to search the space at each temperature.
2382    :param float slope:
2383        Parameter for log schedule
2384    :param bool ranStart=False:
2385        True for set 10% of ranges about x
2386
2387    :returns: (xmin, Jmin, T, feval, iters, accept, retval) where
2388
2389     * xmin (ndarray): Point giving smallest value found.
2390     * Jmin (float): Minimum value of function found.
2391     * T (float): Final temperature.
2392     * feval (int): Number of function evaluations.
2393     * iters (int): Number of cooling iterations.
2394     * accept (int): Number of tests accepted.
2395     * retval (int): Flag indicating stopping condition:
2396
2397              *  0: Points no longer changing
2398              *  1: Cooled to final temperature
2399              *  2: Maximum function evaluations
2400              *  3: Maximum cooling iterations reached
2401              *  4: Maximum accepted query locations reached
2402              *  5: Final point not the minimum amongst encountered points
2403
2404    *Notes*:
2405    Simulated annealing is a random algorithm which uses no derivative
2406    information from the function being optimized. In practice it has
2407    been more useful in discrete optimization than continuous
2408    optimization, as there are usually better algorithms for continuous
2409    optimization problems.
2410
2411    Some experimentation by trying the difference temperature
2412    schedules and altering their parameters is likely required to
2413    obtain good performance.
2414
2415    The randomness in the algorithm comes from random sampling in numpy.
2416    To obtain the same results you can call numpy.random.seed with the
2417    same seed immediately before calling scipy.optimize.anneal.
2418
2419    We give a brief description of how the three temperature schedules
2420    generate new points and vary their temperature. Temperatures are
2421    only updated with iterations in the outer loop. The inner loop is
2422    over xrange(dwell), and new points are generated for
2423    every iteration in the inner loop. (Though whether the proposed
2424    new points are accepted is probabilistic.)
2425
2426    For readability, let d denote the dimension of the inputs to func.
2427    Also, let x_old denote the previous state, and k denote the
2428    iteration number of the outer loop. All other variables not
2429    defined below are input variables to scipy.optimize.anneal itself.
2430
2431    In the 'fast' schedule the updates are ::
2432
2433        u ~ Uniform(0, 1, size=d)
2434        y = sgn(u - 0.5) * T * ((1+ 1/T)**abs(2u-1) -1.0)
2435        xc = y * (upper - lower)
2436        x_new = x_old + xc
2437
2438        c = n * exp(-n * quench)
2439        T_new = T0 * exp(-c * k**quench)
2440
2441
2442    In the 'cauchy' schedule the updates are ::
2443
2444        u ~ Uniform(-pi/2, pi/2, size=d)
2445        xc = learn_rate * T * tan(u)
2446        x_new = x_old + xc
2447
2448        T_new = T0 / (1+k)
2449
2450    In the 'boltzmann' schedule the updates are ::
2451
2452        std = minimum( sqrt(T) * ones(d), (upper-lower) / (3*learn_rate) )
2453        y ~ Normal(0, std, size=d)
2454        x_new = x_old + learn_rate * y
2455
2456        T_new = T0 / log(1+k)
2457
2458    """
2459    x0 = asarray(x0)
2460    lower = asarray(lower)
2461    upper = asarray(upper)
2462
2463    schedule = eval(schedule+'_sa()')
2464    #   initialize the schedule
2465    schedule.init(dims=shape(x0),func=func,args=args,boltzmann=boltzmann,T0=T0,
2466                  learn_rate=learn_rate, lower=lower, upper=upper,
2467                  m=m, n=n, quench=quench, dwell=dwell, slope=slope)
2468
2469    current_state, last_state, best_state = _state(), _state(), _state()
2470    if ranStart:
2471        schedule.set_range(x0,ranRange)
2472    if T0 is None:
2473        x0 = schedule.getstart_temp(best_state)
2474    else:
2475        x0 = random.uniform(size=len(x0))*(upper-lower) + lower
2476        best_state.x = None
2477        best_state.cost = numpy.Inf
2478
2479    last_state.x = asarray(x0).copy()
2480    fval = func(x0,*args)
2481    schedule.feval += 1
2482    last_state.cost = fval
2483    if last_state.cost < best_state.cost:
2484        best_state.cost = fval
2485        best_state.x = asarray(x0).copy()
2486    schedule.T = schedule.T0
2487    fqueue = [100, 300, 500, 700]
2488    iters = 1
2489    keepGoing = True
2490    bestn = 0
2491    while keepGoing:
2492        retval = 0
2493        for n in xrange(dwell):
2494            current_state.x = schedule.update_guess(last_state.x)
2495            current_state.cost = func(current_state.x,*args)
2496            schedule.feval += 1
2497
2498            dE = current_state.cost - last_state.cost
2499            if schedule.accept_test(dE):
2500                last_state.x = current_state.x.copy()
2501                last_state.cost = current_state.cost
2502                if last_state.cost < best_state.cost:
2503                    best_state.x = last_state.x.copy()
2504                    best_state.cost = last_state.cost
2505                    bestn = n
2506                    if best_state.cost < 1.0 and autoRan:
2507                        schedule.set_range(x0,best_state.cost/2.)                       
2508        if dlg:
2509            GoOn = dlg.Update(min(100.,best_state.cost*100),
2510                newmsg='%s%8.5f, %s%d\n%s%8.4f%s'%('Temperature =',schedule.T, \
2511                    'Best trial:',bestn,  \
2512                    'MC/SA Residual:',best_state.cost*100,'%', \
2513                    ))[0]
2514            if not GoOn:
2515                best_state.x = last_state.x.copy()
2516                best_state.cost = last_state.cost
2517                retval = 5
2518        schedule.update_temp()
2519        iters += 1
2520        # Stopping conditions
2521        # 0) last saved values of f from each cooling step
2522        #     are all very similar (effectively cooled)
2523        # 1) Tf is set and we are below it
2524        # 2) maxeval is set and we are past it
2525        # 3) maxiter is set and we are past it
2526        # 4) maxaccept is set and we are past it
2527        # 5) user canceled run via progress bar
2528
2529        fqueue.append(squeeze(last_state.cost))
2530        fqueue.pop(0)
2531        af = asarray(fqueue)*1.0
2532        if retval == 5:
2533            print ' User terminated run; incomplete MC/SA'
2534            keepGoing = False
2535            break
2536        if all(abs((af-af[0])/af[0]) < feps):
2537            retval = 0
2538            if abs(af[-1]-best_state.cost) > feps*10:
2539                retval = 5
2540#                print "Warning: Cooled to %f at %s but this is not" \
2541#                      % (squeeze(last_state.cost), str(squeeze(last_state.x))) \
2542#                      + " the smallest point found."
2543            break
2544        if (Tf is not None) and (schedule.T < Tf):
2545            retval = 1
2546            break
2547        if (maxeval is not None) and (schedule.feval > maxeval):
2548            retval = 2
2549            break
2550        if (iters > maxiter):
2551            print "Warning: Maximum number of iterations exceeded."
2552            retval = 3
2553            break
2554        if (maxaccept is not None) and (schedule.accepted > maxaccept):
2555            retval = 4
2556            break
2557
2558    if full_output:
2559        return best_state.x, best_state.cost, schedule.T, \
2560               schedule.feval, iters, schedule.accepted, retval
2561    else:
2562        return best_state.x, retval
2563
2564def worker(iCyc,data,RBdata,reflType,reflData,covData,out_q,nprocess=-1):
2565    outlist = []
2566    random.seed(int(time.time())%100000+nprocess)   #make sure each process has a different random start
2567    for n in range(iCyc):
2568        result = mcsaSearch(data,RBdata,reflType,reflData,covData,None)
2569        outlist.append(result[0])
2570        print ' MC/SA residual: %.3f%% structure factor time: %.3f'%(100*result[0][2],result[1])
2571    out_q.put(outlist)
2572
2573def MPmcsaSearch(nCyc,data,RBdata,reflType,reflData,covData):
2574    import multiprocessing as mp
2575   
2576    nprocs = mp.cpu_count()
2577    out_q = mp.Queue()
2578    procs = []
2579    iCyc = np.zeros(nprocs)
2580    for i in range(nCyc):
2581        iCyc[i%nprocs] += 1
2582    for i in range(nprocs):
2583        p = mp.Process(target=worker,args=(int(iCyc[i]),data,RBdata,reflType,reflData,covData,out_q,i))
2584        procs.append(p)
2585        p.start()
2586    resultlist = []
2587    for i in range(nprocs):
2588        resultlist += out_q.get()
2589    for p in procs:
2590        p.join()
2591    return resultlist
2592
2593def mcsaSearch(data,RBdata,reflType,reflData,covData,pgbar):
2594    '''default doc string
2595   
2596    :param type name: description
2597   
2598    :returns: type name: description
2599   
2600    '''
2601   
2602    twopi = 2.0*np.pi
2603    global tsum
2604    tsum = 0.
2605   
2606    def getMDparms(item,pfx,parmDict,varyList):
2607        parmDict[pfx+'MDaxis'] = item['axis']
2608        parmDict[pfx+'MDval'] = item['Coef'][0]
2609        if item['Coef'][1]:
2610            varyList += [pfx+'MDval',]
2611            limits = item['Coef'][2]
2612            lower.append(limits[0])
2613            upper.append(limits[1])
2614                       
2615    def getAtomparms(item,pfx,aTypes,SGData,parmDict,varyList):
2616        parmDict[pfx+'Atype'] = item['atType']
2617        aTypes |= set([item['atType'],]) 
2618        pstr = ['Ax','Ay','Az']
2619        XYZ = [0,0,0]
2620        for i in range(3):
2621            name = pfx+pstr[i]
2622            parmDict[name] = item['Pos'][0][i]
2623            XYZ[i] = parmDict[name]
2624            if item['Pos'][1][i]:
2625                varyList += [name,]
2626                limits = item['Pos'][2][i]
2627                lower.append(limits[0])
2628                upper.append(limits[1])
2629        parmDict[pfx+'Amul'] = len(G2spc.GenAtom(XYZ,SGData))
2630           
2631    def getRBparms(item,mfx,aTypes,RBdata,SGData,atNo,parmDict,varyList):
2632        parmDict[mfx+'MolCent'] = item['MolCent']
2633        parmDict[mfx+'RBId'] = item['RBId']
2634        pstr = ['Px','Py','Pz']
2635        ostr = ['Qa','Qi','Qj','Qk']    #angle,vector not quaternion
2636        for i in range(3):
2637            name = mfx+pstr[i]
2638            parmDict[name] = item['Pos'][0][i]
2639            if item['Pos'][1][i]:
2640                varyList += [name,]
2641                limits = item['Pos'][2][i]
2642                lower.append(limits[0])
2643                upper.append(limits[1])
2644        AV = item['Ori'][0]
2645        A = AV[0]
2646        V = AV[1:]
2647        for i in range(4):
2648            name = mfx+ostr[i]
2649            if i:
2650                parmDict[name] = V[i-1]
2651            else:
2652                parmDict[name] = A
2653            if item['Ovar'] == 'AV':
2654                varyList += [name,]
2655                limits = item['Ori'][2][i]
2656                lower.append(limits[0])
2657                upper.append(limits[1])
2658            elif item['Ovar'] == 'A' and not i:
2659                varyList += [name,]
2660                limits = item['Ori'][2][i]
2661                lower.append(limits[0])
2662                upper.append(limits[1])
2663        if 'Tor' in item:      #'Tor' not there for 'Vector' RBs
2664            for i in range(len(item['Tor'][0])):
2665                name = mfx+'Tor'+str(i)
2666                parmDict[name] = item['Tor'][0][i]
2667                if item['Tor'][1][i]:
2668                    varyList += [name,]
2669                    limits = item['Tor'][2][i]
2670                    lower.append(limits[0])
2671                    upper.append(limits[1])
2672        atypes = RBdata[item['Type']][item['RBId']]['rbTypes']
2673        aTypes |= set(atypes)
2674        atNo += len(atypes)
2675        return atNo
2676               
2677    def GetAtomM(Xdata,SGData):
2678        Mdata = []
2679        for xyz in Xdata:
2680            Mdata.append(float(len(G2spc.GenAtom(xyz,SGData))))
2681        return np.array(Mdata)
2682       
2683    def GetAtomT(RBdata,parmDict):
2684        'Needs a doc string'
2685        atNo = parmDict['atNo']
2686        nfixAt = parmDict['nfixAt']
2687        Tdata = atNo*[' ',]
2688        for iatm in range(nfixAt):
2689            parm = ':'+str(iatm)+':Atype'
2690            if parm in parmDict:
2691                Tdata[iatm] = aTypes.index(parmDict[parm])
2692        iatm = nfixAt
2693        for iObj in range(parmDict['nObj']):
2694            pfx = str(iObj)+':'
2695            if parmDict[pfx+'Type'] in ['Vector','Residue']:
2696                if parmDict[pfx+'Type'] == 'Vector':
2697                    RBRes = RBdata['Vector'][parmDict[pfx+'RBId']]
2698                    nAtm = len(RBRes['rbVect'][0])
2699                else:       #Residue
2700                    RBRes = RBdata['Residue'][parmDict[pfx+'RBId']]
2701                    nAtm = len(RBRes['rbXYZ'])
2702                for i in range(nAtm):
2703                    Tdata[iatm] = aTypes.index(RBRes['rbTypes'][i])
2704                    iatm += 1
2705            elif parmDict[pfx+'Type'] == 'Atom':
2706                atNo = parmDict[pfx+'atNo']
2707                parm = pfx+'Atype'              #remove extra ':'
2708                if parm in parmDict:
2709                    Tdata[atNo] = aTypes.index(parmDict[parm])
2710                iatm += 1
2711            else:
2712                continue        #skips March Dollase
2713        return Tdata
2714       
2715    def GetAtomX(RBdata,parmDict):
2716        'Needs a doc string'
2717        Bmat = parmDict['Bmat']
2718        atNo = parmDict['atNo']
2719        nfixAt = parmDict['nfixAt']
2720        Xdata = np.zeros((3,atNo))
2721        keys = {':Ax':Xdata[0],':Ay':Xdata[1],':Az':Xdata[2]}
2722        for iatm in range(nfixAt):
2723            for key in keys:
2724                parm = ':'+str(iatm)+key
2725                if parm in parmDict:
2726                    keys[key][iatm] = parmDict[parm]
2727        iatm = nfixAt
2728        for iObj in range(parmDict['nObj']):
2729            pfx = str(iObj)+':'
2730            if parmDict[pfx+'Type'] in ['Vector','Residue']:
2731                if parmDict[pfx+'Type'] == 'Vector':
2732                    RBRes = RBdata['Vector'][parmDict[pfx+'RBId']]
2733                    vecs = RBRes['rbVect']
2734                    mags = RBRes['VectMag']
2735                    Cart = np.zeros_like(vecs[0])
2736                    for vec,mag in zip(vecs,mags):
2737                        Cart += vec*mag
2738                elif parmDict[pfx+'Type'] == 'Residue':
2739                    RBRes = RBdata['Residue'][parmDict[pfx+'RBId']]
2740                    Cart = np.array(RBRes['rbXYZ'])
2741                    for itor,seq in enumerate(RBRes['rbSeq']):
2742                        QuatA = AVdeg2Q(parmDict[pfx+'Tor'+str(itor)],Cart[seq[0]]-Cart[seq[1]])
2743                        Cart[seq[3]] = prodQVQ(QuatA,Cart[seq[3]]-Cart[seq[1]])+Cart[seq[1]]
2744                if parmDict[pfx+'MolCent'][1]:
2745                    Cart -= parmDict[pfx+'MolCent'][0]
2746                Qori = AVdeg2Q(parmDict[pfx+'Qa'],[parmDict[pfx+'Qi'],parmDict[pfx+'Qj'],parmDict[pfx+'Qk']])
2747                Pos = np.array([parmDict[pfx+'Px'],parmDict[pfx+'Py'],parmDict[pfx+'Pz']])
2748                Xdata.T[iatm:iatm+len(Cart)] = np.inner(Bmat,prodQVQ(Qori,Cart)).T+Pos
2749                iatm += len(Cart)
2750            elif parmDict[pfx+'Type'] == 'Atom':
2751                atNo = parmDict[pfx+'atNo']
2752                for key in keys:
2753                    parm = pfx+key[1:]              #remove extra ':'
2754                    if parm in parmDict:
2755                        keys[key][atNo] = parmDict[parm]
2756                iatm += 1
2757            else:
2758                continue        #skips March Dollase
2759        return Xdata.T
2760       
2761    def getAllTX(Tdata,Mdata,Xdata,SGM,SGT):
2762        allX = np.inner(Xdata,SGM)+SGT
2763        allT = np.repeat(Tdata,allX.shape[1])
2764        allM = np.repeat(Mdata,allX.shape[1])
2765        allX = np.reshape(allX,(-1,3))
2766        return allT,allM,allX
2767
2768    def getAllX(Xdata,SGM,SGT):
2769        allX = np.inner(Xdata,SGM)+SGT
2770        allX = np.reshape(allX,(-1,3))
2771        return allX
2772       
2773    def normQuaternions(RBdata,parmDict,varyList,values):
2774        for iObj in range(parmDict['nObj']):
2775            pfx = str(iObj)+':'
2776            if parmDict[pfx+'Type'] in ['Vector','Residue']:
2777                Qori = AVdeg2Q(parmDict[pfx+'Qa'],[parmDict[pfx+'Qi'],parmDict[pfx+'Qj'],parmDict[pfx+'Qk']])
2778                A,V = Q2AVdeg(Qori)
2779                for i,name in enumerate(['Qa','Qi','Qj','Qk']):
2780                    if i:
2781                        parmDict[pfx+name] = V[i-1]
2782                    else:
2783                        parmDict[pfx+name] = A
2784       
2785    def mcsaCalc(values,refList,rcov,cosTable,ifInv,allFF,RBdata,varyList,parmDict):
2786        ''' Compute structure factors for all h,k,l for phase
2787        input:
2788            refList: [ref] where each ref = h,k,l,m,d,...
2789            rcov:   array[nref,nref] covariance terms between Fo^2 values
2790            ifInv:  bool True if centrosymmetric
2791            allFF: array[nref,natoms] each value is mult*FF(H)/max(mult)
2792            RBdata: [dict] rigid body dictionary
2793            varyList: [list] names of varied parameters in MC/SA (not used here)           
2794            ParmDict: [dict] problem parameters
2795        puts result F^2 in each ref[5] in refList
2796        returns:
2797            delt-F*rcov*delt-F/sum(Fo^2)^2
2798        '''       
2799        global tsum
2800        t0 = time.time()
2801        parmDict.update(dict(zip(varyList,values)))             #update parameter tables
2802        Xdata = GetAtomX(RBdata,parmDict)                       #get new atom coords from RB
2803        allX = getAllX(Xdata,SGM,SGT)                           #fill unit cell - dups. OK
2804        MDval = parmDict['0:MDval']                             #get March-Dollase coeff
2805        HX2pi = 2.*np.pi*np.inner(allX,refList[:3].T)           #form 2piHX for every H,X pair
2806        Aterm = refList[3]*np.sum(allFF*np.cos(HX2pi),axis=0)**2    #compute real part for all H
2807        refList[5] = Aterm
2808        if not ifInv:
2809            refList[5] += refList[3]*np.sum(allFF*np.sin(HX2pi),axis=0)**2    #imaginary part for all H
2810        if len(cosTable):        #apply MD correction
2811            refList[5] *= np.sum(np.sqrt((MDval/(cosTable*(MDval**3-1.)+1.))**3),axis=1)/cosTable.shape[1]
2812        sumFcsq = np.sum(refList[5])
2813        scale = parmDict['sumFosq']/sumFcsq
2814        refList[5] *= scale
2815        refList[6] = refList[4]-refList[5]
2816        M = np.inner(refList[6],np.inner(rcov,refList[6]))
2817        tsum += (time.time()-t0)
2818        return M/np.sum(refList[4]**2)
2819
2820    sq8ln2 = np.sqrt(8*np.log(2))
2821    sq2pi = np.sqrt(2*np.pi)
2822    sq4pi = np.sqrt(4*np.pi)
2823    generalData = data['General']
2824    Amat,Bmat = G2lat.cell2AB(generalData['Cell'][1:7])
2825    Gmat,gmat = G2lat.cell2Gmat(generalData['Cell'][1:7])
2826    SGData = generalData['SGData']
2827    SGM = np.array([SGData['SGOps'][i][0] for i in range(len(SGData['SGOps']))])
2828    SGMT = np.array([SGData['SGOps'][i][0].T for i in range(len(SGData['SGOps']))])
2829    SGT = np.array([SGData['SGOps'][i][1] for i in range(len(SGData['SGOps']))])
2830    fixAtoms = data['Atoms']                       #if any
2831    cx,ct,cs = generalData['AtomPtrs'][:3]
2832    aTypes = set([])
2833    parmDict = {'Bmat':Bmat,'Gmat':Gmat}
2834    varyList = []
2835    atNo = 0
2836    for atm in fixAtoms:
2837        pfx = ':'+str(atNo)+':'
2838        parmDict[pfx+'Atype'] = atm[ct]
2839        aTypes |= set([atm[ct],])
2840        pstr = ['Ax','Ay','Az']
2841        parmDict[pfx+'Amul'] = atm[cs+1]
2842        for i in range(3):
2843            name = pfx+pstr[i]
2844            parmDict[name] = atm[cx+i]
2845        atNo += 1
2846    parmDict['nfixAt'] = len(fixAtoms)       
2847    MCSA = generalData['MCSA controls']
2848    reflName = MCSA['Data source']
2849    phaseName = generalData['Name']
2850    MCSAObjs = data['MCSA']['Models']               #list of MCSA models
2851    upper = []
2852    lower = []
2853    MDvec = np.zeros(3)
2854    for i,item in enumerate(MCSAObjs):
2855        mfx = str(i)+':'
2856        parmDict[mfx+'Type'] = item['Type']
2857        if item['Type'] == 'MD':
2858            getMDparms(item,mfx,parmDict,varyList)
2859            MDvec = np.array(item['axis'])
2860        elif item['Type'] == 'Atom':
2861            getAtomparms(item,mfx,aTypes,SGData,parmDict,varyList)
2862            parmDict[mfx+'atNo'] = atNo
2863            atNo += 1
2864        elif item['Type'] in ['Residue','Vector']:
2865            atNo = getRBparms(item,mfx,aTypes,RBdata,SGData,atNo,parmDict,varyList)
2866    parmDict['atNo'] = atNo                 #total no. of atoms
2867    parmDict['nObj'] = len(MCSAObjs)
2868    aTypes = list(aTypes)
2869    Tdata = GetAtomT(RBdata,parmDict)
2870    Xdata = GetAtomX(RBdata,parmDict)
2871    Mdata = GetAtomM(Xdata,SGData)
2872    allT,allM = getAllTX(Tdata,Mdata,Xdata,SGM,SGT)[:2]
2873    FFtables = G2el.GetFFtable(aTypes)
2874    refs = []
2875    allFF = []
2876    cosTable = []
2877    sumFosq = 0
2878    if 'PWDR' in reflName:
2879        for ref in reflData:
2880            h,k,l,m,d,pos,sig,gam,f = ref[:9]
2881            if d >= MCSA['dmin']:
2882                sig = G2pwd.getgamFW(sig,gam)/sq8ln2        #--> sig from FWHM
2883                SQ = 0.25/d**2
2884                allFF.append(allM*[G2el.getFFvalues(FFtables,SQ,True)[i] for i in allT]/np.max(allM))
2885                refs.append([h,k,l,m,f*m,pos,sig])
2886                sumFosq += f*m
2887                Heqv = np.inner(np.array([h,k,l]),SGMT)
2888                cosTable.append(G2lat.CosAngle(Heqv,MDvec,Gmat))
2889        nRef = len(refs)
2890        cosTable = np.array(cosTable)**2
2891        rcov = np.zeros((nRef,nRef))
2892        for iref,refI in enumerate(refs):
2893            rcov[iref][iref] = 1./(sq4pi*refI[6])
2894            for jref,refJ in enumerate(refs[:iref]):
2895                t1 = refI[6]**2+refJ[6]**2
2896                t2 = (refJ[5]-refI[5])**2/(2.*t1)
2897                if t2 > 10.:
2898                    rcov[iref][jref] = 0.
2899                else:
2900                    rcov[iref][jref] = 1./(sq2pi*np.sqrt(t1)*np.exp(t2))
2901        rcov += (rcov.T-np.diagflat(np.diagonal(rcov)))
2902        Rdiag = np.sqrt(np.diag(rcov))
2903        Rnorm = np.outer(Rdiag,Rdiag)
2904        rcov /= Rnorm
2905    elif 'Pawley' in reflName:  #need a bail out if Pawley cov matrix doesn't exist.
2906        vNames = []
2907        pfx = str(data['pId'])+'::PWLref:'
2908        for iref,refI in enumerate(reflData):           #Pawley reflection set
2909            h,k,l,m,d,v,f,s = refI
2910            if d >= MCSA['dmin'] and v:       #skip unrefined ones
2911                vNames.append(pfx+str(iref))
2912                SQ = 0.25/d**2
2913                allFF.append(allM*[G2el.getFFvalues(FFtables,SQ,True)[i] for i in allT]/np.max(allM))
2914                refs.append([h,k,l,m,f*m,iref,0.])
2915                sumFosq += f*m
2916                Heqv = np.inner(np.array([h,k,l]),SGMT)
2917                cosTable.append(G2lat.CosAngle(Heqv,MDvec,Gmat))
2918        cosTable = np.array(cosTable)**2
2919        nRef = len(refs)
2920        if covData['freshCOV'] and generalData['doPawley'] and MCSA.get('newDmin',True):
2921            vList = covData['varyList']
2922            covMatrix = covData['covMatrix']
2923            rcov = getVCov(vNames,vList,covMatrix)
2924            rcov += (rcov.T-np.diagflat(np.diagonal(rcov)))
2925            Rdiag = np.sqrt(np.diag(rcov))
2926            Rdiag = np.where(Rdiag,Rdiag,1.0)
2927            Rnorm = np.outer(Rdiag,Rdiag)
2928            rcov /= Rnorm
2929            MCSA['rcov'] = rcov
2930            covData['freshCOV'] = False
2931            MCSA['newDmin'] = False
2932        else:
2933            rcov = MCSA['rcov']
2934    elif 'HKLF' in reflName:
2935        for ref in reflData:
2936            [h,k,l,m,d],f = ref[:5],ref[6]
2937            if d >= MCSA['dmin']:
2938                SQ = 0.25/d**2
2939                allFF.append(allM*[G2el.getFFvalues(FFtables,SQ,True)[i] for i in allT]/np.max(allM))
2940                refs.append([h,k,l,m,f*m,0.,0.])
2941                sumFosq += f*m
2942        nRef = len(refs)
2943        rcov = np.identity(len(refs))
2944    allFF = np.array(allFF).T
2945    refs = np.array(refs).T
2946    print ' Minimum d-spacing used: %.2f No. reflections used: %d'%(MCSA['dmin'],nRef)
2947    print ' Number of parameters varied: %d'%(len(varyList))
2948    parmDict['sumFosq'] = sumFosq
2949    x0 = [parmDict[val] for val in varyList]
2950    ifInv = SGData['SGInv']
2951    results = anneal(mcsaCalc,x0,args=(refs,rcov,cosTable,ifInv,allFF,RBdata,varyList,parmDict),
2952        schedule=MCSA['Algorithm'], full_output=True,
2953        T0=MCSA['Annealing'][0], Tf=MCSA['Annealing'][1],dwell=MCSA['Annealing'][2],
2954        boltzmann=MCSA['boltzmann'], learn_rate=0.5, 
2955        quench=MCSA['fast parms'][0], m=MCSA['fast parms'][1], n=MCSA['fast parms'][2],
2956        lower=lower, upper=upper, slope=MCSA['log slope'],ranStart=MCSA.get('ranStart',False),
2957        ranRange=MCSA.get('ranRange',0.10),autoRan=MCSA.get('autoRan',False),dlg=pgbar)
2958    M = mcsaCalc(results[0],refs,rcov,cosTable,ifInv,allFF,RBdata,varyList,parmDict)
2959#    for ref in refs.T:
2960#        print ' %4d %4d %4d %10.3f %10.3f %10.3f'%(int(ref[0]),int(ref[1]),int(ref[2]),ref[4],ref[5],ref[6])
2961#    print np.sqrt((np.sum(refs[6]**2)/np.sum(refs[4]**2)))
2962    Result = [False,False,results[1],results[2],]+list(results[0])
2963    Result.append(varyList)
2964    return Result,tsum
2965
2966       
2967################################################################################
2968##### Quaternion stuff
2969################################################################################
2970
2971def prodQQ(QA,QB):
2972    ''' Grassman quaternion product
2973        QA,QB quaternions; q=r+ai+bj+ck
2974    '''
2975    D = np.zeros(4)
2976    D[0] = QA[0]*QB[0]-QA[1]*QB[1]-QA[2]*QB[2]-QA[3]*QB[3]
2977    D[1] = QA[0]*QB[1]+QA[1]*QB[0]+QA[2]*QB[3]-QA[3]*QB[2]
2978    D[2] = QA[0]*QB[2]-QA[1]*QB[3]+QA[2]*QB[0]+QA[3]*QB[1]
2979    D[3] = QA[0]*QB[3]+QA[1]*QB[2]-QA[2]*QB[1]+QA[3]*QB[0]
2980   
2981#    D[0] = QA[0]*QB[0]-np.dot(QA[1:],QB[1:])
2982#    D[1:] = QA[0]*QB[1:]+QB[0]*QA[1:]+np.cross(QA[1:],QB[1:])
2983   
2984    return D
2985   
2986def normQ(QA):
2987    ''' get length of quaternion & normalize it
2988        q=r+ai+bj+ck
2989    '''
2990    n = np.sqrt(np.sum(np.array(QA)**2))
2991    return QA/n
2992   
2993def invQ(Q):
2994    '''
2995        get inverse of quaternion
2996        q=r+ai+bj+ck; q* = r-ai-bj-ck
2997    '''
2998    return Q*np.array([1,-1,-1,-1])
2999   
3000def prodQVQ(Q,V):
3001    """
3002    compute the quaternion vector rotation qvq-1 = v'
3003    q=r+ai+bj+ck
3004    """
3005    T2 = Q[0]*Q[1]
3006    T3 = Q[0]*Q[2]
3007    T4 = Q[0]*Q[3]
3008    T5 = -Q[1]*Q[1]
3009    T6 = Q[1]*Q[2]
3010    T7 = Q[1]*Q[3]
3011    T8 = -Q[2]*Q[2]
3012    T9 = Q[2]*Q[3]
3013    T10 = -Q[3]*Q[3]
3014    M = np.array([[T8+T10,T6-T4,T3+T7],[T4+T6,T5+T10,T9-T2],[T7-T3,T2+T9,T5+T8]])
3015    VP = 2.*np.inner(V,M)
3016    return VP+V
3017   
3018def Q2Mat(Q):
3019    ''' make rotation matrix from quaternion
3020        q=r+ai+bj+ck
3021    '''
3022    QN = normQ(Q)
3023    aa = QN[0]**2
3024    ab = QN[0]*QN[1]
3025    ac = QN[0]*QN[2]
3026    ad = QN[0]*QN[3]
3027    bb = QN[1]**2
3028    bc = QN[1]*QN[2]
3029    bd = QN[1]*QN[3]
3030    cc = QN[2]**2
3031    cd = QN[2]*QN[3]
3032    dd = QN[3]**2
3033    M = [[aa+bb-cc-dd, 2.*(bc-ad),  2.*(ac+bd)],
3034        [2*(ad+bc),   aa-bb+cc-dd,  2.*(cd-ab)],
3035        [2*(bd-ac),    2.*(ab+cd), aa-bb-cc+dd]]
3036    return np.array(M)
3037   
3038def AV2Q(A,V):
3039    ''' convert angle (radians) & vector to quaternion
3040        q=r+ai+bj+ck
3041    '''
3042    Q = np.zeros(4)
3043    d = nl.norm(np.array(V))
3044    if d:
3045        V /= d
3046        if not A:       #==0.
3047            A = 2.*np.pi
3048        p = A/2.
3049        Q[0] = np.cos(p)
3050        Q[1:4] = V*np.sin(p)
3051    else:
3052        Q[3] = 1.
3053    return Q
3054   
3055def AVdeg2Q(A,V):
3056    ''' convert angle (degrees) & vector to quaternion
3057        q=r+ai+bj+ck
3058    '''
3059    Q = np.zeros(4)
3060    d = nl.norm(np.array(V))
3061    if not A:       #== 0.!
3062        A = 360.
3063    if d:
3064        V /= d
3065        p = A/2.
3066        Q[0] = cosd(p)
3067        Q[1:4] = V*sind(p)
3068    else:
3069        Q[3] = 1.
3070    return Q
3071   
3072def Q2AVdeg(Q):
3073    ''' convert quaternion to angle (degrees 0-360) & normalized vector
3074        q=r+ai+bj+ck
3075    '''
3076    A = 2.*acosd(Q[0])
3077    V = np.array(Q[1:])
3078    V /= sind(A/2.)
3079    return A,V
3080   
3081def Q2AV(Q):
3082    ''' convert quaternion to angle (radians 0-2pi) & normalized vector
3083        q=r+ai+bj+ck
3084    '''
3085    A = 2.*np.arccos(Q[0])
3086    V = np.array(Q[1:])
3087    V /= np.sin(A/2.)
3088    return A,V
3089   
3090def randomQ(r0,r1,r2,r3):
3091    ''' create random quaternion from 4 random numbers in range (-1,1)
3092    '''
3093    sum = 0
3094    Q = np.array(4)
3095    Q[0] = r0
3096    sum += Q[0]**2
3097    Q[1] = np.sqrt(1.-sum)*r1
3098    sum += Q[1]**2
3099    Q[2] = np.sqrt(1.-sum)*r2
3100    sum += Q[2]**2
3101    Q[3] = np.sqrt(1.-sum)*np.where(r3<0.,-1.,1.)
3102    return Q
3103   
3104def randomAVdeg(r0,r1,r2,r3):
3105    ''' create random angle (deg),vector from 4 random number in range (-1,1)
3106    '''
3107    return Q2AVdeg(randomQ(r0,r1,r2,r3))
3108   
3109def makeQuat(A,B,C):
3110    ''' Make quaternion from rotation of A vector to B vector about C axis
3111
3112        :param np.array A,B,C: Cartesian 3-vectors
3113        :returns: quaternion & rotation angle in radians q=r+ai+bj+ck
3114    '''
3115
3116    V1 = np.cross(A,C)
3117    V2 = np.cross(B,C)
3118    if nl.norm(V1)*nl.norm(V2):
3119        V1 /= nl.norm(V1)
3120        V2 /= nl.norm(V2)
3121        V3 = np.cross(V1,V2)
3122    else:
3123        V3 = np.zeros(3)
3124    Q = np.array([0.,0.,0.,1.])
3125    D = 0.
3126    if nl.norm(V3):
3127        V3 /= nl.norm(V3)
3128        D1 = min(1.0,max(-1.0,np.vdot(V1,V2)))
3129        D = np.arccos(D1)/2.0
3130        V1 = C-V3
3131        V2 = C+V3
3132        DM = nl.norm(V1)
3133        DP = nl.norm(V2)
3134        S = np.sin(D)
3135        Q[0] = np.cos(D)
3136        Q[1:] = V3*S
3137        D *= 2.
3138        if DM > DP:
3139            D *= -1.
3140    return Q,D
3141   
3142def annealtests():
3143    from numpy import cos
3144    # minimum expected at ~-0.195
3145    func = lambda x: cos(14.5*x-0.3) + (x+0.2)*x
3146    print anneal(func,1.0,full_output=1,upper=3.0,lower=-3.0,feps=1e-4,maxiter=2000,schedule='cauchy')
3147    print anneal(func,1.0,full_output=1,upper=3.0,lower=-3.0,feps=1e-4,maxiter=2000,schedule='fast')
3148    print anneal(func,1.0,full_output=1,upper=3.0,lower=-3.0,feps=1e-4,maxiter=2000,schedule='boltzmann')
3149
3150    # minimum expected at ~[-0.195, -0.1]
3151    func = lambda x: cos(14.5*x[0]-0.3) + (x[1]+0.2)*x[1] + (x[0]+0.2)*x[0]
3152    print anneal(func,[1.0, 1.0],full_output=1,upper=[3.0, 3.0],lower=[-3.0, -3.0],feps=1e-4,maxiter=2000,schedule='cauchy')
3153    print anneal(func,[1.0, 1.0],full_output=1,upper=[3.0, 3.0],lower=[-3.0, -3.0],feps=1e-4,maxiter=2000,schedule='fast')
3154    print anneal(func,[1.0, 1.0],full_output=1,upper=[3.0, 3.0],lower=[-3.0, -3.0],feps=1e-4,maxiter=2000,schedule='boltzmann')
3155
3156
3157if __name__ == '__main__':
3158    annealtests()
Note: See TracBrowser for help on using the repository browser.