Slabmaker¶
This module contains builder methods for generating grain boundary structures. Currently supports builders for generating symmetric tilt and twist grain boundaries in BCC crystals.
Standard quaternion algebra routines are provided by the transformations.py
library
developed by Christoph Gohlke.
-
class
imeall.slabmaker.gengb_from_quat.
QuaternionGB
[source]¶ imeall.slabmaker.gengb_from_quat.QuaternionGB
contains routines to generate lists of valid grain boundary orientations (i.e. the boundary planes of grains for different orientation axes and misorientation angles) using quaternion algebra. Currently supports routines for symmetric tilt grain boundaries and twist boundaries.-
gen_sym_tilt
(orientation_axis=[0.0, 0.0, 1.0])[source]¶ Generate symmetric tilt grain boundary list. Generate a list of the misorientation angle for the given orientation axis and the relevant boundary plane to generate the symmetric tilt grain boundary bicrystal with
imeall.slabmaker.build_tilt_sym_gb()
.Parameters: orientation_axis (list[int]) – 3d array of floats defining orientation axis. Returns: Dictionary of ‘x’ axis co-ordinates keyed by angle rounded to two floating point numbers. Return type: dict
-
gen_sym_twist
(orientation_axis=[0, 0, 1])[source]¶ Generate symmetric twist grain boundary list. Generate a list of the misorientation angle for the given orientation axis and the relevant boundary plane to generate the symmetric twist grain boundary bicrystal with
imeall.slabmaker.build_twist_sym_gb()
.Parameters: orientation_axis (list[ints]) – orientation axis in ints. Returns: list of [angle, ‘x’-axis-co-ordinates] Return type: list
-
or_axis_to_angle
(m, n, orientation_axis)[source]¶ To generate an approximate spanning set of angles iterate over integers m and n and determine the angle from a combination of these integers according to the chosen orientation axis.
Parameters: Returns: misorientation angle from integers and orientation axis.
Return type:
-
or_axis_to_quat
(orientation_axis)[source]¶ Given an orientation axis converts the vector to quaternion representation defining the orientation axis and a second quaternion giving the 0 of the angle. :
Parameters: orientation_axis (list) – orientation axis. Returns: planequat_1, planequat_2 Return type: quaternion, quaternion
-
-
imeall.slabmaker.slabmaker.
bcc_csl_nn0
(m, n, grain)[source]¶ Rotate grain atom positions according to a rotation matrix generated from a quaternion.
Example for integer m, and orientation axis defined by vector [n n 0], the angle of rotation is given by: arccos((m^{2} - 2n^{2})(m^{2} + 2n^{2})).
Parameters:
-
imeall.slabmaker.slabmaker.
build_tilt_sym_gb
(gbid='', bp=[3, 3, 2], v=[1, -1, 0], min_spacing=12.0, c_space=None, target_dir=None, rbt=None, chem_symbol='Fe', alat=2.83)[source]¶ Generate symmetric tilt grain boundary with appropriate configurations boundary plane (bp) oriented along z axis and orthogonal directions in the the other two planes given the orientation axis (v) and an orthogonal vector bpxv so we have a proper cube. If rbt is not None then rigid body translations are present, this is passed as a list of two numbers abs(rbt[0]) < 1. The cell is then displaced as a function of these numbers. This routine should work for any BCC crystal with the appropriate modification to alat and chem_symbol.
Parameters: - gbid (str) – id label for the grain boundary.
- bp (list, int) – boundary plane normal.
- v (list,int) – orientation axis.
- min_spacing (float) – minimum spacing between grain boundaries.
- c_space (float,optional) – inter-planar spacing if not passed this is calculated automatically.
- target_dir (str) – string specifying grain boundary target directory.
- rbt (list,float) – rigid body translations.
- alat (float) – Lattice constant of the unit cell.
Returns: if target_dir==None returns [z_planes, len(dups), n_grain_unit, grain_c] i.e the location of the interfacial planes, the number of duplicate atoms when the grain is built, the number of atoms in the canonical unit cell, and the bicrystal
ase.Atoms
object. else returnsase.Atoms
.
-
imeall.slabmaker.slabmaker.
build_twist_sym_gb
(gbid='', bp=[0, 0, 1], v=[3, 5, 0], min_spacing=14.0, c_space=None, target_dir=None, rbt=None, chem_symbol='Fe', alat=2.83)[source]¶ Builder for twist boundary structures. For the twist boundaries the orientation axis (or_axis) also defines the boundary plane normal.
Parameters: - gbid (str) –
imeall.gb_models.SubGrainBoundary
id. - bp (list[int]) – Boundary plane.
- v (list[int]) – Defines the ‘x’ axis of the orthorhombic bicrystal.
- min_spacing (float) – Minimum separation between interfaces.
- c_space (float) – Interplanar spacing of boundary planes.
- target_dir (str) – Path str of where to deposit grain boundary.
- rbt (list[float]) – Rigid body translations.
- chem_symbol (str) – chemical symbol of bcc element.
- alat (float) – lattice constanct BCC crystal.
Returns: if target_dir==None returns list [z_planes, sigma_csl, n_grain_unit, grain_c] else returns
ase.Atoms
.- gbid (str) –
-
imeall.slabmaker.slabmaker.
calc_twist_sigma_csl
(bp, v, chem_symbol='Fe', alat=2.83, n=3)[source]¶ Calculate sigma number for a twist symmetric structure.
Parameters: - bp (:numpy:`array`) – Boundary plane.
- v (:numpy:`array`) – x-axis of unit cell.
- chem_symbol (str) – Chemical element abbreviation.
- n (int) – size of unit cell.
Returns: number of atoms, duplicate pairs, Sigma
Return type: list
-
imeall.slabmaker.slabmaker.
csl_lattice_vecs
(m, n)[source]¶ Returns lattice vectors defining a coincident site lattice from quaternion following Zeiner(05).
-
imeall.slabmaker.slabmaker.
csl_tilt_factory
(orientation_axis, boundary_plane, m, n, gbid, grain_a, grain_b, theta=0, target_dir='./', gb_type='tilt')[source]¶ Builds a plot of the coincident site lattice oriented along a suitable plane for viewing purposes. In the gnuplot plane we think of x,y,z orientation. so that the orientation axis = N = z is oriented ‘into’ the page. Vector specifying 0 Angle (for the case of tilt) = v = x y = Nxv.
rotate_plane_z()
takes a grain and rotates it so that the orientation axis of the grain boundary with respect to grain a is orthogonal to the x-y plane. The quaternion to accomplish this rotation is stored in plane_quaternion_z.Parameters: - orientation_axis (list) – orientation axis
- boundary_plane (list) – boundary plane
- m (int) – slope of boundary plane
- n (int) – integer of orientation axis.
- gbid (str) – grain boundary identifying label.
- grain_a (
ase.Atoms
) – blue crystal. - grain_b (
ase.Atoms
) – red crystal. - theta (float) – angle.
- target_dir (str) – target directory to write plotting scripts.
- gb_type (str) – grain boundary type.
-
imeall.slabmaker.slabmaker.
csl_twist_factory
(bp, v, gbid, target_dir, crystal_type='BCC', gb_type='twist')[source]¶ Generate coincident site lattice points for arbitraty twist boundary.
Parameters:
-
imeall.slabmaker.slabmaker.
find_densest_plane
(grain_dict)[source]¶ Return the indices for a set of adjacent planes with the largest number of points. Used in
gen_csl()
.Parameters: grain_dict (dict) – Grain dictionary. Returns: key1 and key2 floating point values representing the z-coordinates of the two adjacent planes along z-axis with the most Atoms on them. Return type: float, float
-
imeall.slabmaker.slabmaker.
gen_canonical_grain_dir
(angle, orientation_axis, boundary_plane, material='alphaFe', target_dir='./', gb_type='tilt')[source]¶ Generate a canonical grain boundary directory. :param angle: misorientation angle in radians. :type angle: float :param orienation_axis: orientation axis. :type orienation_axis: list :param boundary_plane: boundary plane :type boundary_plane: list :param material: material :type material: str, optional :param target_dir: target directory to deposit canonical grain structure. :type target_dir: str, optional :param gb_type: grain boundary type options are ‘tilt’, ‘twist’. :type gb_type: str
-
imeall.slabmaker.slabmaker.
gen_csl
(theta, orientation_axis, boundary_plane, target_dir='./', gbid='0000000000', alat=2.83, chem_symbol='Fe', gb_type='tilt')[source]¶ Generate the coincident site lattice representations.
Parameters: - theta (float) – misorientation angle in radians.
- orientation_axis (list) – orientation axis.
- boundary_plane (list[int]) – boundary plane normal vector.
- target_dir (str) – location to store coincident site lattice files.
- gbid (str) – Grain boundary id.
- alat (float) – Magnitude lattice vector.
- chem_symbol (str) – Chemical symbol.
- gb_type (str) – boundary type [tilt, twist] (Default: ‘tilt’).
-
imeall.slabmaker.slabmaker.
gnu_plot_gb
(boundary_plane, m, invm, gbid, mb=0.0, invmb=0.0, target_dir='./')[source]¶ Generate coincident site lattice plots on the fly in svg format for web based visualization.
Parameters: - boundary_plane (list) – vector defining interfacial plane.
- m (float) – slope of line separating red and blue crystal.
- invm (float) – slope of line perpendicular to m.
- gbid (str) – grain boundary id.
- mb (float,optional) – deprecated.
- invmb (float,optional) – deprecated.
- target_dir (str, optional) – location to write plotting script.
-
imeall.slabmaker.slabmaker.
print_points
(atoms, f, top_grain=True, gb_type='tilt')[source]¶ Print positions of atoms to file.
Parameters:
-
imeall.slabmaker.slabmaker.
rotate_grain
(grain, theta=0.0, x=[0.0, 0.0, 0.0], q=[])[source]¶ Rotate the grain according to quaternion=[Theta, u,v,w]=[w,x,y,z]. Standard routine is passed angle and vector this generates the quaternion to do the rotation however if a quaternion, q!=None, is passed the plane is rotated using q as: qvq^{-1}.
Parameters: Returns: Rotated grain.
-
imeall.slabmaker.slabmaker.
rotate_plane_y
(grain, miller)[source]¶ Rotates atoms in grain so that the miller plane in grain is parallel to the y axis. Returns the quaternion characterizing the rotation which achieves this.
Parameters: - grain (
Atoms
) – grain atoms object. - miller (list,int) – Vector specifying grain boundary interfacial plane.
Returns: quaternion which acheives this rotation.
- grain (
-
imeall.slabmaker.slabmaker.
rotate_plane_z
(grain, miller)[source]¶ Rotates atoms in grain so that planes parallel to the plane defined by the vector miller are parallel to the xy plotting plane.
Parameters: - grain (
Atoms
) – grain structure to be rotated. - miller (numpy array) – 3x1 vector defining plane.
Returns: rotation_quaternion to achieve this rotation.
- grain (
-
imeall.slabmaker.slabmaker.
rotate_vec
(q, vec)[source]¶ Rotate 3 vector with quaternion: qvq^{-1}.
Parameters: - q (quaternion) – rotation quaternion.
- vec (3x1 list) – vector to rotate.
Returns: Rotated 3 vector.