source: trunk/GSASIIimage.py @ 1185

Last change on this file since 1185 was 1185, checked in by vondreele, 10 years ago

revamp strain fitting - now does rings one at a time giving each one strain coeff.

  • Property svn:eol-style set to native
  • Property svn:keywords set to Date Author Revision URL Id
File size: 35.8 KB
Line 
1# -*- coding: utf-8 -*-
2#GSASII image calculations: ellipse fitting & image integration       
3########### SVN repository information ###################
4# $Date: 2014-01-08 21:26:10 +0000 (Wed, 08 Jan 2014) $
5# $Author: vondreele $
6# $Revision: 1185 $
7# $URL: trunk/GSASIIimage.py $
8# $Id: GSASIIimage.py 1185 2014-01-08 21:26:10Z vondreele $
9########### SVN repository information ###################
10'''
11*GSASIIimage: Image calc module*
12================================
13
14Ellipse fitting & image integration
15
16'''
17
18import math
19import time
20import numpy as np
21import numpy.linalg as nl
22import numpy.ma as ma
23import polymask as pm
24from scipy.optimize import leastsq
25import copy
26import GSASIIpath
27GSASIIpath.SetVersionNumber("$Revision: 1185 $")
28import GSASIIplot as G2plt
29import GSASIIlattice as G2lat
30import GSASIIpwd as G2pwd
31import fellipse as fel
32
33# trig functions in degrees
34sind = lambda x: math.sin(x*math.pi/180.)
35asind = lambda x: 180.*math.asin(x)/math.pi
36tand = lambda x: math.tan(x*math.pi/180.)
37atand = lambda x: 180.*math.atan(x)/math.pi
38atan2d = lambda y,x: 180.*math.atan2(y,x)/math.pi
39cosd = lambda x: math.cos(x*math.pi/180.)
40acosd = lambda x: 180.*math.acos(x)/math.pi
41rdsq2d = lambda x,p: round(1.0/math.sqrt(x),p)
42#numpy versions
43npsind = lambda x: np.sin(x*np.pi/180.)
44npasind = lambda x: 180.*np.arcsin(x)/np.pi
45npcosd = lambda x: np.cos(x*np.pi/180.)
46npacosd = lambda x: 180.*np.arccos(x)/np.pi
47nptand = lambda x: np.tan(x*np.pi/180.)
48npatand = lambda x: 180.*np.arctan(x)/np.pi
49npatan2d = lambda y,x: 180.*np.arctan2(y,x)/np.pi
50   
51def pointInPolygon(pXY,xy):
52    'Needs a doc string'
53    #pXY - assumed closed 1st & last points are duplicates
54    Inside = False
55    N = len(pXY)
56    p1x,p1y = pXY[0]
57    for i in range(N+1):
58        p2x,p2y = pXY[i%N]
59        if (max(p1y,p2y) >= xy[1] > min(p1y,p2y)) and (xy[0] <= max(p1x,p2x)):
60            if p1y != p2y:
61                xinters = (xy[1]-p1y)*(p2x-p1x)/(p2y-p1y)+p1x
62            if p1x == p2x or xy[0] <= xinters:
63                Inside = not Inside
64        p1x,p1y = p2x,p2y
65    return Inside
66   
67def peneCorr(tth,dep,tilt=0.,azm=0.):
68    'Needs a doc string'
69    return dep*(1.-npcosd(abs(tilt*npsind(azm))-tth*npcosd(azm)))
70#    return dep*(1.-npcosd(tth))         #best one
71#    return dep*npsind(tth)             #not as good as 1-cos2Q
72       
73def makeMat(Angle,Axis):
74    '''Make rotation matrix from Angle and Axis
75
76    :param float Angle: in degrees
77    :param int Axis: 0 for rotation about x, 1 for about y, etc.
78    '''
79    cs = cosd(Angle)
80    ss = sind(Angle)
81    M = np.array(([1.,0.,0.],[0.,cs,-ss],[0.,ss,cs]),dtype=np.float32)
82    return np.roll(np.roll(M,Axis,axis=0),Axis,axis=1)
83                   
84def FitEllipse(xy):
85   
86    def ellipse_center(p):
87        ''' gives ellipse center coordinates
88        '''
89        b,c,d,f,g,a = p[1]/2, p[2], p[3]/2, p[4]/2, p[5], p[0]
90        num = b*b-a*c
91        x0=(c*d-b*f)/num
92        y0=(a*f-b*d)/num
93        return np.array([x0,y0])
94   
95    def ellipse_angle_of_rotation( p ):
96        ''' gives rotation of ellipse major axis from x-axis
97        range will be -90 to 90 deg
98        '''
99        b,c,d,f,g,a = p[1]/2, p[2], p[3]/2, p[4]/2, p[5], p[0]
100        return 0.5*npatand(2*b/(a-c))
101   
102    def ellipse_axis_length( p ):
103        ''' gives ellipse radii in [minor,major] order
104        '''
105        b,c,d,f,g,a = p[1]/2, p[2], p[3]/2, p[4]/2, p[5], p[0]
106        up = 2*(a*f*f+c*d*d+g*b*b-2*b*d*f-a*c*g)
107        down1=(b*b-a*c)*( (c-a)*np.sqrt(1+4*b*b/((a-c)*(a-c)))-(c+a))
108        down2=(b*b-a*c)*( (a-c)*np.sqrt(1+4*b*b/((a-c)*(a-c)))-(c+a))
109        res1=np.sqrt(up/down1)
110        res2=np.sqrt(up/down2)
111        return np.array([ res2,res1])
112   
113    xy = np.array(xy)
114    x = np.asarray(xy.T[0])[:,np.newaxis]
115    y = np.asarray(xy.T[1])[:,np.newaxis]
116    D =  np.hstack((x*x, x*y, y*y, x, y, np.ones_like(x)))
117    S = np.dot(D.T,D)
118    C = np.zeros([6,6])
119    C[0,2] = C[2,0] = 2; C[1,1] = -1
120    E, V =  nl.eig(np.dot(nl.inv(S), C))
121    n = np.argmax(np.abs(E))
122    a = V[:,n]
123    cent = ellipse_center(a)
124    phi = ellipse_angle_of_rotation(a)
125    radii = ellipse_axis_length(a)
126    phi += 90.
127    if radii[0] > radii[1]:
128        radii = [radii[1],radii[0]]
129        phi -= 90.
130    return cent,phi,radii
131
132def FitDetector(rings,varyList,parmDict,Print=True):
133    'Needs a doc string'
134       
135    def CalibPrint(ValSig):
136        print 'Image Parameters:'
137        ptlbls = 'names :'
138        ptstr =  'values:'
139        sigstr = 'esds  :'
140        for name,fmt,value,sig in ValSig:
141            ptlbls += "%s" % (name.rjust(12))
142            if name == 'phi':
143                ptstr += fmt % (value%360.)
144            else:
145                ptstr += fmt % (value)
146            if sig:
147                sigstr += fmt % (sig)
148            else:
149                sigstr += 12*' '
150        print ptlbls
151        print ptstr
152        print sigstr       
153       
154    def ellipseCalcD(B,xyd,varyList,parmDict):
155        x,y,dsp = xyd
156        wave = parmDict['wave']
157        if 'dep' in varyList:
158            dist,x0,y0,tilt,chi,dep = B[:6]
159        else:
160            dist,x0,y0,tilt,chi = B[:5]
161            dep = parmDict['dep']
162        if 'wave' in varyList:
163            wave = B[-1]
164        phi = chi-90.               #get rotation of major axis from tilt axis
165        tth = 2.0*npasind(wave/(2.*dsp))
166        phi0 = npatan2d(y-y0,x-x0)
167        dxy = peneCorr(tth,dep,tilt,phi0)
168        ttth = nptand(tth)
169        stth = npsind(tth)
170        cosb = npcosd(tilt)
171        tanb = nptand(tilt)       
172        tbm = nptand((tth-tilt)/2.)
173        tbp = nptand((tth+tilt)/2.)
174        sinb = npsind(tilt)
175        d = dist+dxy
176        fplus = d*tanb*stth/(cosb+stth)
177        fminus = d*tanb*stth/(cosb-stth)
178        vplus = d*(tanb+(1+tbm)/(1-tbm))*stth/(cosb+stth)
179        vminus = d*(tanb+(1-tbp)/(1+tbp))*stth/(cosb-stth)
180        R0 = np.sqrt((vplus+vminus)**2-(fplus+fminus)**2)/2.      #+minor axis
181        R1 = (vplus+vminus)/2.                                    #major axis
182        zdis = (fplus-fminus)/2.
183        Robs = np.sqrt((x-x0)**2+(y-y0)**2)
184        rsqplus = R0**2+R1**2
185        rsqminus = R0**2-R1**2
186        R = rsqminus*npcosd(2.*phi0-2.*phi)+rsqplus
187        Q = np.sqrt(2.)*R0*R1*np.sqrt(R-2.*zdis**2*npsind(phi0-phi)**2)
188        P = 2.*R0**2*zdis*npcosd(phi0-phi)
189        Rcalc = (P+Q)/R
190        M = Robs-Rcalc
191        return M
192       
193    names = ['dist','det-X','det-Y','tilt','phi','dep','wave']
194    fmt = ['%12.3f','%12.3f','%12.3f','%12.3f','%12.3f','%12.3f','%12.5f']
195    p0 = [parmDict[key] for key in varyList]
196    result = leastsq(ellipseCalcD,p0,args=(rings.T,varyList,parmDict),full_output=True,ftol=1.e-8)
197    chisq = np.sum(result[2]['fvec']**2)
198    parmDict.update(zip(varyList,result[0]))
199    vals = list(result[0])
200    sig = list(np.sqrt(chisq*np.diag(result[1])))
201    sigList = np.zeros(7)
202    for i,name in enumerate(names):
203        if name in varyList:
204            sigList[i] = sig[varyList.index(name)]
205    ValSig = zip(names,fmt,vals,sig)
206    if Print:
207        CalibPrint(ValSig)
208    return chisq
209                   
210def ImageLocalMax(image,w,Xpix,Ypix):
211    'Needs a doc string'
212    w2 = w*2
213    sizey,sizex = image.shape
214    xpix = int(Xpix)            #get reference corner of pixel chosen
215    ypix = int(Ypix)
216    if (w < xpix < sizex-w) and (w < ypix < sizey-w) and image[ypix,xpix]:
217        Z = image[ypix-w:ypix+w,xpix-w:xpix+w]
218        Zmax = np.argmax(Z)
219        Zmin = np.argmin(Z)
220        xpix += Zmax%w2-w
221        ypix += Zmax/w2-w
222        return xpix,ypix,np.ravel(Z)[Zmax],max(0.0001,np.ravel(Z)[Zmin])   #avoid neg/zero minimum
223    else:
224        return 0,0,0,0     
225   
226def makeRing(dsp,ellipse,pix,reject,scalex,scaley,image):
227    'Needs a doc string'
228    def ellipseC():
229        'compute estimate of ellipse circumference'
230        if radii[0] < 0:        #hyperbola
231            theta = npacosd(1./np.sqrt(1.+(radii[0]/radii[1])**2))
232            print theta
233            return 0
234        apb = radii[1]+radii[0]
235        amb = radii[1]-radii[0]
236        return np.pi*apb*(1+3*(amb/apb)**2/(10+np.sqrt(4-3*(amb/apb)**2)))
237       
238    cent,phi,radii = ellipse
239    cphi = cosd(phi-90.)        #convert to major axis rotation
240    sphi = sind(phi-90.)
241    ring = []
242    C = int(ellipseC())         #ring circumference
243    for i in range(0,C,1):      #step around ring in 1mm increments
244        a = 360.*i/C
245        x = radii[1]*cosd(a)        #major axis
246        y = radii[0]*sind(a)
247        X = (cphi*x-sphi*y+cent[0])*scalex      #convert mm to pixels
248        Y = (sphi*x+cphi*y+cent[1])*scaley
249        X,Y,I,J = ImageLocalMax(image,pix,X,Y)
250        if I and J and I/J > reject:
251            X += .5                             #set to center of pixel
252            Y += .5
253            X /= scalex                         #convert back to mm
254            Y /= scaley
255            if [X,Y,dsp] not in ring:           #no duplicates!
256                ring.append([X,Y,dsp])
257    if len(ring) < 10:
258        ring = []
259    return ring
260   
261def GetEllipse2(tth,dxy,dist,cent,tilt,phi):
262    '''uses Dandelin spheres to find ellipse or hyperbola parameters from detector geometry
263    on output
264    radii[0] (b-minor axis) set < 0. for hyperbola
265   
266    '''
267    radii = [0,0]
268    ttth = tand(tth)
269    stth = sind(tth)
270    ctth = cosd(tth)
271    cosb = cosd(tilt)
272    tanb = tand(tilt)
273    tbm = tand((tth-tilt)/2.)
274    tbp = tand((tth+tilt)/2.)
275    sinb = sind(tilt)
276    d = dist+dxy
277    if tth+abs(tilt) < 90.:      #ellipse
278        fplus = d*tanb*stth/(cosb+stth)
279        fminus = d*tanb*stth/(cosb-stth)
280        vplus = d*(tanb+(1+tbm)/(1-tbm))*stth/(cosb+stth)
281        vminus = d*(tanb+(1-tbp)/(1+tbp))*stth/(cosb-stth)
282        radii[0] = np.sqrt((vplus+vminus)**2-(fplus+fminus)**2)/2.      #+minor axis
283        radii[1] = (vplus+vminus)/2.                                    #major axis
284        zdis = (fplus-fminus)/2.
285    else:   #hyperbola!
286        f = d*abs(tanb)*stth/(cosb+stth)
287        v = d*(abs(tanb)+tand(tth-abs(tilt)))
288        delt = d*stth*(1.+stth*cosb)/(abs(sinb)*cosb*(stth+cosb))
289        eps = (v-f)/(delt-v)
290        radii[0] = -eps*(delt-f)/np.sqrt(eps**2-1.)                     #-minor axis
291        radii[1] = eps*(delt-f)/(eps**2-1.)                             #major axis
292        if tilt > 0:
293            zdis = f+radii[1]*eps
294        else:
295            zdis = -f
296#NB: zdis is || to major axis & phi is rotation of minor axis
297#thus shift from beam to ellipse center is [Z*sin(phi),-Z*cos(phi)]
298    elcent = [cent[0]+zdis*sind(phi),cent[1]-zdis*cosd(phi)]
299    return elcent,phi,radii
300   
301def GetEllipse(dsp,data):
302    '''uses Dandelin spheres to find ellipse or hyperbola parameters from detector geometry
303    as given in image controls dictionary (data) and a d-spacing (dsp)
304    '''
305    cent = data['center']
306    tilt = data['tilt']
307    phi = data['rotation']
308    dep = data['DetDepth']
309    tth = 2.0*asind(data['wavelength']/(2.*dsp))
310    dxy = peneCorr(tth,dep,tilt)
311    dist = data['distance']
312    return GetEllipse2(tth,dxy,dist,cent,tilt,phi)
313       
314def GetDetectorXY(dsp,azm,data):
315    'Needs a doc string'
316   
317    elcent,phi,radii = GetEllipse(dsp,data)
318    phi = data['rotation']-90.          #to give rotation of major axis
319    tilt = data['tilt']
320    dist = data['distance']
321    cent = data['center']
322    tth = 2.0*asind(data['wavelength']/(2.*dsp))
323    ttth = tand(tth)
324    stth = sind(tth)
325    ctth = cosd(tth)
326    cosb = cosd(tilt)
327    if radii[0] > 0.:
328        sinb = sind(tilt)
329        tanb = tand(tilt)
330        fplus = dist*tanb*stth/(cosb+stth)
331        fminus = dist*tanb*stth/(cosb-stth)
332        zdis = (fplus-fminus)/2.
333        rsqplus = radii[0]**2+radii[1]**2
334        rsqminus = radii[0]**2-radii[1]**2
335        R = rsqminus*cosd(2.*azm-2.*phi)+rsqplus
336        Q = np.sqrt(2.)*radii[0]*radii[1]*np.sqrt(R-2.*zdis**2*sind(azm-phi)**2)
337        P = 2.*radii[0]**2*zdis*cosd(azm-phi)
338        radius = (P+Q)/R
339        xy = np.array([radius*cosd(azm),radius*sind(azm)])
340        xy += cent
341    else:   #hyperbola - both branches (one is way off screen!)
342        sinb = abs(sind(tilt))
343        tanb = abs(tand(tilt))
344        f = dist*tanb*stth/(cosb+stth)
345        v = dist*(tanb+tand(tth-abs(tilt)))
346        delt = dist*stth*(1+stth*cosb)/(sinb*cosb*(stth+cosb))
347        ecc = (v-f)/(delt-v)
348        R = radii[1]*(ecc**2-1)/(1-ecc*cosd(azm))
349        if tilt > 0.:
350            offset = 2.*radii[1]*ecc+f      #select other branch
351            xy = [-R*cosd(azm)-offset,R*sind(azm)]
352        else:
353            offset = -f
354            xy = [-R*cosd(azm)-offset,-R*sind(azm)]
355        xy = -np.array([xy[0]*cosd(phi)+xy[1]*sind(phi),xy[0]*sind(phi)-xy[1]*cosd(phi)])
356        xy += cent
357    return xy
358   
359def GetDetXYfromThAzm(Th,Azm,data):
360    'Needs a doc string'
361    dsp = data['wavelength']/(2.0*npsind(Th))   
362    return GetDetectorXY(dsp,azm,data)
363                   
364def GetTthAzmDsp(x,y,data):
365    'Needs a doc string - checked OK for ellipses & hyperbola'
366    wave = data['wavelength']
367    cent = data['center']
368    tilt = data['tilt']
369    dist = data['distance']/cosd(tilt)
370    x0 = data['distance']*tand(tilt)
371    phi = data['rotation']
372    dep = data['DetDepth']
373    LRazim = data['LRazimuth']
374    azmthoff = data['azmthOff']
375    dx = np.array(x-cent[0],dtype=np.float32)
376    dy = np.array(y-cent[1],dtype=np.float32)
377    D = ((dx-x0)**2+dy**2+data['distance']**2)      #sample to pixel distance
378    X = np.array(([dx,dy,np.zeros_like(dx)]),dtype=np.float32).T
379    X = np.dot(X,makeMat(phi,2))
380    Z = np.dot(X,makeMat(tilt,0)).T[2]
381    tth = npatand(np.sqrt(dx**2+dy**2-Z**2)/(dist-Z))
382    dxy = peneCorr(tth,dep,tilt,npatan2d(dy,dx))
383    DX = dist-Z+dxy
384    DY = np.sqrt(dx**2+dy**2-Z**2)
385    tth = npatan2d(DY,DX) 
386    dsp = wave/(2.*npsind(tth/2.))
387    azm = (npatan2d(dy,dx)+azmthoff+720.)%360.
388    G = D/data['distance']**2       #for geometric correction = 1/cos(2theta)^2 if tilt=0.
389    return np.array([tth,azm,G,dsp])
390   
391def GetTth(x,y,data):
392    'Give 2-theta value for detector x,y position; calibration info in data'
393    return GetTthAzmDsp(x,y,data)[0]
394   
395def GetTthAzm(x,y,data):
396    'Give 2-theta, azimuth values for detector x,y position; calibration info in data'
397    return GetTthAzmDsp(x,y,data)[0:2]
398   
399def GetTthAzmD(x,y,data):
400    '''Give 2-theta, azimuth & d-spacing values for detector x,y position;
401     calibration info in data
402    '''
403    return GetTthAzmDsp(x,y,data)[0:3]
404
405def GetDsp(x,y,data):
406    'Give d-spacing value for detector x,y position; calibration info in data'
407    return GetTthAzmDsp(x,y,data)[3]
408       
409def GetAzm(x,y,data):
410    'Give azimuth value for detector x,y position; calibration info in data'
411    return GetTthAzmDsp(x,y,data)[1]
412       
413def ImageCompress(image,scale):
414    'Needs a doc string'
415    if scale == 1:
416        return image
417    else:
418        return image[::scale,::scale]
419       
420def checkEllipse(Zsum,distSum,xSum,ySum,dist,x,y):
421    'Needs a doc string'
422    avg = np.array([distSum/Zsum,xSum/Zsum,ySum/Zsum])
423    curr = np.array([dist,x,y])
424    return abs(avg-curr)/avg < .02
425
426def EdgeFinder(image,data):
427    '''this makes list of all x,y where I>edgeMin suitable for an ellipse search?
428    Not currently used but might be useful in future?
429    '''
430    import numpy.ma as ma
431    Nx,Ny = data['size']
432    pixelSize = data['pixelSize']
433    edgemin = data['edgemin']
434    scalex = pixelSize[0]/1000.
435    scaley = pixelSize[1]/1000.   
436    tay,tax = np.mgrid[0:Nx,0:Ny]
437    tax = np.asfarray(tax*scalex,dtype=np.float32)
438    tay = np.asfarray(tay*scaley,dtype=np.float32)
439    tam = ma.getmask(ma.masked_less(image.flatten(),edgemin))
440    tax = ma.compressed(ma.array(tax.flatten(),mask=tam))
441    tay = ma.compressed(ma.array(tay.flatten(),mask=tam))
442    return zip(tax,tay)
443   
444def MakeFrameMask(data,frame):
445    pixelSize = data['pixelSize']
446    scalex = pixelSize[0]/1000.
447    scaley = pixelSize[1]/1000.
448    blkSize = 512
449    Nx,Ny = data['size']
450    nXBlks = (Nx-1)/blkSize+1
451    nYBlks = (Ny-1)/blkSize+1
452    tam = ma.make_mask_none(data['size'])
453    for iBlk in range(nXBlks):
454        iBeg = iBlk*blkSize
455        iFin = min(iBeg+blkSize,Nx)
456        for jBlk in range(nYBlks):
457            jBeg = jBlk*blkSize
458            jFin = min(jBeg+blkSize,Ny)               
459            nI = iFin-iBeg
460            nJ = jFin-jBeg
461            tax,tay = np.mgrid[iBeg+0.5:iFin+.5,jBeg+.5:jFin+.5]         #bin centers not corners
462            tax = np.asfarray(tax*scalex,dtype=np.float32)
463            tay = np.asfarray(tay*scaley,dtype=np.float32)
464            tamp = ma.make_mask_none((1024*1024))
465            tamp = ma.make_mask(pm.polymask(nI*nJ,tax.flatten(),
466                tay.flatten(),len(frame),frame,tamp)[:nI*nJ])-True  #switch to exclude around frame
467            if tamp.shape:
468                tamp = np.reshape(tamp[:nI*nJ],(nI,nJ))
469                tam[iBeg:iFin,jBeg:jFin] = ma.mask_or(tamp[0:nI,0:nJ],tam[iBeg:iFin,jBeg:jFin])
470            else:
471                tam[iBeg:iFin,jBeg:jFin] = True
472    return tam.T
473   
474def ImageRecalibrate(self,data,masks):
475    'Needs a doc string'
476    import ImageCalibrants as calFile
477    print 'Image recalibration:'
478    time0 = time.time()
479    pixelSize = data['pixelSize']
480    scalex = 1000./pixelSize[0]
481    scaley = 1000./pixelSize[1]
482    pixLimit = data['pixLimit']
483    cutoff = data['cutoff']
484    data['rings'] = []
485    data['ellipses'] = []
486    if not data['calibrant']:
487        print 'no calibration material selected'
488        return True   
489    skip = data['calibskip']
490    dmin = data['calibdmin']
491    Bravais,Cells = calFile.Calibrants[data['calibrant']][:2]
492    HKL = []
493    for bravais,cell in zip(Bravais,Cells):
494        A = G2lat.cell2A(cell)
495        hkl = G2lat.GenHBravais(dmin,bravais,A)
496        HKL += hkl
497    HKL = G2lat.sortHKLd(HKL,True,False)
498    varyList = ['dist','det-X','det-Y','tilt','phi']
499    parmDict = {'dist':data['distance'],'det-X':data['center'][0],'det-Y':data['center'][1],
500        'tilt':data['tilt'],'phi':data['rotation'],'wave':data['wavelength'],'dep':data['DetDepth']}
501    Found = False
502    wave = data['wavelength']
503    frame = masks['Frames']
504    tam = ma.make_mask_none(self.ImageZ.shape)
505    if frame:
506        tam = ma.mask_or(tam,MakeFrameMask(data,frame))
507    for iH,H in enumerate(HKL): 
508        dsp = H[3]
509        tth = 2.0*asind(wave/(2.*dsp))
510        if tth+abs(data['tilt']) > 90.:
511            print 'next line is a hyperbola - search stopped'
512            break
513        ellipse = GetEllipse(dsp,data)
514        Ring = makeRing(dsp,ellipse,pixLimit,cutoff,scalex,scaley,ma.array(self.ImageZ,mask=tam))
515        if Ring:
516            if iH >= skip:
517                data['rings'].append(np.array(Ring))
518            data['ellipses'].append(copy.deepcopy(ellipse+('r',)))
519            Found = True
520        elif not Found:         #skipping inner rings, keep looking until ring found
521            continue
522        else:                   #no more rings beyond edge of detector
523            break
524    rings = np.concatenate((data['rings']),axis=0)
525    if data['DetDepthRef']:
526        varyList.append('dep')
527    FitDetector(rings,varyList,parmDict)
528    data['distance'] = parmDict['dist']
529    data['center'] = [parmDict['det-X'],parmDict['det-Y']]
530    data['rotation'] = np.mod(parmDict['phi'],360.0)
531    data['tilt'] = parmDict['tilt']
532    data['DetDepth'] = parmDict['dep']
533    N = len(data['ellipses'])
534    data['ellipses'] = []           #clear away individual ellipse fits
535    for H in HKL[:N]:
536        ellipse = GetEllipse(H[3],data)
537        data['ellipses'].append(copy.deepcopy(ellipse+('b',)))
538   
539    print 'calibration time = ',time.time()-time0
540    G2plt.PlotImage(self,newImage=True)       
541    return True
542           
543def ImageCalibrate(self,data):
544    'Needs a doc string'
545    import copy
546    import ImageCalibrants as calFile
547    print 'Image calibration:'
548    time0 = time.time()
549    ring = data['ring']
550    pixelSize = data['pixelSize']
551    scalex = 1000./pixelSize[0]
552    scaley = 1000./pixelSize[1]
553    pixLimit = data['pixLimit']
554    cutoff = data['cutoff']
555    if len(ring) < 5:
556        print 'not enough inner ring points for ellipse'
557        return False
558       
559    #fit start points on inner ring
560    data['ellipses'] = []
561    data['rings'] = []
562    outE = FitEllipse(ring)
563    fmt  = '%s X: %.3f, Y: %.3f, phi: %.3f, R1: %.3f, R2: %.3f'
564    fmt2 = '%s X: %.3f, Y: %.3f, phi: %.3f, R1: %.3f, R2: %.3f, chi^2: %.3f, N: %d'
565    if outE:
566        print fmt%('start ellipse: ',outE[0][0],outE[0][1],outE[1],outE[2][0],outE[2][1])
567        ellipse = outE
568    else:
569        return False
570       
571    #setup 360 points on that ring for "good" fit
572    data['ellipses'].append(ellipse[:]+('g',))
573    Ring = makeRing(1.0,ellipse,pixLimit,cutoff,scalex,scaley,self.ImageZ)
574    if Ring:
575        ellipse = FitEllipse(Ring)
576        Ring = makeRing(1.0,ellipse,pixLimit,cutoff,scalex,scaley,self.ImageZ)    #do again
577        ellipse = FitEllipse(Ring)
578    else:
579        print '1st ring not sufficiently complete to proceed'
580        return False
581    print fmt2%('inner ring:    ',ellipse[0][0],ellipse[0][1],ellipse[1],ellipse[2][0],ellipse[2][1],0.,len(Ring))     #cent,phi,radii
582    data['ellipses'].append(ellipse[:]+('r',))
583    data['rings'].append(np.array(Ring))
584    G2plt.PlotImage(self,newImage=True)
585   
586#setup for calibration
587    data['rings'] = []
588    if not data['calibrant']:
589        print 'no calibration material selected'
590        return True
591   
592    skip = data['calibskip']
593    dmin = data['calibdmin']
594#generate reflection set
595    Bravais,Cells = calFile.Calibrants[data['calibrant']][:2]
596    HKL = []
597    for bravais,cell in zip(Bravais,Cells):
598        A = G2lat.cell2A(cell)
599        hkl = G2lat.GenHBravais(dmin,bravais,A)[skip:]
600        HKL += hkl
601    HKL = G2lat.sortHKLd(HKL,True,False)
602    wave = data['wavelength']
603#set up 1st ring
604    elcent,phi,radii = ellipse              #from fit of 1st ring
605    dsp = HKL[0][3]
606    tth = 2.0*asind(wave/(2.*dsp))
607    Ring0 = makeRing(dsp,ellipse,3,cutoff,scalex,scaley,self.ImageZ)
608    ttth = nptand(tth)
609    stth = npsind(tth)
610    ctth = npcosd(tth)
611#1st estimate of tilt; assume ellipse - don't know sign though
612    tilt = npasind(np.sqrt(max(0.,1.-(radii[0]/radii[1])**2))*ctth)
613    if not tilt:
614        print 'WARNING - selected ring was fitted as a circle'
615        print ' - if detector was tilted we suggest you skip this ring - WARNING'
616#1st estimate of dist: sample to detector normal to plane
617    data['distance'] = dist = radii[0]**2/(ttth*radii[1])
618#ellipse to cone axis (x-ray beam); 2 choices depending on sign of tilt
619    zdisp = radii[1]*ttth*tand(tilt)
620    zdism = radii[1]*ttth*tand(-tilt)
621#cone axis position; 2 choices. Which is right?     
622#NB: zdisp is || to major axis & phi is rotation of minor axis
623#thus shift from beam to ellipse center is [Z*sin(phi),-Z*cos(phi)]
624    centp = [elcent[0]+zdisp*sind(phi),elcent[1]-zdisp*cosd(phi)]
625    centm = [elcent[0]+zdism*sind(phi),elcent[1]-zdism*cosd(phi)]
626#check get same ellipse parms either way
627#now do next ring; estimate either way & do a FitDetector each way; best fit is correct one
628    dsp = HKL[1][3]
629    tth = 2.0*asind(wave/(2.*dsp))
630    ellipsep = GetEllipse2(tth,0.,dist,centp,tilt,phi)
631    print fmt%('plus ellipse :',ellipsep[0][0],ellipsep[0][1],ellipsep[1],ellipsep[2][0],ellipsep[2][1])
632    Ringp = makeRing(dsp,ellipsep,3,cutoff,scalex,scaley,self.ImageZ)
633    parmDict = {'dist':dist,'det-X':centp[0],'det-Y':centp[1],
634        'tilt':tilt,'phi':phi,'wave':wave,'dep':0.0}
635    varyList = ['dist','det-X','det-Y','tilt','phi']
636    if len(Ringp) > 10:
637        chip = FitDetector(np.array(Ring0+Ringp),varyList,parmDict,True)
638        tiltp = parmDict['tilt']
639        phip = parmDict['phi']
640        centp = [parmDict['det-X'],parmDict['det-Y']]
641    else:
642        chip = 1e6
643    ellipsem = GetEllipse2(tth,0.,dist,centm,-tilt,phi)
644    print fmt%('minus ellipse:',ellipsem[0][0],ellipsem[0][1],ellipsem[1],ellipsem[2][0],ellipsem[2][1])
645    Ringm = makeRing(dsp,ellipsem,3,cutoff,scalex,scaley,self.ImageZ)
646    if len(Ringm) > 10:
647        parmDict['tilt'] *= -1
648        chim = FitDetector(np.array(Ring0+Ringm),varyList,parmDict,True)
649        tiltm = parmDict['tilt']
650        phim = parmDict['phi']
651        centm = [parmDict['det-X'],parmDict['det-Y']]
652    else:
653        chim = 1e6
654    if chip < chim:
655        data['tilt'] = tiltp
656        data['center'] = centp
657        data['rotation'] = phip
658    else:
659        data['tilt'] = tiltm
660        data['center'] = centm
661        data['rotation'] = phim
662    data['ellipses'].append(ellipsep[:]+('b',))
663    data['rings'].append(np.array(Ringp))
664    data['ellipses'].append(ellipsem[:]+('r',))
665    data['rings'].append(np.array(Ringm))
666    G2plt.PlotImage(self,newImage=True)
667    parmDict = {'dist':data['distance'],'det-X':data['center'][0],'det-Y':data['center'][1],
668        'tilt':data['tilt'],'phi':data['rotation'],'wave':data['wavelength'],'dep':data['DetDepth']}
669    varyList = ['dist','det-X','det-Y','tilt','phi']
670    if data['DetDepthRef']:
671        varyList.append('dep')
672    data['rings'] = []
673    data['ellipses'] = []
674    for i,H in enumerate(HKL):
675        dsp = H[3]
676        tth = 2.0*asind(wave/(2.*dsp))
677        if tth+abs(data['tilt']) > 90.:
678            print 'next line is a hyperbola - search stopped'
679            break
680        print 'HKLD:',H[:4],'2-theta: %.4f'%(tth)
681        elcent,phi,radii = ellipse = GetEllipse(dsp,data)
682        data['ellipses'].append(copy.deepcopy(ellipse+('g',)))
683        print fmt%('predicted ellipse:',elcent[0],elcent[1],phi,radii[0],radii[1])
684        Ring = makeRing(dsp,ellipse,pixLimit,cutoff,scalex,scaley,self.ImageZ)
685        if Ring:
686            data['rings'].append(np.array(Ring))
687            rings = np.concatenate((data['rings']),axis=0)
688            if i:
689                chisq = FitDetector(rings,varyList,parmDict,False)
690                data['distance'] = parmDict['dist']
691                data['center'] = [parmDict['det-X'],parmDict['det-Y']]
692                data['rotation'] = parmDict['phi']
693                data['tilt'] = parmDict['tilt']
694                data['DetDepth'] = parmDict['dep']
695                elcent,phi,radii = ellipse = GetEllipse(dsp,data)
696                print fmt2%('fitted ellipse:   ',elcent[0],elcent[1],phi,radii[0],radii[1],chisq,len(rings))
697            data['ellipses'].append(copy.deepcopy(ellipse+('r',)))
698#            G2plt.PlotImage(self,newImage=True)
699        else:
700            print 'insufficient number of points in this ellipse to fit'
701            break
702    G2plt.PlotImage(self,newImage=True)
703    fullSize = len(self.ImageZ)/scalex
704    if 2*radii[1] < .9*fullSize:
705        print 'Are all usable rings (>25% visible) used? Try reducing Min ring I/Ib'
706    N = len(data['ellipses'])
707    if N > 2:
708        FitDetector(rings,varyList,parmDict)
709        data['distance'] = parmDict['dist']
710        data['center'] = [parmDict['det-X'],parmDict['det-Y']]
711        data['rotation'] = parmDict['phi']
712        data['tilt'] = parmDict['tilt']
713        data['DetDepth'] = parmDict['dep']
714    for H in HKL[:N]:
715        ellipse = GetEllipse(H[3],data)
716        data['ellipses'].append(copy.deepcopy(ellipse+('b',)))
717    print 'calibration time = ',time.time()-time0
718    G2plt.PlotImage(self,newImage=True)       
719    return True
720   
721def Make2ThetaAzimuthMap(data,masks,iLim,jLim):
722    'Needs a doc string'
723    #transforms 2D image from x,y space to 2-theta,azimuth space based on detector orientation
724    pixelSize = data['pixelSize']
725    scalex = pixelSize[0]/1000.
726    scaley = pixelSize[1]/1000.
727   
728    tay,tax = np.mgrid[iLim[0]+0.5:iLim[1]+.5,jLim[0]+.5:jLim[1]+.5]         #bin centers not corners
729    tax = np.asfarray(tax*scalex,dtype=np.float32)
730    tay = np.asfarray(tay*scaley,dtype=np.float32)
731    nI = iLim[1]-iLim[0]
732    nJ = jLim[1]-jLim[0]
733    #make position masks here
734    frame = masks['Frames']
735    tam = ma.make_mask_none((nI,nJ))
736    if frame:
737        tamp = ma.make_mask_none((1024*1024))
738        tamp = ma.make_mask(pm.polymask(nI*nJ,tax.flatten(),
739            tay.flatten(),len(frame),frame,tamp)[:nI*nJ])-True  #switch to exclude around frame
740        tam = ma.mask_or(tam.flatten(),tamp)
741    polygons = masks['Polygons']
742    for polygon in polygons:
743        if polygon:
744            tamp = ma.make_mask_none((1024*1024))
745            tamp = ma.make_mask(pm.polymask(nI*nJ,tax.flatten(),
746                tay.flatten(),len(polygon),polygon,tamp)[:nI*nJ])
747            tam = ma.mask_or(tam.flatten(),tamp)
748    if tam.shape: tam = np.reshape(tam,(nI,nJ))
749    spots = masks['Points']
750    for X,Y,diam in spots:
751        tamp = ma.getmask(ma.masked_less((tax-X)**2+(tay-Y)**2,(diam/2.)**2))
752        tam = ma.mask_or(tam,tamp)
753    TA = np.array(GetTthAzmD(tax,tay,data))     #includes geom. corr. as dist**2/d0**2
754    TA[1] = np.where(TA[1]<0,TA[1]+360,TA[1])
755    return np.array(TA),tam           #2-theta, azimuth & geom. corr. arrays & position mask
756
757def Fill2ThetaAzimuthMap(masks,TA,tam,image):
758    'Needs a doc string'
759    Zlim = masks['Thresholds'][1]
760    rings = masks['Rings']
761    arcs = masks['Arcs']
762    TA = np.dstack((ma.getdata(TA[1]),ma.getdata(TA[0]),ma.getdata(TA[2])))    #azimuth, 2-theta, dist
763    tax,tay,tad = np.dsplit(TA,3)    #azimuth, 2-theta, dist**2/d0**2
764    for tth,thick in rings:
765        tam = ma.mask_or(tam.flatten(),ma.getmask(ma.masked_inside(tay.flatten(),max(0.01,tth-thick/2.),tth+thick/2.)))
766    for tth,azm,thick in arcs:
767        tamt = ma.getmask(ma.masked_inside(tay.flatten(),max(0.01,tth-thick/2.),tth+thick/2.))
768        tama = ma.getmask(ma.masked_inside(tax.flatten(),azm[0],azm[1]))
769        tam = ma.mask_or(tam.flatten(),tamt*tama)
770    taz = ma.masked_outside(image.flatten(),int(Zlim[0]),Zlim[1])
771    tabs = np.ones_like(taz)
772    tam = ma.mask_or(tam.flatten(),ma.getmask(taz))
773    tax = ma.compressed(ma.array(tax.flatten(),mask=tam))
774    tay = ma.compressed(ma.array(tay.flatten(),mask=tam))
775    taz = ma.compressed(ma.array(taz.flatten(),mask=tam))
776    tad = ma.compressed(ma.array(tad.flatten(),mask=tam))
777    tabs = ma.compressed(ma.array(tabs.flatten(),mask=tam))
778    return tax,tay,taz,tad,tabs
779   
780def ImageIntegrate(image,data,masks,blkSize=128,dlg=None):
781    'Needs a doc string'
782    import histogram2d as h2d
783    print 'Begin image integration'
784    LUtth = data['IOtth']
785    LRazm = np.array(data['LRazimuth'],dtype=np.float64)
786    numAzms = data['outAzimuths']
787    numChans = data['outChannels']
788    Dtth = (LUtth[1]-LUtth[0])/numChans
789    Dazm = (LRazm[1]-LRazm[0])/numAzms
790    if data['centerAzm']:
791        LRazm += Dazm/2.
792    NST = np.zeros(shape=(numAzms,numChans),order='F',dtype=np.float32)
793    H0 = np.zeros(shape=(numAzms,numChans),order='F',dtype=np.float32)
794    imageN = len(image)
795    Nx,Ny = data['size']
796    nXBlks = (Nx-1)/blkSize+1
797    nYBlks = (Ny-1)/blkSize+1
798    Nup = nXBlks*nYBlks*3+3
799    t0 = time.time()
800    Nup = 0
801    if dlg:
802        dlg.Update(Nup)
803    for iBlk in range(nYBlks):
804        iBeg = iBlk*blkSize
805        iFin = min(iBeg+blkSize,Ny)
806        for jBlk in range(nXBlks):
807            jBeg = jBlk*blkSize
808            jFin = min(jBeg+blkSize,Nx)               
809            TA,tam = Make2ThetaAzimuthMap(data,masks,(iBeg,iFin),(jBeg,jFin))           #2-theta & azimuth arrays & create position mask
810           
811            Nup += 1
812            if dlg:
813                dlg.Update(Nup)
814            Block = image[iBeg:iFin,jBeg:jFin]
815            tax,tay,taz,tad,tabs = Fill2ThetaAzimuthMap(masks,TA,tam,Block)    #and apply masks
816            Nup += 1
817            if dlg:
818                dlg.Update(Nup)
819            tax = np.where(tax > LRazm[1],tax-360.,tax)                 #put azm inside limits if possible
820            tax = np.where(tax < LRazm[0],tax+360.,tax)
821            if data['SampleAbs'][1]:
822                muR = data['SampleAbs'][0]*(1.+npsind(tax)**2/2.)/(npcosd(tay))
823                tabs = G2pwd.Absorb(data['SampleShape'],muR,tay)
824            if any([tax.shape[0],tay.shape[0],taz.shape[0]]):
825                NST,H0 = h2d.histogram2d(len(tax),tax,tay,taz*tad*tabs,numAzms,numChans,LRazm,LUtth,Dazm,Dtth,NST,H0)
826            Nup += 1
827            if dlg:
828                dlg.Update(Nup)
829    NST = np.array(NST)
830    H0 = np.divide(H0,NST)
831    H0 = np.nan_to_num(H0)
832    del NST
833    if Dtth:
834        H2 = np.array([tth for tth in np.linspace(LUtth[0],LUtth[1],numChans+1)])
835    else:
836        H2 = np.array(LUtth)
837    if Dazm:       
838        H1 = np.array([azm for azm in np.linspace(LRazm[0],LRazm[1],numAzms+1)])
839    else:
840        H1 = LRazm
841    H0 /= npcosd(H2[:-1])**2
842    if data['Oblique'][1]:
843        H0 /= G2pwd.Oblique(data['Oblique'][0],H2[:-1])
844    if 'SASD' in data['type'] and data['PolaVal'][1]:
845        H0 /= np.array([G2pwd.Polarization(data['PolaVal'][0],H2[:-1],Azm=azm)[0] for azm in H1[:-1]+np.diff(H1)/2.])
846    Nup += 1
847    if dlg:
848        dlg.Update(Nup)
849    t1 = time.time()
850    print 'Integration complete'
851    print "Elapsed time:","%8.3f"%(t1-t0), "s"
852    return H0,H1,H2
853   
854def MakeStrStaRing(ring,Image,Controls):
855    ellipse = GetEllipse(ring['Dset'],Controls)
856    pixSize = Controls['pixelSize']
857    scalex = 1000./pixSize[0]
858    scaley = 1000./pixSize[1]
859    Ring = np.array(makeRing(ring['Dset'],ellipse,ring['pixLimit'],ring['cutoff'],scalex,scaley,Image)).T   #returns x,y,dsp for each point in ring
860    if len(Ring):
861        ring['ImxyObs'] = copy.copy(Ring[:2])
862        TA = GetTthAzm(Ring[0],Ring[1],Controls)       #convert x,y to tth,azm
863        TA[0] = Controls['wavelength']/(2.*npsind(TA[0]/2.))      #convert 2th to d
864        ring['ImtaObs'] = TA
865        ring['ImtaCalc'] = np.zeros_like(ring['ImtaObs'])
866        Ring[0] = TA[0]
867        Ring[1] = TA[1]
868        return Ring,ring
869    else:
870        return [],[]    #bad ring; no points found
871   
872def FitStrSta(Image,StrSta,Controls):
873    'Needs a doc string'
874   
875    StaControls = copy.deepcopy(Controls)
876    phi = StrSta['Sample phi']
877    wave = Controls['wavelength']
878    StaControls['distance'] += StrSta['Sample z']*cosd(phi)
879
880    for ring in StrSta['d-zero']:       #get observed x,y,d points for the d-zeros
881        dset = ring['Dset']
882        Ring,R = MakeStrStaRing(ring,Image,StaControls)
883        if len(Ring):
884            ring.update(R)
885            p0 = ring['Emat']
886            ring['Emat'] = FitStrain(Ring,p0,dset,wave,phi)
887    CalcStrSta(StrSta,Controls)
888   
889def CalcStrSta(StrSta,Controls):
890
891    wave = Controls['wavelength']
892    phi = StrSta['Sample phi']
893    E = StrSta['strain']
894    for ring in StrSta['d-zero']:
895        Eij = ring['Emat']
896        E = [[Eij[0],Eij[1],0],[Eij[1],Eij[2],0],[0,0,0]]
897        th,azm = ring['ImtaObs']
898        th0 = np.ones_like(azm)*npasind(wave/(2.*ring['Dset']))
899        V = 1.+np.sum(np.sum(E*calcFij(90.,phi,azm,th0).T/1.e6,axis=2),axis=1)
900        ring['ImtaCalc'] = np.array([V*ring['Dset'],azm])
901
902def calcFij(omg,phi,azm,th):
903    '''Does something...
904
905    Uses parameters as defined by Bob He & Kingsley Smith, Adv. in X-Ray Anal. 41, 501 (1997)
906
907    :param omg: his omega = sample omega rotation; 0 when incident beam || sample surface,
908        90 when perp. to sample surface
909    :param phi: his phi = sample phi rotation; usually = 0, axis rotates with omg.
910    :param azm: his chi = azimuth around incident beam
911    :param th:  his theta = theta
912    '''
913    a = npsind(th)*npcosd(omg)+npsind(azm)*npcosd(th)*npsind(omg)
914    b = -npcosd(azm)*npcosd(th)
915    c = npsind(th)*npsind(omg)-npsind(azm)*npcosd(th)*npcosd(omg)
916    d = a*npcosd(phi)-b*npsind(phi)
917    e = a*npsind(phi)+b*npcosd(phi)
918    Fij = np.array([
919        [d**2,d*e,c*d],
920        [d*e,e**2,c*e],
921        [c*d,c*e,c**2]])
922    return -Fij*nptand(th)
923
924def FitStrain(rings,p0,dset,wave,phi):
925    'Needs a doc string'
926    def StrainPrint(ValSig,dset):
927        print 'Strain tensor for Dset: %.6f'%(dset)
928        ptlbls = 'names :'
929        ptstr =  'values:'
930        sigstr = 'esds  :'
931        for name,fmt,value,sig in ValSig:
932            ptlbls += "%s" % (name.rjust(12))
933            ptstr += fmt % (value)
934            if sig:
935                sigstr += fmt % (sig)
936            else:
937                sigstr += 12*' '
938        print ptlbls
939        print ptstr
940        print sigstr
941       
942    def strainCalc(p,xyd,dset,wave,phi):
943        E = np.array([[p[0],p[1],0],[p[1],p[2],0],[0,0,0]])
944        dspo,azm,dsp = xyd
945        th = npasind(wave/(2.0*dspo))
946        V = 1.+np.sum(np.sum(E*calcFij(90.,phi,azm,th).T/1.e6,axis=2),axis=1)
947        dspc = dset*V
948        return dspo**2-dspc**2
949       
950    names = ['e11','e12','e22']
951    fmt = ['%12.2f','%12.2f','%12.2f']
952    result = leastsq(strainCalc,p0,args=(rings,dset,wave,phi),full_output=True)
953    vals = list(result[0])
954    chisq = np.sum(result[2]['fvec']**2)
955    sig = list(np.sqrt(chisq*np.diag(result[1])))
956    ValSig = zip(names,fmt,vals,sig)
957    StrainPrint(ValSig,dset)
958    return vals
959   
960       
Note: See TracBrowser for help on using the repository browser.