ImageD11 Package

ImageD11 Package

ImageD11 - image analysis software for beamline ID11 at the ESRF

Copyright (C) 2005-2011 Jon Wright

This program is free software; you can redistribute it and/or modify it under the terms of the GNU General Public License as published by the Free Software Foundation; either version 2 of the License, or (at your option) any later version.

This program is distributed in the hope that it will be useful, but WITHOUT ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License for more details.

You should have received a copy of the GNU General Public License along with this program; if not, write to the Free Software Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA

ImageD11_file_series Module

To be moved to fabio sometime

ImageD11.ImageD11_file_series.get_options(parser)
ImageD11.ImageD11_file_series.get_series_from_hdf(hdf_file, dark=None, flood=None)
ImageD11.ImageD11_file_series.get_series_from_options(options, args)

Returns a file series thing - not a fabio one

This gives back a fabioimage object with dark and flood corrected data

ImageD11.ImageD11_file_series.get_series_from_stemnum(options, args, dark=None, flood=None)

Returns a file series thing - not a fabio one

ImageD11.ImageD11_file_series.series_from_fabioseries(fabioseries, dark, flood, options)

ImageD11_thread Module

class ImageD11.ImageD11_thread.ImageD11_thread(myname='ImageD11_thread')

Bases: threading.Thread

Add a stopping mechanism for unhandled exceptions

ImageD11_stop_now()
run()

bisplev Module

Interface to the bispev bit of fitpack.

fitpack (dierckx in netlib) — A Python-C wrapper to FITPACK (by P. Dierckx).
FITPACK is a collection of FORTRAN programs for CURVE and SURFACE FITTING with SPLINES and TENSOR PRODUCT SPLINES.
See
http://www.cs.kuleuven.ac.be/cwis/research/nalag/research/topics/fitpack.html
or
http://www.netlib.org/dierckx/index.html

(I think) This code was copied from fitpack.py and modified by J. P. Wright so as to avoid having to install all of scipy. Hence the following notice:

Copyright 2002 Pearu Peterson all rights reserved, Pearu Peterson <pearu@cens.ioc.ee> Permission to use, modify, and distribute this software is given under the terms of the SciPy (BSD style) license. See LICENSE.txt that came with scipy for specifics.

NO WARRANTY IS EXPRESSED OR IMPLIED. USE AT YOUR OWN RISK.

Pearu Peterson

ImageD11.bisplev.bisplev(x, y, tck, dx=0, dy=0)

Evaluate a bivariate B-spline and its derivatives. Description:

Return a rank-2 array of spline function values (or spline derivative values) at points given by the cross-product of the rank-1 arrays x and y. In special cases, return an array or just a float if either x or y or both are floats.
Inputs:
x, y – Rank-1 arrays specifying the domain over which to evaluate the
spline or its derivative.
tck – A sequence of length 5 returned by bisplrep containing the knot
locations, the coefficients, and the degree of the spline: [tx, ty, c, kx, ky].

dx, dy – The orders of the partial derivatives in x and y respectively.

Outputs: (vals, )
vals – The B-pline or its derivative evaluated over the set formed by
the cross-product of x and y.
ImageD11.bisplev.bisplrep(x, y, z, w=None, xb=None, xe=None, yb=None, ye=None, kx=3, ky=3, task=0, s=None, eps=1e-16, tx=None, ty=None, full_output=0, nxest=None, nyest=None, quiet=1)

Find a bivariate B-spline representation of a surface.

Description:

Given a set of data points (x[i], y[i], z[i]) representing a surface z=f(x,y), compute a B-spline representation of the surface.

Inputs:

x, y, z – Rank-1 arrays of data points. w – Rank-1 array of weights. By default w=ones(len(x)). xb, xe – End points of approximation interval in x. yb, ye – End points of approximation interval in y.

By default xb, xe, yb, ye = x[0], x[-1], y[0], y[-1]
kx, ky – The degrees of the spline (1 <= kx, ky <= 5). Third order
(kx=ky=3) is recommended.
task – If task=0, find knots in x and y and coefficients for a given
smoothing factor, s.
If task=1, find knots and coefficients for another value of the
smoothing factor, s. bisplrep must have been previously called with task=0 or task=1.

If task=-1, find coefficients for a given set of knots tx, ty.

s – A non-negative smoothing factor. If weights correspond to the inverse
of the standard-deviation of the errors in z, then a good s-value should be found in the range (m-sqrt(2*m),m+sqrt(2*m)) where m=len(x)
eps – A threshold for determining the effective rank of an over-determined
linear system of equations (0 < eps < 1) — not likely to need changing.

tx, ty – Rank-1 arrays of the knots of the spline for task=-1 full_output – Non-zero to return optional outputs. nxest, nyest – Over-estimates of the total number of knots. If None then

nxest = max(kx+sqrt(m/2),2*kx+3), nyest = max(ky+sqrt(m/2),2*ky+3)

quiet – Non-zero to suppress printing of messages.

Outputs: (tck, {fp, ier, msg})

tck – A list [tx, ty, c, kx, ky] containing the knots (tx, ty) and
coefficients (c) of the bivariate B-spline representation of the surface along with the degree of the spline.

fp – The weighted sum of squared residuals of the spline approximation. ier – An integer flag about splrep success. Success is indicated if

ier<=0. If ier in [1,2,3] an error occurred but was not raised. Otherwise an error is raised.

msg – A message corresponding to the integer flag, ier.

Remarks:

SEE bisplev to evaluate the value of the B-spline given its tck representation.
ImageD11.bisplev.myasarray(a)

blobcorrector Module

Defines a class for correcting peak positions for spatial distortion via a fit2d spline file

To think about doing - valid regions? What if someone uses a 1K spline file for a 2K image etc?

class ImageD11.blobcorrector.correctorclass(argsplinefile, orientation='edf')

Applies a spatial distortion to a peak position using a fit2d splinefile

correct(xin, yin)

Transform x,y in raw image coordinates into x,y of an idealised image. Returns a tuple (x,y), expects a pair of floats as arguments

distort(xin, yin)

Distort a pair of points xnew, ynew to find where they would be in a raw image

Iterative algorithm...

make_pixel_lut(dims)

Generate an x and y image which maps the array indices into floating point array indices (to be corrected for pixel size later)

returns FIXME - check they are the right way around

add some sort of known splinefile testcase
make_pos_lut(dims)

Generate a look up table of pixel positions in microns # Cache the value in case of multiple calls returns ...

readfit2dspline(name)

Reads a fit2d spline file into a scipy/fitpack tuple, tck A fairly long and dull routine...

test(xin, yin)

Checks that the correct and distort functions are indeed inversely related to each other

class ImageD11.blobcorrector.perfect

Bases: ImageD11.blobcorrector.correctorclass

To use on previously corrected when there is no splinefile Allows pixel size etc to be set

correct(xin, yin)

Do nothing - just return the same values

make_pixel_lut(dims)

Generate an x and y image which maps the array indices into floating point array indices (to be corrected for pixel size later)

returns FIXME - check they are the right way around

add some sort of known splinefile testcase
splinefile = 'NO_CORRECTION_APPLIED'
xsize = 'UNKNOWN'
ysize = 'UNKNOWN'
ImageD11.blobcorrector.readfit2dfloats(filep, nfl)

Interprets a 5E14.7 formatted fortran line filep = file object (has readline method) nfl = the number of floats to read

columnfile Module

columnfile represents an ascii file with titles begining “#” and multiple lines of data

An equals sign “=” on a “#” line implies a parameter = value pair

ImageD11.columnfile.bench()

Compares the timing for reading with columfile versus numpy.loadtxt

ImageD11.columnfile.clean(str_lst)

trim whitespace from titles

ImageD11.columnfile.colfile2db(colfilename, dbname)

Read the columnfile into a database Ignores parameter metadata (not used yet)

ImageD11.columnfile.colfile_from_hdf(hdffile, name=None, obj=None)

Read a columnfile from a hdf file FIXME TODO - add the parameters somewhere (attributes??)

ImageD11.columnfile.colfile_to_hdf(colfile, hdffile, name=None, compression=None, compression_opts=None)

Copy a columnfile into hdf file FIXME TODO - add the parameters somewhere (attributes??)

ImageD11.columnfile.colfileobj_to_hdf(cf, hdffile, name=None)

Save a columnfile into hdf file format FIXME TODO - add the parameters somewhere (attributes??)

class ImageD11.columnfile.columnfile(filename=None, new=False)

Class to represent an ascii file containing multiple named columns

addcolumn(col, name)

Add a new column col to the object with name “name” Overwrites in the case that the column already exists

chkarray()

Ensure self.bigarray holds our attributes

copy()

Returns a (deep) copy of the columnfile

filter(mask)

mask is an nrows long array of true/false

getcolumn(name)

Gets data, if column exists

readfile(filename)

Reads in an ascii columned file

removerows(column_name, values, tol=0)

removes rows where self.column_name == values values is a list of values to remove column name should be in self.titles tol is for floating point (fuzzy) comparisons versus integer

set_attributes()

Set object vars to point into the big array

setcolumn(col, name)

Add a new column col to the object with name “name” Overwrites in the case that the column already exists

writefile(filename)

write an ascii columned file

class ImageD11.columnfile.newcolumnfile(titles)

Bases: ImageD11.columnfile.columnfile

Just like a columnfile, but for creating new files

correct Module

ImageD11.correct.correct(data_object, dark=None, flood=None, do_median=False, monitorval=None, monitorcol=None, filterlist=[])

Does the dark and flood corrections Also PIL filters

fft_index_refac Module

Classes for finding lattice vectors

We assume a very large unit cell Compute the Patterson Function using this cell Origin peaks show up for our crystals These will be real space vectors

ImageD11.fft_index_refac.get_options(parser)
class ImageD11.fft_index_refac.grid(np=128, mr=1.0, nsig=5)
fft()

Compute the Patterson

gv_to_grid_new(gv)

Put gvectors into our grid gv - ImageD11.indexing gvectors

peaksearch(peaksfile)

Peaksearch in the Patterson

props()

Print some properties of the Patterson

pv(v)
read_peaks(peaksfile)

Read in the peaks from a peaksearch

reduce(vecs)
slow_score()
ImageD11.fft_index_refac.test(options)

grain Module

class ImageD11.grain.grain(ubi, translation=None, **kwds)
set_ubi(ubi)
ImageD11.grain.read_grain_file(filename)

read ubifile and return a list of ubi arrays

ImageD11.grain.write_grain_file(filename, list_of_grains)

guicommand Module

Interface between Tkinter gui and the actual useful code.

There should be no scientific algorithms (eventually) on the gui side of this class.

This class will eventually offer macro recording capability.

class ImageD11.guicommand.guicommand

Keeps a log of all commands issued - separates gui code from algorithmical code

execute(obj, command, *args, **kwds)

Pass in object as string [peakmerger|transformer|indexer] Pass in command as string, getattr(command) will be used Returns the return value of the function....

TODO
: change this interface???
eg
: works - returns True
you look for self.lastreturned
fails - returns False
you look for self.lasttraceback
getdata(obj, name)

Allows access to “live” data in the objects wrapped

By passing references back you can circumvent the cleanliness of the interface. Please dont.

Returns object.name

gethistory()

Returns the history of commands run by the gui commander

guiindexer Module

Tkinter wrapper for ImageD11 indexing object

All communication should be via parent guicommander object

Owner of the plot3d window

class ImageD11.guiindexer.guiindexer(parent)
assignpeaks()

see indexing.assigntorings

editparameters()
Has the indexing object update its parameter object
eg : savepars(None)

gets a copy of the parameter object Allows user to edit parameters Has the indexing object update itself from the parameter object

eg : loadpars(None)
find()

see indexing.find

histogram_drlv_fit()

Calls indexer.histogram_drlv_fit Plots indexer.bins versus indexer.histogram

loadfileparameters()

see indexing.loadpars and parameters.loadpars

loadgv()

see indexing.readgvfile

makefriedel()

see indexing.friedelpairs

plotxyz()

Gets gv from indexing object Plots the x,y,z (gv) array in a 3D opengl window

saveindexing()

see indexing.saveindexing

saveparameters()

see indexing.savepars and parameters.savepars

saveubis()

see indexing.saveubis

scorethem()

see indexing.scorethem

guimaker Module

Script to put the menus together and build an appli with each bit being relatively clutterfree

Each of guiindexer, guipeaksearch and guitransformer offer members menuitems to put in menus. The overall gui (imaged11_gui) inherits from this I think

# The code was copied from a book somewhere # probably programming python or maybe cookbook. # Most probably programming python by Mark Lutz

class ImageD11.guimaker.GuiMaker(parent=None)

Bases: Tkinter.Frame

You must inherit from this class and implement the start and makeWidgets methods

addMenuItems(menu, items)
makeMenuBar()
menuBar = []

guipeaksearch Module

class ImageD11.guipeaksearch.guipeaksearcher(parent, quiet='No')

Tkinter wrapper for ImageD11 guipeakmerge object

Note that peaksearching is done via the command line

All communication should be via parent guicommander object

filter()

calls peakmerger.filter (does very little) plots x and y final peak positions

TODO implement filters!!!

harvestpeaks()

gets range from 2d plot (image numbers and omegas) calls peakmerger.harvestpeaks(image_number_range,omega_range)

mergepeaks()

calls peakmerger.mergepeaks and reports number of peaks to user

readpeaks(filename=None)

Runs peakmerger.readpeaks gets names for first and last image for plot title gets omega angles and image numbers for plot plots image number versus omega angle to see if you did a scan

savepeaks(filename=None)

see peakmerger.savepeaks

searchraw()

Explains to user about the command line script

guitransformer Module

class ImageD11.guitransformer.guitransformer(parent, quiet='No')
addcellpeaks()
chooseyz()

choose the columns to use for x / y on detector

computegv()
editparameters()

Gets a copy of the parameter object Allows user to edit parameters

filterhisto()

Call plot histo, then filter on it

fit()
loadfileparameters()
loadfiltered()
plotcols()
plothisto(nbins=None)
plotreta()
plotyz()

Plots the x,y arrays being used

savecolfile()
savegv()
saveparameters(filename=None)
write_graindex_gv()
write_pyFAI()

gv_general Module

Find a generalised way to compute g-vectors backwards and forwards

ImageD11.gv_general.angmod(a)
ImageD11.gv_general.axis_from_matrix(m)

Factory method to create a rotation_axis object from a direction cosine matrix http://en.wikipedia.org/wiki/Rotation_representation_%28mathematics%29

ImageD11.gv_general.chimat(c)
ImageD11.gv_general.chiwedge(chi=0.0, wedge=0.0)

in transform: g = omega . chi . wedge .k

ImageD11.gv_general.g_to_k(g, wavelength, axis=array([ 0., 0., 1.]), pre=None, post=None)
Get the k and angles given the g in
g = pre . rot(axis, angle) . post . k g = omega . chi . wedge . k

...where k will satisfy the Laue equations

The recipe is:
Find the components of g along and perpendicular to the rotation axis
co-ords are a0 = axis,
a1 = axis x g, a2 = a1 x axis = ( axis x g ) x axis

Rotated vector will be a0 + a1 sin(angle) + a2 cos(angle) Laue condition says that [incident beam].[k] = k.sin(theta)

= 2 sin^2(theta) / lambda = sin(t)/d = |g|.sin(t) = |g|*|g|*lambda / 2

Apply any post rotations to the incident beam Apply any pre-rotations to the g vector |g| = [a0 + a1 sin(angle) + a2 cos(angle) ] . [incident beam] => solve for angle

http://www.ping.be/~ping1339/gonio.htm a sin(u) + b cos(u) = c let: tan(u’) = -b/a and: A = a / cos(u’) Finally you’ll find:

A.sin(u-u’) = c A.cos(pi/2 - u - u’) = c
ImageD11.gv_general.k_to_g(k, angles, axis=None, pre=None, post=None)

Computes g = pre . rot(axis, angle) . post . k Typically in ImageD11 post = [wedge][chi] and pre = identity

since g = Omega.Chi.Wedge.k or k = Wedge-1.Chi-1.Omega-1.g that is to say omega is directly under the sample and wedge is at the bottom of the stack
class ImageD11.gv_general.rotation_axis(direction, angle=0.0)
Represents an axis - so far includes:
axis/angle matrix (direction cosines)

potentially can add quaternions and euler etc later?

rotate_vectors(vectors, angles=None)

Given a list of vectors, rotate them to new ones Angle is either self.angle, or 1 per vector, in degrees

http://www.ecse.rpi.edu/Homepages/wrf/Research/Short_Notes/rotation.html p’ = a . p a + (p - a . p a) cos q + a * p sin q

= p cos q + a . p a (1-cos q) + a * p sin q point p, axis a, angle q

http://mathworld.wolfram.com/RotationFormula.html r’ = r cos(t) + n(n.r)(1-cos(t)) + rxn sin(t)

rotate_vectors_inverse(vectors, angles=None)

Same as rotate vectors, but opposing sense

to_matrix()

Returns a 3x3 rotation matrix R = I3cos(t) + ww^T (1-cos(t)) - W sin(t) t = angle W = vector with hat operator = [ 0 -wz wy

wz 0 -wx
-wy wx 0 ]
ImageD11.gv_general.wedgechi(wedge=0.0, chi=0.0)

in transform: g = omega . chi . wedge .k

ImageD11.gv_general.wedgemat(w)

hist Module

Wrapper to ImageD11 histogram generators (c functions)

Eventually to be used to generating sensible image color maps for plotedf

Also try to do a fast median filter based on pixel histograms

ImageD11.hist.hist(data, verbose=0)

Wrapper to ImageD11 histogram generators (c functions) returns histogram if data is a 2D numpy.uint16 array

ImageD11.hist.test_dvhist()

indexing Module

ImageD11.indexing.calc_drlv2(UBI, gv)

Get the difference from being integer hkls UBI: ubi matrix gv: list of g-vectors returns drlv2 = (h_calc - h_int)^2

class ImageD11.indexing.indexer(unitcell=None, gv=None, cosine_tol=0.002, minpks=10, hkl_tol=0.01, ring_1=1, ring_2=2, ds_tol=0.005, wavelength=-1, uniqueness=0.5, eta_range=0.0, max_grains=100)

A class for searching for orientation matrices

assigntorings()

Assign the g-vectors to hkl rings

coverage()

Compute the expected coverage of reciprocal space use the min/max obs values of xp/yp/omega to work out what was measured in the scan? No lambda or

fight_over_peaks()

Get the best ubis from those proposed

find()

Dig out the potential hits

friedelpairs(filename)

Attempt to identify Freidel pairs

Peaks must be assigned to the same powder ring Peaks will be the closest thing to being 180 degrees apart

getind(UBI, tol=None)

Returns the indices of peaks in self.gv indexed by matrix UBI

histogram_drlv_fit(UBI=None, bins=None)

Generate a histogram of |drlv| for a ubi matrix For use in validation of grains

loadpars(filename=None)
out_of_eta_range(eta)

decide if an eta is going to be kept

readgvfile(filename, quiet=False)
refine(UBI)

Refine an orientation matrix and rescore it.

From Paciorek et al Acta A55 543 (1999) UB = R H-1

where: R = sum_n r_n h_n^t H = sum_n h_n h_n^t r = g-vectors h = hkl indices
saveindexing(filename, tol=None)

Save orientation matrices

FIXME
: refactor this into something to do
grain by grain } peak by peak }
savepars(filename=None)
saveubis(filename)

Save the generated ubi matrices into a text file

score(UBI, tol=None)

Decide which are the best orientation matrices

scorethem()
updateparameters()
ImageD11.indexing.mod_360(theta, target)

Find multiple of 360 to add to theta to be closest to target

ImageD11.indexing.myhistogram(data, bins)

The numpy histogram api was changed So here is an api that will not change It is based on that from the old Numeric manual

ImageD11.indexing.readubis(ubifile)

read ubifile and return a list of ubi arrays

ImageD11.indexing.refine(UBI, gv, tol, quiet=True)

Refine an orientation matrix and rescore it.

From Paciorek et al Acta A55 543 (1999)
UB = R H-1
where:
R = sum_n r_n h_n^t H = sum_n h_n h_n^t r = g-vectors h = hkl indices
ImageD11.indexing.ubitoB(ubi)

give the B matrix from ubi

ImageD11.indexing.ubitoRod(ubi)

TODO Testcases!!!

ImageD11.indexing.ubitoU(ubi)

convert ubi to Grainspotter style U The convention is B as being triangular, hopefully as Busing and Levy TODO - make some testcases please!!

ImageD11.indexing.ubitocellpars(ubi)

convert ubi matrix to unit cell

ImageD11.indexing.write_ubi_file(filename, ubilist)

save 3x3 matrices into file

labelimage Module

Class to wrap the connectedpixels c extensions for labelling blobs in images.

ImageD11.labelimage.flip1(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip2(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip3(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip4(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip5(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip6(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip7(x, y)

fast, slow to dety, detz

ImageD11.labelimage.flip8(x, y)

fast, slow to dety, detz

class ImageD11.labelimage.labelimage(shape, fileout=<colorama.ansitowin32.StreamWrapper object at 0x000000000220F0F0>, spatial=<ImageD11.blobcorrector.perfect instance at 0x00000000084DD348>, flipper=<function flip2 at 0x00000000083E0CF8>, sptfile=<colorama.ansitowin32.StreamWrapper object at 0x000000000220F0F0>)

For labelling spots in diffraction images

finalise()

Write out the last frame

format = ' %.4f %.4f %.4f %.0f %.4f %.4f %.4f %.4f %.4f %.4f %.4f %.4f %.4f %.4f %.4f %.4f %.0f %.0f %.4f %.0f %.0f %.0f %.0f %.4f %.4f %.4f %.4f %d %d %d\n'
mergelast()

Merge the last two images searches

output2dpeaks(file_obj)

Write something compatible with the old ImageD11 format which fabian is reading. This is called before mergelast, so we write self.npk/self.res

outputpeaks(peaks)

Peaks are in Numeric arrays nowadays

peaksearch(data, threshold, omega)

# Call the c extensions to do the peaksearch, on entry: # # data = 2D Numeric array (of your data) # threshold = float - pixels above this number are put into objects

titles = '# sc fc omega Number_of_pixels avg_intensity s_raw f_raw sigs sigf covsf sigo covso covfo sum_intensity sum_intensity^2 IMax_int IMax_s IMax_f IMax_o Min_s Max_s Min_f Max_f Min_o Max_o dety detz onfirst onlast spot3d_id\n'

lattice_reduction Module

exception ImageD11.lattice_reduction.BadVectors

Bases: exceptions.Exception

Raised by lattice class when vectors are coplanar or colinear

ImageD11.lattice_reduction.checkvol(v1, v2, v3, min_vec2=1e-18)

Check whether vectors are singular

ImageD11.lattice_reduction.cosangle_vec(ubi, v)

Angle between v in real and reciprocal space eg, is a* parallel to a or not?

ImageD11.lattice_reduction.find_lattice(vecs, min_vec2=1, n_try=None, test_vecs=None, tol=0.1, fraction_indexed=0.9, noisy=False)

vecs - vectors to use to generate the lattice min_vec2 - squared length of min_vec (units as vec) n_try - max number of vecs to test (clips vecs for generating) test_vecs - vectors to index with the unit cell (else use vecs) fraction_indexed - whether or not to return cells

ImageD11.lattice_reduction.fparl(x, y)

fraction of x parallel to y

ImageD11.lattice_reduction.get_options(parser)
ImageD11.lattice_reduction.iter3d(n)

Generate all possible unordered combinations of vectors i,j,k for i,j,k < n

This looping was rewritten thanks to: TAOCP V4 fascicle 3, section 7.2.1.3.

It gives much nicer ordering than previously, as is gradually expands down the list, instead of hammering the start

ImageD11.lattice_reduction.iter3d_old(n)

Generate all possible unordered combinations of vectors i,j,k for i,j,k < n

class ImageD11.lattice_reduction.lattice(v1, v2, v3, direction=None, min_vec2=1e-18)

Bases: object

Represents a 3D crystal lattice built from 3 vectors

flip(v)

See also __init__.__doc__

matrix(direction)
nearest(vecs)

Give back the nearest lattice point indices, in the same direction

remainders(vecs)

Difference between vecs and closest lattice points

score(vecs, tol=0.1, debug=False)

How many peaks have integer error less than tol?

withvec(x, direction='col')

Try to fit x into the lattice Make the remainder together with current vectors Index it as hkl indices whichever vector has the biggest projection is replaced remake the lattice with these 3 vectors

ImageD11.lattice_reduction.mod(x, y)

Returns the part of x with integer multiples of y removed assert that ||mod(x,y)|| <= ||x|| for all y

ImageD11.lattice_reduction.reduce(v1, v2, v3, min_vec2=1e-18)

Try to find a reduced lattice basis of v1, v2, v3

ImageD11.lattice_reduction.rsweep(vl)

One sweep subtracting each from other two This idea comes from Knuth TAOCP sec 3.3.4C

ImageD11.lattice_reduction.search_2folds(ubi)

Inspired by the Yvon Lepage’s method for finding lattice symmetry Check for 2 fold axes by measuring the directions between real and reciprocal vectors with the same indices. In the case of 2 fold axes they should be parallel

ImageD11.lattice_reduction.sortvec_len(vl)

Sorts according to length (d-s-u)

ImageD11.lattice_reduction.sortvec_xyz(vl)

For each choose between x, -x Sorts according to x then y then z components

license Module

listdialog Module

class ImageD11.listdialog.columnchooser(parent, items, title='Choose two columns')

Bases: ImageD11.listdialog.listdialog

Dialog box for setting detector parameters Takes a list of strings and numbers

class ImageD11.listdialog.listdialog(parent, title=None, items=None, logic=None)

Bases: Tkinter.Toplevel

Dialog box for setting detector parameters Takes a list of strings and numbers

apply()
body(master, items, logic=None)
buttonbox()
cancel(event=None)
ok(event=None)
validate()

parameters Module

Class to handle groups of parameters to be saved in and out of files and edited in guis with the fixed/varied info etc

class ImageD11.parameters.par(name, value, helpstring=None, vary=False, can_vary=False, stepsize=None)

Represents a thing which can vary

fromstringlist(sl)

to send to Java

tostringlist()

to catch from Java

class ImageD11.parameters.parameters(**kwds)

Class to hold a set of named parameters

addpar(par)

add a parameter object

dumbtypecheck()

Eventually parameter types (and units and fixed/varied to be specifieable For now it just tries to coerce to float, then does nothing

get(name)
get_parameters()

Returns a dictionary of parameters

get_variable_list()
get_variable_stepsizes()

stepsizes for optimisers

get_variable_values()

values of the parameters

loadparameters(filename)

Load parameters from a file

saveparameters(filename)

Write parameters to a file

set(name, value)
set_parameters(d)

Updates the values of parameters

set_variable_values(values)

set values of the parameters

set_varylist(vl)
update_other(other)

Synchronise an object with the values in this object

update_yourself(other)

Sychronise this parameter objects list of values with another object

ImageD11.parameters.read_par_file(filename)

peakmerge Module

Class for reading peaksearch files, merging peaks on adjacent frames and writing out the result...

This is mostly replaced with the 3D blobsearcher - this algorithm uses a lot of memory

class ImageD11.peakmerge.peak(line, omega, threshold, num, tolerance)

Represents a peak

combine(other)

Combine this peak with another peak (eg merge!) # # If the thresholds are different, we see the same peaks # twice. We will just pick the one with the lower threshold # # Otherwise we try to average them #

class ImageD11.peakmerge.peakmerger(quiet='No')

The useful class - called by the gui to process peaksearch output into spots

filter()

Does nothing really, just makes a self.finalpeaks array

TODO: pass functions to filter/compress finalpeaks

getheaderinfo(key)

try to find “key” in the headers and return it

getheaderkeys()

things in the headers

harvestpeaks(numlim=None, omlim=None, thresholds=None)

Harvests the peaks from the images within a range of imagenumbers and/or omega eg: it actually reads the numbers now

mergepeaks()

Merges peaks together - if they overlap and are on adjacent frames

More complex than initially planned

readpeaks(filename, startom=0.0, omstep=1.0)

Read in the output from the peaksearching script Filename is the output file optionally startom and omstep fill in omega ONLY if missing in file

savepeaks(filename)

# Write out minimal information # list of xcorr,ycorr,omega, try to keep intensity now

setpixeltolerance(tolerance=2)

TODO ah... if only

class ImageD11.peakmerge.pkimage(name)

Contains the header information from an image Also pointers to positions of lines in peaksearch output file

addtoheader(headerline)

save header info in dict headerline is line starting with # in peaksearch output

otherheaderstuff(headerline)

Hard to parse header info - called by addtoheader

ImageD11.peakmerge.roundfloat(x, tol)

Return the float nearest to x stepsize tol

peaksearcher Module

Script for peaksearching images from the command line

Uses the connectedpixels extension for finding blobs above a threshold and the blobcorrector(+splines) for correcting them for spatial distortion

Defines one function (peaksearch) which might be reused

ImageD11.peaksearcher.get_help(usage=True)

return the help string for online help

ImageD11.peaksearcher.get_options(parser)

Add our options to a parser object

ImageD11.peaksearcher.peaksearch(filename, data_object, corrector, thresholds, labims)

filename : The name of the image file for progress info data_object : Fabio object containing data and header corrector : spatial and dark/flood , linearity etc

thresholds : [ float[1], float[2] etc]

labims : label image objects, one for each threshold

ImageD11.peaksearcher.peaksearch_driver(options, args)

To be called with options from command line

class ImageD11.peaksearcher.timer
msg(msg)
tick(msg='')
tock(msg='')

plot3d Module

from example by Tarn Weisner Burton <twburton@users.sourceforge.net> in pyopengl

class ImageD11.plot3d.myOpengl(master=None, cnf={}, **kw)

Bases: OpenGL.Tk.Opengl

tkRedraw(*dummy)

Cause the opengl widget to redraw itself.

class ImageD11.plot3d.plot3d(parent, data=None, lines=None, ubis=None, image=None, pars=None, spline=None)

Bases: Tkinter.Toplevel

changedata(xyz)
go()

Allow the toplevel to return a handle for changing data

goaway()
readimage(image)
readprms(prms)
readspline(spline)
readubis(ubis)
redraw(o)
scorecolor(i=0)
setps()
setupTexture()

rc_array Module

Row/column array type

Inherits from numpy.ndarray

To be used for lattice reduction to apply to g-vectors and patterson peaks equally well and in a coherent way.

class ImageD11.rc_array.rc_array

Bases: numpy.ndarray

Row/Column array Represent a list of row or column vectors

check()

Ensure we have an rc_array which is well behaved Pattern assert(check(v)) should disappear in optimiser

flip(mat)

Flip row to col or col.row

inv()

Inverse matrix of self

nb_vector_axis()

The axis which has the n on it

norm2()

sum(v*v,axis=? for row or col

nvectors()
other_direction()

The one which is not self.direction

vector_axis()

The axis which has the 3 on it

refinegrains Module

ImageD11.refinegrains.compute_lp_factor(colfile, **kwds)

We’ll assume for simplicity that wedge, chi, tilts are all zero Of course we should put that in sometime in the near future

Try using:

Table 6.2.1.1 in Volume C of international tables : (h) Rotation photograph of small crystal, volume V

1. Beam normal to axis rho =

rac{ N^2 e^4 lambda^3 V} { 4 pi m^2 c^4} .

rac{ 1 + cos^2{ heta} } { sin{2 heta } .

rac{ cos{ heta} { sqrt{ cos^2{psi} - sin^2{ heta} } }
p’ |F|^2

rho is Integrated reflection power ratio from a crystal element e, m Electronic charge and mass c Speed of light lambda Wavelength of radiation 2 heta Angle between incident and diffracted beams psi in (h), latitude of reciprocal-lattice point relative to axis of rotation V Volume of crystal, or of irradiated part of powder sample N Number of unit cells per unit volume F Structure factor of hkl reflection p’ Multiplicity factor for single-crystal methods

nb - what is this ? Hopefully would be 1 ?

From this we will take the (1+cos^2(2t))cos(2t)/sqrt(cos^2p-sin^2t)/sin2t

ImageD11.refinegrains.compute_total_intensity(colfile, indices, tth_range, ntrim=2, quiet=False)

Add up the intensity in colfile for peaks given in indices

ImageD11.refinegrains.cosd(x)
ImageD11.refinegrains.cubic(cp)

a=b=c, alpha=beta=gamma=90

ImageD11.refinegrains.hexagonal(cp)

a=b,c, alpha=beta=90,gamma=120

ImageD11.refinegrains.lf(tth, eta)
ImageD11.refinegrains.monoclinic_a(cp)
ImageD11.refinegrains.monoclinic_b(cp)
ImageD11.refinegrains.monoclinic_c(cp)
ImageD11.refinegrains.orthorhombic(cp)

a=b, c, 90,90,90

class ImageD11.refinegrains.refinegrains(tolerance=0.01, intensity_tth_range=(6.1, 6.3), latticesymmetry=<function triclinic at 0x000000000C1F4C88>, OmFloat=True, OmSlop=0.25)

Class to refine the position and orientation of grains and also the detector parameters

applyargs(args)
assignlabels(quiet=False)

Fill out the appropriate labels for the spots

compute_gv(thisgrain, update_columns=False)

Makes self.gv refer be g-vectors computed for this grain in this scan

estimate_steps(gof, guess, steps)
fit(maxiters=100)

Fit the global parameters

generate_grains()
getgrains()
gof(args)

<drlv> for all of the grains in all of the scans

loadfiltered(filename)

Read in file containing filtered and merged peak positions

loadparameters(filename)
makeuniq(symmetry)

Flip orientations to a particular choice

Be aware that if the unit cell is distorted and you do this, you might have problems...

pars = {'wedge': 0.0, 'o12': 0, 'o11': 1, 'fit_tolerance': 0.5, 'cell_alpha': 90.0, 'cell_gamma': 90.0, 'wavelength': 0.2646, 'y_size': 4.6, 'chi': 0.0, 't_z': 0.0, 't_x': 33.0198146824, 't_y': 14.6384893741, 'tilt_z': 0.00995011461696, 'tilt_x': 0.0, 'tilt_y': 0.0623952920101, 'z_size': 4.6, 'cell_lattice_[P,A,B,C,I,F,R]': 'P', 'z_center': 517.007049626, 'omegasign': 1, 'o22': 1, 'o21': 0, 'y_center': 732.950204632, 'cell_beta': 90.0, 'distance': 7367.8452302, 'cell__a': 4.1569162, 'cell__c': 4.1569162, 'cell__b': 4.1569162}
printresult(arg)
readubis(filename)

Read ubi matrices from a text file

refine(ubi, quiet=True)

Fit the matrix without changing the peak assignments

refinepositions(quiet=True, maxiters=100)
refineubis(quiet=True, scoreonly=False)
reset_labels(scanname)
savegrains(filename, sort_npks=True)

Save the refined grains

saveparameters(filename)
set_translation(gr, sc)
stepsizes = {'wedge': 0.0017453292519943296, 'z_center': 0.2, 'wavelength': 0.001, 'y_center': 0.2, 'distance': 200.0, 'chi': 0.0017453292519943296, 't_z': 0.2, 't_x': 0.2, 't_y': 0.2, 'tilt_z': 0.0017453292519943296, 'tilt_x': 0.0017453292519943296, 'tilt_y': 0.0017453292519943296}
ImageD11.refinegrains.sind(x)
ImageD11.refinegrains.test_benoit()
ImageD11.refinegrains.test_nac()
ImageD11.refinegrains.tetragonal(cp)

a=b, c, 90,90,90

ImageD11.refinegrains.triclinic(cp)
ImageD11.refinegrains.trigonalH(cp)

a=b,c, alpha=beta=90,gamma=120

ImageD11.refinegrains.trigonalP(cp)

a=b=c, alpha=beta=gamma

rsv Module

Reciprocal Space Volume

This class is just for holding a volume and loading/saving it to disk Another class (rsv_mapper) will add images into it.

ImageD11.rsv.getbounds(vol, plane)

Returns the extent argument to use for pylab.imshow when plotting a plane

ImageD11.rsv.mem()

debug the memory usage

ImageD11.rsv.readvol(filename, savespace=False)

Read volume from a file returns an rsv object Take care to allocate and read to avoid temporaries

class ImageD11.rsv.rsv(dimensions, bounds, np, **kwds)

Bases: object

A reciprocal space volume

allocate_vol()

Allocates memory for a volume data

normalise(savespace=True)

Return the normalised but avoid divide by zero

plnames = {0: 0, 1: 1, 2: 2, 'H': 0, 'h': 0, 'k': 1, 'K': 1, 'L': 2, 'l': 2}
slice(plane, num)

return signal on plane index num

ImageD11.rsv.writenormedvol(vol, filename)

Write volume in vol to filename - save only the normalised to avoid using so much memory

Compress -1 is for the zeros, which there might be a lot of

ImageD11.rsv.writevol(vol, filename)

Write volume in vol to filename

Compress -1 is for the zeros, which there might be a lot of

rsv_mapper Module

Reciprocal space volume mapper Transfers images into reciprocal space by pixel mapping

ImageD11.rsv_mapper.get_options(parser)

Command line interface for making a mapping Add our options to a parser object

ImageD11.rsv_mapper.main()

A user interface

class ImageD11.rsv_mapper.rsv_mapper(dims, pars, ubi, splinefile=None, np=16, border=10, omegarange=[0, 1, 2, 3, 4, 5, 6, 7, 8, 9, 10, 11, 12, 13, 14, 15, 16, 17, 18, 19, 20, 21, 22, 23, 24, 25, 26, 27, 28, 29, 30, 31, 32, 33, 34, 35, 36, 37, 38, 39, 40, 41, 42, 43, 44, 45, 46, 47, 48, 49, 50, 51, 52, 53, 54, 55, 56, 57, 58, 59, 60, 61, 62, 63, 64, 65, 66, 67, 68, 69, 70, 71, 72, 73, 74, 75, 76, 77, 78, 79, 80, 81, 82, 83, 84, 85, 86, 87, 88, 89, 90, 91, 92, 93, 94, 95, 96, 97, 98, 99, 100, 101, 102, 103, 104, 105, 106, 107, 108, 109, 110, 111, 112, 113, 114, 115, 116, 117, 118, 119, 120, 121, 122, 123, 124, 125, 126, 127, 128, 129, 130, 131, 132, 133, 134, 135, 136, 137, 138, 139, 140, 141, 142, 143, 144, 145, 146, 147, 148, 149, 150, 151, 152, 153, 154, 155, 156, 157, 158, 159, 160, 161, 162, 163, 164, 165, 166, 167, 168, 169, 170, 171, 172, 173, 174, 175, 176, 177, 178, 179, 180, 181, 182, 183, 184, 185, 186, 187, 188, 189, 190, 191, 192, 193, 194, 195, 196, 197, 198, 199, 200, 201, 202, 203, 204, 205, 206, 207, 208, 209, 210, 211, 212, 213, 214, 215, 216, 217, 218, 219, 220, 221, 222, 223, 224, 225, 226, 227, 228, 229, 230, 231, 232, 233, 234, 235, 236, 237, 238, 239, 240, 241, 242, 243, 244, 245, 246, 247, 248, 249, 250, 251, 252, 253, 254, 255, 256, 257, 258, 259, 260, 261, 262, 263, 264, 265, 266, 267, 268, 269, 270, 271, 272, 273, 274, 275, 276, 277, 278, 279, 280, 281, 282, 283, 284, 285, 286, 287, 288, 289, 290, 291, 292, 293, 294, 295, 296, 297, 298, 299, 300, 301, 302, 303, 304, 305, 306, 307, 308, 309, 310, 311, 312, 313, 314, 315, 316, 317, 318, 319, 320, 321, 322, 323, 324, 325, 326, 327, 328, 329, 330, 331, 332, 333, 334, 335, 336, 337, 338, 339, 340, 341, 342, 343, 344, 345, 346, 347, 348, 349, 350, 351, 352, 353, 354, 355, 356, 357, 358, 359], maxpix=None, mask=None)

Bases: object

Maps images into a reciprocal space volume

Basic idea is for a fixed detector, scattering vectors in lab frame are fixed. When rotating sample, just rotate these into crystal frame and use as array indices into volume

add_image(om, data)

RSV = bounds of reciprocal space vol NR = dims of RSV k = scattering vector for image at om == 0 data = 2D image (dims == k.dims) SIG = signal MON = monitor

find_vol(border, omegarange)

find limiting volume The four image corners over 360 degrees

np is the number of pixels per hkl

returns (RSV, NR)
RSV min/max in reciprocal space (np*ubi).gv NR number of points == 1+RSV[i][1]-RSV[i][0]
make_k_vecs()

Generate the k vectors from the experiment parameters given in constructor

writevol(filename)

Save the volume in a hdf file

saintraw Module

Module for read saint raw files from bruker integration

class ImageD11.saintraw.saintraw(filename=None)

Bases: object

condition_filter(name, func)

Remove the peaks according to condition

doc = '\nIHKL 3I4 HKL indices of the reflection \n\n #IMNP 3I4 PRESENT ONLY FOR MODULATED-STRUCTURE DATA\n MNP indices of the reflection. Unused\n indices (in 4- or 5-dimensional cases)\n are written as zero\n\nFI F8.x Integrated intensity in photons, background-\n subtracted, normalized to 1 min/deg, and\n Lp corrected. Number of places to right\n of the decimal point is adjusted to\n balance precision vs. overflow; exponential\n format is used for very large values.\n\nSI F8.x Estimated standard deviation of the\n intensity. Number of places to right\n of the decimal point is adjusted to\n balance precision vs. overflow; exponential\n format is used for very large values.\n\nIBATNO I4 Batch number (for scaling) assigned to\n the reflection\n\nCOSINES 6F8.5 Direction cosines of the incident and\n diffracted beams, relative to the \n unit cell axes. Order is XI,XD,YI,YD,\n ZI,ZD where XI is the X-component of\n the incident beam, XD is the X-component\n of the diffracted beam, etc. Used\n for absorption correction.\n\n\nMSTATUS I3 1X,Z2 Status mask reserved for flagging\n abnormal conditions. There are\n currently none defined.\n\nXO F7.2 Observed X-pixel coordinate of the\n intensity-weighted reflection centroid,\n in reduced pixels (scale of 512x512)\n\nYO F7.2 Observed Y-pixel coordinate of the\n intensity-weighted reflection centroid,\n in reduced pixels (scale of 512x512)\n\nZO F8.2 Observed frame number of the\n intensity-weighted reflection centroid\n\nXP F7.2 Predicted X-pixel of the reflection\n in reduced pixels (scale of 512x512)\n\nYP F7.2 Predicted Y-pixel of the reflection\n in reduced pixels (scale of 512x512)\n\nZP F8.2 Predicted frame number of the reflection\n\nCORLPAF F6.3 Multiplicative correction for Lorentz\n effect, polarization, air absorption,\n and detector faceplate absorption,\n already applied to integrated intensity,\n FI\n\nCORREL F5.2 The correlation coefficient between the\n 3-D profile observed for this reflection\n and the corresponding model 3-D profile\n\nACCUMTIME F7.2 Accumulated hours of exposure\n\nSWING F7.2 Detector swing angle in degrees\n\nROT F7.2 Scan axis (omega or phi) setting in\n degrees at which the reflection was\n observed\n\nIAXIS I2 Scan axis number (2=omega, 3=phi)\n\nISTL I5 SIN(THETA)/LAMBDA times 10,000 \n\nPKSUM I9 Total raw peak counts in photons, with no\n correction for background, normalization, \n or Lp\n\nBGAVG I7 Normally, average BG per pixel in \n photons * 1000. In data sets where the \n background was reported to be very large \n ( > 1024 photons after 1 min/deg \n normalization), the program issues a warning \n message during integration and the scale of\n 1000 is omitted.\n\nROTREL F7.2 Rotation of the scan axis (omega or phi)\n in degrees relative to the start of the\n run\n\nICRYST I4 Crystal number (for scaling)\n\nPKFRAC F6.3 Fraction of the profile volume nearest\n this HKL. 1-PKFRAC is the fraction of\n the intensity which had to be estimated\n from the model profiles because of\n overlap with neighboring spots\n\nIKEY I11 The unique sort key for the group of\n equivalents to which this HKL belongs.\n IKEY = 1048576 * (H+511) + 1024 *\n (K+511) + (L+511), where H,K and L are\n the HKL indices transformed to the base\n asymmetric unit.\n\nIMULT I3 The point-group multiplicity of this HKL\n\nCORLORENTZ F6.3 Lorentz correction\n\nXGEO F8.2 Spatially corrected X relative to beam\n center, in pixels\n\nYGEO F8.2 Spatially corrected Y relative to beam \n center, in pixels\n\nCHI F8.3 Chi setting angle, degrees\n\nOTHER_ANGLE F8.3 The remaining setting angle, degrees. This\n will be phi if scans were in omega, or omega\n if scans were in phi\n\nICOMPONENT I4 Twin component number in SHELXTL HKLF 5\n convention. In a multi-component overlap,\n ICOMPONENT is negated for all but the last\n record of the overlap\n\n'
formats = {}
helps = {}
parsedocs()

Parse the saint documentation for the Bruker format

read(filename)

Read an ascii formatted saint reflection file

sort(name)

Sort according to a column in self.data

take(order)

Put the peaks in the order given in order (indices)

titles = []
tocolumnfile()

Return a columnfile

write(filename)

Write an ascii formatted saint reflection file

scale Module

Class for scaling images with respect to each other due to moco style instabilities.

Depends on there being a large background contribution to work with

Mathematically:

y = a * x + b dy/da = x dy/db = 1 Least squares problem with: matrix = ( dy/da*dy/da , dy/da*dy/db)

( dy/da*dy/db , dy/db*dy/db)
rhs = ( dy/da*y )
( dy/db*y )

Has the option to use only pixels in the image to scale to which are above a threshold (eg the central circle on the bruker)

class ImageD11.scale.scale(im1, threshold=None)
scale(im2)

Fill out RHS and solve returns the scale to apply to the image stored in the class

You probably want the scale to apply to the image you supply ...use scale image for that

scaleimage(im2)

Return a copy of the image scaled to match the class

ImageD11.scale.scaleseries(target, stem, first, last, thresh=None, writeim=None)

Scale a series of [bruker] images to the target TODO - make it work with fabio file series

scalem Module

Module(s) for scaling peaks

  • Read in a series of columnfile type objects
    collection of columnfiles
  • Read in some symmetry information:
    there are 32 crystallographic point groups 11 contain a centre of symmetry and are also laue groups
  • For each peak
    assign it a unique “key” in a list of reflections typically h>=k>=l of the symmetry equivalents Give it a scale factor (or two) Compute the corrected intensity and derivative w.r.t scale factors
  • For each uniq key compute the average group intensity
    with or without median style filtering
  • Compute the overall R(int) merging statistic and minimise that w.r.t scaling
    variables
  • Optimise scale factor variables
class ImageD11.scalem.data_table

Bases: object

Wraps a columnfile/hklf4/brukerraw file Gives h,k,l

  • variables [scales per frame, scan, grain etc]
  • derivatives [of all variables]

To constrain grain scale factors across different scans imply we need to have grains in same datafile.

calc()
get_hkls()

Return the hkl indices of the peaks

set_keys(keys)

Store a list of uniq keys which go with each peak These are used to group peaks into symmetry equivalents

ImageD11.scalem.generate_key(hkl, symmops)

To be redone in C/fortran etc Expect identity to be in symmops All symmetry matrices, not including inversion

class ImageD11.scalem.scalem

Bases: object

This class holds a series of tables (assumes everything fits in memory) ... and does some fitting

add_data(filename, format='ImageD11')

Read a columnfile object. Allowed to have format == ImageD11|bruker|hklf4

calc()
Compute for each data object:
Columns: I_corr, sigma_corr Foreach scaling variable: dI_corr/dvar, dsigma_corr/dvar
set_sym(sym)
Set the symmetry to sym
  • list of matrices (3x3)
  • list of operators (4x4)
  • string (eg m3m etc)

simplex Module

Simplex - a regression method for arbitrary nonlinear function minimization

Simplex minimizes an arbitrary nonlinear function of N variables by the Nelder-Mead Simplex method as described in:

Nelder, J.A. and Mead, R., “A Simplex Method for Function Minimization”,
Computer Journal, Vol. 7, 1965, pp. 308-313.

It makes no assumptions about the smoothness of the function being minimized. It converges to a local minimum which may or may not be the global minimum depending on the initial guess used as a starting point.

# Copyright (c) 2001 Vivake Gupta (v@omniscia.org). All rights reserved. # # This program is free software; you can redistribute it and/or # modify it under the terms of the GNU General Public License as # published by the Free Software Foundation; either version 2 of the # License, or (at your option) any later version. # # This program is distributed in the hope that it will be useful, # but WITHOUT ANY WARRANTY; without even the implied warranty of # MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the # GNU General Public License for more details. # # You should have received a copy of the GNU General Public License # along with this program; if not, write to the Free Software # Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 # USA # # This software is maintained by Vivake (v@omniscia.org) and is available at: # http://www.omniscia.org/~vivake/python/Simplex.py

# Modified (debugged?) 7/16/2004 Michele Vallisneri (vallis@vallis.org)

# Copied into ImageD11 by Jon Wright - This works really really well!!!

class ImageD11.simplex.Simplex(testfunc, guess, increments, kR=-1, kE=2, kC=0.5)
accept_contracted_point()
accept_expanded_point()
accept_reflected_point()
calculate_errors_at_vertices()
contract_simplex()
expand_simplex()
minimize(epsilon=0.0001, maxiters=250, monitor=1)

Walks to the simplex down to a local minima. INPUTS —— epsilon convergence requirement maxiters maximum number of iterations monitor if non-zero, progress info is output to stdout

an array containing the final values lowest value of the error function number of iterations taken to get here

multiple_contract_simplex()
reflect_simplex()
ImageD11.simplex.main()
ImageD11.simplex.myfunc(args)

sym_u Module

ImageD11.sym_u.cubic()
ImageD11.sym_u.find_uniq_hkls(hkls, grp, func=<function hklmax at 0x000000000CAC09E8>)
ImageD11.sym_u.find_uniq_u(u, grp, debug=0, func=<function trace at 0x0000000006573AC8>)
ImageD11.sym_u.fmt(c)
ImageD11.sym_u.generate_group(*args)
ImageD11.sym_u.getgroup(s)

convert a user supplied string to a group ... a little vague still

class ImageD11.sym_u.group(tol=1e-05)

An abstract mathematical finite(?) point rotation groups

additem(item)

add a new member

comp(x, y)

Compare two things for equality

isMember(item)

Decide if item is already in the group

makegroup()

ensure all items = op(x,y) are in group

op(x, y)

Normally multiplication ? Means of generating new thing from two others

ImageD11.sym_u.hexagonal()

P6/mmm 191

ImageD11.sym_u.hklmax(h, hmax=1000)
ImageD11.sym_u.m_from_string(s)

Creates a symmetry operator from a string

ImageD11.sym_u.m_to_string(m)

Creates a symmetry operator string from a matrix

ImageD11.sym_u.monoclinic_a()
ImageD11.sym_u.monoclinic_b()
ImageD11.sym_u.monoclinic_c()
ImageD11.sym_u.orthorhombic()

P222 16

ImageD11.sym_u.rhombohedralP()

R3 primitive

ImageD11.sym_u.test()
ImageD11.sym_u.tetragonal()

P4 75

class ImageD11.sym_u.trans_group(tol=1e-05)

Bases: ImageD11.sym_u.group

Translation group (eg crystal lattice)

FIXME - this is mostly done in lattice_reduction.py instead now

additem(x)

Do lattice reduction before adding as infinite group

isMember(x)
mod(x, y)

Remove y from x to give smallest possible result Find component of x || to y and remove it

op(x, y)

Means of generating new thing from two others In this case add them and mod by group members

reduce(v)

Perform lattice reduction

ImageD11.sym_u.triclinic()
ImageD11.sym_u.trigonal()

P3 143

symops Module

ImageD11.symops.checkop(h, k, l, op, axis)
ImageD11.symops.glide_plane(h, k, l, gtype, axis)
ImageD11.symops.glidetest(gtype, axis)
ImageD11.symops.lattice_centre(h, k, l, ctype)
ImageD11.symops.mirror_plane(h, k, l, axis)
ImageD11.symops.rotation_axis(h, k, l, rtype, axis)
ImageD11.symops.screw_axis(h, k, l, stype, axis)
ImageD11.symops.test_absence(h, k, l, sg)

test_transform Module

class ImageD11.test_transform.testtransform(methodName='runTest')

Bases: unittest.case.TestCase

setUp()
test_compute_tth_eta1()

transform Module

Functions for transforming peaks

ImageD11.transform.compute_g_from_k(k, omega, wedge=0, chi=0)

Compute g-vectors with cached k-vectors

ImageD11.transform.compute_g_vectors(tth, eta, omega, wvln, wedge=0.0, chi=0.0)
Generates spot positions in reciprocal space from
twotheta, wavelength, omega and eta

Assumes single axis vertical ... unless a wedge angle is specified

ImageD11.transform.compute_grain_origins(omega, wedge=0.0, chi=0.0, t_x=0.0, t_y=0.0, t_z=0.0)

# print “Using translations t_x %f t_y %f t_z %f”%(t_x,t_y,t_z) # Compute positions of grains # expecting tx, ty, tz for each diffraction spot # # g = R . W . k # g - is g-vector w.r.t crystal # k is scattering vector in lab # so we want displacement in lab from displacement in sample # shift = W-1 R-1 crystal_translation # # R = ( cos(omega) , sin(omega), 0 ) # (-sin(omega) , cos(omega), 0 ) # ( 0 , 0 , 1 ) # # W = ( cos(wedge) , 0 , sin(wedge) ) # ( 0 , 1 , 0 ) # (-sin(wedge) , 0 , cos(wedge) ) # # C = ( 1 , 0 , 0 ) ??? Use eta0 instead # ( 0 , cos(chi) , sin(chi) ) ??? Use eta0 instead # ( 0 , -sin(chi) , cos(chi) ) ??? Use eta0 instead

ImageD11.transform.compute_k_vectors(tth, eta, wvln)

generate k vectors - scattering vectors in laboratory frame

ImageD11.transform.compute_lorentz_factors(tth, eta, omega, wavelength, wedge=0.0, chi=0.0)

From Kabsch 1988 J. Appl. Cryst. 21 619

Multiply the intensities by: Lorentz = | S.(u x So)| / |S|.|So| S = scattered vector So = incident vector u = unit vector along rotation axis

ImageD11.transform.compute_polarisation_factors(args)

From Kabsch 1988 J. Appl. Cryst. 21 619

DIVIDE the intensities by: <sin2 psi> = (1 - 2p) [ 1 - (n.S/|S|^2) ] + p { 1 + [S.S_0/(|S||S_0|)^2]^2}

p = degree of polarisation (sync = 1, tube = 0.5 , mono + tube in between)
or “probability of finding electric field vector in plane having
normal, n”

S = scattered vector S_0 = incident vector n = normal to polarisation plane, typically perpendicular to S_0

In ImageD11 we normally expect to find: x axis along the beam z axis being up, and parallel to the normal n mentioned above

ImageD11.transform.compute_tth_eta(peaks, y_center=0.0, y_size=0.0, tilt_y=0.0, z_center=0.0, z_size=0.0, tilt_z=0.0, tilt_x=0.0, distance=0.0, o11=1.0, o12=0.0, o21=0.0, o22=-1.0, t_x=0.0, t_y=0.0, t_z=0.0, omega=None, wedge=0.0, chi=0.0, **kwds)

0/10 for style

ImageD11.transform.compute_tth_eta_from_xyz(peaks_xyz, omega, t_x=0.0, t_y=0.0, t_z=0.0, wedge=0.0, chi=0.0, **kwds)

Peaks is a 3 d array of x,y,z peak co-ordinates crystal_translation is the position of the grain giving rise to a diffraction spot

in x,y,z ImageD11 co-ordinates x,y with respect to axis of rotation and or beam centre ?? z with respect to beam height, z centre

omega data needed if crystal translations used

ImageD11.transform.compute_tth_histo(tth, no_bins=100, **kwds)

Compute a histogram of tth values

Returns a normalised histogram (should make this a probability and

For each datapoint, the number of other points in the same bin
ImageD11.transform.compute_xyz_from_tth_eta(tth, eta, omega, t_x=0.0, t_y=0.0, t_z=0.0, wedge=0.0, chi=0.0, **kwds)

Given the tth, eta and omega, compute the xyz on the detector

crystal_translation is the position of the grain giving rise to a diffraction spot

in x,y,z ImageD11 co-ordinates x,y with respect to axis of rotation and or beam centre ?? z with respect to beam height, z centre

omega data needed if crystal translations used

ImageD11.transform.compute_xyz_lab(peaks, y_center=0.0, y_size=0.0, tilt_y=0.0, z_center=0.0, z_size=0.0, tilt_z=0.0, tilt_x=0.0, distance=0.0, o11=1.0, o12=0.0, o21=0.0, o22=-1.0, **kwds)

Peaks is a 2 d array of x,y yc is the centre in y ys is the y pixel size ty is the tilt around y zc is the centre in z zs is the z pixel size tz is the tilt around z dist is the sample - detector distance detector_orientation is a matrix to apply to peaks arg to get ImageD11 convention

(( 0, 1),( 1, 0)) for ( y, x) ((-1, 0),( 0, 1)) for (-x, y) (( 0,-1),(-1, 0)) for (-y,-x)

etc...

ImageD11.transform.cross_product_2x2(a, b)

returns axb for two len(3) vectors a,b

ImageD11.transform.detector_rotation_matrix(tilt_x, tilt_y, tilt_z)

Return the tilt matrix to apply to peaks tilts in radians typically applied to peaks rotating around beam center

ImageD11.transform.uncompute_g_vectors(g, wavelength, wedge=0.0, chi=0.0)

Given g-vectors compute tth,eta,omega assert uncompute_g_vectors(compute_g_vector(tth,eta,omega))==tth,eta,omega

ImageD11.transform.uncompute_one_g_vector(gv, wavelength, wedge=0.0)

Given g-vectors compute tth,eta,omega assert uncompute_g_vectors(compute_g_vector(tth,eta,omega))==tth,eta,omega

transformer Module

class ImageD11.transformer.transformer(parfile=None, fltfile=None)

Handles the algorithmic, fitting and state information for fitting parameters to give experimental calibrations

addcellpeaks(limit=None)

Adds unit cell predicted peaks for fitting against Optional arg limit gives highest angle to use

Depends on parameters:
‘cell__a’,’cell__b’,’cell__c’, ‘wavelength’ # in angstrom ‘cell_alpha’,’cell_beta’,’cell_gamma’ # in degrees ‘cell_lattice_[P,A,B,C,I,F,R]’
addcolumn(col, name)

Return the data

applyargs(args)

for use with simplex/gof function, alter parameters

compute_histo(colname)

Compute the histogram over twotheta for peaks previous read in Filtering is moved to a separate function

compute_tth_eta()

Compute the twotheta and eta for peaks previous read in

compute_tth_histo()

Give hardwire access to tth

computegv()

Using tth, eta and omega angles, compute x,y,z of spot in reciprocal space

filter_min(col, minval)

and filter peaks out falling in bins with less than min_prob

fit(tthmin=0, tthmax=180)

Apply simplex to improve fit of obs/calc tth

get_variable_list()
getaxis()

Compute the rotation axis This handles omegasign in a more elegant way

getcols()
getcolumn(name)

Return the data

getvars()

decide what is refinable

gof(args)

Compute how good is the fit of obs/calc peak positions in tth

loadfileparameters(filename)

Read in beam center etc from file

loadfiltered(filename)

Read in 3D peaks from peaksearch

savegv(filename)

Save g-vectors into a file Use crappy .ass format from previous for now (testing)

saveparameters(filename)

Save beam center etc to file

setvars(varlist)

set the things to refine

setxyomcols(xname, yname, omeganame)
tth_entropy()

Compute the entropy of the two theta values ... this may depend on the number of tth bins (perhaps?) ... to be used for cell parameter free line straightening S = sum -p log p

updateparameters()
write_colfile(filename)

Save out the column file (all info stored we hope)

write_graindex_gv(filename)
write_pyFAI(filename, tthmin=0, tthmax=180)

Write file for Jerome Kieffer’s pyFAI fitting routine to use and run the refinment...

twodplot Module

From the matplotlib examples - modified for mouse

class ImageD11.twodplot.data(x, y, d={})
class ImageD11.twodplot.twodplot(parent=None, data=None, quiet='No')

Bases: Tkinter.Frame

adddata(data)

Takes a tuple of name, data object

autoscale()
autoscaley(e)
bindkeys()
clear()
hideall()
keypress(*arg)
logx()
logy()
on_2(event)
on_3(event)
on_down(event)
on_move(event)
on_up(event)
printplot()
removedata(name)
replot()

ubitool Module

Tools for processing UBI matrices together with data

Uses indexing objects to do refinement

Mostly unimplemented dummies - see also refinegrains

class ImageD11.ubitool.ubitool(ubilist=None, ubifile=None, obsdata=None)
get_orientation(ubi=None)

Give orientation matrix independant of unit cell

get_unit_cell(ubi=None)

Convert UBI representation to give unit cell

read_ubi_file(filename)

Get ubi matrices from a file

refine_translations(ubi=None)

Compute an offset in x/y/z for the origin of the grain with respect to the centre of rotation

validate_peak_assignements()

Make sure each hkl is only assigned to one peak. Offer to merge if there are duplicates

validate_ubi(ubi=None)

Find the number of peaks versus hkl_tol when refining This should plateau when the tolerance is correct and the data are good enough.

validate_ubi_collection()

Check for duplicate orientations

unitcell Module

ImageD11.unitcell.A(h, k, l)
ImageD11.unitcell.B(h, k, l)
ImageD11.unitcell.C(h, k, l)
ImageD11.unitcell.F(h, k, l)
ImageD11.unitcell.I(h, k, l)
ImageD11.unitcell.P(h, k, l)
ImageD11.unitcell.R(h, k, l)
ImageD11.unitcell.cellfromstring(s)
ImageD11.unitcell.cross(a, b)

a x b has length |a||b|sin(theta)

ImageD11.unitcell.degrees(x)
ImageD11.unitcell.norm2(a)

Compute the unit 2 norm

ImageD11.unitcell.radians(x)
ImageD11.unitcell.unit(a)

Normalise vector a to unit length

class ImageD11.unitcell.unitcell(lattice_parameters, symmetry='P', verbose=0)
anglehkls(h1, h2)

Compute the angle between reciprocal lattice vectors h1, h2

ds(h)

computes 1/d for this hkl = hgh

getanglehkls(ring1, ring2)

Cache the last pair called for

gethkls(dsmax)

Generate hkl list Argument dsmax is the d* limit (eg 1/d) Default of zero gives only the (000) reflection

assumes [h|k|l] < 200

gethkls_xfab(dsmax, spg)

Generate hkl list Argument dsmax is the d* limit (eg 1/d) Argument spg is the space group name, e.g. ‘R3-c’

makerings(limit, tol=0.001)

Makes a list of computed powder rings The tolerance is the difference in d* to decide if two peaks overlap

orient(ring1, g1, ring2, g2, verbose=0, all=True)

Compute an orientation matrix using cell parameters and the indexing of two reflections

Orientation matrix: A matrix such that h = A^1 g define 3 vectors t1,t2,t3 t1 is parallel to first peak (unit vector along g1) t2 is in the plane of both (unit vector along g1x(g1xg2)) t3 is perpendicular to both (unit vector along g1xg2)

tostring()

Write out a line containing unit cell information

ImageD11.unitcell.unitcell_from_parameters(pars)

write_graindex_gv Module

Conversion of Soren’s matlab script for making the file final.log which is suitable for input into graindex

ImageD11.write_graindex_gv.get_ds_string(g, ds_list)

Attempt to emulate the graindex (hkl) syntax. Obviously you would not want a (10,0,0) peak!

ImageD11.write_graindex_gv.make_ds_list(cell, limit=2.0)

Generates a list of d-spacings

ImageD11.write_graindex_gv.write_graindex_gv(outfilename, gv, tth, eta, omega, intensity, unitcell)

call with array of gvectors, tth, eta and omega vals

Using an indexing object to get the ring assignments