API Reference#

interfacemaster.cellcalc#

cellcalc.py

class interfacemaster.cellcalc.DSCcalc[source]#

Bases: object

core class computing DSC basis

compute_CNID(hkl, tol=1e-08)[source]#

compute 2D DSC

Parameters:
  • hkl – miller indices of the plane to compute

  • tol (float) – tolerance to judge coincidence, default 1e-8

compute_CSL(tol=1e-08)[source]#

compute CSL

Parameters:

tol (float) – tolerance to judge coincidence, default 1e-8

compute_DSC(to_LLL=True, tol=1e-08)[source]#

solve the DSC equation

Parameters:
  • to_LLL (bool) – whether do LLL algorithm for the found basis, default True

  • tol (float) – tolerance to judge coincidence, default 1e-8

parse_int_U(ai1, ai2, sigma, tol=1e-08)[source]#

initiate computation

Parameters:
  • ai1 (numpy array) – lattice 1 basis set

  • ai2 (numpy array) – lattice 2 basis set

  • sigma (int) – sigma value

  • tol (float) – tolerance to judge coincidence, default 1e-8

interfacemaster.cellcalc.Gram_Schmidt(B0)[source]#

Gram–Schmidt process

Parameters:

B0 (numpy array) – matrix

Returns:

Bstar – the orthogonalized basis for B0

Return type:

numpy array

interfacemaster.cellcalc.LLL(B)[source]#

LLL lattice reduction algorithm

https://en.wikipedia.org/wiki/Lenstra%E2%80%93Lenstra%E2%80%93Lov%C3%A1sz_lattice_basis_reduction_algorithm

Parameters:

B (numpy array) – matrix

Returns:

Bhere – the basis obtained for B

Return type:

numpy array

interfacemaster.cellcalc.MID(lattice, n, tol=1e-08)[source]#

get the Miller indices with a normal n for the lattice

Parameters:
  • lattice (numpy array) – lattice basis set

  • n (numpy array) – normal vector

  • tol (float) – tolerance to judge integer, default 1e-8

Returns:

hkl – the Miller indices of desired plane

Return type:

numpy array

interfacemaster.cellcalc.ang(v1, v2)[source]#

compute the cos of angle between v1 & v2

Parameters:
  • v1 (numpy array) – vectors

  • v2 (numpy array) – vectors

Returns:

cos_12 – the cos of angle between v1 and v2

Return type:

float

interfacemaster.cellcalc.dia_sym_mtx(U)[source]#

return the second-diagonal symmetry transformed matrix of U

Parameters:

U (numpy array) – matrix

Returns:

Ux – the second-diagonal symmetry transformed matrix of U

Return type:

numpy array

interfacemaster.cellcalc.ext_euclid(a, b)[source]#

extended euclidean algorithm

from https://en.wikipedia.org/wiki/Extended_Euclidean_algorithm

Parameters:
  • a (int) – two input integers

  • b (int) – two input integers

Returns:

old_s, old_t, old_r – the output integers

Return type:

int

interfacemaster.cellcalc.find_integer_vectors(v, sigma, tol=1e-08)[source]#

A function for finding the coefficients N so that Nv contains only integer elements and gcd(v1,v2,v3,N) = 1

Parameters:
  • v (numpy array) – a 3D vector

  • sigma (int) – maximum limt of N for searching

  • tol (float) – tolerancea to judge integer, default 1e-8

Returns:

  • Nv (numpy array) – vector of integer elements Nv

  • N (numpy array) – integer coefficients N

interfacemaster.cellcalc.find_integer_vectors_nn(v, sigma, tol=1e-08)[source]#

A function for finding the coefficients N so that Nv contains only integer elements, ignore gcd

Parameters:
  • v (numpy array) – a 3D array

  • sigma (int) – maximum limt of N for searching

  • tol (float) – tolerance, default 1e-8

Returns:

  • Nv (numpy array) – vector of integer elements Nv

  • N (numpy array) – integer coefficients N

interfacemaster.cellcalc.get_indices_from_n_Pc1(n, lattice, Pc1)[source]#

get the Miller indices of certain plane with normal n and one in-plane point Pc1 for certain lattice

Parameters:
  • n (numpy array) – normal vector

  • lattice (numpy array) – lattice basis set

  • Pc1 (numpy array) – an in-plane point in the (hkl) plane

Returns:

hkl – the Miller indices of the desired plane

Return type:

numpy array

interfacemaster.cellcalc.get_normal_from_MI(lattice, hkl)[source]#

get the normal vector of a lattice plane with known miller indices

Parameters:
  • lattice (numpy array) – lattice basis set

  • hkl (numpy array) – the Miller indices

Returns:

n – a normal vector to the (hkl) plane

Return type:

numpy array

interfacemaster.cellcalc.get_normal_index(hkl, lattice)[source]#

get the coordinates in the lattice of a normal vector to the plane (hkl)

Parameters:
  • hkl (array_like) – input index

  • lattice (array_like) – the basis of the lattice

Returns:

v_n – the coordinates in the lattice of a normal vector to the plane (hkl)

Return type:

numpy array

interfacemaster.cellcalc.get_ortho_two_v(B, lim, tol, align_rotation_axis=False, rotation_axis=None)[source]#

get orthogonal cell of a 2D basis

Parameters:
  • B (numpy array) – matrix of 2D basis vectors

  • lim (int) – maximum coefficient

  • tol (float) – tolerance

  • align_rotation_axis (bool) – whether to specify a vector, default False

  • rotation_axis – specified vector, default [1, 1, 1]

Returns:

cell_ortho – a perpendicular 2D basis

Return type:

numpy array

interfacemaster.cellcalc.get_plane(hkl, lattice)[source]#

get the normal vector and one in-plane point for the (hkl) plane of the lattice

Parameters:
  • hkl (array_like) – input index

  • lattice (numpy array) – lattice basis set

Returns:

  • n (numpy array) – the normal vector

  • p (numpy array) – an in-plane point in the (hkl) plane

interfacemaster.cellcalc.get_pri_vec_inplane(hkl, lattice)[source]#

get two primitive lattice vector for the (hkl) plane from Banadaki A D, Patala S. An efficient algorithm for computing the primitive bases of a general lattice plane[J]. Journal of Applied Crystallography, 2015, 48(2): 585-588.

Parameters:
  • hkl (numpy array) – the Miller indices

  • lattice (numpy array) – lattice basis set

Returns:

v1, v2 – the two primitive lattice vectors

Return type:

numpy array

interfacemaster.cellcalc.get_primitive_hkl(hkl, C_lattice, P_lattice, tol=1e-08, n_max=100000)[source]#

convert the miller indices from conventional cell to primitive cell

Parameters:
  • hkl (array_like) – input index

  • C_lattice (numpy array) – conventional cell basis set

  • P_lattice (numpy array) – primitive cell basis set

  • tol (float) – tolerance to judge integer, default 1e-8

  • n_max (int) – maximum limit of N for search

Returns:

hkl_p – the Miller indices in primitive cell

Return type:

numpy array

interfacemaster.cellcalc.get_right_hand(B)[source]#

get right handed lattice of B

Parameters:

B (numpy array) – lattice basis set

Returns:

B_r – the right_handed basis set of lattice B

Return type:

numpy array

interfacemaster.cellcalc.match_rot(deft_rot, axis, tol, exact_rot, av_perpendicular)[source]#

find the intermediate rotation

given an exact rotation matrix and a default rotation matrix, find the intermediate rotation matrix so that the operation combining the default rotation + the intermediate rotation is close to the excat rotation

Parameters:
  • deft_rot (numpy array) – initial orientation matrix

  • axis (numpy array) – the rotation axix

  • tol (float) – tolerance to judge coincidence

  • exact_rot (numpy array) – target orientation matrix

  • av_perpendicular (numpy array) – auxilliary vector

Returns:

inter_rot – intermediate rotation matrix

Return type:

numpy array

interfacemaster.cellcalc.projection(u1, u2)[source]#

get the projection of u1 on u2

Parameters:
  • u1 (numpy array) – vectors

  • u2 (numpy array) – vectors

Returns:

p12 – the projection of u1 on u2

Return type:

numpy array

interfacemaster.cellcalc.rot(a, theta)[source]#

produce a rotation matrix

Parameters:
  • a (numpy array) – rotation axis

  • theta (float) – rotation angle

Returns:

rot_mat – a rotation matrix

Return type:

numpy array

interfacemaster.cellcalc.search_MI_n(lattice, n, tol, lim)[source]#

get the miller indices of planes whose normal vector is close to the input vector n

Parameters:
  • lattice (numpy array) – lattice basis set

  • n (numpy array) – the normal vector

  • tol (float) – tolerance to judge perpendicular

  • lim (int) – upper limit of indices serching

Returns:

hkl – the Miller indices

Return type:

numpy array

interfacemaster.cellcalc.solve_DSC_equations(u, v, w, L, B)[source]#

a function for solving the integer equations for the DSC basis

Parameters:
  • u (int) – integer numbers

  • v (int) – integer numbers

  • w (int) – integer numbers

  • L (int) – integer numbers

  • B (numpy array) – basis set

  • tol (float) – tolerance

Returns:

DSC_basis1, DSC_basis2 – DSC_basis sets

Return type:

numpy array

interfacemaster.interface_generator#

interface_generator.py

interfacemaster.interface_generator.POSCAR_to_cif(Poscar_name, Cif_name)[source]#

generate a cif file for the structure in a POSCAR file

Parameters:
  • Poscar_name (str) – name of a POSCAR file

  • Cif_name (str) – name of a cif file

interfacemaster.interface_generator.Xs_Ys_cell(lattice)[source]#

get the Xs and Ys arrays to draw a cell for a given lattice

interfacemaster.interface_generator.adjust_orientation(lattice)[source]#

adjust the orientation of a lattice so that its first axis is along x-direction and the second axis is in the x-y plane

Parameters:

lattice (numpy array) – a matrix with column lattice vectors

Returns:

lattice, R – rotated lattice, rotation matrix

Return type:

numpy array

interfacemaster.interface_generator.cell_expands(lattice, atoms, elements, xyz)[source]#

expand certain supercell

Parameters:
  • lattice (numpy array) – lattice matrix

  • atoms (numpy array) – list of atom fractional coordinates

  • elements (numpy array) – list of element name of the atoms

  • xyz (list of int) – list of expansion factor for the x, y, z directions

Returns:

  • lattice (numpy array) – lattice matrix

  • atoms (numpy array) – atom coordinates

  • elements (numpy array) – list of element name of the atoms

interfacemaster.interface_generator.clean_fig(num1, num2, axes)[source]#

hide the frames

interfacemaster.interface_generator.convert_vector_index(lattice_0, lattice_f, v_0)[source]#

convert the index of a vector into a different basis

class interfacemaster.interface_generator.core(file_1, file_2, prim_1=True, prim_2=True, verbose=True)[source]#

Bases: object

Core class for dealing with an interface

compute_bicrystal(hkl, lim=20, normal_ortho=False, plane_ortho=False, tol_ortho=1e-10, tol_integer=1e-08, align_rotation_axis=False, rotation_axis=None, inclination_tol=0.7071067811865476)[source]#

compute the transformation to obtain the supercell of the two slabs forming a interface

Parameters:
  • hkl (numpy array) – miller indices of the plane expressed in lattice_1

  • lim (int) – the limit searching for a CSL vector cross the plane, default 20

  • normal_ortho (bool) – whether limit the vector crossing the GB to be normal to the GB, default False

  • plane_ortho (bool) – whether limit the two vectors in the GB plane to be orthogonal, default False

  • tol_ortho (float) – tolerance judging whether orthogonal, default 1e-10

  • tol_integer (float) – tolerance judging integer, default 1e-8

  • align_rotation_axis (bool) – whether to align to rotation axis, default False

  • rotation_axis (list) – rotation axis, defalt [1, 1, 1]

  • inclination_tol (float) – control the angle between the interface and the cross vector, default sqrt(2)/2

compute_bicrystal_two_D(hkl_1, hkl_2, lim=20, normal_ortho=False, tol_ortho=1e-10, tol_integer=1e-08, inclination_tol=0.7071067811865476)[source]#

compute the transformation to obtain the supercell of the two slabs forming a interface (only two_D periodicity)

Parameters:
  • hkl1 (numpy array) – miller indices of the plane expressed in lattice_1 and lattice_2

  • hkl2 (numpy array) – miller indices of the plane expressed in lattice_1 and lattice_2

  • lim (int) – the limit searching for a CSL vector cross the plane, default 20

  • normal_ortho (bool) – whether limit the vector crossing the GB to be normal to the GB, default False

  • tol_ortho (float) – tolerance judging whether orthogonal, default 1e-10

  • tol_integer (float) – tolerance judging integer, default 1e-8

  • inclination_tol (float) – control the angle between the interface and the cross vector, default sqrt(2)/2

define_lammps_regions(region_names, region_los, region_his, ortho=False)[source]#

generate a file defining some regions in the LAMMPS and define the atoms inside these regions into some groups.

Parameters:
  • region_names (list) – list of name of regions

  • region_los (list) – list of the low bounds

  • region_his (list) – list of the hi bounds

  • ortho (bool) – whether the cell is orthogonal

get_bicrystal(dydz=None, dx=0, dp1=0, dp2=0, xyz_1=None, xyz_2=None, vx=0, filename='POSCAR', two_D=False, filetype='VASP', mirror=False, KTI=False)[source]#

generate a cif file for the bicrystal structure

Parameters:
  • dydz (numpy array) – translation vector in the interface

  • dx (float) – translation normal to the interface

  • dp1 (float) – termination of slab 1, 2

  • dp2 (float) – termination of slab 1, 2

  • xyz1 (list) – expansion of slab 1, 2

  • xyz2 (list) – expansion of slab 1, 2

  • vx (float) – vacuum spacing, default 0

  • filename (str) – filename, default ‘POSCAR’

  • two_D (bool) – whether a two CSL

  • filetype (str) – filetype, ‘VASP’ or ‘LAMMPS’, default ‘VASP’

  • mirror (bool) – mirror, default False

  • KTI (bool) – KTI, default False

parse_limit(du, S, sgm1, sgm2, dd)[source]#

set the limitation to accept an appx CSL

Parameters:
  • du (parameters) – see the paper

  • S (parameters) – see the paper

  • sgm1 (parameters) – see the paper

  • sgm2 (parameters) – see the paper

  • dd (parameters) – see the paper

sample_CNID(grid, dx=0, dp1=0, dp2=0, xyz_1=None, xyz_2=None, vx=0, two_D=False, filename='POSCAR', filetype='VASP')[source]#

sampling the CNID and generate POSCARs

Parameters:

grid (numpy array) – 2D grid of sampling

sample_lattice_planes(dx=0, xyz_1=None, xyz_2=None, vx=0, two_D=False, filename='POSCAR', filetype='VASP')[source]#

sampling non-identical lattice planes terminating at the GB

scale(factor_1, factor_2)[source]#

scale the lattice 1 and 2 with specific factors (both the lattice and conventional lattice)

Parameters:
  • factor_1 (float) – scale factor for lattice 1

  • factor_2 (float) – scale factor for lattice 2

search_all_position(axis, theta, theta_range, dtheta, two_D=False)[source]#

main loop finding all the CSL lattices satisfying the limit

Parameters:
  • axis (numpy array) – rotation axis

  • theta (float) – initial rotation angle, in degree

  • theta_range (list) – range varying theta, in degree

  • dtheta (float) – step varying theta, in degree

Notes

log.all_position :

all the searching information

results :

information of the found approximate CSL

search_fixed(R, exact=False, tol=1e-08)[source]#

main loop finding the appx CSL

Parameters:
  • axis (numpy array) – rotation axis

  • theta (float) – initial rotation angle, in degree

  • theta_range (float) – range varying theta, in degree

  • dtheta (float) – step varying theta, in degree

search_one_position(axis, theta, theta_range, dtheta, two_D=False)[source]#

main loop finding the appx CSL

Parameters:
  • axis (numpy array) – rotation axis

  • theta (float) – initial rotation angle, in degree

  • theta_range (float) – range varying theta, in degree

  • dtheta (float) – step varying theta, in degree

search_one_position_2D(hkl_1, hkl_2, theta_range, dtheta, pre_dt=False, pre_R=array([[1., 0., 0.], [0., 1., 0.], [0., 0., 1.]]), integer_tol=1e-08, start=0, exact=False)[source]#

main loop finding the appx CSL

Parameters:
  • axis (numpy array) – rotation axis

  • theta (float) – initial rotation angle, in degree

  • theta_range (float) – range varying theta, in degree

  • dtheta (float) – step varying theta, in degree

set_orientation_axis(axis_1, axis_2, tol=1e-10)[source]#

rotate lattice_2 so that its axis_2 coincident with the axis_1 of lattice_1

interfacemaster.interface_generator.cross_plane(lattice, n, lim, orthogonal, tol, inclination_tol=0.7071067811865476)[source]#

get a primitive lattice vector cross a plane

Parameters:
  • lattice (numpy array) – lattice matrix

  • n (numpy array) – a normal vector

  • lim (int) – control how many vectors to be generated

  • tol (float) – tolerance judging orthogonal

Returns:

n_p – a primitve vector normal to the plane

Return type:

numpy array

interfacemaster.interface_generator.d_hkl(lattice, hkl)[source]#

return the lattice plane spacing of (hkl) of a lattice

interfacemaster.interface_generator.draw_cell(subfig, xs, ys, alpha=0.3, color='k', width=0.5)[source]#

draw the two-dimensional cell

interfacemaster.interface_generator.draw_slab(xs, ys, axes, num, plane_list, lattice_to_screen, elements_list, column, colors, all_elements, l_r, titlesize, legendsize)[source]#

draw the terminating planes in the plane list for a slab

Parameters:
  • xs (numpy array) – x,y arrays to draw the two-dimensional cell

  • ys (numpy array) – x,y arrays to draw the two-dimensional cell

  • axes (list) – list of figures

  • num (int) – number of terminations

  • plane_list (list) – list of planes of atoms

  • elements_list (list) – list of names of elments for these atoms

  • lattice_to_screen (numpy array) – rotated lattice with the interface lying in the screen

  • l_r (str) – left or right slab

  • titlesize (int) – fontsize of title

  • legendsize (int) – fontsize of legend

interfacemaster.interface_generator.draw_slab_dich(xs, ys, c_xs, c_ys, axes, num1, plane_list_1, lattice_to_screen_1, elements_list_1, colors, all_elements, num2, plane_list_2, lattice_to_screen_2, elements_list_2, titlesize)[source]#

draw the dichromatic patterns

Parameters:
  • xs (numpy array) – x,y arrays to draw the interface cell

  • ys (numpy array) – x,y arrays to draw the interface cell

  • c_xs (numpy array) – x,y arrays to draw the CNID cell

  • c_ys (numpy array) – x,y arrays to draw the CNID cell

  • axes (list) – list of figures

  • num1 (int) – number of terminations

  • num2 (int) – number of terminations

  • plane_list_1 (numpy array) – list of atom coordinates of the terminating planes

  • plane_list_2 (numpy array) – list of atom coordinates of the terminating planes

  • elements_list_1 (list) – list of corresponding element names

  • elements_list_2 (list) – list of corresponding element names

  • lattice_to_screen_1 (numpy array) – rotated lattices facing its interface orientation to the screen

  • lattice_to_screen_2 (numpy array) – rotated lattices facing its interface orientation to the screen

  • colors (str) – colors to classify different elements

  • all_elements (list) – list of all the elements for all the atoms

  • titlesize (int) – fontsize of the title

interfacemaster.interface_generator.excess_volume(lattice_1, lattice_bi, atoms_1, atoms_2, dx)[source]#

introduce vacuum between the interfaces

Parameters:
  • lattice_1 (numpy array) – lattice matrix of the first slab

  • lattice_bi (numpy array) – lattice matrix of the bicrystal

  • atoms_1 (numpy array) – atom fractional coordinates of slab 1, slab 2

  • atoms_2 (numpy array) – atom fractional coordinates of slab 1, slab 2

  • dx (float) – length of expands normal to the interface with the same units as lattice para

Returns:

lattice_bi, atoms_1, atoms2 – lattice matrix of bicrystal supercell, atom coordinates of left slab, and atom coordinates of right slab

Return type:

numpyarray

interfacemaster.interface_generator.get_R_to_screen(lattice)[source]#

get a rotation matrix to make the interface plane of the slab located in the screen

interfacemaster.interface_generator.get_ang_list(m1, n)[source]#

return a list of ang cos between one list of vecor and one vector

Parameters:
  • m1 (numpy array) – list of vectors

  • n (numpy array) – a vector

Returns:

c – list of cos

Return type:

numpy array

interfacemaster.interface_generator.get_array_bounds(U)[source]#

get the meshgrid formed by three sets of lower & upper bounds.

the bounds cooresponds to the 8 vertices of the cell consisting of the three column vectors of a matrix

Parameters:

U (numpy array) – integer matrix

Returns:

indice – meshgrid made by the lower and upper bounds of the indices

Return type:

numpy array

interfacemaster.interface_generator.get_disorientation(L1, L2, v1, hkl1, v2, hkl2)[source]#

produce a rotation matrix so that the hkl1 plane overlap with the hkl2 plane and the v1 colinear with v2

Parameters:
  • L1 (numpy array) – lattice basis sets

  • L2 (numpy array) – lattice basis sets

  • v1 (numpy array) – vectors

  • v2 (numpy array) – vectors

  • hkl1 (numpy array) – Miller indices

  • hkl2 (numpy array) – Miller indices

Returns:

rot_mat – a rotation matrix

Return type:

numpy array

interfacemaster.interface_generator.get_height(lattice)[source]#

get the distance of the two surfaces (crossing the first vector) of a cell

interfacemaster.interface_generator.get_nearest_pair(lattice, atoms, indices)[source]#

a function return the indices of two nearest atoms in a periodic block inspired from oekosheri/GB_code

Parameters:
  • lattice (numpy array) – lattice matrix

  • atoms (numpy array) – fractional coordinates

  • indices (numpy array) – indices of atoms

interfacemaster.interface_generator.get_plane_vectors(lattice, n)[source]#

a function get the two vectors normal to a vector

Parameters:
  • lattice (numpy array) – lattice matrix with column lattice vecotrs

  • n (numpy array) – a vector

Returns:

B, indices – B two plane vectors

Return type:

numpy array

interfacemaster.interface_generator.get_sites_elements(structure)[source]#

get the coordinates of atoms and the elements

Parameters:

structure (pymatgen structure class) – input structure

Returns:

  • atoms (numpy array) – fractional coordinates of atoms in the primitive cell

  • elements (numpy aray) – list of element name of the atoms

  • return

interfacemaster.interface_generator.get_surface_slab(structure, hkl, replica=None, inclination_tol=0.7071067811865476, termi_shift=0, vacuum_height=0, plane_normal=False, normal_perp=False, normal_tol=0.001, lim=20, filename='POSCAR', filetype='VASP')[source]#

get a superlattice of a slab containing a desired surface as crystal plane hkl

Parameters:
  • structure (Structure object) – Structure object of the unit cell structure

  • hkl (numpy array) – miller indices of the surface plane

  • replica (list) – expansion of the primitive slab cell, default [1, 1, 1]

  • inclination_tol (float) – required minimum cos value of the angle between the basic crossing vector and the surface, default sqrt(2)/2

  • termi_shift (float) – shift the termination of the surface, default 0

  • vacuum_height (float) – height of vaccum

  • plane_normal (bool) – whether requiring the two vectors in the surface plane to be perpendicular, default False

  • normal_perp (bool) – whether requiring the crossing vector to be perpendicular to the plane, default False

  • normal_tol – tolerance to judge whether perpendicular, default 1e-3

  • lim (int) – control the number of generated vectors to search for the crossing vectors and perpendicular vectors, default 20

  • filename (str) – name of the generated atom files

  • filetype (str, "VASP" or "LAMMPS") – type of the generated atom files (for VASP or LAMMPS)

Returns:

slab_structure – slab structure

Return type:

structure object

interfacemaster.interface_generator.get_unit_mtx(lattice)[source]#

return a unit lattice so that the length of every column vectors is 1

Parameters:

lattice (numpy array) – lattice basis set

Returns:

lattice_return – basis set of a unit lattice

Return type:

numpy array

interfacemaster.interface_generator.reciprocal_lattice(B)[source]#

return the reciprocal lattice of B

interfacemaster.interface_generator.rot(a, theta)[source]#

produce a rotation matrix

Parameters:
  • a (numpy array) – rotation axis

  • Theta (float) – rotation angle

Returns:

rot_mat – a rotation matrix

Return type:

numpy array

interfacemaster.interface_generator.shift_none_copy(lattice, dp, atoms)[source]#

changing terminate involves without requiring to cut or extend the cell for identical interfaces

Parameters:
  • lattice (numpy array) – a matrix with column lattice vectors

  • dp (foat) – height of termination shift

  • atoms (numpy array) – fractional coordinates of the atoms

Returns:

atoms – shifted atom coordinates

Return type:

numpy array

interfacemaster.interface_generator.shift_termi_left(lattice, dp, atoms, elements)[source]#

changing terminate involves requiring to cut or extend the cell for identical interfaces

Parameters:
  • lattice (numpy array) – a matrix with column lattice vectors

  • dp (float) – height of termination shift

  • atoms (numpy array) – fractional coordinates of the atoms

  • elements (list) – list of element names

Returns:

atoms, elements – shifted atom coordinates, corresponding list of elements

Return type:

tuple of numpy array and list

interfacemaster.interface_generator.shift_termi_right(lattice, dp, atoms, elements)[source]#

changing terminate involves requiring to cut or extend the cell for identical interfaces

Parameters:
  • lattice (numpy array) – a matrix with column lattice vectors

  • dp (float) – height of termination shift

  • atoms (numpy array) – fractional coordinates of the atoms

  • elements (list) – list of element names

Returns:

atoms, elements – shifted atom coordinates, corresponding list of elements

Return type:

numpy array

interfacemaster.interface_generator.super_cell(U, lattice, Atoms, elements)[source]#

make supercell

Parameters:
  • U (numpy array) – coefficients of the LC of three vectors from the basis

  • atoms (numpy array) – fractional coordinates of atoms

  • elements (numpy array) – list of the element names of these atoms

Returns:

  • lattice (numpy array) – lattice matrix

  • atoms (numpy array) – atom coordinates

  • elements (numpy array) – list of element name of the atoms

interfacemaster.interface_generator.surface_vacuum(lattice_1, lattice_bi, atoms_bi, vx)[source]#

introduce vacuum at one of the tails of the bicrystal cell

Parameters:
  • lattice_1 (numpy array) – lattice matrix of the first slab

  • lattice_bi (numpy array) – lattice matrix of the bicrystal

  • atoms_bi (numpy array) – atom fractional coordinates of the bicrystal

  • vx (float) – length of the vacuum bulk with units as lattice para

interfacemaster.interface_generator.terminates_scanner_left(slab, atoms, elements, d, round_n=5)[source]#

find all different atomic planes within 1 lattice plane displacing (in the interface), for the left slab

Parameters:
  • slab (numpy array) – basic vectors of the slab

  • atoms (numpy array) – fractional coordinates of atoms

  • elements (list) – list of name of the atoms

  • d (float) – 1 lattice plane displacing

  • round_n (int) – num of bits to round the fraction coordinates to judge identical planes

Returns:

  • plane_list (list) – list of planes of atom fraction coordinates

  • element_list (list) – list of elements in each plane

  • indices_list (list) – list of indices of atoms in each plane

  • dp_list (list) – list of dp parameters as input to select corresponding termination

interfacemaster.interface_generator.terminates_scanner_right(slab, atoms, elements, d, round_n=5)[source]#

find all different atomic planes within 1 lattice plane displacing (in the interface), for the right slab

Parameters:
  • slab (numpy array) – basic vectors of the slab

  • atoms (numpy array) – fractional coordinates of atoms

  • elements (list) – list of name of the atoms

  • d (float) – 1 lattice plane displacing

  • round_n (int) – num of bits to round the fraction coordinates to judge identical planes

Returns:

  • plane_list (list) – list of planes of atom fraction coordinates

  • element_list (list) – list of elements in each plane

  • indices_list (list) – list of indices of atoms in each plane

interfacemaster.interface_generator.terminates_scanner_slab_structure(structure, hkl)[source]#
interfacemaster.interface_generator.three_dot(M1, M2, M3)[source]#

compute the three continuous dot product

Parameters:
  • M1 (numpy array) – matrices

  • M2 (numpy array) – matrices

  • M3 (numpy array) – matrices

Returns:

P – dot product

Return type:

numpy array

interfacemaster.interface_generator.unit_cell_axis(axis)[source]#

get an unit orthogonal cell with the x-axis collinear with certain axis

interfacemaster.interface_generator.unit_v(vector)[source]#

get the unit vector of a vector

interfacemaster.interface_generator.write_LAMMPS(lattice, atoms, elements, filename='lmp_atoms_file', orthogonal=False)[source]#

write LAMMPS input atom file file of a supercell

Parameters:
  • lattice (numpy array) – lattice matrix

  • atoms (numpy array) – fractional atoms coordinates

  • elements (numpy array) – list of element name of the atoms

  • filename (str) – filename of LAMMPS input atom file to write, default “lmp_atoms_file”

  • orthogonal (bool) – whether write orthogonal cell, default False

interfacemaster.interface_generator.write_POSCAR(lattice, atoms, elements, filename='POSCAR')[source]#

write Poscar file of a supercell

Parameters:
  • lattice (numpy array) – lattice matrix

  • atoms (numpy array) – fractional atoms coordinates

  • elements (numpy array) – list of element name of the atoms

  • filename (str) – filename of POSCAR filename to wirte, default “POSCAR”

interfacemaster.interface_generator.write_trans_file(v1, v2, n1, n2)[source]#

write a file including translation information for LAMMPS

Parameters:
  • v1 (numpy array) – CNID vectors

  • v2 (numpy array) – CNID vectors

  • n1 (int) – num of grids for v1 & v2

  • n2 (int) – num of grids for v1 & v2

interfacemaster.symmetric_tilt#

Get the [001], [111], [110] symmetric tilt GBs; Get the csl for bilayer twisted graphenes

interfacemaster.symmetric_tilt.compute_sigma(axis, theta, filename, maxsigma=10000, verbose=True)[source]#

compute sigma values for a given disorientation

Parameters:
  • axis (numpy array) – rotation axis

  • theta (float) – rotation angle

  • maxsigma (int) – maximum sigma value for searching

  • verbose (bool, optional) – If True, the function will print additional information. Default is True.

Returns:

sigma

Return type:

int

interfacemaster.symmetric_tilt.generate_arrays_x_y(x_min, y_min, lim)[source]#

generate x, y meshgrids

Parameters:
  • x_min (int) – minimum x,y

  • y_min (int) – minimum x,y

  • lim (int) – maximum x,y

Returns:

meshgrid – x,y meshgrid

Return type:

numpy array

interfacemaster.symmetric_tilt.get_Ps_sigmas_thetas(lim, axis, maxsigma=100000)[source]#

get the geometric information of all the symmetric tilt GBs for a rotation axis within a searching limitation arguments: lim — control the number of generated referring points axis — rotation axis maxsigma — maximum sigma return: list of referring points, sigma list, angle list, original list of the ‘in-plane’ sigmas

interfacemaster.symmetric_tilt.get_csl_twisted_graphenes(lim, filename, maxsigma=100, verbose=True)[source]#

get the geometric information of all the CS twisted graphene within a searching limitation

Parameters:
  • lim (int) – control the number of generated referring points

  • maxsigma (int) – maximum sigma

Returns:

  • list1 (numpy array) – sigma list

  • list2 (numpy array) – angle list

  • list3 (numpy array) – CNID areas

  • list4 (numpy array) – num of atoms in supercell

  • verbose (bool, optional) – If True, the function will print additional information. Default is True.

interfacemaster.symmetric_tilt.get_hkl(P, axis, tol=0.001)[source]#

given a referring point and a rotation matrix, get the miller indices for this symmetric tilt GB

Parameters:
  • P (numpy array) – referring point

  • axis (numpy array) – rotation axis

Returns:

hkl – miller indices

Return type:

numpy array

interfacemaster.symmetric_tilt.sample_STGB(axis, lim, maxsigma, max_index)[source]#

sampling the symmetric tilt GBs for a given rotation axis

Parameters:
  • axis (numpy array) – rotation axis

  • lim (int) – control the number of generated referring points

  • maxsigma (int) – maximum sigma value

  • max_index (int) – maximum value of miller indices

Returns:

  • list1 (numpy array) – angle list

  • list2 (numpy array) – sigma list

  • list3 (numpy array) – hkl list

interfacemaster.brute_force#

brute_force.py

class interfacemaster.brute_force.GB_runner(my_interface)[source]#

Bases: object

a class doing brute-foce searching for elemental GBs

divide_region()[source]#

divide simulation regions

get_RBT_list(grid_size)[source]#

generate RBT operation list

Parameters:

grid_size (float) – size of grid sampling RBT in angstrom

get_terminations(changing_termination=False)[source]#

generate termination selections

Parameters:

changing_termination (bool) – whether to sample different terminating planes

main_run(core_num)[source]#

main loop doing RBT & merging

Parameters:

core_num (int) – number of CPU cores for simulation

main_run_terminations(core_num)[source]#

main loop by hata’s method

Parameters:

core_num (int) – number of cores for computation

interfacemaster.brute_force.find_pairs_with_closest_distances(atoms_here, bicrystal_lattice)[source]#

Find the pairs of atoms with closest distance

Parameters:
  • atoms_here (numpy array) – cartesian coordinates

  • bicrystal_lattice (numpy array) – lattice matrix of the bicrystal cell

Returns:

  • array_id_del (numpy array) – ids of atoms to delete

  • array_id_dsp (numpy array) – ids of atoms to displace

  • dsps (numpy array) – displacements to be done

  • distances_round[1] (float) – cutoff this turn

  • distances_round[2] (float) – cloest distance after merge

interfacemaster.brute_force.get_lowest()[source]#

get the lowest GB energy

interfacemaster.brute_force.merge_operation(count_start, GB_atoms, bicrystal_lattice, clst_atmc_dstc, RBT, bulk_atoms, core_num, dp1, dp2)[source]#

merge loop

interfacemaster.brute_force.merge_operation_no_RBT(count_start, GB_atoms, bicrystal_lattice, clst_atmc_dstc, bulk_atoms, core_num, dp1, dp2)[source]#

merge loop by hata’s method

interfacemaster.brute_force.read_energy()[source]#

LAMMPS run command

interfacemaster.brute_force.run_LAMMPS(core_num, count)[source]#

LAMMPS run command

interfacemaster.brute_force.screen_out_non_repeated_pairs(ids_1, ids_2)[source]#

input two pairs of ids with self-paring output the indices involving non-repeating pairs

Parameters:
  • ids_1 (numpy arrays) – two pairs of ids

  • ids_2 (numpy arrays) – two pairs of ids

Returns:

screened_ids – ids in ids_1 without repeating pairs

Return type:

numpy array

interfacemaster.brute_force.write_data_here(count, energy, dy, dz, dp1, dp2, delete_cutoff)[source]#

write data for each simulation

Parameters:
  • count (int) – index of simulation

  • energy (float) – GB energy

  • dy (float) – RBT

  • dz (float) – RBT

  • dp1 (float) – terminating shift

  • dp2 (float) – terminating shift

  • delete_cutoff (float) – cutoff to merge atoms