jarvis.analysis.interface.zur ============================= .. py:module:: jarvis.analysis.interface.zur .. autoapi-nested-parse:: Design interface using Zur algorithm and Anderson rule. Classes ------- .. autoapisummary:: jarvis.analysis.interface.zur.ZSLGenerator Functions --------- .. autoapisummary:: jarvis.analysis.interface.zur.gen_sl_transform_matricies jarvis.analysis.interface.zur.rel_strain jarvis.analysis.interface.zur.rel_angle jarvis.analysis.interface.zur.fast_norm jarvis.analysis.interface.zur.vec_angle jarvis.analysis.interface.zur.vec_area jarvis.analysis.interface.zur.reduce_vectors jarvis.analysis.interface.zur.get_factors jarvis.analysis.interface.zur.make_interface jarvis.analysis.interface.zur.get_hetero_type jarvis.analysis.interface.zur.mismatch_strts jarvis.analysis.interface.zur.get_hetero Module Contents --------------- .. py:class:: ZSLGenerator(max_area_ratio_tol=0.09, max_area=400, max_length_tol=0.03, max_angle_tol=0.01) Bases: :py:obj:`object` Uses Zur algorithm to find best matched interfaces. This class is modified from pymatgen. .. py:attribute:: max_area_ratio_tol :value: 0.09 .. py:attribute:: max_area :value: 400 .. py:attribute:: max_length_tol :value: 0.03 .. py:attribute:: max_angle_tol :value: 0.01 .. py:method:: is_same_vectors(vec_set1, vec_set2) Check two sets of vectors are the same. Args: vec_set1(array[array]): an array of two vectors vec_set2(array[array]): second array of two vectors .. py:method:: generate_sl_transformation_sets(film_area, substrate_area) Generate transformation sets for film/substrate. The transformation sets map the film and substrate unit cells to super lattices with a maximum area. Args: film_area(int): the unit cell area for the film. substrate_area(int): the unit cell area for the substrate. Returns: transformation_sets: a set of transformation_sets defined as: 1.) the transformation matricies for the film to create a super lattice of area i*film area 2.) the tranformation matricies for the substrate to create a super lattice of area j*film area .. py:method:: get_equiv_transformations(transformation_sets, film_vectors, substrate_vectors) Apply the transformation_sets to the film and substrate vectors. Generate super-lattices and checks if they matches. Returns all matching vectors sets. Args: transformation_sets(array): an array of transformation sets: each transformation set is an array with the (i,j) indicating the area multipes of the film and subtrate it corresponds to, an array with all possible transformations for the film area multiple i and another array for the substrate area multiple j. film_vectors(array): film vectors to generate super lattices. substrate_vectors(array): substrate vectors to generate super lattices .. py:method:: __call__(film_vectors, substrate_vectors, lowest=False) Run the ZSL algorithm to generate all possible matching. .. py:method:: match_as_dict(film_sl_vectors, substrate_sl_vectors, film_vectors, substrate_vectors, match_area, film_transformation, substrate_transformation) Return dict which contains ZSL match. Args: film_miller(array) substrate_miller(array) .. py:function:: gen_sl_transform_matricies(area_multiple) Generate the transformation matricies. Convert a set of 2D vectors into a super lattice of integer area multiple as proven in Cassels: Cassels, John William Scott. An introduction to the geometry of numbers. Springer Science & Business Media, 2012. Args: area_multiple(int): integer multiple of unit cell area for super lattice area. Returns: matrix_list: transformation matricies to covert unit vectors to super lattice vectors. .. py:function:: rel_strain(vec1, vec2) Calculate relative strain between two vectors. .. py:function:: rel_angle(vec_set1, vec_set2) Calculate the relative angle between two vector sets. Args: vec_set1(array[array]): an array of two vectors. vec_set2(array[array]): second array of two vectors. .. py:function:: fast_norm(a) Much faster variant of numpy linalg norm. .. py:function:: vec_angle(a, b) Calculate angle between two vectors. .. py:function:: vec_area(a, b) Area of lattice plane defined by two vectors. .. py:function:: reduce_vectors(a, b) Generate independent and unique basis vectors based on Zur et al. .. py:function:: get_factors(n) Generate all factors of n. .. py:function:: make_interface(film='', subs='', atol=1, ltol=0.05, max_area=500, max_area_ratio_tol=1.0, seperation=3.0, vacuum=8.0, apply_strain=False) Use as main function for making interfaces/heterostructures. Return mismatch and other information as info dict. Args: film: top/film material. subs: substrate/bottom/fixed material. seperation: minimum seperation between two. vacuum: vacuum will be added on both sides. So 2*vacuum will be added. .. py:function:: get_hetero_type(A={}, B={}) Provide heterojunction classification using Anderson rule. .. py:function:: mismatch_strts() Return mismatch and other information as info dict. Deprecated, preserved here so legacy imports work. .. py:function:: get_hetero() Generate heterostructure. Deprecated also.