source: trunk/GSASIIimage.py @ 4593

Last change on this file since 4593 was 4593, checked in by vondreele, 5 years ago

fix image plotting issues with 2-theta limits, arc & ring masks & azimuth limits with new GetDetectorXY code

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 67.4 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASII image calculations: ellipse fitting & image integration       
3########### SVN repository information ###################
4# $Date: 2020-10-13 19:16:15 +0000 (Tue, 13 Oct 2020) $
5# $Author: vondreele $
6# $Revision: 4593 $
7# $URL: trunk/GSASIIimage.py $
8# $Id: GSASIIimage.py 4593 2020-10-13 19:16:15Z vondreele $
9########### SVN repository information ###################
10'''
11*GSASIIimage: Image calc module*
12================================
13
14Ellipse fitting & image integration
15
16'''
17from __future__ import division, print_function
18import math
19import time
20import numpy as np
21import numpy.linalg as nl
22import numpy.ma as ma
23from scipy.optimize import leastsq
24import scipy.interpolate as scint
25import scipy.special as sc
26import copy
27import GSASIIpath
28GSASIIpath.SetVersionNumber("$Revision: 4593 $")
29try:
30    import GSASIIplot as G2plt
31except ImportError: # expected in scriptable w/o matplotlib and/or wx
32    pass
33import GSASIIlattice as G2lat
34import GSASIIpwd as G2pwd
35import GSASIIspc as G2spc
36#import GSASIImath as G2mth
37import GSASIIfiles as G2fil
38
39# trig functions in degrees
40sind = lambda x: math.sin(x*math.pi/180.)
41asind = lambda x: 180.*math.asin(x)/math.pi
42tand = lambda x: math.tan(x*math.pi/180.)
43atand = lambda x: 180.*math.atan(x)/math.pi
44atan2d = lambda y,x: 180.*math.atan2(y,x)/math.pi
45cosd = lambda x: math.cos(x*math.pi/180.)
46acosd = lambda x: 180.*math.acos(x)/math.pi
47rdsq2d = lambda x,p: round(1.0/math.sqrt(x),p)
48#numpy versions
49npsind = lambda x: np.sin(x*np.pi/180.)
50npasind = lambda x: 180.*np.arcsin(x)/np.pi
51npcosd = lambda x: np.cos(x*np.pi/180.)
52npacosd = lambda x: 180.*np.arccos(x)/np.pi
53nptand = lambda x: np.tan(x*np.pi/180.)
54npatand = lambda x: 180.*np.arctan(x)/np.pi
55npatan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
56nxs = np.newaxis
57debug = False
58   
59def pointInPolygon(pXY,xy):
60    'Needs a doc string'
61    #pXY - assumed closed 1st & last points are duplicates
62    Inside = False
63    N = len(pXY)
64    p1x,p1y = pXY[0]
65    for i in range(N+1):
66        p2x,p2y = pXY[i%N]
67        if (max(p1y,p2y) >= xy[1] > min(p1y,p2y)) and (xy[0] <= max(p1x,p2x)):
68            if p1y != p2y:
69                xinters = (xy[1]-p1y)*(p2x-p1x)/(p2y-p1y)+p1x
70            if p1x == p2x or xy[0] <= xinters:
71                Inside = not Inside
72        p1x,p1y = p2x,p2y
73    return Inside
74   
75def peneCorr(tth,dep,dist):
76    'Needs a doc string'
77    return dep*(1.-npcosd(tth))*dist**2/1000.         #best one
78       
79def makeMat(Angle,Axis):
80    '''Make rotation matrix from Angle and Axis
81
82    :param float Angle: in degrees
83    :param int Axis: 0 for rotation about x, 1 for about y, etc.
84    '''
85    cs = npcosd(Angle)
86    ss = npsind(Angle)
87    M = np.array(([1.,0.,0.],[0.,cs,-ss],[0.,ss,cs]),dtype=np.float32)
88    return np.roll(np.roll(M,Axis,axis=0),Axis,axis=1)
89                   
90def FitEllipse(xy):
91   
92    def ellipse_center(p):
93        ''' gives ellipse center coordinates
94        '''
95        b,c,d,f,a = p[1]/2., p[2], p[3]/2., p[4]/2., p[0]
96        num = b*b-a*c
97        x0=(c*d-b*f)/num
98        y0=(a*f-b*d)/num
99        return np.array([x0,y0])
100   
101    def ellipse_angle_of_rotation( p ):
102        ''' gives rotation of ellipse major axis from x-axis
103        range will be -90 to 90 deg
104        '''
105        b,c,a = p[1]/2., p[2], p[0]
106        return 0.5*npatand(2*b/(a-c))
107   
108    def ellipse_axis_length( p ):
109        ''' gives ellipse radii in [minor,major] order
110        '''
111        b,c,d,f,g,a = p[1]/2., p[2], p[3]/2., p[4]/2, p[5], p[0]
112        up = 2*(a*f*f+c*d*d+g*b*b-2*b*d*f-a*c*g)
113        down1=(b*b-a*c)*( (c-a)*np.sqrt(1+4*b*b/((a-c)*(a-c)))-(c+a))
114        down2=(b*b-a*c)*( (a-c)*np.sqrt(1+4*b*b/((a-c)*(a-c)))-(c+a))
115        res1=np.sqrt(up/down1)
116        res2=np.sqrt(up/down2)
117        return np.array([ res2,res1])
118   
119    xy = np.array(xy)
120    x = np.asarray(xy.T[0])[:,np.newaxis]
121    y = np.asarray(xy.T[1])[:,np.newaxis]
122    D =  np.hstack((x*x, x*y, y*y, x, y, np.ones_like(x)))
123    S = np.dot(D.T,D)
124    C = np.zeros([6,6])
125    C[0,2] = C[2,0] = 2; C[1,1] = -1
126    E, V =  nl.eig(np.dot(nl.inv(S), C))
127    n = np.argmax(np.abs(E))
128    a = V[:,n]
129    cent = ellipse_center(a)
130    phi = ellipse_angle_of_rotation(a)
131    radii = ellipse_axis_length(a)
132    phi += 90.
133    if radii[0] > radii[1]:
134        radii = [radii[1],radii[0]]
135        phi -= 90.
136    return cent,phi,radii
137
138def FitDetector(rings,varyList,parmDict,Print=True,covar=False):
139    '''Fit detector calibration parameters
140
141    :param np.array rings: vector of ring positions
142    :param list varyList: calibration parameters to be refined
143    :param dict parmDict: all calibration parameters
144    :param bool Print: set to True (default) to print the results
145    :param bool covar: set to True to return the covariance matrix (default is False)
146    :returns: [chisq,vals,sigList] unless covar is True, then
147        [chisq,vals,sigList,coVarMatrix] is returned
148    '''
149       
150    def CalibPrint(ValSig,chisq,Npts):
151        print ('Image Parameters: chi**2: %12.3g, Np: %d'%(chisq,Npts))
152        ptlbls = 'names :'
153        ptstr =  'values:'
154        sigstr = 'esds  :'
155        for name,value,sig in ValSig:
156            ptlbls += "%s" % (name.rjust(12))
157            if name == 'phi':
158                ptstr += Fmt[name] % (value%360.)
159            else:
160                ptstr += Fmt[name] % (value)
161            if sig:
162                sigstr += Fmt[name] % (sig)
163            else:
164                sigstr += 12*' '
165        print (ptlbls)
166        print (ptstr)
167        print (sigstr)       
168       
169    def ellipseCalcD(B,xyd,varyList,parmDict):
170       
171        x,y,dsp = xyd
172        varyDict = dict(zip(varyList,B))
173        parms = {}
174        for parm in parmDict:
175            if parm in varyList:
176                parms[parm] = varyDict[parm]
177            else:
178                parms[parm] = parmDict[parm]
179        phi = parms['phi']-90.               #get rotation of major axis from tilt axis
180        tth = 2.0*npasind(parms['wave']/(2.*dsp))
181        phi0 = npatan2d(y-parms['det-Y'],x-parms['det-X'])
182        dxy = peneCorr(tth,parms['dep'],parms['dist'])
183        stth = npsind(tth)
184        cosb = npcosd(parms['tilt'])
185        tanb = nptand(parms['tilt'])       
186        tbm = nptand((tth-parms['tilt'])/2.)
187        tbp = nptand((tth+parms['tilt'])/2.)
188        d = parms['dist']+dxy
189        fplus = d*tanb*stth/(cosb+stth)
190        fminus = d*tanb*stth/(cosb-stth)
191        vplus = d*(tanb+(1+tbm)/(1-tbm))*stth/(cosb+stth)
192        vminus = d*(tanb+(1-tbp)/(1+tbp))*stth/(cosb-stth)
193        R0 = np.sqrt((vplus+vminus)**2-(fplus+fminus)**2)/2.      #+minor axis
194        R1 = (vplus+vminus)/2.                                    #major axis
195        zdis = (fplus-fminus)/2.
196        Robs = np.sqrt((x-parms['det-X'])**2+(y-parms['det-Y'])**2)
197        rsqplus = R0**2+R1**2
198        rsqminus = R0**2-R1**2
199        R = rsqminus*npcosd(2.*phi0-2.*phi)+rsqplus
200        Q = np.sqrt(2.)*R0*R1*np.sqrt(R-2.*zdis**2*npsind(phi0-phi)**2)
201        P = 2.*R0**2*zdis*npcosd(phi0-phi)
202        Rcalc = (P+Q)/R
203        M = (Robs-Rcalc)*25.        #why 25? does make "chi**2" more reasonable
204        return M
205       
206    names = ['dist','det-X','det-Y','tilt','phi','dep','wave']
207    fmt = ['%12.3f','%12.3f','%12.3f','%12.3f','%12.3f','%12.4f','%12.6f']
208    Fmt = dict(zip(names,fmt))
209    p0 = [parmDict[key] for key in varyList]
210    result = leastsq(ellipseCalcD,p0,args=(rings.T,varyList,parmDict),full_output=True,ftol=1.e-8)
211    chisq = np.sum(result[2]['fvec']**2)/(rings.shape[0]-len(p0))   #reduced chi^2 = M/(Nobs-Nvar)
212    parmDict.update(zip(varyList,result[0]))
213    vals = list(result[0])
214    if not len(vals):
215        sig = []
216        ValSig = []
217        sigList = []
218    else:   
219        sig = list(np.sqrt(chisq*np.diag(result[1])))
220        sigList = np.zeros(7)
221        for i,name in enumerate(varyList):
222            sigList[i] = sig[varyList.index(name)]
223        ValSig = zip(varyList,vals,sig)
224    if Print:
225        if len(sig):
226            CalibPrint(ValSig,chisq,rings.shape[0])
227        else:
228            print(' Nothing refined')
229    if covar:
230        return [chisq,vals,sigList,result[1]]
231    else:
232        return [chisq,vals,sigList]
233
234def FitMultiDist(rings,varyList,parmDict,Print=True,covar=False):
235    '''Fit detector calibration parameters with multi-distance data
236
237    :param np.array rings: vector of ring positions (x,y,dist,d-space)
238    :param list varyList: calibration parameters to be refined
239    :param dict parmDict: calibration parameters
240    :param bool Print: set to True (default) to print the results
241    :param bool covar: set to True to return the covariance matrix (default is False)
242    :returns: [chisq,vals,sigDict] unless covar is True, then
243        [chisq,vals,sigDict,coVarMatrix] is returned
244    '''
245       
246    def CalibPrint(parmDict,sigDict,chisq,Npts):
247        ptlbls = 'names :'
248        ptstr =  'values:'
249        sigstr = 'esds  :'
250        for d in sorted(set([i[5:] for i in parmDict.keys() if 'det-X' in i]),key=lambda x:int(x)):
251            fmt = '%12.3f'
252            for key in 'det-X','det-Y','delta':
253                name = key+d
254                if name not in parmDict: continue
255                ptlbls += "%12s" % name
256                ptstr += fmt % (parmDict[name])
257                if name in sigDict:
258                    sigstr += fmt % (sigDict[name])
259                else:
260                    sigstr += 12*' '
261                if len(ptlbls) > 68:
262                    print()
263                    print (ptlbls)
264                    print (ptstr)
265                    print (sigstr)
266                    ptlbls = 'names :'
267                    ptstr =  'values:'
268                    sigstr = 'esds  :'
269        if len(ptlbls) > 8:
270            print()
271            print (ptlbls)
272            print (ptstr)
273            print (sigstr)
274        print ('\nImage Parameters: chi**2: %12.3g, Np: %d'%(chisq,Npts))
275        ptlbls = 'names :'
276        ptstr =  'values:'
277        sigstr = 'esds  :'
278        names = ['wavelength', 'dep', 'phi', 'tilt']
279        if 'deltaDist' in parmDict:
280            names += ['deltaDist']
281        for name in names:
282            if name == 'wavelength':
283                fmt = '%12.6f'
284            elif name == 'dep':
285                fmt = '%12.4f'
286            else:
287                fmt = '%12.3f'
288               
289            ptlbls += "%s" % (name.rjust(12))
290            if name == 'phi':
291                ptstr += fmt % (parmDict[name]%360.)
292            else:
293                ptstr += fmt % (parmDict[name])
294            if name in sigDict:
295                sigstr += fmt % (sigDict[name])
296            else:
297                sigstr += 12*' '
298        print (ptlbls)
299        print (ptstr)
300        print (sigstr)
301        print()
302
303    def ellipseCalcD(B,xyd,varyList,parmDict):
304        x,y,dist,dsp = xyd
305        varyDict = dict(zip(varyList,B))
306        parms = {}
307        for parm in parmDict:
308            if parm in varyList:
309                parms[parm] = varyDict[parm]
310            else:
311                parms[parm] = parmDict[parm]
312        # create arrays with detector center values
313        detX = np.array([parms['det-X'+str(int(d))] for d in dist])
314        detY = np.array([parms['det-Y'+str(int(d))] for d in dist])
315        if 'deltaDist' in parms:
316            deltaDist = parms['deltaDist']
317        else:
318            deltaDist = np.array([parms['delta'+str(int(d))] for d in dist])
319               
320        phi = parms['phi']-90.               #get rotation of major axis from tilt axis
321        tth = 2.0*npasind(parms['wavelength']/(2.*dsp))
322        phi0 = npatan2d(y-detY,x-detX)
323        dxy = peneCorr(tth,parms['dep'],dist-deltaDist)
324        stth = npsind(tth)
325        cosb = npcosd(parms['tilt'])
326        tanb = nptand(parms['tilt'])       
327        tbm = nptand((tth-parms['tilt'])/2.)
328        tbp = nptand((tth+parms['tilt'])/2.)
329        d = (dist-deltaDist)+dxy
330        fplus = d*tanb*stth/(cosb+stth)
331        fminus = d*tanb*stth/(cosb-stth)
332        vplus = d*(tanb+(1+tbm)/(1-tbm))*stth/(cosb+stth)
333        vminus = d*(tanb+(1-tbp)/(1+tbp))*stth/(cosb-stth)
334        R0 = np.sqrt((vplus+vminus)**2-(fplus+fminus)**2)/2.      #+minor axis
335        R1 = (vplus+vminus)/2.                                    #major axis
336        zdis = (fplus-fminus)/2.
337        Robs = np.sqrt((x-detX)**2+(y-detY)**2)
338        rsqplus = R0**2+R1**2
339        rsqminus = R0**2-R1**2
340        R = rsqminus*npcosd(2.*phi0-2.*phi)+rsqplus
341        Q = np.sqrt(2.)*R0*R1*np.sqrt(R-2.*zdis**2*npsind(phi0-phi)**2)
342        P = 2.*R0**2*zdis*npcosd(phi0-phi)
343        Rcalc = (P+Q)/R
344        return (Robs-Rcalc)*25.        #why 25? does make "chi**2" more reasonable
345       
346    p0 = [parmDict[key] for key in varyList]
347    result = leastsq(ellipseCalcD,p0,args=(rings.T,varyList,parmDict),full_output=True,ftol=1.e-8)
348    chisq = np.sum(result[2]['fvec']**2)/(rings.shape[0]-len(p0))   #reduced chi^2 = M/(Nobs-Nvar)
349    parmDict.update(zip(varyList,result[0]))
350    vals = list(result[0])
351    if chisq > 1:
352        sig = list(np.sqrt(chisq*np.diag(result[1])))
353    else:
354        sig = list(np.sqrt(np.diag(result[1])))
355    sigDict = {name:s for name,s in zip(varyList,sig)}
356    if Print:
357        CalibPrint(parmDict,sigDict,chisq,rings.shape[0])
358    if covar:
359        return [chisq,vals,sigDict,result[1]]
360    else:
361        return [chisq,vals,sigDict]
362   
363def ImageLocalMax(image,w,Xpix,Ypix):
364    'Needs a doc string'
365    w2 = w*2
366    sizey,sizex = image.shape
367    xpix = int(Xpix)            #get reference corner of pixel chosen
368    ypix = int(Ypix)
369    if not w:
370        ZMax = np.sum(image[ypix-2:ypix+2,xpix-2:xpix+2])
371        return xpix,ypix,ZMax,0.0001
372    if (w2 < xpix < sizex-w2) and (w2 < ypix < sizey-w2) and image[ypix,xpix]:
373        ZMax = image[ypix-w:ypix+w,xpix-w:xpix+w]
374        Zmax = np.argmax(ZMax)
375        ZMin = image[ypix-w2:ypix+w2,xpix-w2:xpix+w2]
376        Zmin = np.argmin(ZMin)
377        xpix += Zmax%w2-w
378        ypix += Zmax//w2-w
379        return xpix,ypix,np.ravel(ZMax)[Zmax],max(0.0001,np.ravel(ZMin)[Zmin])   #avoid neg/zero minimum
380    else:
381        return 0,0,0,0     
382   
383def makeRing(dsp,ellipse,pix,reject,scalex,scaley,image,mul=1):
384    'Needs a doc string'
385    def ellipseC():
386        'compute estimate of ellipse circumference'
387        if radii[0] < 0:        #hyperbola
388#            theta = npacosd(1./np.sqrt(1.+(radii[0]/radii[1])**2))
389#            print (theta)
390            return 0
391        apb = radii[1]+radii[0]
392        amb = radii[1]-radii[0]
393        return np.pi*apb*(1+3*(amb/apb)**2/(10+np.sqrt(4-3*(amb/apb)**2)))
394       
395    cent,phi,radii = ellipse
396    cphi = cosd(phi-90.)        #convert to major axis rotation
397    sphi = sind(phi-90.)
398    ring = []
399    C = int(ellipseC())*mul         #ring circumference in mm
400    azm = []
401    for i in range(0,C,1):      #step around ring in 1mm increments
402        a = 360.*i/C
403        x = radii[1]*cosd(a-phi+90.)        #major axis
404        y = radii[0]*sind(a-phi+90.)
405        X = (cphi*x-sphi*y+cent[0])*scalex      #convert mm to pixels
406        Y = (sphi*x+cphi*y+cent[1])*scaley
407        X,Y,I,J = ImageLocalMax(image,pix,X,Y)
408        if I and J and float(I)/J > reject:
409            X += .5                             #set to center of pixel
410            Y += .5
411            X /= scalex                         #convert back to mm
412            Y /= scaley
413            if [X,Y,dsp] not in ring:           #no duplicates!
414                ring.append([X,Y,dsp])
415                azm.append(a)
416    if len(ring) < 10:
417        ring = []
418        azm = []
419    return ring,azm
420   
421def GetEllipse2(tth,dxy,dist,cent,tilt,phi):
422    '''uses Dandelin spheres to find ellipse or hyperbola parameters from detector geometry
423    on output
424    radii[0] (b-minor axis) set < 0. for hyperbola
425   
426    '''
427    radii = [0,0]
428    stth = sind(tth)
429    cosb = cosd(tilt)
430    tanb = tand(tilt)
431    tbm = tand((tth-tilt)/2.)
432    tbp = tand((tth+tilt)/2.)
433    sinb = sind(tilt)
434    d = dist+dxy
435    if tth+abs(tilt) < 90.:      #ellipse
436        fplus = d*tanb*stth/(cosb+stth)
437        fminus = d*tanb*stth/(cosb-stth)
438        vplus = d*(tanb+(1+tbm)/(1-tbm))*stth/(cosb+stth)
439        vminus = d*(tanb+(1-tbp)/(1+tbp))*stth/(cosb-stth)
440        radii[0] = np.sqrt((vplus+vminus)**2-(fplus+fminus)**2)/2.      #+minor axis
441        radii[1] = (vplus+vminus)/2.                                    #major axis
442        zdis = (fplus-fminus)/2.
443    else:   #hyperbola!
444        f = d*abs(tanb)*stth/(cosb+stth)
445        v = d*(abs(tanb)+tand(tth-abs(tilt)))
446        delt = d*stth*(1.+stth*cosb)/(abs(sinb)*cosb*(stth+cosb))
447        eps = (v-f)/(delt-v)
448        radii[0] = -eps*(delt-f)/np.sqrt(eps**2-1.)                     #-minor axis
449        radii[1] = eps*(delt-f)/(eps**2-1.)                             #major axis
450        if tilt > 0:
451            zdis = f+radii[1]*eps
452        else:
453            zdis = -f
454#NB: zdis is || to major axis & phi is rotation of minor axis
455#thus shift from beam to ellipse center is [Z*sin(phi),-Z*cos(phi)]
456    elcent = [cent[0]+zdis*sind(phi),cent[1]-zdis*cosd(phi)]
457    return elcent,phi,radii
458   
459def GetEllipse(dsp,data):
460    '''uses Dandelin spheres to find ellipse or hyperbola parameters from detector geometry
461    as given in image controls dictionary (data) and a d-spacing (dsp)
462    '''
463    cent = data['center']
464    tilt = data['tilt']
465    phi = data['rotation']
466    dep = data.get('DetDepth',0.0)
467    tth = 2.0*asind(data['wavelength']/(2.*dsp))
468    dist = data['distance']
469    dxy = peneCorr(tth,dep,dist)
470    return GetEllipse2(tth,dxy,dist,cent,tilt,phi)
471
472def GetDetectorXY(dsp,azm,data):
473    '''Get detector x,y position from d-spacing (dsp), azimuth (azm,deg)
474    & image controls dictionary (data) - new version
475    it seems to be only used in plotting
476    '''
477    def LinePlaneCollision(planeNormal, planePoint, rayDirection, rayPoint, epsilon=1e-6):
478     
479        ndotu = planeNormal.dot(rayDirection)
480        if ndotu < epsilon:   
481                return None
482     
483        w = rayPoint - planePoint
484        si = -planeNormal.dot(w) / ndotu
485        Psi = w + si * rayDirection + planePoint
486        return Psi
487   
488   
489    dist = data['distance']
490    cent = data['center']
491    T = makeMat(data['tilt'],0)
492    R = makeMat(data['rotation'],2)
493    MN = np.inner(R,np.inner(R,T))
494    iMN= nl.inv(MN)
495    tth = 2.0*npasind(data['wavelength']/(2.*dsp))
496    vect = np.array([npsind(tth)*npcosd(azm),npsind(tth)*npsind(azm),npcosd(tth)])
497    dxyz0 = np.inner(np.array([0.,0.,1.0]),MN)    #tilt detector normal
498    dxyz0 += np.array([0.,0.,dist])                 #translate to distance
499    dxyz0 = np.inner(dxyz0,makeMat(data['det2theta'],1).T)   #rotate on 2-theta
500    dxyz1 = np.inner(np.array([cent[0],cent[1],0.]),MN)    #tilt detector cent
501    dxyz1 += np.array([0.,0.,dist])                 #translate to distance
502    dxyz1 = np.inner(dxyz1,makeMat(data['det2theta'],1).T)   #rotate on 2-theta
503    xyz = LinePlaneCollision(dxyz0,dxyz1,vect,2.*dist*vect)
504    if xyz is None:
505        return None
506    xyz = np.inner(xyz,makeMat(-data['det2theta'],1).T)
507    xyz -= np.array([0.,0.,dist])                 #translate back
508    xyz = np.inner(xyz,iMN)
509    return np.squeeze(xyz)[:2]+cent
510       
511def GetDetectorXY2(dsp,azm,data):
512    '''Get detector x,y position from d-spacing (dsp), azimuth (azm,deg)
513    & image controls dictionary (data)
514    it seems to be only used in plotting
515    '''   
516    elcent,phi,radii = GetEllipse(dsp,data)
517    phi = data['rotation']-90.          #to give rotation of major axis
518    tilt = data['tilt']
519    dist = data['distance']
520    cent = data['center']
521    tth = 2.0*asind(data['wavelength']/(2.*dsp))
522    stth = sind(tth)
523    cosb = cosd(tilt)
524    if radii[0] > 0.:
525        sinb = sind(tilt)
526        tanb = tand(tilt)
527        fplus = dist*tanb*stth/(cosb+stth)
528        fminus = dist*tanb*stth/(cosb-stth)
529        zdis = (fplus-fminus)/2.
530        rsqplus = radii[0]**2+radii[1]**2
531        rsqminus = radii[0]**2-radii[1]**2
532        R = rsqminus*cosd(2.*azm-2.*phi)+rsqplus
533        Q = np.sqrt(2.)*radii[0]*radii[1]*np.sqrt(R-2.*zdis**2*sind(azm-phi)**2)
534        P = 2.*radii[0]**2*zdis*cosd(azm-phi)
535        radius = (P+Q)/R
536        xy = np.array([radius*cosd(azm),radius*sind(azm)])
537        xy += cent
538    else:   #hyperbola - both branches (one is way off screen!)
539        sinb = abs(sind(tilt))
540        tanb = abs(tand(tilt))
541        f = dist*tanb*stth/(cosb+stth)
542        v = dist*(tanb+tand(tth-abs(tilt)))
543        delt = dist*stth*(1+stth*cosb)/(sinb*cosb*(stth+cosb))
544        ecc = (v-f)/(delt-v)
545        R = radii[1]*(ecc**2-1)/(1-ecc*cosd(azm))
546        if tilt > 0.:
547            offset = 2.*radii[1]*ecc+f      #select other branch
548            xy = [-R*cosd(azm)-offset,-R*sind(azm)]
549        else:
550            offset = -f
551            xy = [-R*cosd(azm)-offset,R*sind(azm)]
552        xy = -np.array([xy[0]*cosd(phi)+xy[1]*sind(phi),xy[0]*sind(phi)-xy[1]*cosd(phi)])
553        xy += cent
554    if data['det2theta']:
555        xy[0] += dist*nptand(data['det2theta']+data['tilt']*npsind(data['rotation']))
556    return xy
557   
558def GetDetXYfromThAzm(Th,Azm,data):
559    '''Computes a detector position from a 2theta angle and an azimultal
560    angle (both in degrees) - apparently not used!
561    '''
562    dsp = data['wavelength']/(2.0*npsind(Th))   
563    return GetDetectorXY(dsp,Azm,data)
564                   
565def GetTthAzmDsp2(x,y,data): #expensive
566    '''Computes a 2theta, etc. from a detector position and calibration constants - checked
567    OK for ellipses & hyperbola.
568
569    :returns: np.array(tth,azm,G,dsp) where tth is 2theta, azm is the azimutal angle,
570       G is ? and dsp is the d-space
571    '''
572    wave = data['wavelength']
573    cent = data['center']
574    tilt = data['tilt']
575    dist = data['distance']/cosd(tilt)
576    x0 = dist*tand(tilt)
577    phi = data['rotation']
578    dep = data.get('DetDepth',0.)
579    azmthoff = data['azmthOff']
580    dx = np.array(x-cent[0],dtype=np.float32)
581    dy = np.array(y-cent[1],dtype=np.float32)
582    D = ((dx-x0)**2+dy**2+dist**2)      #sample to pixel distance
583    X = np.array(([dx,dy,np.zeros_like(dx)]),dtype=np.float32).T
584    X = np.dot(X,makeMat(phi,2))
585    Z = np.dot(X,makeMat(tilt,0)).T[2]
586    tth = npatand(np.sqrt(dx**2+dy**2-Z**2)/(dist-Z))
587    dxy = peneCorr(tth,dep,dist)
588    DX = dist-Z+dxy
589    DY = np.sqrt(dx**2+dy**2-Z**2)
590    tth = npatan2d(DY,DX) 
591    dsp = wave/(2.*npsind(tth/2.))
592    azm = (npatan2d(dy,dx)+azmthoff+720.)%360.
593    G = D/dist**2       #for geometric correction = 1/cos(2theta)^2 if tilt=0.
594    return np.array([tth,azm,G,dsp])
595   
596def GetTthAzmDsp(x,y,data): #expensive
597    '''Computes a 2theta, etc. from a detector position and calibration constants - checked
598    OK for ellipses & hyperbola.
599
600    :returns: np.array(tth,azm,G,dsp) where tth is 2theta, azm is the azimutal angle,
601       G is ? and dsp is the d-space
602    '''
603   
604    def costth(xyz):
605        u = xyz/nl.norm(xyz,axis=-1)[:,:,nxs]
606        return np.dot(u,np.array([0.,0.,1.]))
607       
608#zero detector 2-theta: tested with tilted images - perfect integrations
609    wave = data['wavelength']
610    dx = x-data['center'][0]
611    dy = y-data['center'][1]
612    tilt = data['tilt']
613    dist = data['distance']/npcosd(tilt)    #sample-beam intersection point
614    T = makeMat(tilt,0)
615    R = makeMat(data['rotation'],2)
616    MN = np.inner(R,np.inner(R,T))
617    dxyz0 = np.inner(np.dstack([dx,dy,np.zeros_like(dx)]),MN)    #correct for 45 deg tilt
618    dxyz0 += np.array([0.,0.,dist])
619    if data['DetDepth']:
620        ctth0 = costth(dxyz0)
621        tth0 = npacosd(ctth0)
622        dzp = peneCorr(tth0,data['DetDepth'],dist)
623        dxyz0[:,:,2] += dzp
624#non zero detector 2-theta:
625    if data['det2theta']:       
626        tthMat = makeMat(data['det2theta'],1)
627        dxyz = np.inner(dxyz0,tthMat.T)
628    else:
629        dxyz = dxyz0
630    ctth = costth(dxyz)
631    tth = npacosd(ctth)
632    dsp = wave/(2.*npsind(tth/2.))
633    azm = (npatan2d(dxyz[:,:,1],dxyz[:,:,0])+data['azmthOff']+720.)%360.       
634# G-calculation       
635    x0 = data['distance']*nptand(tilt)
636    x0x = x0*npcosd(data['rotation'])
637    x0y = x0*npsind(data['rotation'])
638    distsq = data['distance']**2
639    G = ((dx-x0x)**2+(dy-x0y)**2+distsq)/distsq       #for geometric correction = 1/cos(2theta)^2 if tilt=0.
640    return tth,azm,G,dsp
641   
642def GetTth(x,y,data):
643    'Give 2-theta value for detector x,y position; calibration info in data'
644    return GetTthAzmDsp(x,y,data)[0]
645   
646def GetTthAzm(x,y,data):
647    'Give 2-theta, azimuth values for detector x,y position; calibration info in data'
648    return GetTthAzmDsp(x,y,data)[0:2]
649   
650def GetTthAzmG2(x,y,data):
651    '''Give 2-theta, azimuth & geometric corr. values for detector x,y position;
652     calibration info in data - only used in integration - old version
653    '''
654    'Needs a doc string - checked OK for ellipses & hyperbola'
655    tilt = data['tilt']
656    dist = data['distance']/npcosd(tilt)
657    MN = -np.inner(makeMat(data['rotation'],2),makeMat(tilt,0))
658    dx = x-data['center'][0]
659    dy = y-data['center'][1]
660    dz = np.dot(np.dstack([dx.T,dy.T,np.zeros_like(dx.T)]),MN).T[2]   
661    xyZ = dx**2+dy**2-dz**2   
662    tth0 = npatand(np.sqrt(xyZ)/(dist-dz))
663    dzp = peneCorr(tth0,data['DetDepth'],dist)
664    tth = npatan2d(np.sqrt(xyZ),dist-dz+dzp) 
665    azm = (npatan2d(dy,dx)+data['azmthOff']+720.)%360.
666   
667    distsq = data['distance']**2
668    x0 = data['distance']*nptand(tilt)
669    x0x = x0*npcosd(data['rotation'])
670    x0y = x0*npsind(data['rotation'])
671    G = ((dx-x0x)**2+(dy-x0y)**2+distsq)/distsq       #for geometric correction = 1/cos(2theta)^2 if tilt=0.
672    return tth,azm,G
673
674def GetTthAzmG(x,y,data):
675    '''Give 2-theta, azimuth & geometric corr. values for detector x,y position;
676     calibration info in data - only used in integration
677     checked OK for ellipses & hyperbola
678     This is the slow step in image integration
679     '''     
680    def costth(xyz):
681        u = xyz/nl.norm(xyz,axis=-1)[:,:,nxs]
682        return np.dot(u,np.array([0.,0.,1.]))
683       
684#zero detector 2-theta: tested with tilted images - perfect integrations
685    dx = x-data['center'][0]
686    dy = y-data['center'][1]
687    tilt = data['tilt']
688    dist = data['distance']/npcosd(tilt)    #sample-beam intersection point
689    T = makeMat(tilt,0)
690    R = makeMat(data['rotation'],2)
691    MN = np.inner(R,np.inner(R,T))
692    dxyz0 = np.inner(np.dstack([dx,dy,np.zeros_like(dx)]),MN)    #correct for 45 deg tilt
693    dxyz0 += np.array([0.,0.,dist])
694    if data['DetDepth']:
695        ctth0 = costth(dxyz0)
696        tth0 = npacosd(ctth0)
697        dzp = peneCorr(tth0,data['DetDepth'],dist)
698        dxyz0[:,:,2] += dzp
699#non zero detector 2-theta:
700    if data['det2theta']:       
701        tthMat = makeMat(data['det2theta'],1)
702        dxyz = np.inner(dxyz0,tthMat.T)
703    else:
704        dxyz = dxyz0
705    ctth = costth(dxyz)
706    tth = npacosd(ctth)
707    azm = (npatan2d(dxyz[:,:,1],dxyz[:,:,0])+data['azmthOff']+720.)%360.       
708# G-calculation       
709    x0 = data['distance']*nptand(tilt)
710    x0x = x0*npcosd(data['rotation'])
711    x0y = x0*npsind(data['rotation'])
712    distsq = data['distance']**2
713    G = ((dx-x0x)**2+(dy-x0y)**2+distsq)/distsq       #for geometric correction = 1/cos(2theta)^2 if tilt=0.
714    return tth,azm,G
715
716def GetDsp(x,y,data):
717    'Give d-spacing value for detector x,y position; calibration info in data'
718    return GetTthAzmDsp(x,y,data)[3]
719       
720def GetAzm(x,y,data):
721    'Give azimuth value for detector x,y position; calibration info in data'
722    return GetTthAzmDsp(x,y,data)[1]
723   
724def meanAzm(a,b):
725    AZM = lambda a,b: npacosd(0.5*(npsind(2.*b)-npsind(2.*a))/(np.pi*(b-a)/180.))/2.
726    azm = AZM(a,b)
727#    quad = int((a+b)/180.)
728#    if quad == 1:
729#        azm = 180.-azm
730#    elif quad == 2:
731#        azm += 180.
732#    elif quad == 3:
733#        azm = 360-azm
734    return azm     
735       
736def ImageCompress(image,scale):
737    ''' Reduces size of image by selecting every n'th point
738    param: image array: original image
739    param: scale int: intervsl between selected points
740    returns: array: reduced size image
741    '''
742    if scale == 1:
743        return image
744    else:
745        return image[::scale,::scale]
746       
747def checkEllipse(Zsum,distSum,xSum,ySum,dist,x,y):
748    'Needs a doc string'
749    avg = np.array([distSum/Zsum,xSum/Zsum,ySum/Zsum])
750    curr = np.array([dist,x,y])
751    return abs(avg-curr)/avg < .02
752
753def GetLineScan(image,data):
754    Nx,Ny = data['size']
755    pixelSize = data['pixelSize']
756    scalex = 1000./pixelSize[0]         #microns --> 1/mm
757    scaley = 1000./pixelSize[1]
758    wave = data['wavelength']
759    numChans = data['outChannels']
760    LUtth = np.array(data['IOtth'],dtype=np.float)
761    azm = data['linescan'][1]-data['azmthOff']
762    Tx = np.array([tth for tth in np.linspace(LUtth[0],LUtth[1],numChans+1)])
763    Ty = np.zeros_like(Tx)
764    dsp = wave/(2.0*npsind(Tx/2.0))
765    xy = np.array([GetDetectorXY(d,azm,data) for d in dsp]).T
766    xy[1] *= scalex
767    xy[0] *= scaley
768    xy = np.array(xy,dtype=int)
769    Xpix = ma.masked_outside(xy[1],0,Ny-1)
770    Ypix = ma.masked_outside(xy[0],0,Nx-1)
771    xpix = Xpix[~(Xpix.mask+Ypix.mask)].compressed()
772    ypix = Ypix[~(Xpix.mask+Ypix.mask)].compressed()
773    Ty = image[xpix,ypix]
774    Tx = ma.array(Tx,mask=Xpix.mask+Ypix.mask).compressed()
775    return [Tx,Ty]
776
777def EdgeFinder(image,data):
778    '''this makes list of all x,y where I>edgeMin suitable for an ellipse search?
779    Not currently used but might be useful in future?
780    '''
781    import numpy.ma as ma
782    Nx,Ny = data['size']
783    pixelSize = data['pixelSize']
784    edgemin = data['edgemin']
785    scalex = pixelSize[0]/1000.
786    scaley = pixelSize[1]/1000.   
787    tay,tax = np.mgrid[0:Nx,0:Ny]
788    tax = np.asfarray(tax*scalex,dtype=np.float32)
789    tay = np.asfarray(tay*scaley,dtype=np.float32)
790    tam = ma.getmask(ma.masked_less(image.flatten(),edgemin))
791    tax = ma.compressed(ma.array(tax.flatten(),mask=tam))
792    tay = ma.compressed(ma.array(tay.flatten(),mask=tam))
793    return zip(tax,tay)
794   
795def MakeFrameMask(data,frame):
796    import polymask as pm
797    pixelSize = data['pixelSize']
798    scalex = pixelSize[0]/1000.
799    scaley = pixelSize[1]/1000.
800    blkSize = 512
801    Nx,Ny = data['size']
802    nXBlks = (Nx-1)//blkSize+1
803    nYBlks = (Ny-1)//blkSize+1
804    tam = ma.make_mask_none(data['size'])
805    for iBlk in range(nXBlks):
806        iBeg = iBlk*blkSize
807        iFin = min(iBeg+blkSize,Nx)
808        for jBlk in range(nYBlks):
809            jBeg = jBlk*blkSize
810            jFin = min(jBeg+blkSize,Ny)               
811            nI = iFin-iBeg
812            nJ = jFin-jBeg
813            tax,tay = np.mgrid[iBeg+0.5:iFin+.5,jBeg+.5:jFin+.5]         #bin centers not corners
814            tax = np.asfarray(tax*scalex,dtype=np.float32)
815            tay = np.asfarray(tay*scaley,dtype=np.float32)
816            tamp = ma.make_mask_none((1024*1024))
817            tamp = ma.make_mask(pm.polymask(nI*nJ,tax.flatten(),
818                tay.flatten(),len(frame),frame,tamp)[:nI*nJ])^True  #switch to exclude around frame
819            if tamp.shape:
820                tamp = np.reshape(tamp[:nI*nJ],(nI,nJ))
821                tam[iBeg:iFin,jBeg:jFin] = ma.mask_or(tamp[0:nI,0:nJ],tam[iBeg:iFin,jBeg:jFin])
822            else:
823                tam[iBeg:iFin,jBeg:jFin] = True
824    return tam.T
825   
826def ImageRecalibrate(G2frame,ImageZ,data,masks,getRingsOnly=False):
827    '''Called to repeat the calibration on an image, usually called after
828    calibration is done initially to improve the fit.
829
830    :param G2frame: The top-level GSAS-II frame or None, to skip plotting
831    :param np.Array ImageZ: the image to calibrate
832    :param dict data: the Controls dict for the image
833    :param dict masks: a dict with masks
834    :returns: a list containing vals,varyList,sigList,parmDict,covar or rings
835      (with an array of x, y, and d-space values) if getRingsOnly is True
836      or an empty list, in case of an error
837    '''
838    import ImageCalibrants as calFile
839    if not getRingsOnly:
840        G2fil.G2Print ('Image recalibration:')
841    time0 = time.time()
842    pixelSize = data['pixelSize']
843    scalex = 1000./pixelSize[0]
844    scaley = 1000./pixelSize[1]
845    pixLimit = data['pixLimit']
846    cutoff = data['cutoff']
847    data['rings'] = []
848    data['ellipses'] = []
849    if data['DetDepth'] > 0.5:          #patch - redefine DetDepth
850        data['DetDepth'] /= data['distance']
851    if not data['calibrant']:
852        G2fil.G2Print ('warning: no calibration material selected')
853        return []   
854    skip = data['calibskip']
855    dmin = data['calibdmin']
856    if data['calibrant'] not in calFile.Calibrants:
857        G2fil.G2Print('Warning: %s not in local copy of image calibrants file'%data['calibrant'])
858        return []
859    calibrant = calFile.Calibrants[data['calibrant']]
860    Bravais,SGs,Cells = calibrant[:3]
861    HKL = []
862    for bravais,sg,cell in zip(Bravais,SGs,Cells):
863        A = G2lat.cell2A(cell)
864        if sg:
865            SGData = G2spc.SpcGroup(sg)[1]
866            hkl = G2pwd.getHKLpeak(dmin,SGData,A,Inst=None,nodup=True)
867            HKL += list(hkl)
868        else:
869            hkl = G2lat.GenHBravais(dmin,bravais,A)
870            HKL += list(hkl)
871    if len(calibrant) > 5:
872        absent = calibrant[5]
873    else:
874        absent = ()
875    HKL = G2lat.sortHKLd(HKL,True,False)
876    varyList = [item for item in data['varyList'] if data['varyList'][item]]
877    parmDict = {'dist':data['distance'],'det-X':data['center'][0],'det-Y':data['center'][1],
878        'setdist':data.get('setdist',data['distance']),
879        'tilt':data['tilt'],'phi':data['rotation'],'wave':data['wavelength'],'dep':data['DetDepth']}
880    Found = False
881    wave = data['wavelength']
882    frame = masks['Frames']
883    tam = ma.make_mask_none(ImageZ.shape)
884    if frame:
885        tam = ma.mask_or(tam,MakeFrameMask(data,frame))
886    for iH,H in enumerate(HKL):
887        if debug:   print (H) 
888        dsp = H[3]
889        tth = 2.0*asind(wave/(2.*dsp))
890        if tth+abs(data['tilt']) > 90.:
891            G2fil.G2Print ('next line is a hyperbola - search stopped')
892            break
893        ellipse = GetEllipse(dsp,data)
894        if iH not in absent and iH >= skip:
895            Ring = makeRing(dsp,ellipse,pixLimit,cutoff,scalex,scaley,ma.array(ImageZ,mask=tam))[0]
896        else:
897            Ring = makeRing(dsp,ellipse,pixLimit,1000.0,scalex,scaley,ma.array(ImageZ,mask=tam))[0]
898        if Ring:
899            if iH not in absent and iH >= skip:
900                data['rings'].append(np.array(Ring))
901            data['ellipses'].append(copy.deepcopy(ellipse+('r',)))
902            Found = True
903        elif not Found:         #skipping inner rings, keep looking until ring found
904            continue
905        else:                   #no more rings beyond edge of detector
906            data['ellipses'].append([])
907            continue
908    if not data['rings']:
909        G2fil.G2Print ('no rings found; try lower Min ring I/Ib',mode='warn')
910        return []   
911       
912    rings = np.concatenate((data['rings']),axis=0)
913    if getRingsOnly:
914        return rings,HKL
915    [chisq,vals,sigList,covar] = FitDetector(rings,varyList,parmDict,True,True)
916    data['wavelength'] = parmDict['wave']
917    data['distance'] = parmDict['dist']
918    data['center'] = [parmDict['det-X'],parmDict['det-Y']]
919    data['rotation'] = np.mod(parmDict['phi'],360.0)
920    data['tilt'] = parmDict['tilt']
921    data['DetDepth'] = parmDict['dep']
922    data['chisq'] = chisq
923    N = len(data['ellipses'])
924    data['ellipses'] = []           #clear away individual ellipse fits
925    for H in HKL[:N]:
926        ellipse = GetEllipse(H[3],data)
927        data['ellipses'].append(copy.deepcopy(ellipse+('b',)))   
928    G2fil.G2Print ('calibration time = %.3f'%(time.time()-time0))
929    if G2frame:
930        G2plt.PlotImage(G2frame,newImage=True)       
931    return [vals,varyList,sigList,parmDict,covar]
932
933def ImageCalibrate(G2frame,data):
934    '''Called to perform an initial image calibration after points have been
935    selected for the inner ring.
936    '''
937    import ImageCalibrants as calFile
938    G2fil.G2Print ('Image calibration:')
939    time0 = time.time()
940    ring = data['ring']
941    pixelSize = data['pixelSize']
942    scalex = 1000./pixelSize[0]
943    scaley = 1000./pixelSize[1]
944    pixLimit = data['pixLimit']
945    cutoff = data['cutoff']
946    varyDict = data['varyList']
947    if varyDict['dist'] and varyDict['wave']:
948        G2fil.G2Print ('ERROR - you can not simultaneously calibrate distance and wavelength')
949        return False
950    if len(ring) < 5:
951        G2fil.G2Print ('ERROR - not enough inner ring points for ellipse')
952        return False
953       
954    #fit start points on inner ring
955    data['ellipses'] = []
956    data['rings'] = []
957    outE = FitEllipse(ring)
958    fmt  = '%s X: %.3f, Y: %.3f, phi: %.3f, R1: %.3f, R2: %.3f'
959    fmt2 = '%s X: %.3f, Y: %.3f, phi: %.3f, R1: %.3f, R2: %.3f, chi**2: %.3f, Np: %d'
960    if outE:
961        G2fil.G2Print (fmt%('start ellipse: ',outE[0][0],outE[0][1],outE[1],outE[2][0],outE[2][1]))
962        ellipse = outE
963    else:
964        return False
965       
966    #setup 360 points on that ring for "good" fit
967    data['ellipses'].append(ellipse[:]+('g',))
968    Ring = makeRing(1.0,ellipse,pixLimit,cutoff,scalex,scaley,G2frame.ImageZ)[0]
969    if Ring:
970        ellipse = FitEllipse(Ring)
971        Ring = makeRing(1.0,ellipse,pixLimit,cutoff,scalex,scaley,G2frame.ImageZ)[0]    #do again
972        ellipse = FitEllipse(Ring)
973    else:
974        G2fil.G2Print ('1st ring not sufficiently complete to proceed',mode='warn')
975        return False
976    if debug:
977        G2fil.G2Print (fmt2%('inner ring:    ',ellipse[0][0],ellipse[0][1],ellipse[1],
978            ellipse[2][0],ellipse[2][1],0.,len(Ring)))     #cent,phi,radii
979    data['ellipses'].append(ellipse[:]+('r',))
980    data['rings'].append(np.array(Ring))
981    G2plt.PlotImage(G2frame,newImage=True)
982   
983#setup for calibration
984    data['rings'] = []
985    if not data['calibrant']:
986        G2fil.G2Print ('Warning: no calibration material selected')
987        return True
988   
989    skip = data['calibskip']
990    dmin = data['calibdmin']
991#generate reflection set
992    calibrant = calFile.Calibrants[data['calibrant']]
993    Bravais,SGs,Cells = calibrant[:3]
994    HKL = []
995    for bravais,sg,cell in zip(Bravais,SGs,Cells):
996        A = G2lat.cell2A(cell)
997        if sg:
998            SGData = G2spc.SpcGroup(sg)[1]
999            hkl = G2pwd.getHKLpeak(dmin,SGData,A,Inst=None,nodup=True)
1000            HKL += list(hkl)
1001        else:
1002            hkl = G2lat.GenHBravais(dmin,bravais,A)
1003            HKL += list(hkl)
1004    HKL = G2lat.sortHKLd(HKL,True,False)[skip:]
1005#set up 1st ring
1006    elcent,phi,radii = ellipse              #from fit of 1st ring
1007    dsp = HKL[0][3]
1008    G2fil.G2Print ('1st ring: try %.4f'%(dsp))
1009    if varyDict['dist']:
1010        wave = data['wavelength']
1011        tth = 2.0*asind(wave/(2.*dsp))
1012    else:   #varyDict['wave']!
1013        dist = data['distance']
1014        tth = npatan2d(radii[0],dist)
1015        data['wavelength'] = wave =  2.0*dsp*sind(tth/2.0)
1016    Ring0 = makeRing(dsp,ellipse,3,cutoff,scalex,scaley,G2frame.ImageZ)[0]
1017    ttth = nptand(tth)
1018    ctth = npcosd(tth)
1019#1st estimate of tilt; assume ellipse - don't know sign though
1020    if varyDict['tilt']:
1021        tilt = npasind(np.sqrt(max(0.,1.-(radii[0]/radii[1])**2))*ctth)
1022        if not tilt:
1023            G2fil.G2Print ('WARNING - selected ring was fitted as a circle')
1024            G2fil.G2Print (' - if detector was tilted we suggest you skip this ring - WARNING')
1025    else:
1026        tilt = data['tilt']
1027#1st estimate of dist: sample to detector normal to plane
1028    if varyDict['dist']:
1029        data['distance'] = dist = radii[0]**2/(ttth*radii[1])
1030    else:
1031        dist = data['distance']
1032    if varyDict['tilt']:
1033#ellipse to cone axis (x-ray beam); 2 choices depending on sign of tilt
1034        zdisp = radii[1]*ttth*tand(tilt)
1035        zdism = radii[1]*ttth*tand(-tilt)
1036#cone axis position; 2 choices. Which is right?     
1037#NB: zdisp is || to major axis & phi is rotation of minor axis
1038#thus shift from beam to ellipse center is [Z*sin(phi),-Z*cos(phi)]
1039        centp = [elcent[0]+zdisp*sind(phi),elcent[1]-zdisp*cosd(phi)]
1040        centm = [elcent[0]+zdism*sind(phi),elcent[1]-zdism*cosd(phi)]
1041#check get same ellipse parms either way
1042#now do next ring; estimate either way & do a FitDetector each way; best fit is correct one
1043        fail = True
1044        i2 = 1
1045        while fail:
1046            dsp = HKL[i2][3]
1047            G2fil.G2Print ('2nd ring: try %.4f'%(dsp))
1048            tth = 2.0*asind(wave/(2.*dsp))
1049            ellipsep = GetEllipse2(tth,0.,dist,centp,tilt,phi)
1050            G2fil.G2Print (fmt%('plus ellipse :',ellipsep[0][0],ellipsep[0][1],ellipsep[1],ellipsep[2][0],ellipsep[2][1]))
1051            Ringp = makeRing(dsp,ellipsep,3,cutoff,scalex,scaley,G2frame.ImageZ)[0]
1052            parmDict = {'dist':dist,'det-X':centp[0],'det-Y':centp[1],
1053                'tilt':tilt,'phi':phi,'wave':wave,'dep':0.0}       
1054            varyList = [item for item in varyDict if varyDict[item]]
1055            if len(Ringp) > 10:
1056                chip = FitDetector(np.array(Ring0+Ringp),varyList,parmDict,True)[0]
1057                tiltp = parmDict['tilt']
1058                phip = parmDict['phi']
1059                centp = [parmDict['det-X'],parmDict['det-Y']]
1060                fail = False
1061            else:
1062                chip = 1e6
1063            ellipsem = GetEllipse2(tth,0.,dist,centm,-tilt,phi)
1064            G2fil.G2Print (fmt%('minus ellipse:',ellipsem[0][0],ellipsem[0][1],ellipsem[1],ellipsem[2][0],ellipsem[2][1]))
1065            Ringm = makeRing(dsp,ellipsem,3,cutoff,scalex,scaley,G2frame.ImageZ)[0]
1066            if len(Ringm) > 10:
1067                parmDict['tilt'] *= -1
1068                chim = FitDetector(np.array(Ring0+Ringm),varyList,parmDict,True)[0]
1069                tiltm = parmDict['tilt']
1070                phim = parmDict['phi']
1071                centm = [parmDict['det-X'],parmDict['det-Y']]
1072                fail = False
1073            else:
1074                chim = 1e6
1075            if fail:
1076                i2 += 1
1077        if chip < chim:
1078            data['tilt'] = tiltp
1079            data['center'] = centp
1080            data['rotation'] = phip
1081        else:
1082            data['tilt'] = tiltm
1083            data['center'] = centm
1084            data['rotation'] = phim
1085        data['ellipses'].append(ellipsep[:]+('b',))
1086        data['rings'].append(np.array(Ringp))
1087        data['ellipses'].append(ellipsem[:]+('r',))
1088        data['rings'].append(np.array(Ringm))
1089        G2plt.PlotImage(G2frame,newImage=True)
1090    if data['DetDepth'] > 0.5:          #patch - redefine DetDepth
1091        data['DetDepth'] /= data['distance']
1092    parmDict = {'dist':data['distance'],'det-X':data['center'][0],'det-Y':data['center'][1],
1093        'tilt':data['tilt'],'phi':data['rotation'],'wave':data['wavelength'],'dep':data['DetDepth']}
1094    varyList = [item for item in varyDict if varyDict[item]]
1095    data['rings'] = []
1096    data['ellipses'] = []
1097    for i,H in enumerate(HKL):
1098        dsp = H[3]
1099        tth = 2.0*asind(wave/(2.*dsp))
1100        if tth+abs(data['tilt']) > 90.:
1101            G2fil.G2Print ('next line is a hyperbola - search stopped')
1102            break
1103        if debug:   print ('HKLD:'+str(H[:4])+'2-theta: %.4f'%(tth))
1104        elcent,phi,radii = ellipse = GetEllipse(dsp,data)
1105        data['ellipses'].append(copy.deepcopy(ellipse+('g',)))
1106        if debug:   print (fmt%('predicted ellipse:',elcent[0],elcent[1],phi,radii[0],radii[1]))
1107        Ring = makeRing(dsp,ellipse,pixLimit,cutoff,scalex,scaley,G2frame.ImageZ)[0]
1108        if Ring:
1109            data['rings'].append(np.array(Ring))
1110            rings = np.concatenate((data['rings']),axis=0)
1111            if i:
1112                chisq = FitDetector(rings,varyList,parmDict,False)[0]
1113                data['distance'] = parmDict['dist']
1114                data['center'] = [parmDict['det-X'],parmDict['det-Y']]
1115                data['rotation'] = parmDict['phi']
1116                data['tilt'] = parmDict['tilt']
1117                data['DetDepth'] = parmDict['dep']
1118                data['chisq'] = chisq
1119                elcent,phi,radii = ellipse = GetEllipse(dsp,data)
1120                if debug:   print (fmt2%('fitted ellipse:   ',elcent[0],elcent[1],phi,radii[0],radii[1],chisq,len(rings)))
1121            data['ellipses'].append(copy.deepcopy(ellipse+('r',)))
1122#            G2plt.PlotImage(G2frame,newImage=True)
1123        else:
1124            if debug:   print ('insufficient number of points in this ellipse to fit')
1125#            break
1126    G2plt.PlotImage(G2frame,newImage=True)
1127    fullSize = len(G2frame.ImageZ)/scalex
1128    if 2*radii[1] < .9*fullSize:
1129        G2fil.G2Print ('Are all usable rings (>25% visible) used? Try reducing Min ring I/Ib')
1130    N = len(data['ellipses'])
1131    if N > 2:
1132        FitDetector(rings,varyList,parmDict)[0]
1133        data['wavelength'] = parmDict['wave']
1134        data['distance'] = parmDict['dist']
1135        data['center'] = [parmDict['det-X'],parmDict['det-Y']]
1136        data['rotation'] = parmDict['phi']
1137        data['tilt'] = parmDict['tilt']
1138        data['DetDepth'] = parmDict['dep']
1139    for H in HKL[:N]:
1140        ellipse = GetEllipse(H[3],data)
1141        data['ellipses'].append(copy.deepcopy(ellipse+('b',)))
1142    G2fil.G2Print ('calibration time = %.3f'%(time.time()-time0))
1143    G2plt.PlotImage(G2frame,newImage=True)       
1144    return True
1145   
1146def Make2ThetaAzimuthMap(data,iLim,jLim): #most expensive part of integration!
1147    'Needs a doc string'
1148    #transforms 2D image from x,y space to 2-theta,azimuth space based on detector orientation
1149    pixelSize = data['pixelSize']
1150    scalex = pixelSize[0]/1000.
1151    scaley = pixelSize[1]/1000.
1152    tay,tax = np.mgrid[iLim[0]+0.5:iLim[1]+.5,jLim[0]+.5:jLim[1]+.5]         #bin centers not corners
1153    tax = np.asfarray(tax*scalex,dtype=np.float32).flatten()
1154    tay = np.asfarray(tay*scaley,dtype=np.float32).flatten()
1155    nI = iLim[1]-iLim[0]
1156    nJ = jLim[1]-jLim[0]
1157    TA = np.empty((4,nI,nJ))
1158    TA[:3] = np.array(GetTthAzmG(np.reshape(tax,(nI,nJ)),np.reshape(tay,(nI,nJ)),data))     #includes geom. corr. as dist**2/d0**2 - most expensive step
1159    TA[1] = np.where(TA[1]<0,TA[1]+360,TA[1])
1160    TA[3] = G2pwd.Polarization(data['PolaVal'][0],TA[0],TA[1]-90.)[0]
1161    return TA           #2-theta, azimuth & geom. corr. arrays
1162
1163def MakeMaskMap(data,masks,iLim,jLim,tamp):
1164    import polymask as pm
1165    pixelSize = data['pixelSize']
1166    scalex = pixelSize[0]/1000.
1167    scaley = pixelSize[1]/1000.
1168   
1169    tay,tax = np.mgrid[iLim[0]+0.5:iLim[1]+.5,jLim[0]+.5:jLim[1]+.5]         #bin centers not corners
1170    tax = np.asfarray(tax*scalex,dtype=np.float32).flatten()
1171    tay = np.asfarray(tay*scaley,dtype=np.float32).flatten()
1172    nI = iLim[1]-iLim[0]
1173    nJ = jLim[1]-jLim[0]
1174    #make position masks here
1175    frame = masks['Frames']
1176    tam = ma.make_mask_none((nI*nJ))
1177    if frame:
1178        tam = ma.mask_or(tam,ma.make_mask(pm.polymask(nI*nJ,tax,
1179            tay,len(frame),frame,tamp)[:nI*nJ])^True)
1180    polygons = masks['Polygons']
1181    for polygon in polygons:
1182        if polygon:
1183            tam = ma.mask_or(tam,ma.make_mask(pm.polymask(nI*nJ,tax,
1184                tay,len(polygon),polygon,tamp)[:nI*nJ]))
1185    for X,Y,rsq in masks['Points'].T:
1186        tam = ma.mask_or(tam,ma.getmask(ma.masked_less((tax-X)**2+(tay-Y)**2,rsq)))
1187    if tam.shape: 
1188        tam = np.reshape(tam,(nI,nJ))
1189    else:
1190        tam = ma.make_mask_none((nI,nJ))
1191    for xline in masks.get('Xlines',[]):    #a y pixel position
1192        if iLim[0] <= xline <= iLim[1]:
1193            tam[xline-iLim[0],:] = True
1194    for yline in masks.get('Ylines',[]):    #a x pixel position
1195        if jLim[0] <= yline <= jLim[1]:
1196            tam[:,yline-jLim[0]] = True           
1197    return tam           #position mask
1198
1199def Fill2ThetaAzimuthMap(masks,TA,tam,image):
1200    'Needs a doc string'
1201    Zlim = masks['Thresholds'][1]
1202    rings = masks['Rings']
1203    arcs = masks['Arcs']
1204    TA = np.dstack((ma.getdata(TA[1]),ma.getdata(TA[0]),ma.getdata(TA[2]),ma.getdata(TA[3])))    #azimuth, 2-theta, dist, pol
1205    tax,tay,tad,pol = np.dsplit(TA,4)    #azimuth, 2-theta, dist**2/d0**2, pol
1206    for tth,thick in rings:
1207        tam = ma.mask_or(tam.flatten(),ma.getmask(ma.masked_inside(tay.flatten(),max(0.01,tth-thick/2.),tth+thick/2.)))
1208    for tth,azm,thick in arcs:
1209        tamt = ma.getmask(ma.masked_inside(tay.flatten(),max(0.01,tth-thick/2.),tth+thick/2.))
1210        tama = ma.getmask(ma.masked_inside(tax.flatten(),azm[0],azm[1]))
1211        tam = ma.mask_or(tam.flatten(),tamt*tama)
1212    taz = ma.masked_outside(image.flatten(),int(Zlim[0]),Zlim[1])
1213    tabs = np.ones_like(taz)
1214    tam = ma.mask_or(tam.flatten(),ma.getmask(taz))
1215    tax = ma.compressed(ma.array(tax.flatten(),mask=tam))   #azimuth
1216    tay = ma.compressed(ma.array(tay.flatten(),mask=tam))   #2-theta
1217    taz = ma.compressed(ma.array(taz.flatten(),mask=tam))   #intensity
1218    tad = ma.compressed(ma.array(tad.flatten(),mask=tam))   #dist**2/d0**2
1219    tabs = ma.compressed(ma.array(tabs.flatten(),mask=tam)) #ones - later used for absorption corr.
1220    pol = ma.compressed(ma.array(pol.flatten(),mask=tam))   #polarization
1221    return tax,tay,taz/pol,tad,tabs
1222
1223def MakeUseTA(data,blkSize=128):
1224    Nx,Ny = data['size']
1225    nXBlks = (Nx-1)//blkSize+1
1226    nYBlks = (Ny-1)//blkSize+1
1227    useTA = []
1228    for iBlk in range(nYBlks):
1229        iBeg = iBlk*blkSize
1230        iFin = min(iBeg+blkSize,Ny)
1231        useTAj = []
1232        for jBlk in range(nXBlks):
1233            jBeg = jBlk*blkSize
1234            jFin = min(jBeg+blkSize,Nx)
1235            TA = Make2ThetaAzimuthMap(data,(iBeg,iFin),(jBeg,jFin))          #2-theta & azimuth arrays & create position mask
1236            useTAj.append(TA)
1237        useTA.append(useTAj)
1238    return useTA
1239
1240def MakeUseMask(data,masks,blkSize=128):
1241    Masks = copy.deepcopy(masks)
1242    Masks['Points'] = np.array(Masks['Points']).T           #get spots as X,Y,R arrays
1243    if np.any(masks['Points']):
1244        Masks['Points'][2] = np.square(Masks['Points'][2]/2.)
1245    Nx,Ny = data['size']
1246    nXBlks = (Nx-1)//blkSize+1
1247    nYBlks = (Ny-1)//blkSize+1
1248    useMask = []
1249    tamp = ma.make_mask_none((1024*1024))       #NB: this array size used in the fortran histogram2d
1250    for iBlk in range(nYBlks):
1251        iBeg = iBlk*blkSize
1252        iFin = min(iBeg+blkSize,Ny)
1253        useMaskj = []
1254        for jBlk in range(nXBlks):
1255            jBeg = jBlk*blkSize
1256            jFin = min(jBeg+blkSize,Nx)
1257            mask = MakeMaskMap(data,Masks,(iBeg,iFin),(jBeg,jFin),tamp)          #2-theta & azimuth arrays & create position mask
1258            useMaskj.append(mask)
1259        useMask.append(useMaskj)
1260    return useMask
1261
1262def ImageIntegrate(image,data,masks,blkSize=128,returnN=False,useTA=None,useMask=None):
1263    'Integrate an image; called from OnIntegrateAll and OnIntegrate in G2imgGUI'    #for q, log(q) bins need data['binType']
1264    import histogram2d as h2d
1265    G2fil.G2Print ('Begin image integration; image range: %d %d'%(np.min(image),np.max(image)))
1266    CancelPressed = False
1267    LUtth = np.array(data['IOtth'])
1268    LRazm = np.array(data['LRazimuth'],dtype=np.float64)
1269    numAzms = data['outAzimuths']
1270    numChans = (data['outChannels']//4)*4
1271    Dazm = (LRazm[1]-LRazm[0])/numAzms
1272    if '2-theta' in data.get('binType','2-theta'):
1273        lutth = LUtth               
1274    elif 'log(q)' in data['binType']:
1275        lutth = np.log(4.*np.pi*npsind(LUtth/2.)/data['wavelength'])
1276    elif 'q' == data['binType'].lower():
1277        lutth = 4.*np.pi*npsind(LUtth/2.)/data['wavelength']
1278    dtth = (lutth[1]-lutth[0])/numChans
1279    muT = data.get('SampleAbs',[0.0,''])[0]
1280    if data['DetDepth'] > 0.5:          #patch - redefine DetDepth
1281        data['DetDepth'] /= data['distance']
1282    if 'SASD' in data['type']:
1283        muT = -np.log(muT)/2.       #Transmission to 1/2 thickness muT
1284    Masks = copy.deepcopy(masks)
1285    Masks['Points'] = np.array(Masks['Points']).T           #get spots as X,Y,R arrays
1286    if np.any(masks['Points']):
1287        Masks['Points'][2] = np.square(Masks['Points'][2]/2.)
1288    NST = np.zeros(shape=(numAzms,numChans),order='F',dtype=np.float32)
1289    H0 = np.zeros(shape=(numAzms,numChans),order='F',dtype=np.float32)
1290    H2 = np.linspace(lutth[0],lutth[1],numChans+1)
1291    Nx,Ny = data['size']
1292    nXBlks = (Nx-1)//blkSize+1
1293    nYBlks = (Ny-1)//blkSize+1
1294    tbeg = time.time()
1295    times = [0,0,0,0,0]
1296    tamp = ma.make_mask_none((1024*1024))       #NB: this array size used in the fortran histogram2d
1297    for iBlk in range(nYBlks):
1298        iBeg = iBlk*blkSize
1299        iFin = min(iBeg+blkSize,Ny)
1300        for jBlk in range(nXBlks):
1301            jBeg = jBlk*blkSize
1302            jFin = min(jBeg+blkSize,Nx)
1303            # next is most expensive step!
1304
1305            t0 = time.time()
1306            if useTA:
1307                TA = useTA[iBlk][jBlk]
1308            else:
1309                TA = Make2ThetaAzimuthMap(data,(iBeg,iFin),(jBeg,jFin))           #2-theta & azimuth arrays & create position mask
1310            times[1] += time.time()-t0      #xy->th,azm
1311
1312            t0 = time.time()
1313            if useMask:
1314                tam = useMask[iBlk][jBlk]
1315            else:
1316                tam = MakeMaskMap(data,Masks,(iBeg,iFin),(jBeg,jFin),tamp)
1317            times[0] += time.time()-t0      #apply masks
1318
1319            t0 = time.time()
1320            Block = image[iBeg:iFin,jBeg:jFin]          #apply image spotmask here
1321            tax,tay,taz,tad,tabs = Fill2ThetaAzimuthMap(Masks,TA,tam,Block)    #and apply masks
1322            tax = np.where(tax > LRazm[1],tax-360.,tax)                 #put azm inside limits if possible
1323            tax = np.where(tax < LRazm[0],tax+360.,tax)
1324            if data.get('SampleAbs',[0.0,''])[1]:
1325                if 'Cylind' in data['SampleShape']:
1326                    muR = muT*(1.+npsind(tax)**2/2.)/(npcosd(tay))      #adjust for additional thickness off sample normal
1327                    tabs = G2pwd.Absorb(data['SampleShape'],muR,tay)
1328                elif 'Fixed' in data['SampleShape']:    #assumes flat plate sample normal to beam
1329                    tabs = G2pwd.Absorb('Fixed',muT,tay)
1330            if 'log(q)' in data.get('binType',''):
1331                tay = np.log(4.*np.pi*npsind(tay/2.)/data['wavelength'])
1332            elif 'q' == data.get('binType','').lower():
1333                tay = 4.*np.pi*npsind(tay/2.)/data['wavelength']
1334            times[2] += time.time()-t0          #fill map
1335
1336            t0 = time.time()
1337            taz = np.array((taz*tad/tabs),dtype='float32')
1338            if any([tax.shape[0],tay.shape[0],taz.shape[0]]):
1339                NST,H0 = h2d.histogram2d(len(tax),tax,tay,taz,
1340                    numAzms,numChans,LRazm,lutth,Dazm,dtth,NST,H0)
1341            times[3] += time.time()-t0          #binning
1342    G2fil.G2Print('End integration loops')
1343
1344    t0 = time.time()
1345    #prepare masked arrays of bins with pixels for interpolation setup
1346    H2msk = [ma.array(H2[:-1],mask=np.logical_not(nst)) for nst in NST]
1347#    H0msk = [ma.array(np.divide(h0,nst),mask=np.logical_not(nst)) for nst,h0 in zip(NST,H0)]
1348    H0msk = [ma.array(h0/nst,mask=np.logical_not(nst)) for nst,h0 in zip(NST,H0)]
1349    #make linear interpolators; outside limits give NaN
1350    H0 = []
1351    problemEntries = []
1352    for i,(h0msk,h2msk) in enumerate(zip(H0msk,H2msk)):
1353        try:
1354            h0int = scint.interp1d(h2msk.compressed(),h0msk.compressed(),bounds_error=False)
1355            #do interpolation on all points - fills in the empty bins; leaves others the same
1356            h0 = h0int(H2[:-1])
1357            H0.append(h0)
1358        except ValueError:
1359            problemEntries.append(i)
1360            H0.append(np.zeros(numChans,order='F',dtype=np.float32))
1361    H0 = np.nan_to_num(np.array(H0))       
1362    if 'log(q)' in data.get('binType',''):
1363        H2 = 2.*npasind(np.exp(H2)*data['wavelength']/(4.*np.pi))
1364    elif 'q' == data.get('binType','').lower():
1365        H2 = 2.*npasind(H2*data['wavelength']/(4.*np.pi))
1366    if Dazm:       
1367        H1 = np.array([azm for azm in np.linspace(LRazm[0],LRazm[1],numAzms+1)])
1368    else:
1369        H1 = LRazm
1370    if 'SASD' not in data['type']:
1371        H0 *= np.array(G2pwd.Polarization(data['PolaVal'][0],H2[:-1],0.)[0])
1372    H0 /= np.abs(npcosd(H2[:-1]-np.abs(data['det2theta'])))           #parallax correction
1373    if 'SASD' in data['type']:
1374        H0 /= npcosd(H2[:-1])           #one more for small angle scattering data?
1375    if data['Oblique'][1]:
1376        H0 /= G2pwd.Oblique(data['Oblique'][0],H2[:-1])
1377    times[4] += time.time()-t0          #cleanup
1378
1379    G2fil.G2Print ('Step times: \n apply masks  %8.3fs xy->th,azm   %8.3fs fill map     %8.3fs \
1380        \n binning      %8.3fs cleanup      %8.3fs'%(times[0],times[1],times[2],times[3],times[4]))
1381    G2fil.G2Print ("Elapsed time:","%8.3fs"%(time.time()-tbeg))
1382    G2fil.G2Print ('Integration complete')
1383    if problemEntries:
1384        msg = ""
1385        for i in problemEntries:
1386            if msg: msg += ', '
1387            msg += "{:.2f}".format((H1[i+1]+H1[i])/2.)
1388        print('\nWarning, integrations have no pixels at these azimuth(s)\n\t',msg,'\n')
1389    if returnN:     #As requested by Steven Weigand
1390        return H0,H1,H2,NST,CancelPressed
1391    else:
1392        return H0,H1,H2,CancelPressed
1393   
1394def MakeStrStaRing(ring,Image,Controls):
1395    ellipse = GetEllipse(ring['Dset'],Controls)
1396    pixSize = Controls['pixelSize']
1397    scalex = 1000./pixSize[0]
1398    scaley = 1000./pixSize[1]
1399    Ring = np.array(makeRing(ring['Dset'],ellipse,ring['pixLimit'],ring['cutoff'],scalex,scaley,Image)[0]).T   #returns x,y,dsp for each point in ring
1400    if len(Ring):
1401        ring['ImxyObs'] = copy.copy(Ring[:2])
1402        TA = GetTthAzm(Ring[0],Ring[1],Controls)       #convert x,y to tth,azm
1403        TA[0] = Controls['wavelength']/(2.*npsind(TA[0]/2.))      #convert 2th to d
1404        ring['ImtaObs'] = TA
1405        ring['ImtaCalc'] = np.zeros_like(ring['ImtaObs'])
1406        Ring[0] = TA[0]
1407        Ring[1] = TA[1]
1408        return Ring,ring
1409    else:
1410        ring['ImxyObs'] = [[],[]]
1411        ring['ImtaObs'] = [[],[]]
1412        ring['ImtaCalc'] = [[],[]]
1413        return [],[]    #bad ring; no points found
1414   
1415def FitStrSta(Image,StrSta,Controls):
1416    'Needs a doc string'
1417   
1418    StaControls = copy.deepcopy(Controls)
1419    phi = StrSta['Sample phi']
1420    wave = Controls['wavelength']
1421    pixelSize = Controls['pixelSize']
1422    scalex = 1000./pixelSize[0]
1423    scaley = 1000./pixelSize[1]
1424    StaType = StrSta['Type']
1425    StaControls['distance'] += StrSta['Sample z']*cosd(phi)
1426
1427    for ring in StrSta['d-zero']:       #get observed x,y,d points for the d-zeros
1428        dset = ring['Dset']
1429        Ring,R = MakeStrStaRing(ring,Image,StaControls)
1430        if len(Ring):
1431            ring.update(R)
1432            p0 = ring['Emat']
1433            val,esd,covMat = FitStrain(Ring,p0,dset,wave,phi,StaType)
1434            ring['Emat'] = val
1435            ring['Esig'] = esd
1436            ellipse = FitEllipse(R['ImxyObs'].T)
1437            ringxy,ringazm = makeRing(ring['Dcalc'],ellipse,0,0.,scalex,scaley,Image)
1438            ring['ImxyCalc'] = np.array(ringxy).T[:2]
1439            ringint = np.array([float(Image[int(x*scalex),int(y*scaley)]) for y,x in np.array(ringxy)[:,:2]])
1440            ringint /= np.mean(ringint)
1441            ring['Ivar'] = np.var(ringint)
1442            ring['covMat'] = covMat
1443            G2fil.G2Print ('Variance in normalized ring intensity: %.3f'%(ring['Ivar']))
1444    CalcStrSta(StrSta,Controls)
1445   
1446def IntStrSta(Image,StrSta,Controls):
1447    StaControls = copy.deepcopy(Controls)
1448    pixelSize = Controls['pixelSize']
1449    scalex = 1000./pixelSize[0]
1450    scaley = 1000./pixelSize[1]
1451    phi = StrSta['Sample phi']
1452    StaControls['distance'] += StrSta['Sample z']*cosd(phi)
1453    RingsAI = []
1454    for ring in StrSta['d-zero']:       #get observed x,y,d points for the d-zeros
1455        Ring,R = MakeStrStaRing(ring,Image,StaControls)
1456        if len(Ring):
1457            ellipse = FitEllipse(R['ImxyObs'].T)
1458            ringxy,ringazm = makeRing(ring['Dcalc'],ellipse,0,0.,scalex,scaley,Image,5)
1459            XY = np.array(ringxy).T
1460            Th,Azm = GetTthAzm(XY[0],XY[1],Controls)
1461            pola = G2pwd.Polarization(Controls['PolaVal'][0],Th,Azm-90.)[0]
1462            ring['ImxyCalc'] = np.array(ringxy).T[:2]
1463            ringint = np.array([float(Image[int(x*scalex),int(y*scaley)]) for y,x in np.array(ringxy)[:,:2]])
1464            ringint /= np.mean(ringint)
1465            ringint /= pola
1466            G2fil.G2Print (' %s %.3f %s %.3f %s %d'%('d-spacing',ring['Dcalc'],'sig(MRD):',np.sqrt(np.var(ringint)),'# points:',len(ringint)))
1467            RingsAI.append(np.array(list(zip(ringazm,ringint))).T)
1468    return RingsAI
1469   
1470def CalcStrSta(StrSta,Controls):
1471
1472    wave = Controls['wavelength']
1473    phi = StrSta['Sample phi']
1474    StaType = StrSta['Type']
1475    for ring in StrSta['d-zero']:
1476        Eij = ring['Emat']
1477        E = [[Eij[0],Eij[1],0],[Eij[1],Eij[2],0],[0,0,0]]
1478        th,azm = ring['ImtaObs']
1479        th0 = np.ones_like(azm)*npasind(wave/(2.*ring['Dset']))
1480        V = -np.sum(np.sum(E*calcFij(90.,phi,azm,th0).T/1.e6,axis=2),axis=1)
1481        if StaType == 'True':
1482            ring['ImtaCalc'] = np.array([np.exp(V)*ring['Dset'],azm])
1483        else:
1484            ring['ImtaCalc'] = np.array([(V+1.)*ring['Dset'],azm])
1485        dmin = np.min(ring['ImtaCalc'][0])
1486        dmax = np.max(ring['ImtaCalc'][0])
1487        if ring.get('fixDset',True):
1488            if abs(Eij[0]) < abs(Eij[2]):         #tension
1489                ring['Dcalc'] = dmin+(dmax-dmin)/4.
1490            else:                       #compression
1491                ring['Dcalc'] = dmin+3.*(dmax-dmin)/4.
1492        else:
1493            ring['Dcalc'] = np.mean(ring['ImtaCalc'][0])
1494
1495def calcFij(omg,phi,azm,th):
1496    '''    Uses parameters as defined by Bob He & Kingsley Smith, Adv. in X-Ray Anal. 41, 501 (1997)
1497
1498    :param omg: his omega = sample omega rotation; 0 when incident beam || sample surface,
1499        90 when perp. to sample surface
1500    :param phi: his phi = sample phi rotation; usually = 0, axis rotates with omg.
1501    :param azm: his chi = azimuth around incident beam
1502    :param th:  his theta = theta
1503    '''
1504    a = npsind(th)*npcosd(omg)+npsind(azm)*npcosd(th)*npsind(omg)
1505    b = -npcosd(azm)*npcosd(th)
1506    c = npsind(th)*npsind(omg)-npsind(azm)*npcosd(th)*npcosd(omg)
1507    d = a*npsind(phi)+b*npcosd(phi)
1508    e = a*npcosd(phi)-b*npsind(phi)
1509    Fij = np.array([
1510        [d**2,d*e,c*d],
1511        [d*e,e**2,c*e],
1512        [c*d,c*e,c**2]])
1513    return -Fij
1514
1515def FitStrain(rings,p0,dset,wave,phi,StaType):
1516    'Needs a doc string'
1517    def StrainPrint(ValSig,dset):
1518        print ('Strain tensor for Dset: %.6f'%(dset))
1519        ptlbls = 'names :'
1520        ptstr =  'values:'
1521        sigstr = 'esds  :'
1522        for name,fmt,value,sig in ValSig:
1523            ptlbls += "%s" % (name.rjust(12))
1524            ptstr += fmt % (value)
1525            if sig:
1526                sigstr += fmt % (sig)
1527            else:
1528                sigstr += 12*' '
1529        print (ptlbls)
1530        print (ptstr)
1531        print (sigstr)
1532       
1533    def strainCalc(p,xyd,dset,wave,phi,StaType):
1534        E = np.array([[p[0],p[1],0],[p[1],p[2],0],[0,0,0]])
1535        dspo,azm,dsp = xyd
1536        th = npasind(wave/(2.0*dspo))
1537        V = -np.sum(np.sum(E*calcFij(90.,phi,azm,th).T/1.e6,axis=2),axis=1)
1538        if StaType == 'True':
1539            dspc = dset*np.exp(V)
1540        else:
1541            dspc = dset*(V+1.)
1542        return dspo-dspc
1543       
1544    names = ['e11','e12','e22']
1545    fmt = ['%12.2f','%12.2f','%12.2f']
1546    result = leastsq(strainCalc,p0,args=(rings,dset,wave,phi,StaType),full_output=True)
1547    vals = list(result[0])
1548    chisq = np.sum(result[2]['fvec']**2)/(rings.shape[1]-3)     #reduced chi^2 = M/(Nobs-Nvar)
1549    covM = result[1]
1550    covMatrix = covM*chisq
1551    sig = list(np.sqrt(chisq*np.diag(result[1])))
1552    ValSig = zip(names,fmt,vals,sig)
1553    StrainPrint(ValSig,dset)
1554    return vals,sig,covMatrix
1555
1556def FitImageSpots(Image,ImMax,ind,pixSize,nxy,spotSize=1.0):
1557   
1558    def calcMean(nxy,pixSize,img):
1559        gdx,gdy = np.mgrid[0:nxy,0:nxy]
1560        gdx = ma.array((gdx-nxy//2)*pixSize[0]/1000.,mask=~ma.getmaskarray(ImBox))
1561        gdy = ma.array((gdy-nxy//2)*pixSize[1]/1000.,mask=~ma.getmaskarray(ImBox))
1562        posx = ma.sum(gdx)/ma.count(gdx)
1563        posy = ma.sum(gdy)/ma.count(gdy)
1564        return posx,posy
1565   
1566    def calcPeak(values,nxy,pixSize,img):
1567        back,mag,px,py,sig = values
1568        peak = np.zeros([nxy,nxy])+back
1569        nor = 1./(2.*np.pi*sig**2)
1570        gdx,gdy = np.mgrid[0:nxy,0:nxy]
1571        gdx = (gdx-nxy//2)*pixSize[0]/1000.
1572        gdy = (gdy-nxy//2)*pixSize[1]/1000.
1573        arg = (gdx-px)**2+(gdy-py)**2       
1574        peak += mag*nor*np.exp(-arg/(2.*sig**2))
1575        return ma.compressed(img-peak)/np.sqrt(ma.compressed(img))
1576   
1577    def calc2Peak(values,nxy,pixSize,img):
1578        back,mag,px,py,sigx,sigy,rho = values
1579        peak = np.zeros([nxy,nxy])+back
1580        nor = 1./(2.*np.pi*sigx*sigy*np.sqrt(1.-rho**2))
1581        gdx,gdy = np.mgrid[0:nxy,0:nxy]
1582        gdx = (gdx-nxy//2)*pixSize[0]/1000.
1583        gdy = (gdy-nxy//2)*pixSize[1]/1000.
1584        argnor = -1./(2.*(1.-rho**2))
1585        arg = (gdx-px)**2/sigx**2+(gdy-py)**2/sigy**2-2.*rho*(gdx-px)*(gdy-py)/(sigx*sigy)       
1586        peak += mag*nor*np.exp(argnor*arg)
1587        return ma.compressed(img-peak)/np.sqrt(ma.compressed(img))       
1588   
1589    nxy2 = nxy//2
1590    ImBox = Image[ind[1]-nxy2:ind[1]+nxy2+1,ind[0]-nxy2:ind[0]+nxy2+1]
1591    back = np.min(ImBox)
1592    mag = np.sum(ImBox-back)
1593    vals = [back,mag,0.,0.,0.2,0.2,0.]
1594    ImBox = ma.array(ImBox,dtype=float,mask=ImBox>0.75*ImMax)
1595    px = (ind[0]+.5)*pixSize[0]/1000.
1596    py = (ind[1]+.5)*pixSize[1]/1000.
1597    if ma.any(ma.getmaskarray(ImBox)):
1598        vals = calcMean(nxy,pixSize,ImBox)
1599        posx,posy = [px+vals[0],py+vals[1]]
1600        return [posx,posy,spotSize]
1601    else:
1602        result = leastsq(calc2Peak,vals,args=(nxy,pixSize,ImBox),full_output=True)
1603        vals = result[0]
1604        ratio = vals[4]/vals[5]
1605        if 0.5 < ratio < 2.0 and vals[2] < 2. and vals[3] < 2.:
1606            posx,posy = [px+vals[2],py+vals[3]]
1607            return [posx,posy,min(6.*vals[4],spotSize)]
1608        else:
1609            return None
1610   
1611def AutoSpotMask(Image,Masks,Controls,numChans,dlg=None):
1612   
1613    frame = Masks['Frames']
1614    tam = ma.make_mask_none(Image.shape)
1615    if frame:
1616        tam = ma.mask_or(tam,MakeFrameMask(Controls,frame))
1617    LUtth = np.array(Controls['IOtth'])
1618    dtth = (LUtth[1]-LUtth[0])/numChans
1619    esdMul = Masks['SpotMask']['esdMul']
1620    prob = 100.*sc.erf(esdMul/np.sqrt(2.))
1621    print(' Spots greater than %.2f of band intensity are masked'%prob)
1622    mask = ma.make_mask_none(Image.shape)
1623    band = ma.array(Image,mask=ma.nomask)
1624    TA = Make2ThetaAzimuthMap(Controls,(0,Image.shape[0]),(0,Image.shape[1]))[0]    #2-theta array
1625    TThs = np.linspace(LUtth[0],LUtth[1],numChans,False)
1626    for it,TTh in enumerate(TThs):
1627        band.mask = ma.masked_outside(TA,TTh,TTh+dtth).mask+tam
1628        pcmax = np.percentile(band.compressed(),[prob,50.])
1629        mband = ma.masked_greater(band,pcmax[0])
1630        std = ma.std(mband)
1631        anom = ma.masked_greater((band-pcmax[1])/std,esdMul)
1632        mask ^= (anom.mask^band.mask)
1633        if not dlg is None:
1634            GoOn = dlg.Update(it,newmsg='Processed 2-theta rings = %d'%(it))
1635            if not GoOn[0]:
1636                break
1637    return mask
1638
1639def DoPolaCalib(ImageZ,imageData,arcTth):
1640    ''' Determine image polarization by successive integrations with & without preset arc mask.
1641    After initial search, does a set of five with offset azimuth to get mean(std) result.
1642    '''
1643    from scipy.optimize import minimize_scalar
1644    print(' Polarization fit for mask at %.1f deg 2-theta'%arcTth)
1645    data = copy.deepcopy(imageData)
1646    data['IOtth'] = [arcTth-2.,arcTth+2.]
1647    data['fullIntegrate'] = True
1648    data['LRazimuth'] = [0.,360.]
1649    data['outChannels'] = 200
1650    data['outAzimuths'] = 1
1651    Arc = [arcTth,[85.,95.],2.0]
1652    print(' Integration 2-theta test range %.1f - %.1f in 200 steps'%(data['IOtth'][0],data['IOtth'][1]))
1653    print(' Mask azimuth range: %.1f - %.1f'%(Arc[1][0],Arc[1][1]))
1654    Masks = {'Points':[],'Rings':[],'Arcs':[],'Polygons':[],'Frames':[],
1655        'Thresholds':imageData['range'],'SpotMask':{'esdMul':3.,'spotMask':None}}
1656    AMasks = {'Points':[],'Rings':[],'Arcs':[Arc,],'Polygons':[],'Frames':[],
1657        'Thresholds':imageData['range'],'SpotMask':{'esdMul':3.,'spotMask':None}}
1658   
1659    def func(p):
1660        p = min(1.0,max(p,0.0))
1661        data['PolaVal'][0] = p
1662        H0 = ImageIntegrate(ImageZ,data,Masks)[0]
1663        H0A = ImageIntegrate(ImageZ,data,AMasks)[0]
1664        M = np.sum(H0)-np.sum(H0A)
1665        print(' Polarization %.4f, fxn: %.1f'%(p,M))
1666        return M**2
1667    time0 = time.time()
1668    res = minimize_scalar(func,bracket=[1.,.999],tol=.0001)
1669    print(res)
1670    pola = min(1.0,max(res.x,.0))
1671   
1672    Pola = []
1673    for arc in [75,80,85,90,95]:
1674        Arc = [arcTth,[arc,arc+10.],2.0]
1675        AMasks = {'Points':[],'Rings':[],'Arcs':[Arc,],'Polygons':[],'Frames':[],
1676            'Thresholds':imageData['range'],'SpotMask':{'esdMul':3.,'spotMask':None}}
1677        res = minimize_scalar(func,bracket=[pola-.001,pola],tol=.0001)
1678        print(' Mask azimuth range: %.1f - %.1f'%(Arc[1][0],Arc[1][1]))
1679        print(' pola: %.5f'%res.x)
1680        Pola.append(res.x)
1681    Pola = np.array(Pola)
1682    mean = np.mean(Pola)
1683    std = int(10000*np.std(Pola))
1684    print(' Polarization: %.4f(%d)'%(mean,std))
1685    print(' time: %.2fs'%(time.time()-time0))
1686    imageData['PolaVal'][0] = mean
1687   
Note: See TracBrowser for help on using the repository browser.