source: trunk/GSASIIimage.py @ 1169

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

fix space group ops. cif output

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