API¶
Latticegeneration module¶
Code for generating (quasi-)lattices and the corresponding k-vectors
-
latticegen.latticegeneration.
anylattice_gen
(r_k, theta, order, symmetry=6, size=500, kappa=1.0, psi=0.0, shift=array([0, 0]), normalize=False, chunks=(- 1, - 1))¶ Generate a regular lattice of any symmetry. The lattice is generated from the symmetry 360/symmetry degree rotated k-vectors of length r_k, further rotated by theta degrees. Size either an int for a size*size lattice, or tuple (N,M) for a rectangle N*M. Optionally, the lattice can be strained by a factor kappa in direction psi[1].
With higher order frequency components upto order order
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants. 1/r_k is the line spacing in pixels in the resulting image.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
symmetry (int) – symmetry of the lattice.
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
kappa (float, default: 1) – strain/deformation magnitude. 1 corresponds to no strain.
psi (float, default: 0) – Principal strain direction with respect to horizontal.
shift (iterable or array, optional) – shift of the lattice in pixels. Either a pair (x,y) global shift, or an (2xNxM) array where (NxM) corresponds to size.
normalize (bool, default: False) – if true, normalize the output values to the interval [0,1].
chunks (int or pair of int, optional) – dask chunks in which to divide the returned lattice.
- Returns
lattice – The generated lattice.
- Return type
Dask array
See also
References
[1] T. Benschop et al., 2020, https://doi.org/10.1103/PhysRevResearch.3.013153
-
latticegen.latticegeneration.
anylattice_gen_np
(r_k, theta, order=1, symmetry=6, size=50, kappa=1.0, psi=0.0, shift=array([0, 0]))¶ Generate a regular lattice of any symmetry in pure numpy.
The lattice is generated from the symmetry 360/symmetry degree rotated k-vectors of length r_k, further rotated by theta degrees. Size either an int for a size*size lattice, or tuple (N,M) for a rectangle N*M. Optionally, the lattice can be strained by a factor kappa in direction psi[1].
With higher order frequency components upto order order
The generated lattice gets returned as a numpy array. Only usable for small values of size and order.
See also
-
latticegen.latticegeneration.
combine_ks
(kvecs, order=1, return_counts=False)¶ Generate all possible different sums of kvecs upto order.
- Parameters
kvecs (array-like 2xN) – k vectors to combine.
order (int, default: 1) – Number of different kvecs to combine for each resulting tks.
return_counts (bool, default False) – if True, also return the number of possible combinations for each different combination.
- Returns
tks ((2xM) array of float) – All possible unique combinations of kvecs upto order
counts ((1xM) array of int, optional) – Number of combinations for each vectors in tks.
-
latticegen.latticegeneration.
generate_ks
(r_k, theta, kappa=1.0, psi=0.0, sym=6)¶ Generate k-vectors from given parameters.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants. 1/r_k is the line spacing in pixels in the resulting image.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
kappa (float, default: 1) – strain/deformation magnitude. 1 corresponds to no strain. Larger values corresponds to stretching along the psi direction in real space, so compression along the same direction in k-space.
psi (float, default: 0) – Principal strain direction with respect to horizontal in degrees.
sym (int, default 6) – Rotational symmetry of the unstrained lattice.
- Returns
ks
- Return type
np.array (2x(sym + 1))
-
latticegen.latticegeneration.
hexlattice_gen
(r_k, theta, order, size=500, kappa=1.0, psi=0.0, shift=array([0, 0]), **kwargs)¶ Generate a regular hexagonal lattice. The lattice is generated from the six 60 degree rotated k-vectors of length r_k, further rotated by theta degrees. Optionally, the lattice can be strained by a factor kappa in direction psi [1]. Rendered with higher order frequency components upto order order.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants. 1/r_k is the line spacing in pixels in the resulting image.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
kappa (float, default: 1) – strain/deformation magnitude. 1 corresponds to no strain.
psi (float, default: 0) – Principal strain direction with respect to horizontal in degrees.
shift (iterable or array, optional) – shift of the lattice in pixels. Either a pair (x,y) global shift, or an (2xNxM) array where (NxM) corresponds to size.
**kwargs (dict) – Keyword arguments to be passed to anylattice_gen
- Returns
lattice – The generated lattice.
- Return type
Dask array
See also
-
latticegen.latticegeneration.
hexlattice_gen_fast
(r_k, theta, order, size=500, kappa=1.0, psi=0.0, shift=array([0, 0]))¶ Generate a regular hexagonal lattice.
Speed optimized version, losing some mathematical precision. Tested to be accurate down to max(1e-10*r_k, 1e-10) w.r.t. the regular function. The lattice is generated from the six 60 degree rotated k-vectors of length r_k, further rotated by theta degrees. Size either an int for a size*size lattice, or tuple (N,M) for a rectangle N*M. Optionally, the lattice can be strained by a factor kappa in direction psi[1].
With higher order frequency components upto order order
The generated lattice gets returned as a dask array.
See also
-
latticegen.latticegeneration.
physical_lattice_gen
(a_0, theta, order, pixelspernm=10, symmetry='hexagonal', size=500, epsilon=None, delta=0.16, **kwargs)¶ Generate a physical lattice
Wraps anylattice_gen. Using a lattice constant a_0 in nm and a resolution in pixels per nm, generate a rendering of a lattice of size pixels. Optionally, the lattice can be strained by a factor kappa in direction psi[1].
With higher order frequency components upto order order
- Parameters
a_0 (float) – lattice constant in nm
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
pixelspernm (float, default: 10) – number of pixels per nanometer.
symmetry ({'hexagonal', 'trigonal', 'square'}) – symmetry of the lattice.
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
epsilon (float) – Lattice strain
delta (float, default=0.16) – Poisson ratio for lattice strain. Default value corresponds to graphene.
**kwargs (dict) – Keyword arguments to be passed to anylattice_gen
- Returns
lattice – The generated lattice.
- Return type
Dask array
See also
-
latticegen.latticegeneration.
squarelattice_gen
(r_k, theta, order, size=500, kappa=1.0, psi=0.0, shift=array([0, 0]), **kwargs)¶ Generate a regular square lattice.
The lattice is generated from the four 90 degree rotated k-vectors of length r_k, further rotated by theta degrees. Optionally, the lattice can be strained by a factor kappa in direction psi`[1]. Rendered with higher order frequency components upto order `order.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
kappa (float, default: 1) – strain/deformation magnitude. 1 corresponds to no strain.
psi (float, default: 0) – Principal strain direction with respect to horizontal.
shift (iterable or array, optional) – shift of the lattice in pixels. Either a pair (x,y) global shift, or an (2xNxM) array where (NxM) corresponds to size.
**kwargs (dict) – Keyword arguments to be passed to anylattice_gen
- Returns
lattice – The generated lattice.
- Return type
Dask array
See also
-
latticegen.latticegeneration.
trilattice_gen
(r_k, theta, order, size=500, kappa=1.0, psi=0.0, shift=array([0, 0]), **kwargs)¶ Generate a regular trigonal lattice.
The lattice is generated from the six 60 degree rotated k-vectors of length r_k, further rotated by theta degrees. Optionally, the lattice can be strained by a factor kappa in direction psi`[1]. Rendered with higher order frequency components upto order `order.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
kappa (float, default: 1) – strain/deformation magnitude. 1 corresponds to no strain.
psi (float, default: 0) – Principal strain direction with respect to horizontal.
shift (iterable or array, optional) – shift of the lattice in pixels. Either a pair (x,y) global shift, or an (2xNxM) array where (NxM) corresponds to size.
**kwargs (dict) – Keyword arguments to be passed to anylattice_gen
- Returns
lattice – The generated lattice.
- Return type
Dask array
See also
Singularities module¶
Code to generate edge dislocations/singularities in lattices
-
latticegen.singularities.
gen_dists_image
(peaks, imageshape, rmax=55)¶ Generate squared distance to peaks
- Parameters
peaks ((N,2) array_like) – peak locations in image
imageshape (pair of int) – shape of the original image
rmax (float, default 55) – maximum radius to use around each peak
- Returns
res – Array of shape imageshape, with the squared distance to the nearest peak in float, with a maximum value of rmax.
- Return type
array_like
-
latticegen.singularities.
hexlattice_gen_singularity
(r_k, theta, order, size=250, position=[0, 0], shift=array([0, 0]), **kwargs)¶ Generate a hexagonal lattice with a singularity.
Singularity is shifted position from the center. Not yet equivalent to hexlattice_gen_singularity_legacy.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
kappa (float, default: 1) – strain/deformation magnitude. 1 corresponds to no strain.
position (iterable, default [0, 0]) – [x, y] position of the singularity in pixels with respect to the center.
shift (iterable or array, optional) – shift of the lattice. Either a pair (x,y) global shift, or an (2xNxM) array where (NxM) corresponds to size.
**kwargs (dict) – Keyword arguments to be passed to hexlattice_gen
- Returns
lattice – The generated lattice with singularity
- Return type
Dask array
See also
hexlattice_gen
,singularity_shift
-
latticegen.singularities.
hexlattice_gen_singularity_legacy
(r_k, theta, order, size=250)¶ Generate a regular hexagonal lattice of 2*size times 2*size. The lattice is generated from the six 60 degree rotations of k0, further rotated by theta degrees. With higher order frequency components upto order order and containing a singularity in the center. theta != 0 does not yield a correct lattice.
The generated lattice gets returned as a dask array.
-
latticegen.singularities.
refine_peaks
(image, peaks)¶ Refine peak locations using a bivariate spline
Uses scipy.interpolate.RectBivariateSpline to interpolate image, and scipy.optimize.minimize with bounds to find the interpolated maximum.
- Parameters
image (2D array) – image to interpolate
peaks ((2,N) array of ints) – list of peaks to refine.
- Returns
interppeaks – The refined peak locations.
- Return type
(2,N) array of floats
-
latticegen.singularities.
refined_singularity
(r_k, theta=0, order=3, S=500, remove_center_atom=True, position=array([0., 0.]), **kwargs)¶ Created a refined singularity without distorted atoms in the center
Generate a lattice with a dislocation, then extract atom positions using subpixel_peak_local_max. Optionally remove the atom at the position of the dislocation, add a gaussian blob at each atom position.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
S (int, default: 500) – Size of the resulting lattice in pixels. The returned lattice will be square.
position (iterable, default [0, 0]) – [x, y] position of the singularity in pixels with respect to the center.
remove_center_atom (bool, default=True) – whether to remove the atom at the center of the dislocation
**kwargs (dict) – Keyword arguments to be passed to hexlattice_gen
- Returns
lattice – The generated refined lattice with singularity
- Return type
numpy array
See also
-
latticegen.singularities.
singularity_shift
(r_k, theta, size=500, position=[0, 0], alpha=0.0, symmetry=6)¶ Generate the shift of an edge dislocation.
Generate the shift / displacement field for a hexagonal lattice of size, with the Burgers vector at angle alpha (radians) w.r.t. theta (in degrees). Burgers vector has length np.sin(2*np.pi/symmetry) / r_k, to create a first order edge dislocation.
- Parameters
r_k (float) – length of lattice vectors in k-space. Larger r_k correspond to smaller real space lattice constants.
theta (float) – Angle of the first lattice vector with respect to positive horizontal.
order (int) – Order upto which to generate higher frequency components by combining lattice vectors
size (int, or pair of int, default: 500) – Size of the resulting lattice in pixels. if int, the returned lattice will be square.
position (iterable, default: [0, 0]) – [x, y] position of the singularity in pixels with respect to the center.
alpha (float, default 0.0) – angle in radians of the Burgers vector with respect to the first lattice vector.
symmetry (int, default: 6) – symmetry of the lattice.
- Returns
shift – The generated shift corresponding to a singularity, which can be passed to anylattice_gen
- Return type
ndarray
See also
anylattice_gen
-
latticegen.singularities.
subpixel_peak_local_max
(image, **kwargs)¶ A subpixel accurate local peak find
Wraps skimage.feature.peak_local_max to return subpixel accurate peak locations from a BivariateSpline interpolation.
- Parameters
image (2D array) – image to find maxima in
**kwargs (dict) – keyword arguments passed to skimage.feature.peak_local_max. indices=False is not supported.
- Returns
Containing the coordinates of the peaks
- Return type
(N,2) array
Transformations module¶
Common code for geometrical transformations
-
latticegen.transformations.
a_0_to_r_k
(a_0, symmetry=6)¶ Transform realspace lattice constant to r_k
Where $r_k = (a_0sin(2pi/symmetry))^{-1}$. i.e. r_k is 1/(2pi) the reciprocal lattice constant.
- Parameters
a_0 (float) – realspace lattice constant
symmetry (integer) – symmetry of the lattice
- Returns
r_k – length of k-vector in frequency space
- Return type
float
-
latticegen.transformations.
apply_transformation_matrix
(vecs, matrix)¶ Apply transformation matrix to a list of vectors.
Apply transformation matrix matrix to a list of vectors vecs. vecs can either be a list of vectors, or a NxM array, where N is the number of M-dimensional vectors.
- Parameters
vecs (2D array or iterable) – list of vectors to be transformed
matrix (2D array) – Array corresponding to the transformation matrix
- Returns
Transformed vectors
- Return type
2D array
-
latticegen.transformations.
epsilon_to_kappa
(r_k, epsilon, delta=0.16)¶ Convert frequency r_k and strain epsilon to corresponding r_k and kappa as consumed by functions in latticegeneration.
- Returns
r_k2 (float)
kappa (float)
See also
latticegeneration.generate_ks
-
latticegen.transformations.
r_k_to_a_0
(r_k, symmetry=6)¶ Transform r_k to a realspace lattice constant a_0
- Parameters
r_k (float) – length of k-vector in frequency space as used by lattigeneration functions
symmetry (integer) – symmetry of the lattice
- Returns
a_0 – realspace lattice constant
- Return type
float
-
latticegen.transformations.
rotate
(vecs, angle)¶ Rotate 2D vectors vecs by angle radians around the origin.
-
latticegen.transformations.
rotation_matrix
(angle)¶ Create a rotation matrix
an array corresponding to the 2D transformation matrix of a rotation over angle.
- Parameters
angle (float) – rotation angle in radians
- Returns
2D transformation matrix corresponding to the rotation
- Return type
ndarray (2x2)
-
latticegen.transformations.
scaling_matrix
(kappa, dims=2)¶ Create a scaling matrix
Creates a numpy array containing the dims-dimensional scaling matrix scaling the first dimension by kappa.
- Parameters
kappa (float) – scaling factor of first dimension
dims (int, default: 2) – number of dimensions
- Returns
scaling matrix corresponding to scaling the first dimension by a factor kappa
- Return type
ndarray (dims x dims)
-
latticegen.transformations.
strain_matrix
(epsilon, delta=0.16, axis=0, diff=True)¶ Create a scaling matrix corresponding to uniaxial strain
Only works for the 2D case.
- Parameters
epsilon (float) – applied strain
delta (float, default 0.16) – Poisson ratio. default value corresponds to graphene
axis ({0, 1}) – Axis along which to apply the strain.
diff (bool, default True) – Whether to apply in diffraction or real space.
- Returns
scaling matrix corresponding to epsilon strain along axis
- Return type
ndarray (2 x 2)
-
latticegen.transformations.
wrapToPi
(x)¶ Wrap all values of x to the interval [-pi,pi)