hexrd.transforms.xf module

hexrd.transforms.xf.anglesToGVec(angs, bHat_l, eHat_l, rMat_s=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), rMat_c=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]))[source]

from ‘eta’ frame out to lab (with handy kwargs to go to crystal or sample)

hexrd.transforms.xf.angles_in_range(angles, starts, stops, degrees=True)[source]

Determine whether angles lie in or out of specified ranges

angles - a list/array of angles starts - a list of range starts stops - a list of range stops

OPTIONAL ARGS: degrees - [True] angles & ranges in degrees (or radians)

hexrd.transforms.xf.angularDifference(angList0, angList1, units='radians')[source]

Calculate the proper (acute) angular difference in the context of a branch cut.

Parameters

angList0array_like

An (n, ) array on input angles in radians.

angList1array_like

An (n, ) array on input angles in radians.

unitsstr, optional

Either ‘degrees’ or ‘radians’. The default is ‘radians’.

Returns

numpy.ndarray

The (n, ) array of acute angular differences bwetween the inputs on the circle [-pi, pi].

hexrd.transforms.xf.arccosSafe(temp)[source]

Protect against numbers slightly larger than 1 in magnitude due to round-off

hexrd.transforms.xf.columnNorm(a)[source]

normalize array of column vectors (hstacked, axis = 0)

hexrd.transforms.xf.detectorXYToGvec(xy_det, rMat_d, rMat_s, tVec_d, tVec_s, tVec_c, distortion=None, beamVec=array([[-0.], [-0.], [-1.]]), etaVec=array([[1.], [0.], [0.]]), output_ref=False)[source]

Takes a list cartesian (x, y) pairs in the detector coordinates and calculates the associated reciprocal lattice (G) vectors and (bragg angle, azimuth) pairs with respect to the specified beam and azimth (eta) reference directions.

Parameters

xy_detnumpy.ndarray

DESCRIPTION.

rMat_dnumpy.ndarray

(3, 3) array, the COB taking DETECTOR FRAME components to LAB FRAME.

rMat_snumpy.ndarray

(3, 3) array, the COB taking SAMPLE FRAME components to LAB FRAME.

tVec_dnumpy.ndarray

(3, 1) array, the translation vector connecting LAB to DETECTOR.

tVec_snumpy.ndarray

(3, 1) array, the translation vector connecting LAB to SAMPLE.

tVec_cnumpy.ndarray

(3, 1) array, the translation vector connecting SAMPLE to CRYSTAL.

distortionclass, optional

the distortion operator, if applicable. The default is None.

beamVecnumpy.ndarray, optional

(3, 1) array, the vector pointing to the X-ray source in the LAB FRAME. The default is bVec_ref.

etaVecnumpy.ndarray, optional

(3, 1) array, the vector defining eta=0 in the LAB FRAME. The default is eta_ref.

output_refbool, optional

if true, will output calculated angles in the LAB FRAME. The default is False.

Returns

tth_eta_reftuple, optional

if output_ref is True, (n, 2) ndarray containing the (tTh, eta) pairs in the LAB REFERENCE FRAME associated with each (x, y)

tth_etatuple

(n, 2) ndarray containing the (tTh, eta) pairs associated with each (x, y).

gVec_lnumpy.ndarray

(3, n) ndarray containing the associated G vector directions in the LAB FRAME associated with gVecs

hexrd.transforms.xf.gvecToDetectorXY(gVec_c, rMat_d, rMat_s, rMat_c, tVec_d, tVec_s, tVec_c, beamVec=array([[-0.], [-0.], [-1.]]))[source]

Takes a list of unit reciprocal lattice vectors in crystal frame to the specified detector-relative frame.

The following conditions must be satisfied:

  1. the reciprocal lattice vector must be able to satisfy a bragg condition

  2. the associated diffracted beam must intersect the detector plane

Parameters

gVec_cnumpy.ndarray

(3, n) array of n reciprocal lattice vectors in the CRYSTAL FRAME.

rMat_dnumpy.ndarray

(3, 3) array, the COB taking DETECTOR FRAME components to LAB FRAME.

rMat_snumpy.ndarray

(3, 3) array, the COB taking SAMPLE FRAME components to LAB FRAME.

rMat_cnumpy.ndarray

(3, 3) array, the COB taking CRYSTAL FRAME components to SAMPLE FRAME.

tVec_dnumpy.ndarray

(3, 1) array, the translation vector connecting LAB to DETECTOR.

tVec_snumpy.ndarray

(3, 1) array, the translation vector connecting LAB to SAMPLE.

tVec_cnumpy.ndarray

(3, 1) array, the translation vector connecting SAMPLE to CRYSTAL.

beamVecnumpy.ndarray, optional

(3, 1) array, the vector pointing to the X-ray source in the LAB FRAME. The default is bVec_ref.

Returns

numpy.ndarray

(3, m)darray containing the intersections of m <= n diffracted beams associated with gVecs.

hexrd.transforms.xf.makeBinaryRotMat(axis)[source]
hexrd.transforms.xf.makeDetectorRotMat(tiltAngles)[source]

Form the (3, 3) tilt rotations from the tilt angle list:

tiltAngles = [gamma_Xl, gamma_Yl, gamma_Zl] in radians

hexrd.transforms.xf.makeEtaFrameRotMat(bHat_l, eHat_l)[source]

make eta basis COB matrix with beam antiparallel with Z

takes components from ETA frame to LAB

**NO EXCEPTION HANDLING FOR COLINEAR ARGS IN NUMBA VERSION!

…put checks for non-zero magnitudes and non-colinearity in wrapper?

hexrd.transforms.xf.makeGVector(hkl, bMat)[source]

take a CRYSTAL RELATIVE B matrix onto a list of hkls to output unit reciprocal lattice vectors (a.k.a. lattice plane normals)

Required Arguments: hkls – (3, n) ndarray of n hstacked reciprocal lattice vector component

triplets

bMat – (3, 3) ndarray representing the matirix taking reciprocal lattice

vectors to the crystal reference frame

Output: gVecs – (3, n) ndarray of n unit reciprocal lattice vectors

(a.k.a. lattice plane normals)

To Do: * might benefit from some assert statements to catch improperly shaped

input.

hexrd.transforms.xf.makeOscillRotMat(oscillAngles)[source]

oscillAngles = [chi, ome]

hexrd.transforms.xf.makeRotMatOfExpMap(expMap)[source]

Calculate the rotation matrix of an exponential map.

Parameters

expMaparray_like

A 3-element sequence representing the exponential map n*phi.

Returns

rMatnp.ndarray

The (3, 3) array representing the rotation matrix of the input exponential map parameters.

hexrd.transforms.xf.mapAngle(ang, *args, **kwargs)[source]

Utility routine to map an angle into a specified period

actual function is mapAngle(ang[, range], units=angularUnits). range is optional and defaults to the appropriate angle for the unit centered on 0.

hexrd.transforms.xf.oscillAnglesOfHKLs(hkls, chi, rMat_c, bMat, wavelength, vInv=array([[1.], [1.], [1.], [0.], [0.], [0.]]), beamVec=array([[-0.], [-0.], [-1.]]), etaVec=array([[1.], [0.], [0.]]))[source]

Parameters

hklsnumpy.ndarray

(3, n_) array of reciprocal lattice vector components.

chiscalar

the inclination angle of the SAMPLE FRAME about the LAB X.

rMat_cnumpy.ndarray

(3, 3) array, the COB taking CRYSTAL FRAME components to SAMPLE FRAME.

bMatnumpy.ndarray

(3, 3) array, the COB matrix taking RECIPROCAL FRAME components to the CRYSTAL FRAME.

wavelengthscalar

the X-ray wavelength in Ångstroms.

vInvnumpy.ndarray, optional

(6, 1) array representing the components of the inverse stretch tensor in the SAMPLE FRAME. The default is vInv_ref.

beamVecnumpy.ndarray, optional

(3, 1) array, the vector pointing to the X-ray source in the LAB FRAME. The default is bVec_ref.

etaVecnumpy.ndarray, optional

(3, 1) array, the vector defining eta=0 in the LAB FRAME. The default is eta_ref.

Returns

retvaltuple

Two (3, n) numpy.ndarrays representing the pair of feasible (tTh, eta, ome) solutions for the n input hkls.

Notes

The reciprocal lattice vector, G, will satisfy the the Bragg condition when:

b.T * G / ||G|| = -sin(theta)

where b is the incident beam direction (k_i) and theta is the Bragg angle consistent with G and the specified wavelength. The components of G in the lab frame in this case are obtained using the crystal orientation, Rc, and the single-parameter oscillation matrix, Rs(ome):

Rs(ome) * Rc * G / ||G||

The equation above can be rearranged to yeild an expression of the form:

a*sin(ome) + b*cos(ome) = c

which is solved using the relation:

a*sin(x) + b*cos(x) = sqrt(a**2 + b**2) * sin(x + alpha)

—> sin(x + alpha) = c / sqrt(a**2 + b**2)

where:

alpha = arctan2(b, a)

The solutions are:

/ | arcsin(c / sqrt(a**2 + b**2)) - alpha

x = <
pi - arcsin(c / sqrt(a**2 + b**2)) - alpha

There is a double root in the case the reflection is tangent to the Debye-Scherrer cone (c**2 = a**2 + b**2), and no solution if the Laue condition cannot be satisfied (filled with NaNs in the results array here)

hexrd.transforms.xf.polarRebin(thisFrame, npdiv=2, mmPerPixel=(0.2, 0.2), convertToTTh=False, rMat_d=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), tVec_d=array([0., 0., -1000.]), beamVec=array([[-0.], [-0.], [-1.]]), etaVec=array([[1.], [0.], [0.]]), rhoRange=array([20, 200]), numRho=1000, etaRange=array([-0.08726646, 6.19591884]), numEta=36, verbose=True, log=None)[source]

Performs polar rebinning of an input image.

Parameters

thisFramearray_like

An (n, m) array representing an image.

npdivint, optional

The number of pixel subdivision. The default is 2.

mmPerPixelarray_like, optional

A 2-element sequence with the pixel pitch in mm for the row and column dimensions, respectively. The default is (0.2, 0.2).

convertToTThbool, optional

If true, output abcissa is in angular (two-theta) coordinates. The default is False.

rMat_dnumpy.ndarray, optional

(3, 3) array, the COB taking DETECTOR FRAME components to LAB FRAME. The default is I3.

tVec_dnumpy.ndarray, optional

(3, 1) array, the translation vector connecting LAB to DETECTOR. The default is np.r_[0., 0., -1000.].

beamVecTYPE, optional

DESCRIPTION. The default is bVec_ref.

etaVecTYPE, optional

DESCRIPTION. The default is eta_ref.

rhoRangearray_like, optional

the min/max limits of the rebinning region in pixels. The default is np.r_[20, 200].

numRhoint, optional

the number of bins to use in the radial dimension. The default is 1000.

etaRangearray_like, optional

A 2-element sequence describing the min/max azimuthal angles in radians. The default is (-0.08726646259971647, 6.19591884457987).

numEtaint, optional

The number of azimuthal bins. The default is 36.

verbosebool, optional

DESCRIPTION. The default is True.

logTYPE, optional

DESCRIPTION. The default is None.

Raises

RuntimeError

DESCRIPTION.

Returns

polImgTYPE

DESCRIPTION.

hexrd.transforms.xf.quat_distance(q1, q2, qsym)[source]
hexrd.transforms.xf.quat_product_matrix(q, mult='right')[source]

Form 4 x 4 array to perform the quaternion product

USAGE

qmat = quatProductMatrix(q, mult=’right’)

INPUTS
  1. quats is (4,), an iterable representing a unit quaternion horizontally concatenated

  2. mult is a keyword arg, either ‘left’ or ‘right’, denoting the sense of the multiplication:

    / quatProductMatrix(h, mult=’right’) * q

    q * h –> <

    quatProductMatrix(q, mult=’left’) * h

OUTPUTS
  1. qmat is (4, 4), the left or right quaternion product operator

NOTES
*) This function is intended to replace a cross-product based

routine for products of quaternions with large arrays of quaternions (e.g. applying symmetries to a large set of orientations).

hexrd.transforms.xf.reg_grid_indices(edges, points_1d)[source]

get indices in a 1-d regular grid.

edges are just that:

point: x (2.5)

edges: |1 |2 |3 |4 |5
indices: | 0 | 1 | 2 | 3 |

above the deltas are + and the index for the point is 1

point: x (2.5)

edges: |5 |4 |3 |2 |1
indices: | 0 | 1 | 2 | 3 |

here the deltas are - and the index for the point is 2

  • can handle grids with +/- deltas

  • be careful when using with a cyclical angular array! edges and points must be mapped to the same branch cut, and abs(edges[0] - edges[-1]) = 2*pi

hexrd.transforms.xf.rotate_vecs_about_axis(angle, axis, vecs)[source]

Rotate vectors about an axis

INPUTS angle - array of angles (len == 1 or n) axis - array of unit vectors (shape == (3, 1) or (3, n)) vec - array of vectors to be rotated (shape = (3, n))

Quaternion formula: if we split v into parallel and perpedicular components w.r.t. the axis of quaternion q,

v = a + n

then the action of rotating the vector dot(R(q), v) becomes

v_rot = (q0**2 - |q|**2)(a + n) + 2*dot(q, a)*q + 2*q0*cross(q, n)

hexrd.transforms.xf.rowNorm(a)[source]

normalize array of row vectors (vstacked, axis = 1)

hexrd.transforms.xf.unitVector(a)[source]

normalize array of column vectors (hstacked, axis = 0)

hexrd.transforms.xf.validateAngleRanges(angList, startAngs, stopAngs, ccw=True)[source]

A better way to go. find out if an angle is in the range CCW or CW from start to stop

There is, of course, an ambigutiy if the start and stop angle are the same; we treat them as implying 2*pi having been mapped