You can adapt this file completely to your liking, but it should at least contain the root toctree directive.

Welcome to GBpy’s documentation!

Tutorials

For a complete list of tutorials please click here.

Lattice Class

class lattice.Lattice(*args)[source]

This class contains all the crystallographic information required for each atom type. Currently there are only two pre-configured atoms available in class i.e. ‘Al’ and ‘Mg’. Up on need user can create a new instance of this class with the same attributes. The attributes of this class are: ...

Notes

Examples of elem_type elem_type = ‘Mg’;

elem_type = ‘Al’;

elem_type = ‘Cu’;

elem_type = ‘Ni’;

elem_type = ‘cF_Id’;

elem_type = ‘cI_Id’;

elem_type = ‘cP_Id’;

elem_type = ‘hP_Id’;

Attributes

elem_type: string Element of Interest
pearson: string Pearson symbol for the lattice
lat_params: dictionary Lattice parameters (‘a’, ‘b’, ‘c’, ‘alpha’, ‘beta’, ‘gamma’)
l_g_go: numpy array Primitve basis of the lattice
basis_atoms: Location of the basis atoms in the primitive lattice
cryst_ptgrp: string Crystallographic point group of the lattice
burgers_mag: float The smallest burgers vector in the lattice
eam_file: eam_file name for atomistic simulations

Methods

str  
Method for printing the lattice class  

Geometry Tools

class quaternion.quaternion(*args, **kwargs)[source]

The quaternion class is defined to define the quaternion parameterization of rotation and their operations ...

Attributes

quaternion: numpy array 5 x n dimensions

Methods

all([axis, out]) Returns True if all elements evaluate to True.
any([axis, out]) Returns True if any of the elements of a evaluate to True.
argmax([axis, out]) Return indices of the maximum values along the given axis.
argmin([axis, out]) Return indices of the minimum values along the given axis of a.
argpartition(kth[, axis, kind, order]) Returns the indices that would partition this array.
argsort([axis, kind, order]) Returns the indices that would sort this array.
astype(dtype[, order, casting, subok, copy]) Copy of the array, cast to a specified type.
byteswap(inplace) Swap the bytes of the array elements
choose(choices[, out, mode]) Use an index array to construct a new array from a set of choices.
clip(a_min, a_max[, out]) Return an array whose values are limited to [a_min, a_max].
compress(condition[, axis, out]) Return selected slices of this array along given axis.
conj() Complex-conjugate all elements.
conjugate() Return the complex conjugate, element-wise.
copy([order]) Return a copy of the array.
cumprod([axis, dtype, out]) Return the cumulative product of the elements along the given axis.
cumsum([axis, dtype, out]) Return the cumulative sum of the elements along the given axis.
diagonal([offset, axis1, axis2]) Return specified diagonals.
dot(b[, out]) Dot product of two arrays.
dump(file) Dump a pickle of the array to the specified file.
dumps() Returns the pickle of the array as a string.
fill(value) Fill the array with a scalar value.
flatten([order]) Return a copy of the array collapsed into one dimension.
getfield(dtype[, offset]) Returns a field of the given array as a certain type.
item(*args) Copy an element of an array to a standard Python scalar and return it.
itemset(*args) Insert scalar into an array (scalar is cast to array’s dtype, if possible)
max([axis, out]) Return the maximum along a given axis.
mean([axis, dtype, out]) Returns the average of the array elements along given axis.
min([axis, out]) Return the minimum along a given axis.
newbyteorder([new_order]) Return the array with the same data viewed with a different byte order.
nonzero() Return the indices of the elements that are non-zero.
partition(kth[, axis, kind, order]) Rearranges the elements in the array in such a way that value of the element in kth position is in the position it would be in a sorted array.
prod([axis, dtype, out]) Return the product of the array elements over the given axis
ptp([axis, out]) Peak to peak (maximum - minimum) value along a given axis.
put(indices, values[, mode]) Set a.flat[n] = values[n] for all n in indices.
ravel([order]) Return a flattened array.
repeat(repeats[, axis]) Repeat elements of an array.
reshape(shape[, order]) Returns an array containing the same data with a new shape.
resize(new_shape[, refcheck]) Change shape and size of array in-place.
round([decimals, out]) Return a with each element rounded to the given number of decimals.
searchsorted(v[, side, sorter]) Find indices where elements of v should be inserted in a to maintain order.
setfield(val, dtype[, offset]) Put a value into a specified place in a field defined by a data-type.
setflags([write, align, uic]) Set array flags WRITEABLE, ALIGNED, and UPDATEIFCOPY, respectively.
sort([axis, kind, order]) Sort an array, in-place.
squeeze([axis]) Remove single-dimensional entries from the shape of a.
std([axis, dtype, out, ddof]) Returns the standard deviation of the array elements along given axis.
sum([axis, dtype, out]) Return the sum of the array elements over the given axis.
swapaxes(axis1, axis2) Return a view of the array with axis1 and axis2 interchanged.
take(indices[, axis, out, mode]) Return an array formed from the elements of a at the given indices.
tofile(fid[, sep, format]) Write array to a file as text or binary (default).
tolist() Return the array as a (possibly nested) list.
tostring([order]) Construct a Python string containing the raw data bytes in the array.
trace([offset, axis1, axis2, dtype, out]) Return the sum along diagonals of the array.
transpose(*axes) Returns a view of the array with axes transposed.
var([axis, dtype, out, ddof]) Returns the variance of the array elements, along given axis.
view([dtype, type]) New view of array with the same data.

Integer Manipulations

integer_manipulations.gcd_array(input, order='all')[source]

The function computes the GCD of an array of numbers.

Parameters:

input : numpy array or list of intgers

Input n-D array of integers (most suitable for 1D and 2D arrays)

order : {‘rows’, ‘columns’, ‘col’, ‘all’}, optional

Returns:

Agcd: numpy array

An array of greatest common divisors of the input

See also

gcd
from fractions module for computing gcd of two integers

Notes

  • If order = all, the input array is flattened and the GCD is calculated
  • If order = rows, GCD of elements in each row is calculated
  • If order = columns or cols, GCD of elements in each column is calculated
integer_manipulations.lcm_vec(a, b)[source]

The function computes the LCM of two 1D array of integers of length and retruns a 1D array of lcm values

Parameters:

a, b : numpy array

Input 1D arrays of integers

Returns:

lcm_vector: numpy array

Output 1D array of integers

See also

lcm_arry

integer_manipulations.lcm_array(input, order='all')[source]

The function computes the LCM of an array of numbers.

Parameters:

input : numpy array or list of intgers

Input n-D array of integers (most suitable for 1D and 2D arrays)

order : {‘rows’, ‘columns’, ‘col’, ‘all’}, optional

Returns:

Alcm: numpy array

An array of least common multiples of the input

See also

gcd_array

Notes

  • If order = all, the input array is flattened and the LCM is calculated
  • If order = rows, LCM of elements in each row is calculated
  • If order = columns or cols, LCM of elements in each column is calculated
integer_manipulations.int_check(input, precis=6)[source]

Checks whether the input variable (arrays) is an interger or not. A precision value is specified and the integer check is performed up to that decimal point.

Parameters:

input : numpy array or list

Input n-D array of floats.

precis : Integer

Default = 6. A value that specifies the precision to which the number is an integer. precis = 6 implies a precision of \(10^{-6}\).

Returns:

cond: Boolean

True if the element is an integer to a certain precision, False otherwise

integer_manipulations.rat(input, tol=1e-06)[source]

The function returns a rational (p/q) approximation of a given floating point array to a given precision

Parameters:

input : numpy array or list of real numbers

tol : floating point tolerance value

Default = 1e-06

Returns:

N, D: Integer numpy arrays

N and D contain the numerators (p) and denominators (q) of the rational approximations

integer_manipulations.int_finder(input_v, tol=1e-06, order='all', tol1=1e-06)[source]

The function computes the scaling factor required to multiply the given input array to obtain an integer array. The integer array is returned.

Parameters:

input1 : numpy array or list of real numbers

tol : floating point tolerance value

Default = 1e-06

order : {‘rows’, ‘columns’, ‘col’, ‘all’}

Defualt = ‘all’

tol1:

Returns:

output: numpy float array

An array of integers obtained by scaling input

See also

gcd_array

Notes

  • If order = all, the input array is flattened and then scaled
  • If order = rows, elements in each row are scaled
  • If order = columns or cols, elements in each column are scaled
integer_manipulations.int_mult(input, tol=1e-06)[source]

The function computes the scaling factor required to multiply the given input array to obtain an integer array. The integer array is returned.

Parameters:

input : numpy array or list of real numbers

tol : floating point tolerance value

Default = 1e-06

Returns:

N: numpy float array

An array of integers obtained by scaling input

Int_Mat: numpy float array

An array of integers obtained by scaling input

See also

int_finder

Notes

Change this function to accept rows and columns as input

CSL Utility Function

csl_utility_functions.csl_rotations(sigma, sig_type, lat_type)[source]

The function computes the CSL rotation matrices r_g1tog2_g1 corresponding to a give sigma and lattice

Parameters:

sigma : int

Sigma corresponding to the transformation matrix

sig_type: {‘common’, ‘specific’}

If the sigma generating function depends on the lattice type, then sig_type is ‘specific’, otherwise it is ‘common’

lat_type: Lattice class

Attributes of the underlying lattice

Returns:

sig_rots: dictionary

keys: ‘N’, ‘D’

sig_rots[‘N’], sig_rots[‘D’]: Numerator and Integer matrices

The transformation matrix is N/D in the g1 reference frame (i.e. r_g1tog2_g1)

Notes

The following steps are considered to obtain the sigma rotation:

  • compute_inp_params: computes tau and kmax that fixes the range of integer qudruples sampled

  • mesh_muvw: Creates the integer quadruples that depend on sigma, tau, kmax, crystallographic point group

  • eliminate_idrots: Eliminates Identity rotations

  • If specific rotations are desired:
    • mesh_muvw_fz: Restricts quadruples to fundamental zone
    • check_fsig_int: Filters out quadruples that do not meet the condition specified in this function
  • sigtype_muvw: Filters out quadruple combinations depending on the type of sigma rotation

  • eliminate_mults: Eliminates integer quadruples that are same except for a scaling factor

  • check_sigma: Returns integer quadruples that result in the sigma rotation

  • compute_tmat: Computes the transformation matrix from the integer quadruple

  • disorient_sigmarots: Converts all the transformations to the fundamental zone of the corresponding crystallogrphic point group

  • check_sigma_rots: Checks that the transformation matrix is a sigma rotation and returns them as numerator and denominator matrices

csl_utility_functions.proper_ptgrp(cryst_ptgrp)[source]

Returns the proper point group corresponding to a crystallographic point group

Parameters:

cryst_ptgrp: string

Crystallogrphic point group in Schoenflies notation

Returns:

proper_ptgrp: string

Proper point group in Schoenflies notation

csl_utility_functions.largest_odd_factor(var_arr)[source]

Function that computes the larges odd factors of an array of integers

Parameters:

var_arr: numpy array

Array of integers whose largest odd factors needs to be computed

Returns:

odd_d: numpy array

Array of largest odd factors of each integer in var_arr

csl_utility_functions.compute_inp_params(lattice, sig_type)[source]

tau and kmax necessary for possible integer quadruple combinations are computed

Parameters:

lattice: Lattice class

Attributes of the underlying lattice

sig_type: {‘common’, ‘specific’}

Returns:

tau: float

tau is a rational number \(= \frac{\nu}{\mu}\)

kmax: float

kmax is an integer that depends on \(\mu \ , \nu\)

csl_utility_functions.mesh_muvw(cryst_ptgrp, sigma, sig_type, *args)[source]

Compute max allowed values of [m,U,V,W] and generates an array of integer quadruples

Parameters:

cryst_ptgrp: string

Proper point group in Schoenflies notation

sigma: integer

Sigma number

sig_type: {‘common’, ‘specific’}

args[0]: dictionary

keys: ‘nu’, ‘mu’, ‘kmax’

Returns:

Integer quadruple numpy array

csl_utility_functions.mesh_muvw_fz(quad_int, cryst_ptgrp, sig_type, *args)[source]

For given integer quadruples, the set belonging to the corresponding fundamental zone are separated out and retruned.

Parameters:

quad_int: numpy array

Integer quadruples

cryst_ptgrp: string

Proper point group in Schoenflies notation

sig_type: {‘common’, ‘specific’}

args[0]: dictionary

keys: ‘nu’, ‘mu’, ‘kmax’

Returns:

Integer quadruple numpy array belonging to the fundamental zone

of the corresponding crystallographic point group

csl_utility_functions.check_fsig_int(quad_int, cryst_ptgrp, sigma, *args)[source]

For specific sigma rotations, a function of m, U, V, W (fsig) is computed. The ratio of fsig and sigma should be a divisor of kmax. This condition is checked and those integer quadruples that satisfy this condition are returned

Parameters:

quad_int: numpy array

Integer quadruples

cryst_ptgrp: string

Proper point group in Schoenflies notation

sigma: float

sigma number

args[0]: dictionary

keys: ‘nu’, ‘mu’, ‘kmax’

Returns:

quad_int: numpy array

Integer quadruple array that satisfy the above mentioned condition

csl_utility_functions.eliminate_idrots(quad_int)[source]

Eliminate the roations that belong to the identity matrix and return the integer quadruples

csl_utility_functions.sigtype_muvw(quad_int, cryst_ptgrp, sig_type)[source]

The type of integer quadruples are different for common and specific sigma rotations. For example, for D4 point group, common rotations satisfy the condition u = 0 and v = 0 or m = 0 and w = 0. The specific rotations belong to the complimentary set. Depending on the sig_type (common, specific), the appropriate set of the integer quadruples are returned.

Parameters:

quad_int: numpy array

Integer quadruples.

cryst_ptgrp: string

Proper point group in Schoenflies notation.

sig_type: {‘common’, ‘specific’}

Returns:

quad_int: numpy array

Integer quadruple array that satisfy the above mentioned condition.

csl_utility_functions.eliminate_mults(quad_int)[source]

Divide all the integer quadruples by their corresponding least common multiples and return the unique set of integer quadruples

csl_utility_functions.check_sigma(quad_int, sigma, cryst_ptgrp, sig_type, *args)[source]

The integer quadruples that correspond to a sigma rotation satisfy certain conditions. These conditions are checked and all the quadruples that do not meet these requirements are filtered out. These conditions depend on the rotation type (common or specific) and the lattice type (crystallogrphic point group and mu, nu)

Parameters:

quad_int: numpy array

Integer quadruples.

cryst_ptgrp: string

Proper point group in Schoenflies notation.

sig_type: {‘common’, ‘specific’}

args[0]: dictionary

keys: ‘nu’, ‘mu’, ‘kmax’

Returns:

quad_int: numpy array

Integer quadruple array that satisfy the above mentioned condition.

See also

check_fsig_int

csl_utility_functions.gcd1d_arr(arr_tup)[source]

A tuple of one-D arrays are passed with equal size and the gcd of their rows is computed

Parameters:

arr_tup: tuple

one-D arrays of integers of equal size.

Returns:

GCD of rows of 1D arrays of integers

csl_utility_functions.compute_tmat(quad_int, tau, lat_type)[source]

The transformation matrix (r_g1tog2_g1) corresponding to the integer quadruple is computed. The matrix elements depend on m, U, V, W and the crystallographic point group and tau = (nu/mu)

Parameters:

quad_int: numpy array

Integer quadruples.

tau: float

\(\frac{\nu}{\mu}\)

lat_type: Lattice class

Attributes of the underlying lattice

Returns:

g : numpy array

dimension = 3, n x 3 x 3 transformation matrices

csl_utility_functions.disorient_sigmarots(r_g1tog2_g1, l_g_go, cryst_ptgrp)[source]

The disorientation corresponding to each rotation matrix is computed and the unique set is returned

Parameters:

r_g1tog2_g1: numpy array (n x 3 x 3)

Transformation matrices in g1 reference frame

l_g_go: numpy array

The primitive basis vectors of the underlying lattice in the orthogonal reference frame.

cryst_ptgrp: string

Proper point group in Schoenflies notation

Returns:

rots_g1tog2_g1: numpy array (n x 3 x 3)

Transformation matrices in g1 reference frame in the fundamental zone

csl_utility_functions.check_sigma_rots(r_g1tog2_g1, sigma)[source]

The sigma transformation matrix has the property that sigma is the smallest integer such that sigma*T is an integer matrix. This condition is checked and the numerator and denominatr(sigma) matrices are returned

Parameters:

r_g1tog2_g1: numpy array (n x 3 x 3)

Transformation matrices in g1 reference frame

sigma: float

sigma number

Returns:

{‘N’: rots_n, ‘D’: rots_d}: dictionary

rots_n: numpy array

numerator matrices n x 3 x3

rots_d: numpy array

denominator matrices n x 3 x3

CSL/DSC Computation

find_csl_dsc.find_csl_dsc(L_G1_GO1, R_G1ToG2_G1)[source]

This function calls the csl_finder and dsc_finder and returns the CSL and DSC basis vectors in ‘g1’ reference frame.

Parameters:

L_G1_GO1: numpy array

The three basis vectors for the primitive unit cell (as columns) are given with respect to the GO1 reference frame.

R_G1ToG2_G1: 3X3 numpy array

The rotation matrix defining the transformation in ‘G1’ reference frame. The subscript ‘G1’ refers to the primitive unit cell of G lattice.

Returns

l_csl_g1, l_dsc_g1: numpy arrays

The basis vectors of csl and dsc lattices in the g1 reference frame

find_csl_dsc.sigma_calc(t_g1tog2_g1)[source]

Computes the sigma of the transformation matrix

  • if det(T) = det(T^{-1}) then sigma1 = sigma2 is returned
  • if det(T) neq det(T^{-1}) then max(sigma1, sigma2) is returned
find_csl_dsc.reciprocal_mat(l_g_go)[source]

The reciprocal matrix with reciprocal basis vectors is computed for the input matrix with primitve basis vectors

Parameters:

l_g_go: numpy array

The primitive basis vectors b1x, b1y, b1z

Returns:

rl_g_go: numpy array

The primitve reciprocal basis vectors

find_csl_dsc.csl_elem_div_thm_l1(T0, l_g1n_g1)[source]

The csl basis vectors are obtained from the diagonal matrix using the algorithm specified in doi:10.1107/S056773947601231X. There are two algorithms specified based on numerators or denominators of the T0 matrix. The numerators are used in this function.

Parameters:

T0: numpy array

The transformation matrix in G1n reference frame

l_g1n_g1: numpy array

The ‘new’ basis vectors of g1 lattice (g1n) in g1 reference frame

Returns:

l_csl_g1: numpy array

The CSL basis vectors in g1 reference frame

find_csl_dsc.csl_elem_div_thm_l2(t0, l_g2n_g2)[source]

The csl basis vectors are obtained from the diagonal matrix using the algorithm specified in doi:10.1107/S056773947601231X. There are two algorithms specified based on numerators or denominators of the T0 matrix. The denominators are used in this function.

Parameters:

T0: numpy array

The transformation matrix in G1n reference frame

l_g2n_g2: numpy array

The ‘new’ basis vectors of g2 lattice (g2n) in g2 reference frame

Returns:

l_csl_g2: numpy array

The CSL basis vectors in g2 reference frame

find_csl_dsc.csl_finder_smith(r_g1tog2_g1)[source]

This funciton extracts the CSL basis when transformation between the two lattices is given (r_g1tog2_g1). The algorithms used are based on the following article: doi:10.1107/S056773947601231X)

Parameters:

r_g1tog2_g1: numpy array

The 3x3 transformation matrix in g1 reference frame

Returns:

l_csl_g1: numpy array

3 x 3 matrix with the csl basis vectors as columns

Notes

The “Reduced” refer to the use of LLL algorithm to compute a basis that is as close to orthogonal as possible. (Refer to http://en.wikipedia.org/wiki/Lattice_reduction) for further detials on the concept of Lattice Reduction

find_csl_dsc.check_csl_finder_smith(r_g1tog2_g1, Sigma, L_G1_GO1, L_CSL_G1)[source]

This function checks the obtained CSL basis vectors are correct by using the following conditions: * The CSL basis vectors are integer combinations of basis vectors of lattice 1 * The CSL basis vectors are integer combinations of basis vectors of lattice 2 * The volume enclosed by the CSL is sigma times the volume of lattice 1

find_csl_dsc.dsc_finder(L_G2_G1, L_G1_GO1)[source]

The DSC lattice is computed for the bi-crystal, if the transformation matrix l_g2_g1 is given and the basis vectors of the underlying crystal l_g_go (in the orthogonal reference go frame) are known. The following relationship is used: The reciprocal of the coincidence site lattice of the reciprocal lattices is the DSC lattice

Parameters:

l_g2_g1: numpy array

transformation matrix (r_g1tog2_g1)

l_g1_go1: numpy array

basis vectors (as columns) of the underlying lattice expressed in the orthogonal ‘go’ reference frame

Returns:

l_dsc_g1: numpy array

The dsc lattice basis vectors (as columns) expressed in the g1 reference

Notes

The “Reduced” refer to the use of LLL algorithm to compute a basis that is as close to orthogonal as possible. (Refer to http://en.wikipedia.org/wiki/Lattice_reduction) for further detials on the concept of Lattice Reduction

find_csl_dsc.check_dsc_finder(R_G1ToG2_G1, Sigma, L_G1_GO1, L_DSC_G1, L_CSL_G1)[source]

This function checks the obtained DSC basis vectors are correct by using the following conditions: * Lattice 1 basis vectors are integer combinations of basis vectors of the DSC lattice * Lattice 2 basis vectors are integer combinations of basis vectors of the DSC lattice * The volume enclosed by the DSC is 1/sigma times the volume of lattice 1

Boundary Plane Basis

bp_basis.check_int_mats(l1, l2)[source]

The function checks whether or not the elements of $A_{left}^{-1}$*B are integers

Parameters:l1, l2 : numpy arrays of same dimensions
Returns:cond: numpy boolean array
bp_basis.check_2d_csl(l_pl1_g1, l_pl2_g1, l_csl_g1)[source]

The function checks whether or not the CSL basis may be expressed as a linear integer combination of the plane bases of planes 1 and 2

Parameters:

l_pl1_g1, l_pl2_g1: numpy arrays of basis vectors for plane 1 and 2

in the g1 reference frame

bp_basis.lbi_dioph_soln(a, b, c)[source]

Computes the diophantaine solution for the equation ax + by = c

bp_basis.compute_basis_vec(d_eq)[source]

The function computes y1, y2, y3 such that h*y1 + k*y2 + l*y3 = 0 and modulus of y1 is a minimum

Parameters:

d_eq: numpy array or list of size 3 and dimension 1

h = d_eq[0], k = d_eq[1], l = d_eq[2]

Returns:

np.array([y1, y2, y3])

bp_basis.bp_basis(miller_ind)[source]

The function computes the primitve basis of the plane if the boundary plane indices are specified

Parameters:

miller_ind: numpy array

Miller indices of the plane (h k l)

Returns:

l_pl_g1: numpy array

The primitive basis of the plane in ‘g1’ reference frame

bp_basis.pl_density(l_pl_g1, l_g1_go1)[source]

For a given two-dimensional plane basis, the planar density is computed

Parameters:

l_pl_g1: numpy array

l_g1_go1: numpy array

Basis vectors of the underlying lattice with respect to the orthogonal reference frame ‘go1’

Returns:

pd: float

Planar density = (1/area covered by plane basis)

bp_basis.csl_finder_2d(l_pl1_g1, l_pl2_g1)[source]

Given two plane bases, the 2D CSL bases are obtined by utilizing the smith normal form of the transformation between the two bases

Parameters:

l_pl1_g1, l_pl2_g1: numpy array

Basis vectors of planes 1 and 2 expressed in g1 reference frame

Returns:

l_2d_csl_g1: numpy array

The basis vectors of the 2D CSL expressed in g1 reference frame

bp_basis.gb_2d_csl(inds, t_mat, l_g_go, inds_type='miller_index', mat_ref='g1')[source]

For a given boundary plane normal ‘bp1_g1’ and the misorientation matrix ‘t_g1tog2_g1’, the two-dimensional CSL lattice is computed

Parameters:

inds: numpy array

The boundary plane indices

inds_type: string

{‘miller_index’, ‘normal_go’, ‘normal_g’}

t_mat: numpy array

Transformation matrix from g1 to g2 in ‘mat_ref’ reference frame

mat_ref: string

{‘go1’, ‘g1’}

lattice: Lattice class

Attributes of the underlying lattice

Returns:

l_2d_csl_g1, l_pl1_g1, l_pl2_g1: numpy arrays

l_2d_csl_g1 is the 2d CSL in g1 ref frame.

l_pl1_g1 is the plane 1 basis in g1 ref frame.

l_pl2_g1 is the plane 2 basis in g1 ref frame.

bp_basis.bicryst_planar_den(inds, t_mat, l_g_go, inds_type='miller_index', mat_ref='go1')[source]

The function computes the planar densities of the planes 1 and 2 and the two-dimensional CSL

Parameters:

inds: numpy array

The boundary plane indices.

inds_type: string

{‘miller_index’, ‘normal_go’, ‘normal_g’}

t_mat: numpy array

Transformation matrix from g1 to g2 in go1 (or g1) reference frame.

mat_ref: string

{‘go1’, ‘g1’}

lattice: Lattice class

Attributes of the underlying lattice.

Returns:

pl_den_pl1, pl_den_pl2: numpy array

The planar density of planes 1 and 2.

pl_den_csl: numpy array

The planare density of the two-dimensional CSL.

Misorientation Fundamental Zones

misorient_fz.misorient_fz(misquats, cryst_ptgrp, tol=1e-12)[source]

The function takes as input the misorientations and the corresponding crystallographic point group. It converts them using symmetry operations and returns the disorientations

Parameters:

misquats: Quaternion class

Quaternion misorientations

cryst_ptgrp: string

Crystallogrphic point group in Schoenflies notation

tol: float

Tolerance for the disorientation to belong in the fundamental zone

Returns:

disquats: quaternion class

Disorientations for the given misorientations

misorient_fz.check_cond(g, cryst_ptgrp, tol)[source]
Parameters:

g: quaternion object

Misorientation

cryst_ptgrp: string

Crystallogrphic point group in Schoenflies notation

tol: float

Tolerance for the misorientation to belong in the fundamental zone

Returns:

True or False: Boolean

Depending on whether or not the misorientation is a disorientation

Generate Symmetry Operators

generate_symm_ops.generate_symm_mats(cryst_ptgrp, tol=1e-10)[source]

Give crystallographic point group, this function generates all the symmetry operations (as matrices) that belong to the point group using ‘generators’

Parameters:

cryst_ptgrp: string

Crystallogrphic point group in Schoenflies notation

tol: float

The tolerance used to check if two matrices are the same

Returns:

symm_mat: numpy array

Size: n x 3 x3

Symmetry operations as matrices for the corresponding point group

generate_symm_ops.generate_symm_quats(cryst_ptgrp, tol=1e-10)[source]

Give crystallographic point group, this function generates all the symmetry operations (as quaternions) that belong to the point group using ‘generators’

Parameters:

cryst_ptgrp: string

Crystallogrphic point group in Schoenflies notation

tol: float

The tolerance used to check if two matrices are the same

Returns:

symm_quat: quaternion array

Size: n x 5

Symmetry operations as matrices for the corresponding point group

generate_symm_ops.save_symm_pkl(cryst_ptgrp, op_type)[source]

A pkl file with the symmetry operations of op_type (matrices or quaternions) are created and stored in the ‘pkl_files’ directory

Parameters:

cryst_ptgrp: string

Crystallogrphic point group in Schoenflies notation

op_type: {‘matrices’, ‘quats’}

Creates matrices or quaternion symmetry operations depending on op_type

Tools

tools.vrrotmat2vec(mat1, rot_type='proper')[source]

Create an axis-angle np.array from Rotation Matrix:

Parameters:

mat1: nx3x3 numpy array

The nx3x3 rotation matrices to convert

rot_type: string (‘proper’ or ‘improper’)

improper if there is a possibility of having improper matrices in the input, proper otherwise.

Default: proper

Returns:

ax_ang: numpy 5xn array

The 3D rotation axis and angle (ax_ang)

5 entries:

First 3: axis

4: angle

5: 1 for proper and -1 for improper

tools.vrrotvec2mat(ax_ang)[source]

Create a Rotation Matrix from Axis-Angle vector:

Parameters:

``ax_ang``: numpy 5xn array

The 3D rotation axis and angle (ax_ang)

5 entries:

First 3: axis

4: angle

5: 1 for proper and -1 for improper

Returns:

mtx: nx3x3 numpy array

3x3 rotation matrices

tools.unique_rows_tol(data, tol=1e-12, return_index=False, return_inverse=False)[source]

This function returns the unique rows of the input matrix within that are within the specified tolerance.

Parameters:

data: numpy array (m x n)

tol: double

tolerance of comparison for each rows Default: 1e-12

return_index: Boolean

flag to return the index of unique rows based on the indices of the output

return_inverse: Boolean

flag to return the index of unique rows based on the indices of the input

Returns:

unique_rows: numpy array (m’ x n)

ia: numpy array, integer (m’ x 1)

unique rows based on the indices of the output

ic: numpy array, integer (m x 1)

unique rows based on the indices of the input

See also

unique

tools.axang2quat(ax_ang)[source]

Create a quaternion corresponding to the rotation specified by an axis and an angle

Parameters:ax_ang: numpy array or a list of (4 x 1)
Returns:quaternion_rep: numpy array (5 x 1)
tools.mat2quat(mat, rot_type='proper')[source]

Convert Rotation Matrices to Quaternions

Parameters:

mat: numpy array or a list of (3 x 3)

rotation matrix

rot_type: string (‘proper’ or ‘improper’)

improper if there is a possibility of having improper matrices in the input, proper otherwise.

Default: proper

Returns:

quaternion_rep: numpy array (5 x 1)

See also

quat2mat, axang2quat

tools.quat2mat(q)[source]

Convert Quaternion Arrays to Rotation Matrix

Parameters:

q: numpy array (5 x 1)

quaternion

Returns:

g: numpy array (3 x 3)

rotation matrix

See also

mat2quat, axang2quat

tools.lll_reduction(matrix, delta=0.75)[source]

This function has been borrowed from pymatgen project with slight changes in input and output handling. Link to the project:

Link to the source code:
https://github.com/materialsproject/pymatgen/blob/92ee88ab6a6ec6e27b717150931e6d484d37a4e6/pymatgen/core/lattice.py

Performs a Lenstra-Lenstra-Lovasz lattice basis reduction to obtain a c-reduced basis. This method returns a basis which is as “good” as possible, with “good” defined by orthongonality of the lattice vectors.

Parameters:

delta: float

Reduction parameter.

Default of 0.75 is usually fine.

Returns:

a: numpy array

Reduced lattice:

tools.eq(m1, m2, tol)[source]

Check if the two rotation matrices are the same

tools.message_display(CheckMatrix, Checknumber, Message, Precis)[source]

This function displays a Message (passed as input) and gives and error in case the matrix passed to it is not integral.`

tools.ehermite(a, b)[source]

Elementary Hermite tranformation. For integers a and b, E = ehermite(a,b) returns an integer matrix with determinant 1 such that E * [a;b] = [g;0], where g is the gcd of a and b. E = ehermite(a,b) This function is in some ways analogous to GIVENS. John Gilbert, 415-812-4487, December 1993 gilbert@parc.xerox.com Xerox Palo Alto Research Center

Parameters:

a, b: integers

Returns:

E: numpy array 3x3

integer matrix with determinant 1 such that E * [a;b] = [g;0], where g is the gcd of a and b.

tools.smith_nf(matrix)[source]

Smith normal form of an integer matrix. [U,S,V] = smith(A) returns integer matrices U, S, and V such that A = U*S*V’, S is diagonal and nonnegative, S(i,i) divides S(i+1,i+1) for all i, det U =+-1, and det V =+-1. s = smith(A) just returns diag(S). Uses function ehermite. [U,S,V] = smith(A);

This function is in some ways analogous to SVD. Originally implemented by: John Gilbert, 415-812-4487, December 1993 gilbert@parc.xerox.com Xerox Palo Alto Research Center

Parameters:

matrix: numpy array

Returns:

S: numpy array

S is diagonal and nonnegative, S(i,i) divides S(i+1,i+1) for all i

U: numpy array

det(U) =+-1

V: numpy array

det(V) =+-1

tools.extgcd(x, y)[source]

Return a tuple (u, v, d); they are the greatest common divisor d of two integers x and y and u, v such that d = x * u + y * v.

class tools.Col[source]

This class is defined to ouput a word or sentence in a different color to the standard shell. The colors available are: pink, blue, green, dgrn: dark green, yellow, amber

Methods

c_prnt(text, color) Print a string in color,
c_prnt(text, color)[source]

Print a string in color,

Parameters:

Col: Col class instance

an instance of the Col class

text: string

Text to be shown in color.

color: string

Returns:

N/A

Indices and tables