source: trunk/GSASIIimage.py @ 3842

Last change on this file since 3842 was 3842, checked in by vondreele, 4 years ago

force number of integration bins to be divisible by 4 - solves crash problem from fortran histogram routine

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