https://github.com/FABLE-3DXRD/xrd_simulator/blob/main/docs/source/images/logo.png?raw=true https://img.shields.io/pypi/pyversions/xrd-simulator.svg? https://github.com/FABLE-3DXRD/xrd_simulator/actions/workflows/python-package-run-tests-linux-py38.yml/badge.svg? https://github.com/FABLE-3DXRD/xrd_simulator/actions/workflows/python-package-run-tests-macos-py38.yml/badge.svg? https://github.com/FABLE-3DXRD/xrd_simulator/actions/workflows/python-package-conda-linux-py38.yml/badge.svg? https://github.com/FABLE-3DXRD/xrd_simulator/actions/workflows/python-package-conda-macos-py38.yml/badge.svg? https://github.com/FABLE-3DXRD/xrd_simulator/actions/workflows/pages/pages-build-deployment/badge.svg? https://badge.fury.io/py/xrd-simulator.svg? https://anaconda.org/conda-forge/vsc-install/badges/platforms.svg? https://anaconda.org/conda-forge/xrd_simulator/badges/latest_release_relative_date.svg?

Simulate X-ray Diffraction from Polycrystals in 3D.

https://img.shields.io/badge/stability-alpha-f4d03f.svg?

The X-R ay D iffraction SIMULATOR package defines polycrystals as a mesh of tetrahedral single crystals and simulates diffraction as collected by a 2D discretized detector array while the sample is rocked around an arbitrary rotation axis. The full journal paper associated to the release of this code can be found here:

xrd_simulator: 3D X-ray diffraction simulation software supporting 3D polycrystalline microstructure morphology descriptions Henningsson, A. & Hall, S. A. (2023). J. Appl. Cryst. 56, 282-292. https://doi.org/10.1107/S1600576722011001

xrd_simulator was originally developed with the hope to answer questions about measurement optimization in scanning x-ray diffraction experiments. However, xrd_simulator can simulate a wide range of experimental diffraction setups. The essential idea is that the sample and beam topology can be arbitrarily specified, and their interaction simulated as the sample is rocked. This means that standard “non-powder” experiments such as scanning-3dxrd and full-field 3dxrd (or HEDM if you like) can be simulated as well as more advanced measurement sequences such as helical scans for instance. It is also possible to simulate powder like scenarios using orientation density functions as input.

Introduction

Before reading all the boring documentation (which is hosted here) let’s dive into some end to end examples to get us started on a good flavour.

The xrd_simulator is built around four python objects which reflect a diffraction experiment:

  • A beam of x-rays (using the xrd_simulator.beam module)

  • A 2D area detector (using the xrd_simulator.detector module)

  • A 3D polycrystal sample (using the xrd_simulator.polycrystal module)

  • A rigid body sample motion (using the xrd_simulator.motion module)

Once these objects are defined it is possible to let the detector collect scattering of the polycrystal as the sample undergoes the prescribed rigid body motion while being illuminated by the xray beam.

Let’s go ahead and build ourselves some x-rays:

import numpy as np
from xrd_simulator.beam import Beam
# The beam of xrays is represented as a convex polyhedron
# We specify the vertices in a numpy array.
beam_vertices = np.array([
    [-1e6, -500., -500.],
    [-1e6, 500., -500.],
    [-1e6, 500., 500.],
    [-1e6, -500., 500.],
    [1e6, -500., -500.],
    [1e6, 500., -500.],
    [1e6, 500., 500.],
    [1e6, -500., 500.]])

beam = Beam(
    beam_vertices,
    xray_propagation_direction=np.array([1., 0., 0.]),
    wavelength=0.28523,
    polarization_vector=np.array([0., 1., 0.]))

We will also need to define a detector:

from xrd_simulator.detector import Detector
# The detector plane is defined by it's corner coordinates det_corner_0,det_corner_1,det_corner_2
detector = Detector(pixel_size_z=75.0,
                    pixel_size_y=55.0,
                    det_corner_0=np.array([142938.3, -38400., -38400.]),
                    det_corner_1=np.array([142938.3, 38400., -38400.]),
                    det_corner_2=np.array([142938.3, -38400., 38400.]))

Next we go ahead and produce a sample, to do this we need to first define a mesh that describes the topology of the sample, in this example we make the sample shaped as a ball:

from xrd_simulator.mesh import TetraMesh
# xrd_simulator supports several ways to generate a mesh, here we
# generate meshed solid sphere using a level set.
mesh = TetraMesh.generate_mesh_from_levelset(
    level_set=lambda x: np.linalg.norm(x) - 768.0,
    bounding_radius=769.0,
    max_cell_circumradius=450.)

Every element in the sample is composed of some material, or “phase”, we define the present phases in a list of xrd_simulator.phase.Phase objects, in this example only a single phase is present:

from xrd_simulator.phase import Phase
quartz = Phase(unit_cell=[4.926, 4.926, 5.4189, 90., 90., 120.],
               sgname='P3221',  # (Quartz)
               path_to_cif_file=None  # phases can be defined from crystalographic information files
               )

The polycrystal sample can now be created. In this example the crystal elements have random orientations and the strain is uniformly zero in the sample:

from scipy.spatial.transform import Rotation as R
from xrd_simulator.polycrystal import Polycrystal
orientation = R.random(mesh.number_of_elements).as_matrix()
polycrystal = Polycrystal(mesh,
                          orientation,
                          strain=np.zeros((3, 3)),
                          phases=quartz,
                          element_phase_map=None)

We may save the polycrystal to disc by using the builtin save() command as

polycrystal.save('my_polycrystal', save_mesh_as_xdmf=True)

We can visualize the sample by loading the .xdmf file into your favorite 3D rendering program. In paraview the sampled colored by one of its Euler angles looks like this:

https://github.com/FABLE-3DXRD/xrd_simulator/blob/main/docs/source/images/example_polycrystal_readme.png?raw=true

We can now define some motion of the sample over which to integrate the diffraction signal:

from xrd_simulator.motion import RigidBodyMotion
motion = RigidBodyMotion(rotation_axis=np.array([0, 1/np.sqrt(2), -1/np.sqrt(2)]),
                         rotation_angle=np.radians(1.0),
                         translation=np.array([123, -153.3, 3.42]))

Now that we have an experimental setup we may collect diffraction by letting the beam and detector interact with the sample:

polycrystal.diffract(beam, detector, motion)
diffraction_pattern = detector.render(frames_to_render=0,
                                        lorentz=False,
                                        polarization=False,
                                        structure_factor=False,
                                        method="project")

The resulting rendered detector frame will look something like the below. Note that the positions of the diffraction spots may vary as the crystal orientations were randomly generated!:

import matplotlib.pyplot as plt
fig,ax = plt.subplots(1,1)
ax.imshow(diffraction_pattern, cmap='gray')
plt.show()
https://github.com/FABLE-3DXRD/xrd_simulator/blob/main/docs/source/images/diffraction_pattern.png?raw=true

To compute several frames simply change the motion and collect the diffraction again. The sample may be moved before each computation using the same or another motion.

polycrystal.transform(motion, time=1.0)
polycrystal.diffract(beam, detector, motion)

Many more options for experimental setups and intensity rendering exist, have fun experimenting! The above example code can be found as a single .py file here.

Installation

Anaconda installation (Linux and Macos)

xrd_simulator is distributed on the conda-forge channel and the preferred way to install the xrd_simulator package is via Anaconda:

conda install -c conda-forge xrd_simulator
conda create -n xrd_simulator
conda activate xrd_simulator

This is meant to work across OS-systems and requires an Anaconda installation.

(The conda-forge feedstock of xrd_simulator can be found here.)

Anaconda installation (Windows)

To install with anaconda on windows you must make sure that external dependencies of pygalmesh are preinstalled on your system. Documentation on installing these package can be found elsewhere.

Pip Installation

Pip installation is possible, however, external dependencies of pygalmesh must the be preinstalled on your system. Installation of these will be OS dependent and documentation can be found elsewhere.:

pip install xrd-simulator

Source installation

Naturally one may also install from the sources:

git clone https://github.com/FABLE-3DXRD/xrd_simulator.git
cd xrd_simulator
python setup.py install

This will then again require the pygalmesh dependencies to be resolved beforehand.

Credits

xrd_simulator makes good use of xfab and pygalmesh. The source code of these repos can be found here:

Citation

If you feel that xrd_simulator was helpful in your research we would love for you to cite us.

xrd_simulator: 3D X-ray diffraction simulation software supporting 3D polycrystalline microstructure morphology descriptions Henningsson, A. & Hall, S. A. (2023). J. Appl. Cryst. 56, 282-292. https://doi.org/10.1107/S1600576722011001

Documentation

polycrystal

The polycrystal module is used to represent a polycrystalline sample. The xrd_simulator.polycrystal.Polycrystal object holds the function xrd_simulator.polycrystal.Polycrystal.diffract() which may be used to compute diffraction. To move the sample spatially, the function xrd_simulator.polycrystal.Polycrystal.transform() can be used. Here is a minimal example of how to instantiate a polycrystal object and save it to disc:

Examples:
import numpy as np
from scipy.spatial.transform import Rotation as R
from xrd_simulator.mesh import TetraMesh
from xrd_simulator.phase import Phase
from xrd_simulator.polycrystal import Polycrystal

# The toplogy of the polycrystal is described by a tetrahedral mesh,
# xrd_simulator supports several ways to generate a mesh, here we
# generate meshed solid sphere using a level set.
mesh = TetraMesh.generate_mesh_from_levelset(
    level_set=lambda x: np.linalg.norm(x) - 768.0,
    bounding_radius=769.0,
    max_cell_circumradius=450.)

# Each element of the mesh is a single crystal with properties defined
# by an xrd_simulator.phase.Phase object.
quartz = Phase(unit_cell=[4.926, 4.926, 5.4189, 90., 90., 120.],
               sgname='P3221',  # (Quartz)
               path_to_cif_file=None  # phases can be defined from crystalographic information files
               )

# The polycrystal can now map phase(s) (only quartz here), orientations and
# strains to the tetrahedral mesh elements.
orientation = R.random(mesh.number_of_elements).as_matrix()
polycrystal = Polycrystal(mesh,
                          orientation,
                          strain=np.zeros((3, 3)),
                          phases=quartz,
                          element_phase_map=None)

# The polycrystal may be saved to disc for later usage.
polycrystal.save('my_polycrystal', save_mesh_as_xdmf=True)
polycrystal_loaded_from_disc = Polycrystal.load('my_polycrystal.pc')

Below follows a detailed description of the polycrystal class attributes and functions.

class xrd_simulator.polycrystal.Polycrystal(mesh, orientation, strain, phases, element_phase_map=None)[source]

Represents a multi-phase polycrystal as a tetrahedral mesh where each element can be a single crystal

The polycrystal is created in laboratory coordinates. At instantiation it is assumed that the sample and lab coordinate systems are aligned.

Parameters
  • mesh (xrd_simulator.mesh.TetraMesh) – Object representing a tetrahedral mesh which defines the geometry of the sample. (At instantiation it is assumed that the sample and lab coordinate systems are aligned.)

  • orientation (numpy array) – Per element orientation matrices (sometimes known by the capital letter U), (shape=(N,3,3)) or (shape=(3,3)) if the orientation is the same for all elements. The orientation matrix maps from crystal coordinates to sample coordinates.

  • strain (numpy array) – Per element (Green-Lagrange) strain tensor, in lab coordinates, (shape=(N,3,3)) or (shape=(3,3)) if the strain is the same for all elements elements.

  • phases (xrd_simulator.phase.Phase or list of xrd_simulator.phase.Phase) – Phase of the polycrystal, or for multiphase samples, a list of all phases present in the polycrystal.

  • element_phase_map (numpy array) – Index of phase that elements belong to such that phases[element_phase_map[i]] gives the xrd_simulator.phase.Phase object of element number i. None if the sample is composed of a single phase. (Defaults to None)

mesh_lab[source]

Object representing a tetrahedral mesh which defines the geometry of the sample in a fixed lab frame coordinate system. This quantity is updated when the sample transforms.

Type

xrd_simulator.mesh.TetraMesh

mesh_sample[source]

Object representing a tetrahedral mesh which defines the geometry of the sample in a sample coordinate system. This quantity is not updated when the sample transforms.

Type

xrd_simulator.mesh.TetraMesh

orientation_lab[source]

Per element orientation matrices mapping from the crystal to the lab coordinate system, this quantity is updated when the sample transforms. (shape=(N,3,3)).

Type

numpy array

orientation_sample[source]

Per element orientation matrices mapping from the crystal to the sample coordinate system., this quantity is not updated when the sample transforms. (shape=(N,3,3)).

Type

numpy array

strain_lab[source]

Per element (Green-Lagrange) strain tensor in a fixed lab frame coordinate system, this quantity is updated when the sample transforms. (shape=(N,3,3)).

Type

numpy array

strain_sample[source]

Per element (Green-Lagrange) strain tensor in a sample coordinate system., this quantity is not updated when the sample transforms. (shape=(N,3,3)).

Type

numpy array

phases[source]

List of all unique phases present in the polycrystal.

Type

list of xrd_simulator.phase.Phase

element_phase_map[source]

Index of phase that elements belong to such that phases[element_phase_map[i]] gives the xrd_simulator.phase.Phase object of element number i.

Type

numpy array

diffract(beam, detector, rigid_body_motion, min_bragg_angle=0, max_bragg_angle=None, verbose=False, number_of_processes=1, number_of_frames=1)[source]

Compute diffraction from the rotating and translating polycrystal while illuminated by an xray beam.

The xray beam interacts with the polycrystal producing scattering units which are stored in a detector frame. The scattering units may be rendered as pixelated patterns on the detector by using a detector rendering option.

Parameters
  • beam (xrd_simulator.beam.Beam) – Object representing a monochromatic beam of xrays.

  • detector (xrd_simulator.detector.Detector) – Object representing a flat rectangular detector.

  • rigid_body_motion (xrd_simulator.motion.RigidBodyMotion) – Rigid body motion object describing the polycrystal transformation as a function of time on the domain (time=[0,1]) over which diffraction is to be computed.

  • min_bragg_angle (float) – Minimum Bragg angle (radians) below which to not compute diffraction. Defaults to 0.

  • max_bragg_angle (float) – Maximum Bragg angle (radians) after which to not compute diffraction. By default the max_bragg_angle is approximated by wrapping the detector corners in a cone with apex at the sample-beam intersection centroid for time=0 in the rigid_body_motion.

  • verbose (bool) – Prints progress. Defaults to True.

  • number_of_processes (int) – Optional keyword specifying the number of desired processes to use for diffraction computation. Defaults to 1, i.e a single processes.

  • number_of_frames (int) – Optional keyword specifying the number of desired temporally equidistantly spaced frames to be collected. Defaults to 1, which means that the detector reads diffraction during the full rigid body motion and integrates out the signal to a single frame. The number_of_frames keyword primarily allows for single rotation axis full 180 dgrs or 360 dgrs sample rotation data sets to be computed rapidly and convinently.

transform(rigid_body_motion, time)[source]

Transform the polycrystal by performing a rigid body motion (translation + rotation)

This function will update the polycrystal mesh (update in lab frame) with any dependent quantities, such as face normals etc. Likewise, it will update the per element crystal orientation matrices (U) as well as the lab frame description of strain tensors.

Parameters
  • rigid_body_motion (xrd_simulator.motion.RigidBodyMotion) – Rigid body motion object describing the polycrystal transformation as a function of time on the domain time=[0,1].

  • time (float) – Time between [0,1] at which to call the rigid body motion.

save(path, save_mesh_as_xdmf=True)[source]

Save polycrystal to disc (via pickling).

Parameters
  • path (str) – File path at which to save, ending with the desired filename.

  • save_mesh_as_xdmf (bool) – If true, saves the polycrystal mesh with associated strains and crystal orientations as a .xdmf for visualization (sample coordinates). The results can be vizualised with for instance paraview (https://www.paraview.org/). The resulting data fields of the mesh data are the 6 unique components of the strain tensor (in sample coordinates) and the 3 Bunge Euler angles (Bunge, H. J. (1982). Texture Analysis in Materials Science. London: Butterworths.). Additionally a single field specifying the material phases of the sample will be saved.

classmethod load(path)[source]

Load polycrystal from disc (via pickling).

Parameters

path (str) – File path at which to load, ending with the desired filename.

Warning

This function will unpickle data from the provied path. The pickle module is not intended to be secure against erroneous or maliciously constructed data. Never unpickle data received from an untrusted or unauthenticated source.

mesh

The mesh module is used to represent the morphology of a polycrystalline sample. Once created and linked to a polycrystal the mesh can be accessed directly through the xrd_simulator.polycrystal.Polycrystal. Here is a minimal example of how to instantiate a mesh and save it to disc:

Examples:
import numpy as np
from xrd_simulator.mesh import TetraMesh

# Generate mesh with 2 elements from nodal coordinates.
nodal_coordinates = np.array([[0,0,0],[0,1,0],[1,0,0],[0,0,1],[0,0,-1]])
element_node_map  = np.array([[0,1,2,3],[0,1,2,4]])
mesh = TetraMesh.generate_mesh_from_vertices(nodal_coordinates, element_node_map)

# The mesh may be saved to disc for later usage or visualization.
mesh.save('my_mesh')
mesh_loaded_from_disc = TetraMesh.load('my_mesh.xdmf')

This should look somethign like this in a 3D viewer like paraview:

https://github.com/FABLE-3DXRD/xrd_simulator/blob/main/docs/source/images/mesh_example.png?raw=true

Below follows a detailed description of the mesh class attributes and functions.

class xrd_simulator.mesh.TetraMesh[source]

Defines a 3D tetrahedral mesh with associated geometry data such face normals, centroids, etc.

For level-set mesh generation the TetraMesh uses the meshio package: For more meshing tools please see this package directly (which itself is a wrapper of CGAL)

coord[source]

Nodal coordinates, shape=(nenodes, 3). Each row in coord defines the coordinates of a mesh node.

Type

numpy array

enod[source]

Tetra element nodes shape=(nelm, nenodes).e.g enod[i,:] gives the nodal indices of element i.

Type

numpy array

dof[source]

Per node degrees of freedom, i.e dof[i,:] gives the degrees of freedom of node i.

Type

numpy array

efaces[source]

Element faces nodal indices, shape=(nelm, nenodes, 3). e.g efaces[i,j,:] gives the nodal indices of face j of element i.

Type

numpy array

enormals[source]

Element faces outwards normals (nelm, nefaces, 3). e.g enormals[i,j,:] gives the normal of face j of element i.

Type

numpy array

ecentroids[source]

Per element centroids, shape=(nelm, 3).

Type

numpy array

eradius[source]

Per element bounding ball radius, shape=(nelm, 1).

Type

numpy array

espherecentroids[source]

Per element bounding ball centroids, shape=(nelm, 3).

Type

numpy array

evolumes[source]

Per element volume, shape=(nelm,).

Type

numpy array

centroid[source]

Global centroid of the entire mesh, shape=(3,)

Type

numpy array

number_of_elements[source]

Number of tetrahedral elements in the mesh.

Type

int

classmethod generate_mesh_from_vertices(coord, enod)[source]

Generate a mesh from vertices using the meshio package:

Parameters
  • coord (numpy array) – Nodal coordinates, shape=(nenodes, 3). Each row in coord defines the coordinates of a mesh node.

  • enod (numpy array) – Tetra element nodes shape=(nelm, nenodes).e.g enod[i,:] gives the nodal indices of element i.

classmethod generate_mesh_from_levelset(level_set, bounding_radius, max_cell_circumradius)[source]

Generate a mesh from a level set using the pygalmesh package:

Parameters
  • level_set (callable) – Level set, level_set(x) should give a negative output on the exterior of the mesh and positive on the interior.

  • bounding_radius (float) – Bounding radius of mesh.

  • max_cell_circumradius (float) – Bound for element radii.

translate(translation_vector)[source]

Translate the mesh.

Parameters

translation_vector (numpy.array) – [x,y,z] translation vector, shape=(3,)

rotate(rotation_axis, angle)[source]

Rotate the mesh.

Parameters
  • rotation_axis (numpy array) – Rotation axis shape=(3,)

  • rotation_angle (float) – Radians for rotation.

update(rigid_body_motion, time)[source]

Apply a rigid body motion transformation to the mesh.

Parameters
  • rigid_body_motion (xrd_simulator.motion.RigidBodyMotion) – Rigid body motion object describing the polycrystal transformation as a function of time on the domain time=[0,1].

  • time (float) – Time between [0,1] at which to call the rigid body motion.

save(file, element_data=None)[source]

Save the tetra mesh to .xdmf paraview readable format for visualization.

Parameters
  • file (str) – Absolute path to save the mesh in .xdmf format.

  • element_data (dict of list or numpy array) – Data associated to the elements.

classmethod load(path)[source]

Load a mesh from a saved mesh file set using the meshio package:

Parameters

file (str) – Absolute path to the mesh file.

phase

The phase module is used to represent material phase. Each element of the xrd_simulator.polycrystal.Polycrystal mesh is linked to a xrd_simulator.phase.Phase describing things like the lattice reference unit cell and generating Miller indices of scattering planes. By providing a crystalographic information file (.cif) structure factors can be computed.

Here is a minimal example of how to instantiate a Phase object

Examples:
from xrd_simulator.phase import Phase

# The phase can be specified via a crystalographic information file.
quartz = Phase(unit_cell=[4.926, 4.926, 5.4189, 90., 90., 120.],
               sgname='P3221',
               path_to_cif_file='quartz.cif'
               )

The .cif file used in the above example can be found here.

Below follows a detailed description of the Phase class attributes and functions.

class xrd_simulator.phase.Phase(unit_cell, sgname, path_to_cif_file=None)[source]

Defines properties related to a crystal class.

The barebone Phase object holds the space group name and unit cell of the crystal. From this it is possible to compute a set of Miller indices for a given wavelength that can give diffraction. In addition, it is possible to compute unit cell structure factors for the hkls which is usefull to model the scattered intensity. This will however require a CIF file.

Parameters
  • unit_cell (list of float) – Crystal unit cell representation of the form [a,b,c,alpha,beta,gamma], where alpha,beta and gamma are in units of degrees while a,b and c are in units of anstrom.

  • sgname (string) – Name of space group , e.g ‘P3221’ for quartz, SiO2, for instance

  • path_to_cif_file (string) – Path to CIF file. Defaults to None, in which case no structure factors are computed, i.e structure_factors=None.

unit_cell[source]

Crystal unit cell representation of the form [a,b,c,alpha,beta,gamma], where alpha,beta and gamma are in units of degrees while a,b and c are in units of anstrom.

Type

list of float

sgname[source]

Name of space group , e.g ‘P3221’ for quartz, SiO2, for instance

Type

string

miller_indices[source]

Allowable integer Miller indices (h,k,l) of shape=(n,3).

Type

numpy array

structure_factors[source]

Structure factors of allowable Miller indices (`miller_indices`) of shape=(n,2). structure_factors[i,0] gives the real structure factor of hkl=miller_indices[i,:] while structure_factors[i,0] gives the corresponding imaginary part of the structure factor.

Type

numpy array

path_to_cif_file[source]

Path to CIF file.

Type

string

setup_diffracting_planes(wavelength, min_bragg_angle, max_bragg_angle)[source]

Generates all Miller indices (h,k,l) that will diffract given wavelength and Bragg angle bounds.

If self.path_to_cif_file is not None, structure factors are computed in addition to the hkls.

Parameters
  • wavelength (float) – xray wavelength in units of anstrom.

  • min_bragg_angle (float) – Maximum Bragg angle, in radians, allowed to be taken during diffraction.

  • max_bragg_angle (float) – Minimum Bragg angle, in radians, allowed to be taken during diffraction.

NOTE: This function will skip Miller indices that have a zero intensity due to the unit cell structure factor vanishing, i.e forbidden reflections, such as a 100 in an fcc for instance, will not be included.

beam

The beam module is used to represent a beam of xrays. The idea is to create a xrd_simulator.beam.Beam object and pass it along to the xrd_simulator.polycrystal.Polycrystal.diffract() function to compute diffraction from the sample for the specified xray beam geometry. Here is a minimal example of how to instantiate a beam object and save it to disc:

Examples:
import numpy as np
from xrd_simulator.beam import Beam

# The beam of xrays is represented as a convex polyhedron
# We specify the vertices in a numpy array.
beam_vertices = np.array([
    [-5., 0., 0.],
    [-5., 1., 0.],
    [-5., 0., 1.],
    [-5., 1., 1.],
    [5., 0., 0.],
    [5., 1., 0.],
    [5., 0., 1.],
    [5., 1., 1.]])

# The xray beam object is instantiated
beam = Beam(
    beam_vertices,
    xray_propagation_direction=np.array([1., 0., 0.]),
    wavelength=0.28523,
    polarization_vector=np.array([0., 1., 0.]))

# The xray beam may be saved to disc for later usage.
beam.save('my_xray_beam')
beam_loaded_from_disc = Beam.load('my_xray_beam.beam')

Below follows a detailed description of the beam class attributes and functions.

class xrd_simulator.beam.Beam(beam_vertices, xray_propagation_direction, wavelength, polarization_vector)[source]

Represents a monochromatic xray beam as a convex polyhedra with uniform intensity.

The beam is described in the laboratory coordinate system.

Parameters
  • beam_vertices (numpy array) – Vertices of the xray beam in units of microns, shape=(N,3).

  • xray_propagation_direction (numpy array) – Propagation direction of xrays, shape=(3,).

  • wavelength (float) – Xray wavelength in units of angstrom.

  • polarization_vector (numpy array) – Beam linear polarization unit vector shape=(3,). Must be orthogonal to the xray propagation direction.

vertices[source]

Vertices of the xray beam in units of microns, shape=(N,3).

Type

numpy array

wavelength[source]

Xray wavelength in units of angstrom.

Type

float

wave_vector[source]

Beam wavevector with norm 2*pi/wavelength, shape=(3,)

Type

numpy array

polarization_vector[source]

Beam linear polarization unit vector shape=(3,). Must be orthogonal to the xray propagation direction.

Type

numpy array

centroid[source]

Beam convex hull centroid shape=(3,)

Type

numpy array

halfspaces[source]

Beam halfspace equation coefficients shape=(N,3). A point x is on the interior of the halfspace if: halfspaces[i,:-1].dot(x) + halfspaces[i,-1] <= 0.

Type

numpy array

set_beam_vertices(beam_vertices)[source]

Set the beam vertices defining the beam convex hull and update all dependent quantities.

Parameters

beam_vertices (numpy array) – Vertices of the xray beam in units of microns, shape=(N,3).

contains(points)[source]

Check if the beam contains a number of point(s).

Parameters

points (numpy array) – Point(s) to evaluate shape=(3,n) or shape=(3,).

Returns

numpy array with 1 if the point is contained by the beam and 0 otherwise, if single point is passed this returns scalar 1 or 0.

intersect(vertices)[source]

Compute the beam intersection with a convex polyhedra.

Parameters

vertices (numpy array) – Vertices of a convex polyhedra with shape=(N,3).

Returns

A scipy.spatial.ConvexHull object formed from the vertices of the intersection between beam vertices and input vertices.

save(path)[source]

Save the xray beam to disc (via pickling).

Parameters

path (str) – File path at which to save, ending with the desired filename.

classmethod load(path)[source]

Load the xray beam from disc (via pickling).

Parameters

path (str) – File path at which to load, ending with the desired filename.

Warning

This function will unpickle data from the provied path. The pickle module is not intended to be secure against erroneous or maliciously constructed data. Never unpickle data received from an untrusted or unauthenticated source.

detector

The detector module is used to represent a 2D area detector. After diffraction from a xrd_simulator.polycrystal.Polycrystal has been computed, the detector can render the scattering as a pixelated image via the xrd_simulator.detector.Detector.render() function.

Here is a minimal example of how to instantiate a detector object and save it to disc:

Examples:
import numpy as np
from xrd_simulator.detector import Detector

# The detector is defined by it's corner coordinates det_corner_0, det_corner_1, det_corner_2
detector = Detector(pixel_size_z=75.0,
                    pixel_size_y=55.0,
                    det_corner_0=np.array([142938.3, -38400., -38400.]),
                    det_corner_1=np.array([142938.3, 38400., -38400.]),
                    det_corner_2=np.array([142938.3, -38400., 38400.]))

# The detector may be saved to disc for later usage.
detector.save('my_detector')
detector_loaded_from_disc = Detector.load('my_detector.det')

Below follows a detailed description of the detector class attributes and functions.

class xrd_simulator.detector.Detector(pixel_size_z, pixel_size_y, det_corner_0, det_corner_1, det_corner_2)[source]

Represents a rectangular 2D area detector.

The detector collects xrd_simulator.scattering_unit.ScatteringUnit during diffraction from a xrd_simulator.polycrystal.Polycrystal. The detector implements various rendering of the scattering as 2D pixelated images. The detector geometry is described by specifying the locations of three detector corners. The detector is described in the laboratory coordinate system.

Parameters
  • pixel_size_z (float) – Pixel side length along zdhat (rectangular pixels) in units of microns. (zdhat is the unit vector from det_corner_0 towards det_corner_2)

  • pixel_size_y (float) – Pixel side length along ydhat (rectangular pixels) in units of microns. (ydhat is the unit vector from det_corner_0 towards det_corner_1)

  • det_corner_0 (numpy array) – Detector corner 3d coordinates shape=(3,). The origin of the detector is at det_corner_0.

  • det_corner_1 (numpy array) – Detector corner 3d coordinates shape=(3,). The origin of the detector is at det_corner_0.

  • det_corner_2 (numpy array) – Detector corner 3d coordinates shape=(3,). The origin of the detector is at det_corner_0.

pixel_size_z[source]

Pixel side length along zdhat (rectangular pixels) in units of microns.

Type

float

pixel_size_y[source]

Pixel side length along ydhat (rectangular pixels) in units of microns.

Type

float

det_corner_0,det_corner_1,det_corner_2

Detector corner 3d coordinates shape=(3,). The origin of the detector is at det_corner_0.

Type

numpy array

frames[source]

Analytical diffraction patterns which

Type

list of list of scattering_unit.ScatteringUnit

zdhat,ydhat

Detector basis vectors.

Type

numpy array

normal[source]

Detector normal, fromed as the cross product: numpy.cross(self.zdhat, self.ydhat)

Type

numpy array

zmax,ymax

Detector width and height.

Type

numpy array

pixel_coordinates[source]

Real space 3d detector pixel coordinates. shape=(n,3)

Type

numpy array

point_spread_function[source]

Scalar point spread function called as point_spread_function(z, y). The z and y coordinates are assumed to be in local units of pixels. I.e point_spread_function(0, 0) returns the value of the pointspread function at the location of the point being spread. This is meant to model blurring due to detector optics. Defaults to a Gaussian with standard deviation 1.0 and mean at (z,y)=(0,0).

Type

callable

property point_spread_kernel_shape[source]

Number of pixels in zdhat and ydhat over which to apply the pointspread function for each scattering event. I.e the shape of the kernel that will be convolved with the diffraction pattern. The values of the kernel is defined by the point_spread_function. Defaults to shape (5, 5).

NOTE: The point_spread_function is automatically normalised over the point_spread_kernel_shape domain such that the

final convolution over the detector diffraction pattern is intensity preserving.

Type

point_spread_kernel_shape (tuple)

render(frames_to_render, lorentz=True, polarization=True, structure_factor=True, method='centroid', verbose=True, number_of_processes=1)[source]

Render a pixelated diffraction pattern onto the detector plane .

NOTE: The value read out on a pixel in the detector is an approximation of the integrated number of counts over the pixel area.

Parameters
  • frames_to_render (int or iterable of int or str) – Indices of the frame in the frames list to be rendered. Optionally the keyword string ‘all’ can be passed to render all frames of the detector.

  • lorentz (bool) – Weight scattered intensity by Lorentz factor. Defaults to False.

  • polarization (bool) – Weight scattered intensity by Polarization factor. Defaults to False.

  • structure_factor (bool) – Weight scattered intensity by Structure Factor factor. Defaults to False.

  • method (str) – Rendering method, must be one of `project` , `centroid` or `centroid_with_scintillator`. Defaults to `centroid`. The default,```method=centroid```, is a simple deposit of intensity for each scattering_unit onto the detector by tracing a line from the sample scattering region centroid to the detector plane. The intensity is deposited into a single detector pixel regardless of the geometrical shape of the scattering_unit. If instead `method=project` the scattering regions are projected onto the detector depositing a intensity over possibly several pixels as weighted by the optical path lengths of the rays diffracting from the scattering region. If `method=centroid_with_scintillator a centroid type raytracing is used, but the scintillator point spread is applied before deposition unto the pixel grid, revealing sub pixel shifts in the scattering events.

  • verbose (bool) – Prints progress. Defaults to True.

  • number_of_processes (int) – Optional keyword specifying the number of desired processes to use for diffraction computation. Defaults to 1, i.e a single processes.

Returns

A pixelated frame as a (numpy array) with shape inferred form the detector geometry and pixel size.

NOTE: This function can be overwitten to do more advanced models for intensity.

pixel_index_to_theta_eta(incoming_wavevector, pixel_zd_index, pixel_yd_index, scattering_origin=array([0, 0, 0]))[source]

Compute bragg angle and azimuth angle for a detector pixel index.

Parameters
  • pixel_zd_index (float) – Coordinate in microns along detector zd axis.

  • pixel_yd_index (float) – Coordinate in microns along detector yd axis.

  • (obj (scattering_origin) – numpy array): Origin of diffraction in microns. Defaults to np.array([0, 0, 0]).

Returns

(tuple) Bragg angle theta and azimuth angle eta (measured from det_corner_1 - det_corner_0 axis) in radians

pixel_coord_to_theta_eta(incoming_wavevector, pixel_zd_coord, pixel_yd_coord, scattering_origin=array([0, 0, 0]))[source]

Compute bragg angle and azimuth angle for a detector coordinate.

Parameters
  • pixel_zd_coord (float) – Coordinate in microns along detector zd axis.

  • pixel_yd_coord (float) – Coordinate in microns along detector yd axis.

  • (obj (scattering_origin) – numpy array): Origin of diffraction in microns. Defaults to np.array([0, 0, 0]).

Returns

(tuple) Bragg angle theta and azimuth angle eta (measured from det_corner_1 - det_corner_0 axis) in radians

get_intersection(ray_direction, source_point)[source]

Get detector intersection in detector coordinates of a single ray originating from source_point.

Parameters
  • ray_direction (numpy array) – Vector in direction of the xray propagation

  • source_point (numpy array) – Origin of the ray.

Returns

(tuple) zd, yd in detector plane coordinates.

contains(zd, yd)[source]

Determine if the detector coordinate zd,yd lies within the detector bounds.

Parameters
  • zd (float) – Detector z coordinate

  • yd (float) – Detector y coordinate

Returns

(boolean) True if the zd,yd is within the detector bounds.

project(scattering_unit, box)[source]

Compute parametric projection of scattering region unto detector.

Parameters
  • scattering_unit (xrd_simulator.ScatteringUnit) – The scattering region.

  • box (tuple of int) – indices of the detector frame over which to compute the projection. i.e the subgrid of the detector is taken as: array[[box[0]:box[1], box[2]:box[3]].

Returns

(numpy array) clip lengths between scattering_unit polyhedron and rays traced from the detector.

get_wrapping_cone(k, source_point)[source]

Compute the cone around a wavevector such that the cone wraps the detector corners.

Parameters
  • k (numpy array) – Wavevector forming the central axis of cone `shape=(3,)`.

  • source_point (numpy array) – Origin of the wavevector `shape=(3,)`.

Returns

(float) Cone opening angle divided by two (radians), corresponding to a maximum bragg angle after

which scattering will systematically miss the detector.

save(path)[source]

Save the detector object to disc (via pickling).

Parameters

path (str) – File path at which to save, ending with the desired filename.

classmethod load(path)[source]

Load the detector object from disc (via pickling).

Parameters

path (str) – File path at which to load, ending with the desired filename.

Warning

This function will unpickle data from the provied path. The pickle module is not intended to be secure against erroneous or maliciously constructed data. Never unpickle data received from an untrusted or unauthenticated source.

motion

The motion module is used to represent a rigid body motion. During diffraction from a xrd_simulator.polycrystal.Polycrystal the xrd_simulator.motion.RigidBodyMotion object describes how the sample is translating and rotating. The motion can be used to update the polycrystal position via the xrd_simulator.polycrystal.Polycrystal.transform() function.

Here is a minimal example of how to instantiate a rigid body motion object, apply the motion to a pointcloud and save the motion to disc:

Examples:
import numpy as np
from xrd_simulator.motion import RigidBodyMotion

# The motion is described by a rotation axis and a translation vector
motion = RigidBodyMotion(rotation_axis=np.array([0, 1/np.sqrt(2), -1/np.sqrt(2)]),
                         rotation_angle=np.radians(2.0),
                         translation=np.array([123, -153.3, 3.42]))

# The motion can be applied to transform a set of points
points = np.random.rand(3, 10)
transformed_points = motion( points, time=0.421 )

# The motion may be saved to disc for later usage.
motion.save('my_motion')
motion_loaded_from_disc = RigidBodyMotion.load('my_motion.motion')

Below follows a detailed description of the RigidBodyMotion class attributes and functions.

class xrd_simulator.motion.RigidBodyMotion(rotation_axis, rotation_angle, translation, origin=array([0., 0., 0.]))[source]

Rigid body transformation of euclidean points by an euler axis rotation and a translation.

A rigid body motion is defined in the laboratory coordinates system.

The Motion is parametric in the interval time=[0,1] and will perform a rigid body transformation of a point x by linearly uniformly rotating it from [0, rotation_angle] and translating [0, translation]. I.e if called at a time time=t the motion will first rotate the point t*rotation_angle radians around rotation_axis and next translate the point by the vector t*translation.

Parameters
  • rotation_axis (numpy array) – Rotation axis shape=(3,)

  • rotation_angle (float) – Radians for final rotation, when time=1.

  • translation (numpy array) – Translation vector shape=(3,)

  • origin (numpy array) – Point in space about which the rigid body motion is defined Defaults to the origin (0,0,0). All translations are executed in relation to the origin and all rotation are rotations about the point of origin. shape=(3,)

rotation_axis[source]

Rotation axis shape=(3,)

Type

numpy array

rotation_angle[source]

Radians for final rotation, when time=1.

Type

float

translation[source]

Translation vector shape=(3,)

Type

numpy array

origin[source]

Point in space about which the rigid body motion is defined Defaults to the origin (0,0,0). All translations are executed in relation to the origin and all rotation are rotations about the point of origin. shape=(3,)

Type

numpy array

__call__(vectors, time)[source]

Find the transformation of a set of points at a prescribed time.

Calling this method will execute the rigid body motion with respect to the

currently set origin.

Parameters
  • vectors (numpy array) – A set of points to be transformed (shape=(3,N))

  • time (float) – Time to compute for.

Returns

Transformed vectors (numpy array) of shape=(3,N).

rotate(vectors, time)[source]

Find the rotational transformation of a set of vectors at a prescribed time.

NOTE: This function only applies the rigid body rotation and will not respect the

origin of the motion! This function is intended for rotation of diffraction and wavevectors. Use the __call__ method to perform a physical rigidbody motion respecting the origin.

Parameters
  • vectors (numpy array) – A set of points in 3d euclidean space to be rotated (shape=(3,N))

  • time (float) – Time to compute for.

Returns

Transformed vectors (numpy array) of shape=(3,N).

translate(vectors, time)[source]

Find the translational transformation of a set of points at a prescribed time.

NOTE: This function only applies the rigid body translation.

Parameters
  • vectors (numpy array) – A set of points in 3d euclidean space to be rotated (shape=(3,N))

  • time (float) – Time to compute for.

Returns

Transformed vectors (numpy array) of shape=(3,N).

inverse()[source]

Create an instance of the inverse motion, defined by negative translation- and rotation-axis vectors.

Returns

(xrd_simulator.RigidBodyMotion) The inverse motion with a reversed rotation and translation.

save(path)[source]

Save the motion object to disc (via pickling).

Parameters

path (str) – File path at which to save, ending with the desired filename.

classmethod load(path)[source]

Load the motion object from disc (via pickling).

Parameters

path (str) – File path at which to load, ending with the desired filename.

Warning

This function will unpickle data from the provied path. The pickle module is not intended to be secure against erroneous or maliciously constructed data. Never unpickle data received from an untrusted or unauthenticated source.

templates

The templates module allows for fast creation of a few select sample types and diffraction geometries without having to worry about any of the “under the hood” scripting.

xrd_simulator.templates.s3dxrd(parameters)[source]

Construct a scaning-three-dimensional-xray diffraction experiment.

This is a helper/utility function for quickly creating an experiment. For full controll over the diffraction geometry consider custom creation of the primitive quantities: (xrd_simulator.beam.Beam), and (xrd_simulator.detector.Detector) seperately.

Parameters

parameters (dict) –

Dictionary with fields as

"detector_distance"(float) Distance form sample origin to

detector centre in units of microns.

"number_of_detector_pixels_z"(int) Number of detector pixels

along z-axis.

"number_of_detector_pixels_y"(int) Number of detector pixels

along y-axis.

"detector_center_pixel_z"(float) Intersection pixel coordinate between

beam centroid line and detector along z-axis.

"detector_center_pixel_y"(float) Intersection pixel coordinate between

beam centroid line and detector along y-axis.

"pixel_side_length_z"(float) Detector pixel side length in units

of microns along z-axis.

"pixel_side_length_y"(float) Detector pixel side length in units

of microns along y-axis.

"wavelength" : (float) Wavelength in units of Angstrom.

"beam_side_length_z"(float) Beam side length in units

of microns.

"beam_side_length_y"(float) Beam side length in units

of microns.

"rotation_step"(float) Angular frame integration step in

units of radians.

"rotation_axis"(numpy array) Axis around which to

positively rotate the sample by rotation_step radians.

Returns

(xrd_simulator.beam.Beam), (xrd_simulator.detector.Detector), (xrd_simulator.motion.RigidBodyMotion).

Return type

(xrd_simulator) objects defining an experiment

Examples

import numpy as np
from xrd_simulator import templates

parameters = {
    "detector_distance": 191023.9164,
    "detector_center_pixel_z": 256.2345,
    "detector_center_pixel_y": 255.1129,
    "pixel_side_length_z": 181.4234,
    "pixel_side_length_y": 180.2343,
    "number_of_detector_pixels_z": 512,
    "number_of_detector_pixels_y": 512,
    "wavelength": 0.285227,
    "beam_side_length_z": 512 * 200.,
    "beam_side_length_y": 512 * 200.,
    "rotation_step": 1.0,
    "rotation_axis": np.array([0., 0., 1.0])
}

beam, detector, motion = templates.s3dxrd(parameters)
xrd_simulator.templates.polycrystal_from_odf(orientation_density_function, number_of_crystals, sample_bounding_cylinder_height, sample_bounding_cylinder_radius, unit_cell, sgname, path_to_cif_file=None, maximum_sampling_bin_seperation=0.08726646259971647, strain_tensor=<function <lambda>>)[source]

Fill a cylinder with crystals from a given orientation density function.

The orientation_density_function is sampled by discretizing orientation space over the unit quarternions. Each bin is assigned its apropiate probability, assuming the orientation_density_function is approximately constant over a single bin. Each sampled orientation is constructed by first drawing a random bin and next drawing uniformly from within that bin, again assuming that orientation_density_function is approximately constant over a bin.

Parameters
  • orientation_density_function (callable) – orientation_density_function(x, q) -> float where input variable x is a numpy array of shape (3,) representing a spatial coordinate in the cylinder (x,y,z) and q is a numpy array of shape (4,) representing a orientation in so3 by a unit quarternion. The format of the quarternion is “scalar last” (same as in scipy.spatial.transform.Rotation).

  • number_of_crystals (int) – Approximate number of crystal elements to compose the cylinder volume.

  • sample_bounding_cylinder_height (float) – Height of sample cylinder in units of microns.

  • sample_bounding_cylinder_radius (float) – Radius of sample cylinder in units of microns.

  • unit_cell (list of float) – Crystal unit cell representation of the form [a,b,c,alpha,beta,gamma], where alpha,beta and gamma are in units of degrees while a,b and c are in units of anstrom.

  • sgname (string) – Name of space group , e.g ‘P3221’ for quartz, SiO2, for instance

  • path_to_cif_file (string) – Path to CIF file. Defaults to None, in which case no structure factors are computed.

  • maximum_sampling_bin_seperation (float) – Discretization steplength of orientation space using spherical coordinates over the unit quarternions in units of radians. A smaller steplength gives more accurate sampling of the input orientation_density_function but is computationally slower.

  • strain_tensor (callable) – Strain tensor field over sample cylinder. strain_tensor(x) -> numpy array of shape (3,3) where input variable x is a numpy array of shape (3,) representing a spatial coordinate in the cylinder (x,y,z).

Returns

(xrd_simulator.polycrystal.Polycrystal)

Examples

import numpy as np
from xrd_simulator.templates import polycrystal_from_odf


# uniform orientation distribution function.
def ODF(x, q): return 1. / (np.pi**2)


number_of_crystals = 500
bounding_height = 50.0
bounding_radius = 25.0
unit_cell = [4.926, 4.926, 5.4189, 90., 90., 120.]
sgname = 'P3221',  # Quartz
max_bin = np.radians(10.0)
path_to_cif_file = None

def strain_tensor(x): return np.array([[0, 0, 0.02 * x[2] / bounding_height],
                                       [0, 0, 0],
                                       [0, 0, 0]])  # Linear strain gradient along rotation axis.


polycrystal = polycrystal_from_odf(ODF,
                                   number_of_crystals,
                                   bounding_height,
                                   bounding_radius,
                                   unit_cell,
                                   sgname,
                                   path_to_cif_file,
                                   max_bin,
                                   strain_tensor)

More complicated ODFs are also possible to use. Why not use a von-Mises-Fisher distribution for instance:

Examples

from scipy.special import iv
import numpy as np
# Von Mises-Fisher orientation distribution function.
# https://en.wikipedia.org/wiki/Von_Mises%E2%80%93Fisher_distribution
def ODF(x, q):
    kappa = 0.064 * x.dot(x) + 0.001
    mu    = np.array([ 0.,  0.95689793,  0.28706938, -0.0440173 ])
    p     = 4.0
    I     = iv( p/2 - 1, kappa)
    Cp    = kappa**(p/2 - 1) / ( (2*np.pi)**(p/2) * I )
    return Cp * ( np.exp( kappa * np.dot( q, mu ) ) + np.exp( -(kappa * np.dot( q, mu )) ) )
xrd_simulator.templates.get_uniform_powder_sample(sample_bounding_radius, number_of_grains, unit_cell, sgname, strain_tensor=array([[0., 0., 0.], [0., 0., 0.], [0., 0., 0.]]), path_to_cif_file=None)[source]

Generate a polycyrystal with grains overlayed at the origin and orientations drawn uniformly.

Parameters
  • sample_bounding_radius (float) – Bounding radius of sample. All tetrahedral crystal elements will be overlayed within a sphere of sample_bounding_radius radius.

  • number_of_grains (int) – Number of grains composing the polycrystal sample.

  • unit_cell (list of float) – Crystal unit cell representation of the form [a,b,c,alpha,beta,gamma], where alpha,beta and gamma are in units of degrees while a,b and c are in units of anstrom.

  • sgname (string) – Name of space group , e.g ‘P3221’ for quartz, SiO2, for instance

  • strain_tensor (numpy array) – Strain tensor to apply to all tetrahedral crystal elements contained within the sample. shape=(3,3).

  • path_to_cif_file (string) – Path to CIF file. Defaults to None, in which case no structure factors are computed.

Returns

(xrd_simulator.polycrystal) A polycyrystal sample with number_of_grains grains.

Examples

import numpy as np
from xrd_simulator import templates

polycrystal = templates.get_uniform_powder_sample(
    sample_bounding_radius=1.203,
    number_of_grains=15,
    unit_cell=[4.926, 4.926, 5.4189, 90., 90., 120.],
    sgname='P3221',
    strain_tensor=np.array([[0, 0, -0.0012],
                            [0, 0, 0],
                            [0, 0, 0]])
)

scattering_unit

A scattering unit is generated by the xrd_simulator.polycrystal.Polycrystal as the beam interacts with the sample. For advanced users it can be usefull to interact with the scattering units, this is preferably done via the xrd_simulator.detector.Detector which will hold all the scatterign units created during diffraction from the polycrystal.

class xrd_simulator.scattering_unit.ScatteringUnit(convex_hull, scattered_wave_vector, incident_wave_vector, wavelength, incident_polarization_vector, rotation_axis, time, phase, hkl_indx, element_index)[source]

Defines a scattering region in space as a single crystal as a convex polyhedra.

The scattering unit is described in laboratory coordinates.

Parameters
  • convex_hull (scipy.spatial.ConvexHull) – Object describing the convex hull of the scattering unit.

  • scattered_wave_vector (numpy array) – Scattered wavevector `shape=(3,)`

  • incident_wave_vector (numpy array) – Incident wavevector `shape=(3,)`

  • wavelength (float) – Wavelength of xrays in units of angstrom.

  • incident_polarization_vector (numpy array) – Unit vector of linear polarization `shape=(3,)`

  • rotation_axis (numpy array) – Sample motion rotation axis `shape=(3,)`

  • time (numpy array) – Parametric value in range [0,1] determining at which instance during sample motion the scattering occured.

  • phase (xrd_simulator.phase.Phase) – The phase of the scattering unit.

  • hkl_indx (int) – Index of Miller index in the phase.miller_indices list.

convex_hull[source]

Object describing the convex hull of the scattering unit.

Type

scipy.spatial.ConvexHull

scattered_wave_vector[source]

Scattered wavevector `shape=(3,)`

Type

numpy array

incident_wave_vector[source]

Incident wavevector `shape=(3,)`

Type

numpy array

wavelength[source]

Wavelength of xrays in units of angstrom.

Type

float

incident_polarization_vector[source]

Unit vector of linear polarization `shape=(3,)`

Type

numpy array

rotation_axis[source]

Sample motion rotation axis `shape=(3,)`

Type

numpy array

time[source]

Parametric value in range [0,1] determining at which instance during sample motion the scattering occured.

Type

numpy array

phase[source]

The phase of the scattering unit.

Type

xrd_simulator.phase.Phase

hkl_indx[source]

Index of Miller index in the phase.miller_indices list.

Type

int

element_index[source]

Index of mesh tetrahedral element refering to a xrd_simulator.polycrystal.Polycrystal object from which the scattering unit originated.

Type

int

property hkl[source]

Miller indices [h,k,l] shape=(3,).

Type

hkl (numpy array)

property real_structure_factor[source]

Real part of unit cell structure factor

Type

hkl (numpy array)

property imaginary_structure_factor[source]

Imaginary part of unit cell structure factor

Type

hkl (numpy array)

property lorentz_factor[source]

Compute the Lorentz intensity factor for a scattering_unit.

property polarization_factor[source]

Compute the Polarization intensity factor for a scattering_unit.

property centroid[source]

centroid of the scattering region. shape=(3,)

Type

centroid (numpy array)

property volume[source]

volume of the scattering region volume

Type

volume (float)

laue

Collection of functions for solving time dependent Laue equations for arbitrary rigid body motions. This module is mainly used internally by the xrd_simulator.polycrystal.Polycrystal. However, for the advanced user, access to these functions may be of interest.

xrd_simulator.laue.get_G(U, B, G_hkl)[source]

Compute the diffraction vector

\[\boldsymbol{G} = \boldsymbol{U}\boldsymbol{B}\boldsymbol{G}_{hkl}\]
Parameters
  • U (numpy array) Orientation matrix of shape=(3,3) (unitary) –

  • B (numpy array) – Reciprocal to grain coordinate mapping matrix of shape=(3,3).

  • G_hkl (numpy array) – Miller indices, i.e the h,k,l integers (shape=(3,n)).

Returns

Sample coordinate system diffraction vector. (shape=(3,n))

Return type

G (numpy array)

xrd_simulator.laue.get_bragg_angle(G, wavelength)[source]

Compute a Bragg angle given a diffraction (scattering) vector.

Parameters
  • G (numpy array) – Sample coordinate system diffraction vector. (shape=(3,n))

  • wavelength (float) – Photon wavelength in units of angstrom.

Returns

in units of radians. (shape=(n,))

Return type

Bragg angles (float)

xrd_simulator.laue.get_sin_theta_and_norm_G(G, wavelength)[source]

Compute a Bragg angle given a diffraction (scattering) vector.

Parameters
  • G (numpy array) – Sample coordinate system diffraction vector.

  • wavelength (float) – Photon wavelength in units of angstrom.

Returns

in units of radians and ||G||.

Return type

sin(Bragg angle) (float)

xrd_simulator.laue.find_solutions_to_tangens_half_angle_equation(rho_0, rho_1, rho_2, delta_omega)[source]

Find all solutions, t, to the equation (maximum 2 solutions exists)

\[\rho_0 \cos(t \Delta \omega) + \rho_1 \sin(t \Delta \omega) + \rho_2 = 0. \quad\quad (1)\]

by rewriting as

\[(\rho_2 - \rho_0) s^2 + 2 \rho_1 s + (\rho_0 + \rho_2) = 0. \quad\quad (2)\]

where

\[s = \tan(t \Delta \omega / 2). \quad\quad (3)\]

and

\[\Delta \omega\]

is a rotation angle

Parameters
  • rho_0 (float) – Coefficients rho_0,rho_1 and rho_2 of equation (1).

  • rho_1 (float) – Coefficients rho_0,rho_1 and rho_2 of equation (1).

  • rho_2 (float) – Coefficients rho_0,rho_1 and rho_2 of equation (1).

  • delta_omega (float) – Radians of rotation.

Returns

2 Arrays of solutions of matching shape to input. if no solutions exist on an interval the corresponding instances in the solution array holds np.nan values.

Return type

(tuple of numpy.array)