pynx.scattering: scattering calculations from vectors of xyz atomic positions and hkl values

Description

This module aims to help computing scattering (X-ray or neutrons) for atomic structures, especially if they are distorted or disordered.

The library uses GPU computing (although parallel CPU computing is also available as a fall-back), with the following platforms:

  • nVidia’s CUDA toolkit and the pyCUDA library

  • OpenCL language, along with pyOpenCL library.

Using GPU computing, PyNX provides fast parallel computation of scattering from large assemblies of atoms (>>1000 atoms) and 1D, 2D or 3D coordinates (>>1000) in reciprocal space.

Typical computing speeds on GPUs more than 10^11 reflections.atoms/s on nVidia cards (3.5x10^11 on 2xTitan, 2x10^11 on a GTX 690, 5x10^10 on a GTX295, 2.5x10^10 on a GTX465), more than 2 orders of magnitude faster than on a CPU.

Note that the main interest of pynx.scattering is the ability to compute scattering from any assembly of atoms (not regularly-spaced) to any set of points in reciprocal space. While a FFT will always be faster, it is much more restrictive since the FFT imposes a strict relation between the sampling in real (atomic positions) and reciprocal (hkl coordinates) space.

API Reference

This is the reference for scattering calculations. Note that this part of the code is older than the cdi, ptycho and wavefront modules, so that the API is much more basic.

Structure factor calculations

pynx.scattering.fhkl.Fhkl_thread(h, k, l, x, y, z, occ=None, verbose=False, gpu_name='GTX 295', nbCPUthread=None, sz_imag=None, language='OpenCL', cl_platform='')

Compute F(hkl)=SUM_i exp(-2j*pi*(h*x_i + k*y_i + l*z_i)) or equivalently: F(k) =SUM_i exp(-2j*pi*(sx*x_i + sy*y_i + sz*z_i))

nbCPUthread can be used only for CPU computing, when no GPU is available. nbCPUthread can be set to the number of cores available. Using None (the default) makes the program recognize the number of available processors/cores

If sz_imag is not equal to None, then it means that the scattering vector s has an imaginary component (absorption), e.g. because of grazing incidence condition.

Parameters:
  • h – numpy vector of h coordinates (reciprocal lattice units)

  • k – numpy vector of k coordinates (reciprocal lattice units)

  • l – numpy vector of l coordinates (reciprocal lattice units)

  • x – numpy vector of x coordinates (fractional coordinates in crystal space)

  • y – numpy vector of x coordinates (fractional coordinates in crystal space)

  • z – numpy vector of x coordinates (fractional coordinates in crystal space)

  • occ – numpy vector of occupancy. Can be None

  • verbose – if True, will show some computing/threading messages

  • gpu_name – name of the GPu to be used

  • nbCPUthread – if using a CPU, number of threads to use

  • sz_imag – see grazing incidence example

  • language – either ‘cuda’, ‘opencl’ or ‘cpu’

  • cl_platform – the opencl platform to use

Returns:

a tuple with (calculated structure factor as a numpy complex64 array, elapsed time in seconds)

pynx.scattering.fhkl.Fhkl_gold(h, k, l, x, y, z, occ=None, verbose=False, dtype=<class 'numpy.float64'>)
Compute reference value (using CPU) of:

F(hkl)=SUM_i exp(2j*pi*(h*x_i + k*y_i + l*z_i))

Grazing incidence diffraction

This module requires the cctbx library.

class pynx.scattering.gid.Crystal(unitcell, spacegroup, scatterers)

Description of a material content.

GetAtomDensity()

Number of atoms per cubic meter

GetCriticalAngle()

Calculate the critical angle for this material :return: the critical angle, in radians

GetDensity()

kg/m^3

GetF0(scatterer, h, k, l)

Get the thomson scattering factor for a scatterer at given h,k,l reciprocal unit coordinates

Parameters:
  • scatterer – the scatterer

  • h – reciprocal lattice unit

  • k – reciprocal lattice unit

  • l – reciprocal lattice unit

Returns:

a numpy array of the thomson scattering factor with the same shape as h, k, l.

GetLinearAbsorptionCoeff()

Calculate the linear absorption coefficient for this material :return:

GetRefractionIndexDeltaBeta()

Refraction index, n=1-delta-i*beta. returns a tuple delta,beta

class pynx.scattering.gid.DistortedWave(material0, material1, wave, delta=None, beta=None)

Reflected and refracted wave at an interface.

This will compute the reflection and refraction coefficients, as well as the transmitted (complex) wavevector).

Calculation can be done by supplying the Wave object, and either:

  • the materials used as pynx.gid.Crystal objects: - material0: upper (incoming) layer material, can be None (air or vacuum) - material1: lower layer material

  • the delta, beta values for the difference of the refraction index between the lower and the upper layer (delta: real part, beta: imaginary part)

If delta and beta are supplied, then material0 and material1 are ignored, otherwise the refraction index delta and beta will be calculated using cctbx from material1 and material0.

class pynx.scattering.gid.Scatterer(label, site, occup, u, energy, fp=None, fs=None)

Class used to store scatterer properties, including f’ and f’’ values.

Parameters:
  • label – name of the scatterer

  • site – the site to occupy

  • occup – occupancy

  • u – isotropic displacement parameter

  • energy – energy in keV

  • fp – f’

  • fs – f’’

GetF0(stol)

Get the thomson scattering factor for this scatterer

Parameters:

stol – sin(theta)/lambda value, either a scalar or a numpy array

Returns:

pynx.scattering.gid.W2E(x)

nrj->wavelength or wavelength->nrj , nrj in eV, wavelength in Angstroems

pynx.scattering.gid.fhkl_dwba4(x, y, z, h, k, l=None, occ=None, alphai=0.2, alphaf=None, substrate=None, wavelength=1.0, e_par=0.0, e_perp=1.0, gpu_name='CPU', use_fractional=True, language='OpenCL', cl_platform='', separate_paths=False, **kwargs)

Calculate the grazing-incidence X-ray scattered intensity taking into account 4 scattering paths, for a nanostructure object located above a given substrate. The 5th path is the scattering from the substrate, assumed to be below the interface at z=0.

Parameters:
  • x,y,z – coordinates of the atoms in fractionnal coordinates (relative to the substrate unit cell)- if use_fractional==False, these should be given in Angstroems

  • h,k,l – reciprocal space coordinates. If use_fractional==False, these should be given in inverse Angstroems (multiplied by 2pi - physicist ‘k’ convention, abs(k)=4pisin(theta)/lambda, i.e. these correspond to k_x,k_y,k_z).

  • alphaf (alphai,) – incident and outgoing angles, in radians

  • substrate – the substrate material, as a pynx.gid.Crystal object - this will be used to calculate the material refraction index.

  • wavelength – in Angstroems

  • e_par,e_perp – percentage of polarisation parallel and perpendicular to the incident plane

  • use_fractional – if True (the default), then coordinates for atoms and reciprocal space are given relatively to the unit cell, otherwise in Angstroems and 2pi*inverse Angstroems.

Note: Either l or alphaf must be supplied - it is assumed that the lattice coordinates are such that the [001] direction is perpendicular to the surface.

pynx.scattering.gid.fhkl_dwba5(x, y, z, h, k, l=None, occ=None, alphai=0.2, alphaf=None, substrate=None, wavelength=1.0, e_par=0.0, e_perp=1.0, gpu_name='CPU', use_fractional=True, verbose=False, language='OpenCL', cl_platform='', separate_paths=False)

WARNING: this code is still in development, and needs to be checked !

Calculate the grazing-incidence X-ray scattered intensity taking into account 5 scattering paths, for a nanostructure object located above a given substrate. All atoms with z>0 are assumed to be above the surface, and their scattering is computed using the 4 DWBA paths. Atoms with z<=0 are below the surface, and their scattering is computed using a single path, taking into account the refraction and the attenuation length.

Parameters:
  • x,y,z – coordinates of the atoms in fractionnal coordinates (relative to the substrate unit cell)

  • h,k,l – reciprocal space coordinates

  • alphaf (alphai,) – incident and outgoing angles, in radians

  • substrate – the substrate material, as a pynx.gid.Crystal object - this will be used to calculate the material refraction index.

  • wavelength – in Angstroems

  • e_par,e_perp – percentage of polarisation parallel and perpendicular to the incident plane

Note: Either l OR alphaf must be supplied - it is assumed that the lattice coordinates are such that the [001] direction is perpendicular to the surface.

Thomson scattering factor

pynx.scattering.fthomson.f_thomson(s, element)

Calculate the Thomson scattering factor for one atom or ion, for a given value of sin(theta)/lambda.

Parameters:
  • s\(sin(\vartheta)/\lambda\)

  • element – string of the element name (‘H’, ‘Fe’, ‘Cu1+’,…)

Returns: