Welcome to the s3dxrd package

This is a scientific code originally developed to adress scanning-3d-xray-diffraction (s3dxrd) strain measurements of polycrystalline materials.

Intragranular strain is computed based on a series of line integral measurements. The s3dxrd package supports regression either by a simple weighted least squares approach or alternatively by a Gaussian Proccess. The later statistical model uses spatial correlation assumptions and an equlibrium prior to find good fits to data.

If you want to use this code, it is strongly recomended that you have a look at the underlying publication: describing the weighted least squares approach (named “ASR” in the paper)

Reconstructing intragranular strain fields in polycrystalline materials from scanning 3DXRD data, Henningsson, N. A., Hall, S. A., Wright, J. P. & Hektor, J. (2020). J. Appl. Cryst. 53, 314-325.

The paper describing the Gaussian Process regression procedure is also available here:

Intragranular Strain Estimation in Far-Field Scanning X-ray Diffraction using a Gaussian Processes, Henningsson, A. & Hendriks, J. (2021). J. Appl. Cryst. 54.

This paper may also help the user to understand some of the mathematical notation hinted at in the code.

Installation

Installation via pip is available as

pip3 install s3dxrd

For manuall installation, first get the code to your local machine by:

git clone https://github.com/AxelHenningsson/scanning-xray-diffraction.git

Next go to the repository and try to install

cd scanning-xray-diffraction
python setup.py build install

You will now recieve messages about dependecies that need be installed first. Go through these untill the build succeeds.

Documentation

Documentation is hosted seperately at github pages:

measurement

measurement.Id11.peaks_to_vectors(flt_paths, zpos, param_path, ubi_paths, omegastep, ymin, ymax, ystep, hkltol=0.05, nmedian=5, rcut=0.2, save_path=None)

Convert an x-ray diffraction dataset saved in the ImageD11 format to vector format

Based on an inital guess of ubi matrices all data contained in a sereis of Id11 peak files is analysed. Grian shapes are reconstructed and grain orientations refined. The data is converted into a vector format where each measurement is associated to a parametric line trhough the grain as well as an average approximate strain along the diffracting planes.

The proccessed output from this function can be used to further perform regression for strain tensor maps using any of the in package supported regression techinques.

Parameters:
  • flt_paths (list of string) – Absolute file paths to Id11 peaks files. These must contain a column named `dty` that stores the sample y-translations in units of microns. These translations should be centered at the rotation axis (dty=0 means the xray go through the rotation axis.)

  • zpos (list of float) – z-translations of sample corresponding to peaks in flt_paths.

  • param_path (string) – Absolute file path to Id11 parameters file.

  • ubi_paths (list of string) – Absolute file paths to Id11 ubi matrices file.

  • omegastep (float) – Rotation intervall in degrees over wich detector images are where read out.

  • ymin (float) – Minimum y-translation of smaple in units of microns.

  • ymax (float) – Maximum y-translation of smaple in units of microns.

  • ystep (float) – Sample y-translation stepsize.

  • hkltol (float) – Miller index tolerance for assigning peaks to a grain. Defaults to 0.05.

  • nmedian (float) – Number of median deviations of a peak to discard as outlier. Defaults to 5.

  • rcut (float) – Relative threshold for defining grain shape from tomographic reconstrction. Defaults to 0.2.

  • save_path (string) – Path at which to save all computed for arrays. Defaults to None.

Returns:

dictionary with keys and values:

Y (numpy array): Average scalar strains along scanned lines. shape=(N,)

sig_m (numpy array): Per measurement standard deviations. shape=(N,)

entry (numpy array): Entry points for line integral meaurements. shape=(k,N)

exit (numpy array): Exit points for line integral meaurements. shape=(k,N)

nhat (numpy array): X-ray beam direction normal in sample coordinate system. shape=(3,N)

kappa (numpy array): Direction of strain (Y) in sample coordinate system. shape=(3,N)

L (numpy array): Lengths of scanned lines. shape=(N,)

nsegs (numpy array): Number of segments for each line integral. shape=(N,)

polygons (dict of dict of Shapely Polygons): Polygonal boundary representation of grians.

polygon_orientations (dict of dict of numpy array): Crystal orientation of each grain slice in polygons.

orientations (numpy array): Per measurement crystal orientation (uniform along scan-line).

measurement_grain_map (numpy array): dictionary keys as integers corresponding to polygons, such that polygons[str(measurement_grain_map[i])] gives the polygon slices dictionary of the grain that gave rise to measurement number i. shape=(N,)

labeled_volume (numpy array): 3d numpy array of sample with values corresponding to grain index.

labeled_grains (dict of dict of ImageD11 Grain): Id11 grain for each grain slice in polygons.

Return type:

(dict)

regression

regression.bayesian.BeltramiApproxAnisoRayMeas(A, B, C, D, E, F, H, entry, exit, nsegs, L, kappa, subset_idx=None)
regression.bayesian.BeltramiApproxAnisoStrainField(A, B, C, D, E, F, x, y, z, H)

Calculates regressor matrix given hyperparameters for beltrami components and spatial locations.

This function works for aniostropic materials as the compliance H can be supplied.

Parameters:
  • A (dict) – Basis fuctions for Beltrami component A

  • B (dict) – Basis fuctions for Beltrami component B

  • C (dict) – Basis fuctions for Beltrami component C

  • D (dict) – Basis fuctions for Beltrami component D

  • E (dict) – Basis fuctions for Beltrami component E

  • F (dict) – Basis fuctions for Beltrami component F

  • x (torch.tensor) – x coordinate of points to calculate basis functions at

  • y (torch.tensor) – y coordinate of points to calculate basis functions at

  • z (torch.tensor) – z coordinate of points to calculate basis functions at

  • H (torch.tensor) – the compliance matrix, either a `shape=(6,6)` average for the material or `shape=(N,6,6)` individual compliances for each point

Returns:

The regressor matrix

regression.bayesian.BeltramiApproxAnisoStressField(A, B, C, D, E, F, x, y, z)

Calculates regressor matrix given hyperparameters for beltrami components and spatial locations.

This function works for aniostropic materials as the compliance H can be supplied.

Parameters:
  • A (dict) – Basis fuctions for Beltrami component A

  • B (dict) – Basis fuctions for Beltrami component B

  • C (dict) – Basis fuctions for Beltrami component C

  • D (dict) – Basis fuctions for Beltrami component D

  • E (dict) – Basis fuctions for Beltrami component E

  • F (dict) – Basis fuctions for Beltrami component F

  • x (torch.tensor) – x coordinate of points to calculate basis functions at

  • y (torch.tensor) – y coordinate of points to calculate basis functions at

  • z (torch.tensor) – z coordinate of points to calculate basis functions at

Returns:

The regressor matrix

NOTE: Approximating stress is simply done by approximating strain with H=I unity.

regression.bayesian.basisFuncDInts(lambdas, Lx, Ly, Lz, entry, exit, req, nsegs, Q1, Q2, Q3, Q4)

Computes line integrals of the derivatives of basis functions.

Parameters:
  • lambdas (torch.tensor) – placement of each basis function (frequencies)

  • Lx (torch.tensor) – domain size in x,y,z directions

  • Ly (torch.tensor) – domain size in x,y,z directions

  • Lz (torch.tensor) – domain size in x,y,z directions

  • entry (torch.tensor) – entry points of each segment of each measurement ray

  • exit (torch.tensor) – exit points of each segment of each measurement ray

  • req (list of bool) – Components desired to calculate for.

  • nsegs (torch.tensor) – Number of segments for each line integral. shape=(N,)

  • Q1 (torch.tensor) – Integral partial results computed by defineFourierBasis().

  • Q2 (torch.tensor) – Integral partial results computed by defineFourierBasis().

  • Q3 (torch.tensor) – Integral partial results computed by defineFourierBasis().

  • Q4 (torch.tensor) – Integral partial results computed by defineFourierBasis().

Returns:

torch.tensor of `shape=(n, mm_adj, 6)`, the line integrals of the derivatives of each basis function

regression.bayesian.basisFuncDerivs(lambdas, Lx, Ly, Lz, x, y, z, req)

Calculates the basis function derivatives Beltrami basis components.

Parameters:
  • lambdas (torch.tensor) – Basis function frequencies.

  • Lx (torch.tensor) – Basis function phases in x,y,z domain.

  • Ly (torch.tensor) – Basis function phases in x,y,z domain.

  • Lz (torch.tensor) – Basis function phases in x,y,z domain.

  • x (torch.tensor) – Coordinates at which to calculate basis functions at.

  • y (torch.tensor) – Coordinates at which to calculate basis functions at.

  • z (torch.tensor) – Coordinates at which to calculate basis functions at.

  • req (list of bool) – Components desired to calculate for.

Returns:

Partial derivatives of the basis function, `shape=6,n,m)`.

Return type:

(torch.tensor of dict)

regression.bayesian.defineBeltramiBasis(m, theta, nhat=None)

Defines a set of Fourier basis functions used in the bayesian regression.

Parameters:
  • m (torch.tensor) – Number of basis functions in each frequency dimension.

  • theta (torch.tensor) – Hyper parameters defining covariance function.

  • nhat (torch.tensor) – X-ray beam directionin sample coordinate system.

Returns:

tuple of dictionaries dictionaries of type defineFourierBasis(). one dictionary for each beltrami tensor componenent. Also returns a concatenated tensor with all dictionaries contained.

Return type:

(tuple of dict)

regression.bayesian.defineFourierBasis(m, lx, ly, lz, sig_f, nhat=None, subset=None)

Defines a set of Fourier basis functions used in the bayesian regression.

Parameters:
  • m (torch.tensor) – number of basis functions along each axis

  • lx (torch.tensor) – lengthscale along x axis

  • ly (torch.tensor) – lengthscale along y axis

  • lz (torch.tensor) – lengthscale along z axis

  • sig_f (torch.tensor) – prior signal uncertainty

  • nhat (torch.tensor) – the unit vector directions of the measurements

Returns:

dictionary with keys and values:

lambdas (torch.tensor): placement of each basis function (frequencies)

SLambda (torch.tensor): spectral densities for each basis function (i.e. prior variances)

Lx,Ly,Lz (torch.tensor): domain size in x,y,z directions

Q1,Q2,Q3,Q4 (torch.tensor): Integral partial results to be used in basisFuncDInts().

Return type:

(dict)

regression.bayesian.negativeLogLikelihood(Phi, Y, sig_m, v, r)

Computes the negative log likelihood based on gaussian measurement model.

cost = 0.5*error^T Sigma^{-1} error + 0.5 log(det(Sigma)) + C where Sigma is the predictions of the measurement variances Sigma = Variance of the predicted value + variance of the measurements (sig_m^2) This computation is done in square root form for computational efficiency and numerical accuracy

Parameters:
  • Phi (torch.tensor) – the regressor matrix

  • Y (torch.tensor) – the measurement values

  • sig_m (torch.tensor) – the measurement standard deviations, either a scalar or an N,1 vector of individual stds

  • v (torch.tensor) – the fit basis function coefficients

  • r (torch.tensor) – the cholesky factor calculated during regression

Returns:

The negative log likelihood

regression.bayesian.nll_scipy_wrapper_RT_aniso(logtheta, Y_train, entry_train, exit_train, nsegs_train, L_train, kappa_train, Y_test, entry_test, exit_test, nsegs_test, L_test, kappa_test, sig_m_train, sig_m_test, Hray_train, Hray_test, m, nhat, train_idx, test_idx)

Computes the negative log likelihood based on gaussian measurement model.

regression.bayesian.predict(Phi, v, r)

Predicts the point location values of a field given regressor matrix.

Based on the fitted basis function coefficients v, and the cholesky factor r also returns the cholesky factor of the variances of these predictions

Parameters:
  • Phi (torch.tensor) – The regressor matrix

  • v (torch.tensor) – The fit basis function coefficients

  • r (torch.tensor) – The cholesky factor calculated during regression

Returns:

predictions, chol_var

regression.bayesian.regression(Phi, Y, sig_m, SLambda)

Performs linear bayesian regression to fit basis function coefficients.

Parameters:
  • Phi (torch.tensor) – The regressor matrix

  • Y (torch.tensor) – The measurements

  • sig_m (torch.tensor) – The measurement standard deviation, either a scalar or a N,1 array of individual stds

  • SLambda (torch.tensor) – The prior covariances for the basis functions

Returns:

(v,r) where v is the basis function coefficients and r is the cholesky factor required form stuffs

regression.wlsq.trust_constr_solve(mesh, directions, strains, entry, exit, nhat, weights, beam_width, grad_constraint, maxiter, verbose=True, nprocs=1)

Compute a voxelated strain-tensor field weighted least squares fit to a series of line integral strain measures.

Assigns strain tensors over a 2d input mesh on a per element basis.

Parameters:
  • mesh (numpy array) – Coordinate of mesh element nodes, `shape=(N,4)` for quads.

  • directions (numpy array) – unit vectors along which strains apply.

  • strains (numpy array) – Average strains along lines.

  • entry (numpy array) – Xray grain entry points per line.

  • exit (numpy array) – Xray grain exit points per line.

  • nhat (numpy array) – Xray direction per integral in sample coordinate system.

  • weights (numpy array) – Per measurement weights, higher weight gives more impact of equation.

  • beam_width (float) – Width of beam in units of microns.

  • grad_constraint (float) – In x-y plane smoothing constraint of strain reconstruction. Deviations between two neighboring pixels cannot exceed this value for any strain component.

  • maxiter (int) – Maximum number of WLSQ iterations to perform.

  • verbose (bool) – If to print progress. Defaults to True.

  • nprocs (int) – Number of threads to use when assembling the system matrix. Defaults to 1. (Remember to wrap your code in a if __name__="__main__": when running on multiple threads.)

Returns:

List of strain tensor components each with `shape=mesh.shape`. The order of components is “XX”,”YY”,”ZZ”,”YZ”,”XZ”,”XY”.

Return type:

(list of numpy array)

utils

utils.mesher.mesh_from_polygon(w, poly)

Generate a quadratile 2d mesh over a Shapely polygon object.

Parameters:
  • w (float) – Element side length.

  • poly (Shapely Polygon) – Polygon object to be meshed. The output mesh follows the poly coordinate system.

Returns:

(numpy array) The nodal coordinates of each element in the mesh.

utils.mesher.mesh_to_pixels(mesh, element_width, shape, values)

Generate pixelated map from a mesh coordinate matrix.

Parameters:
  • mesh (numpy array) – Element nodal coordinates

  • element_width (float) – Element side length.

  • shape (tuple) – Desired output array shape, (always uses odd shape)

  • values (numpy array) – Per elemetn values to fill the output array with.

Returns:

(numpy array) Pixelated map of input field stored in values.

utils.save.as_vtk_voxel_volume(file, arrays, names)

Save arrays as fields that are paraview compatible (vtk files).

Parameters:
  • file (string) – Asolute path to save the output at ending with desired filename.

  • arrays (list of numpy array) – The arrays to be saved.

  • names (list of string) – Names of the arrays refered to by paraview. These will also be the default titles in paraview.

utils.save.load_vectors(path)

Load vectors produced and saved by s3dxrd.Id11.peaks_to_vectors().

Parameters:

path (string) – Asolute path ending with filename of file to load.

Returns:

dictionary with fields as specified in s3dxrd.Id11.peaks_to_vectors().

Return type:

(dict)

utils.save.mesh2voxels(file, mesh, zpos, values, names)

Save arrays as fields that are paraview compatible (vtk files).

Parameters:
  • file (string) – Asolute path to save the output at ending with desired filename.

  • mesh (list of numpy array) – Series of 2d x-y element mesh coordinate arrays.

  • zpos (list of float) – z-cordinates of each of the meshes in ´´´mesh´´´.

  • values (list of list of numpy array) – For each mesh this list contains a list for each tensor component and each mesh element. I.e values[i][j][k] gives the value associated to the k:th element in the mesh[i] for the tensor component j.

  • names (list of string) – Names of the tensor component arrays refered to by paraview. These will also be the default titles in paraview.

utils.compliance.build_per_point_compliance_and_mask(compliance, Xgrid, Ygrid, Zgrid, polygons, polygon_orientations, zpos, ystep)

Given a material compliance and orientation field, build a per point complance field and density mask.

This function iterates over `polygons` to check if the points in `Xgrid`, `Ygrid`, `Zgrid` are contained by any polygon. If so, the compliance of the given polygon is assigned to the point, as specified by polygon_orientations and compliance.

Parameters:
  • compliance (numpy array) – Compliance matrix, shape=(6,6).

  • Xgrid (numpy array) – x,y,z coordinate arrays, shape=(m,n,l).

  • Ygrid (numpy array) – x,y,z coordinate arrays, shape=(m,n,l).

  • Zgrid (numpy array) – x,y,z coordinate arrays, shape=(m,n,l).

  • polygons (dict) – See s3dxrd.measurements.Id11.peaks_to_vectors()

  • polygon_orientations (dict) – See s3dxrd.measurements.Id11.peaks_to_vectors()

  • zpos (list) – z-coordinates of polygons.

  • ystep (float) – Distance in microns between gridpoints.

Returns:

tuple with values:

Hpoint (numpy array): Per gridpoint compliance of shape=(m,n,l,6,6). Zero outside of mask.

mask (numpy array): Density mask of shape=(m,n,l). np.nan is assigned to empty points.

Return type:

(tuple)

utils.compliance.rotate_compliance(U, H)

Rotate a series of 6x6 “tensor matrices” according to a series of 3x3 rotation matrices.

Parameters:
  • U (numpy array) – Orientation/rotation tensors, `shape=(Nx3x3)`

  • H (numpy array) – Compliance matrices, `shape=(Nx6x6)`

Returns:

(numpy array) of `shape=(Nx6x6)` rotated compliance matrices.

Indices and tables