Welcome to phenum’s documentation!

Contents:

Methods used for generating the symmetry group. All of these methods were adapted to python from their original fortran implementations which can be found at: https://github.com/msg-byu/symlib

phenum.symmetry.bring_into_cell(vec, cart_to_latt, latt_to_cart, eps)

This subroutine translates a point back into the unit cell in lattice coordinates, the coefficients of the point must all be less than one and at least zero.

Parameters:
  • vec (list of int) – The atom’s position vector
  • cart_to_latt (numpy ndarray) – The matrix that transforms from cartesian to lattice coordinates
  • latt_to_cart (numpy ndarray) – The matrix that transforms form lattice to cartesian coordinates
  • eps (float) – Finite precision tolerance.
Returns:

vec – The vector brought back into the cell.

Return type:

list of float

phenum.symmetry.get_concs_for_size(size, nspecies, res_concs, nB, concs)

Gets the concentration ranges for the atoms within the cells of certain sizes given the constraints provided such as concentration restrictions and the presence of arrows. Code rewritten from the get_conetration_list subroutine of: https://github.com/msg-byu/enumlib/blob/master/src/derivative_structure_generator.f90

Parameters:
  • size (int) – The cell size in integer form
  • nspecies (int) – the number of atomic species in the system
  • nB (int) – the number of basis vectors being used
  • res_concs (bool) – a logical that indicates of the concentrations are being restricted
  • concs (list of int) – A 2D integer array that contains the concentration ranges for each atom in the system
Returns:

c_list – The allowed concentration ranges.

Return type:

array-like

phenum.symmetry.get_lattice_pointGroup(a_vecs, eps=1e-10)
This routine returns only the point group of the lattite rather
than the space group of the given crystal structure.
Parameters:
  • a_vecs (array-like) – The 2D array that contains the parent lattice vectors as row vectors.
  • eps (float, optional) – Finite precision tolerance
Returns:

lattpg_op – The point group for the lattice in cartesian coordinates.

Return type:

array-like

phenum.symmetry.get_spaceGroup(par_lat, atomType, bas_vecs, eps=1e-10, lattcoords=False)

This routine takes a crystal structure (basis vectors and basis atoms and positions) and returns the point operators and fractional translations of the space group. The routine assumes that the given crystal structure is already primitive.

Parameters:
  • par_lat (array-like) – A 2D array that contains the parent lattice vectors.
  • atomType (list of int) – Integer array representing the type of each basis atom.
  • bas_vecs (array-like) – A 2D array that contains the basis vectors for the cell.
  • eps (float, optional) – Finite precisions tolerance.
  • lattcoords (bool, optional) – True if vectors are in lattice coordinates rather than cartesian.
Returns:

(sg_ops, sg_fracts) – The rotation and mirror operations of the space group, and the translation operations of the space group.

Return type:

array-like, array-like

phenum.element_data.get_lattice_parameter(elements, concentrations, default_title)

Finds the lattice parameters for the provided atomic species using Vagars law.

Parameters:
  • elements (list of str) – A list of the elements in the system.
  • default_title (str) – The default system title.
  • concentrations (list of int) – The concentrations of each element.
Returns:

lat_param – The lattice parameter of the system.

Return type:

float

Raises:

ValueError – if the number of elements doesn’t match the len of the concentration array.

Methods needed for the enumeration of arrows.

phenum.phonons.add_arrows(col, agroup, dim, translations, accept=None, nested=False, num_wanted=None, small=False, supers=False)

Finds the unique arrangements of arrows for a given configuration.

Parameters:
  • col (list) – A 2D integer array of the initial labeling.
  • dim (int) – The number of directions the arrows can point.
  • translations (list) – The translations of the lattice.
  • accept (float, optional) – acceptance rate of configurations for large enumerations.
  • agroup (list) – The stabilizers for the colors only with the arrow permutations.
  • nested (bool, optionl) – Set to True if this is called from within the context of another progress bar.
  • supers (bool, optional) – True if we want to include super periodic arrangements in the enumeration.
Returns:

arsurvivors – The unique arrow arrangements.

Return type:

list

phenum.phonons.arrow_concs(cList, aconcs)

Uses the concentrations of the atoms and the arrows to make a labeling for the system. It then sorts them to be in the correct order for the code.

Parameters:
  • cListr (list) – An integer array the concentration of the colors.
  • aconcs (list) – An array of floats of the fractional number of arrows for each color
Returns:

decoration_w_arrows – The labeling of both colors and arrows.

Return type:

list

phenum.phonons.enum_sys(groupfile, concs, a_concs, num_wanted, HNF, params, supers, accept=None)

Enumerates a random subset of the unique structures that have the shape defined by the symmetry group and the specified concentration.

Parameters:
  • groupfile (str) – Path to the file containing the symmetry group.
  • concs (list) – Integer array of the concentrations for each species.
  • a_concs (list) – Integer array of the arrow concentrations for each species.
  • HNF (list) – The HNF matrix for the system we’re currently enumerating
  • params (dict) – The dictionary of parameters read in from lattice.in
  • supers (bool) – True if superperiodic structures are to be kept.
  • accept (float, optional) – The acceptance rate for this enumeration.
  • num_wanted (int) – The number of structures to pick randomly from the enumerated list.
Returns:

configs – A list of the unique labelings in this system.

Return type:

list

Raises:

ValueError – if the number of configurations found doesn’t match the number requested.

phenum.phonons.get_arrow_concs(params)

If the concentrations are being restricted then find the correct arrow for each species included.

Parameters:params (dict) – The parameters read in from lattice.in.
Returns:a_concs – The arrow concentration for each species.
Return type:list
phenum.phonons.how_many_arrows(tcol)

Determines the number of colors that have arrows on them.

Parameters:tcol (list) – A 2D array of the labeling that contains the colors and arrows for each site.
Returns:arrows – The number of arrows in the system. n_species (int): The number of different atomic species with arrows. aconcs (list): The concentration of each arrow species.
Return type:int

Methods for reading and writing enumeration and polya-counting results.

phenum.structures.make_enum_in(infile, distribution, outfile, number=None, sizes=None, save=True, seed=None, restrict=None)

Makes an enum.in file if the distrubiton type is all with the desired number of structures. Otherwise prints the distribution information to the screen for the user.

Parameters:
  • infile (str) – The name of the ‘polya.out’ style input file.
  • n (int) – The number of structures or ‘all’.
  • outfile (str) – The name of the output file for the distribution
  • save (bool) – True if the data is to be saved to file.
  • seed (float) – The seed for the random number generator.
  • distribution (list) – The parameters that the distribution is over (‘shape’, ‘conc’, ‘size’, ‘all’).
  • sizes (list) – when specified, limit the distribution to these integer cell sizes; otherwise, look for all cell sizes we have data for.
  • restrict (list) – A list containing the restriction type (‘shape’,’conc’) and the file of allowed values.
Raises:

ValueError – if the ‘polya.*’ file cannot be found.

Contains routines used to map the enumerated structures to real space.

phenum.vector_utils.cartesian2direct(sLV, aBas, eps)

This routine takes three lattice vectors and a list of atomic basis vector in Cartesian coordinates and converts them to direct (“lattice”) coordinates.

Parameters:
  • sLV (list) – The superlattice vectors in cartesian coordinates as rows of a matrix.
  • aBas (list) – Atomic basis vectors in cartesian coordinates.
  • eps (float) – Finite precision tolerance.
Returns:

aBas – The atomic basis vectors in direct coordinates.

Return type:

list

phenum.vector_utils.map_enumStr_to_real_space(system_data, structure_data, minkowskiReduce)

Maps an enumerated structure back to real space. Returns a dictionary containing the real space data.

Parameters:
  • system_data (dict) –

    A dictionary containing all the information about the sysytem with keys:

    “plattice”: The parrent lattice vectors.

    “dvecs”: The atomic basis vectors.

    “title”: The name of the system.

    “bulksurf”: ‘bulk’ or ‘surface’.

    “nD”: The number of atomic basis vectors.

    “k”: The number of atomic species.

    “eps”: Finite precision tolerance.

  • sturture_data (dict) –

    A dictionary containing the information for this structure with keys:

    “strN”: The structure number.

    “hnfN”: The HNFs number.

    “hnf_degen”: The HNFs degeneracy.

    “lab_degen”: The label degeneracy.

    “tot_degen”: The total degeneracy for the structure.

    “sizeN”: The system size.

    “n”: Number of structure within this size.

    “pgOps”: The number of point group operations.

    “diag”: The diagonal of the SNF.

    “HNF”: The HNF matrix.

    “L”: The left transform matrix.

    “labeling”: The atomic labeling for the structure.

    “directions”: The arrow labeling for the structure.

  • minkowskiReduce (bool) – Logical indicating if basis should be reduced.
Returns:

space_data

A dictionary of the structure mapped to real space with keys:

“sLV”: The super lattice vectors as rows in a matrix.

“aBas”: The attomic basis vectors.

“spin”: A list of the occupanices.

“gIndx”: The integer group index.

“x”: A list of the concentrations of each atom type.

Return type:

dict

The heart of the code is here. These methods build the tree and check the canfigurations for uniqueness.

phenum.tree.brancher(concs, group, colors_w_arrows, dim, supers, cellsize, total=0, subset=None, accept=None, seed=None)

This routine navigates the tree and saves the unique configurations to an array survivors.

Parameters:
  • concs (list) – The concentrations of the colors or atoms.
  • dim (int) – The number of directions the arrows can point.
  • supers (bool) – True if superperiodic structures are to be kept.
  • cellsize (int) – The number of cells in the system.
  • total (int, optional) – The total number predicted by polya. Default is 0.
  • group (list) – The symmetry group for the system, including the effect on the arrows (displacement directions).
  • colors_w_arrows (list) – An integer 2D array that indiciates which atoms are being displaced, i.e., where the arrows are.
  • subset (list, optional) – An integer array of the subset of unique arrangements wanted. Default is None.
  • accept (float, optional) – for large enumerations, how often to accept configurations. Default is None.
  • seed (int, optaion) – The random seed. Default is None.
Returns:

survivors – A 3D array of the unique arrangements of colors and arrows.

Return type:

list

phenum.tree.guess_and_check_brancher(concs, group, colors_w_arrows, dim, supers, cellsize, num_wanted)

When the number of configurations wanted is sufficiently small relative to the total number of unique configurations then it is faster to ‘pick’ them from the list of possible configurations and then verify that none of your selected configurations are equivalent.

Parameters:
  • concs (list) – An integer arrray of the concentrations of the colors or atoms.
  • dim (int) – The number of directions the arrows can point.
  • supers (bool) – Logical that indicates if super periodic structures are to be kept.
  • cellsize (int) – The number of cells in the system.
  • group (list) – The symmetry group for the system, including the effect on the arrows (displacement directions).
  • colors_w_arrows (list) – An integer 2D array that indiciates which atoms are being displaced, i.e., where the arrows are.
Returns:

survivors – A 3D array of the unique arrangements of colors and arrows.

Return type:

list

Methods for taking an HNF and finding the site permutations and the arrow permutations. All of these have been modified for python from their original fortran implementations which can be found at: https://github.com/msg-byu/enumlib/. Of these only get_rotation get_rotation_perms_lists() has been modified in order to produce the arrow group. The class ArrowPerm has also been added and implemented inside the RotPermList class.

class phenum.grouptheory.ArrowPerm(site_perm=None, arrow_perm=None)

ArrowPerm pairs a site permutation with an arrow permutation.

site_perm

list – A 1D array containing a site permutation.

arrow_perm

list – A 1D array containing an arrow permutation.

class phenum.grouptheory.OpList(rot=None, shift=None)

OpList maintains the data structure of the fixing operations form the fortran code.

rot

array-like – A 3D array containing the rotations.

shift

array-like – A 2D array containing the translations of the lattice.

class phenum.grouptheory.RotPermList(nL=None, v=None, perm=None, RotIndx=None, arrows=None)

RotPermList maintains the data structure for the rotations and permutations from the fortran code.

nL

int – An integer indicating the number of operations.

v

array-like – A 2D array containing the lattice vectors.

perm

ArrowPerm – An ArrowPerm object.

RotIndx

list – List of the rotation indices.

phenum.grouptheory.SmithNormalForm(HNF)

Finds the Smith Normal Form (SNF) of an HNF as well as the left and right transforms.

Parameters:HNF (list of list) – An integer matrix for which the SNF is to be found.
Returns:(M, A, B) – The M is the Smith Normal Form matrix of the input matrix, A is the left transform matrix, and B is the right transform matrix. The matrices are such that np.dot(np.dot(A,HNF),B) = M.
Return type:matrix, matrix, matrix
phenum.grouptheory.a_group(trans, rots)

This subroutine combines that translations of the lattice with the rotatians of the lattice.

Parameters:
  • trans (list) – A 2D array where each row is an translation of the lattice.
  • rots (list) – A 3D array. Each row’s first entry is the site permutations and the second entry is the arrow permutations.
Returns:

groupi – The full symmetry group.

Return type:

list

phenum.grouptheory.a_group_gen(trans, rots)

This subroutine combines that translations of the lattice with the rotatians of the lattice for the generators of the system.

Parameters:
  • trans (list) – A 2D array where each row is an translation of the lattice.
  • rots (list) – A 3D array. Each row’s first entry is the site permutations and the second entry is the arrow permutations.
Returns:

groupi – The full symmetry group.

Return type:

list

phenum.grouptheory.get_full_HNF(HNF)

Turns the reduced HNF to a full HNF matrix.

Parameters:HNF (list) – The 1D array that contains all the lower diagonal entries of the HNF matrix.
Returns:full_hnf – A 2D list containing the HNF matrix.
Return type:list of list
phenum.grouptheory.get_sym_group(par_lat, bas_vecs, HNF, LatDim, arrows=True)

Generates the symmetry group for a given lattice and HNF

Parameters:
  • par_lat (numpy ndarray) – a 2D array that contains the parrent lattice vectors
  • bas_vecs (list) – A list of the atomic basis vectors for the lattice.
  • HNF (numpy ndarray) – The HNF matrix that defines the supercell.
  • LatDim (int) – 2 if a surface case case 3 if bulk.
  • arrows (bool, optional) – True if the arrow group is to be found as well.
Returns:

sym_group

The symmetry operations represented as permutations

of a list.

Return type:

RotPermList

Methods for reading and writing enumeration and polya-counting results.

phenum.io_utils.create_labeling(config)

This routine takes a string of atomic locations as a vector and the atomic concentrations and returns a unique labeling.

Parameters:config (list of int) – A list of integers describing the arrangement of atoms on the lattice.
Returns:label – The atomic labeling for the arrangement. arrow (str): The arrow directions for the arrangement.
Return type:str
phenum.io_utils.read_enum(filename='enum.in')

Reads the list of structures and the number of random candidates to enumerate from the ‘enum.in’ styled file. It should contain only integers and have:

HNF concs num_wanted

where,

HNF: the list of 6 independent entries in the HNF matrix; concs: list of integer species concentrations; num_wanted: number of random structures to draw from the enumerated list.

Parameters:filename (str, optional) – The filename and path. Default ‘enum.in’.
Returns:systems – A list containing the lists of [HNF, concs, num_wanted].
Return type:list of lists
phenum.io_utils.read_enum_out(args)

Reads the enum.out file and builds a dictionary with the needed information to construct a POSCAR.

Parameters:args (dict) – The makeStr.py input arguments.
Returns:system

A dictionary with the keys:

“plattice”: The parrent lattice vectors.

“dvecs”: The atomic basis vectors.

“title”: The name of the system.

“bulksurf”: ‘bulk’ or ‘surface’.

“nD”: The number of atomic basis vectors.

“k”: The number of atomic species.

“eps”: Finite precision tolerance.

structure_data (list of dict): A list of dictionaries for each desired structure. Each
dictionary contains:
“strN”: The structure number.

“hnfN”: The HNFs number.

“hnf_degen”: The HNFs degeneracy.

“lab_degen”: The label degeneracy.

“tot_degen”: The total degeneracy for the structure.

“sizeN”: The system size.

“n”: Number of structure within this size.

“pgOps”: The number of point group operations.

“diag”: The diagonal of the SNF.

“HNF”: The HNF matrix.

“L”: The left transform matrix.

“labeling”: The atomic labeling for the structure.

“directions”: The arrow labeling for the structure.

Return type:dict
phenum.io_utils.read_group(fname)

Reads the symmetry group in from the ‘rot_perms’ styled group output by enum.x.

Parameters:fname (str) – Path to the file to read the group from.
Returns:group – The symmetry group.
Return type:numpy ndarray
phenum.io_utils.read_lattice(filename='lattice.in')

Reads the lattice.in file.

Parameters:filename (str, optional) – The filename to be read in. Default is ‘lattice.in’.
Returns:result

A dictionary with the following fields:

“sizes”: the range of cell sizes,

“lat_vecs”: lattice vectors of the parent cell,

“nspecies”: the number of atomic species in the enumeration,

“basis_vecs”: basis vectors of the parent cell,

“is_crestricted”: logical that indicates if the concentrations will be restricted,

“arrows”: logical that indicates if arrows are present,

“concs”: array of the concentrations in format [1,3].

Return type:dict
phenum.io_utils.which(program)

Determines if and where an executable exists on the users path. This code was contributed by Jay at http://stackoverflow.com/a/377028

Parameters:program (str) – The name, or path for the program.
Returns:The program or executable.
phenum.io_utils.write_POSCAR(system_data, space_data, structure_data, args)

Writes a vasp POSCAR style file for the input structure and system data.

Parameters:
  • system_data (dict) – A dictionary of the system_data.
  • space_data (dict) – A dictionary containing the spacial data
  • structure_data (dict) – : a dictionary of the data for this structure
  • args (dict) – Dictionary of user supplied input.
phenum.io_utils.write_enum(params, outfile='enum.out')

Writes a ‘struct_enum.out’ style preamble to the specified output file using the enumeration parameters gleaned from a ‘lattice.in’ styled file.

Parameters:
  • params (dict) – values returned from method:read_lattice().
  • outfile (str, optional) – path to desired output file. Default is ‘enum.out’.
phenum.io_utils.write_struct_enum(params)

Writes a ‘struct_enum.in’ file for the executable enum.x.

Parameters:params (dict) – values returned from method:read_lattice().

Numerical methods for mathematical functions neede for the program

phenum.numerics.binomial_coefficient(n, k)

Finds the binomial coefficient n choose k. See https://en.wikipedia.org/wiki/Binomial_coefficient for details.

Parameters:
  • n (int) – An integer.
  • k (int) – An integer.
Returns:

t – The binomial coefficient, n choose k.

Return type:

int

phenum.numerics.factorial(num)

Finds the factorial of the input integer.

Parameters:num (int) – The integer to find the factorial of.
Returns:fact – The factorial of num.
Return type:int

AUTHORS: Conrad W. Rosenbrock, Wiley S. Morgan (June 2015)

Classes to support the calculation of coefficients for specific terms in a product of multinomials. Construct a product class by specifying the exponent and target term and then add multinomials using the Product instance’s append(). The coefficient is then available from the coeff().

class phenum.polyaburnside.Multinomial(power, coeff, arrowings, exponent=1)

Represents a multinomial expansion.

power

list – The power on each of the unexpanded variables in the multinomial; of the form (x^2+y^2+z^2) => 2.

coef

list – The coefficient on each of the variables in the multinomial; e.g. (x^2 + a y^2) => a. For more than two variables, each color that has arrowing gets the coefficient while the others get 1.

arrowings

list – The exponents of the arrows terms.

exponent

int – The exponent of the entire multinomial, default value is 1.

powersum

list – The integer value that all term exponents in the multinomial should sum to (or be less than).

possible_powers

list – The possible powers based on the exponent in the multinomial.

static nchoosek(n, k)

This implementation was taken from “Binomial Coefficient Computation: Recursion or Iteration?” by Yannis Manolopoulos, ACM SIGCSE Bulletin InRoads, Vol.34, No.4, December 2002. http://delab.csd.auth.gr/papers/SBI02m.pdf It is supposed to be robust against large, intermediate values and to have optimal complexity.

Parameters:
  • n (int) – The value of n in n choose k.
  • k (int) – The value of k in n choose k.
Returns:

t – The value of n choose k.

Return type:

int

nchoosekm(sequence)

Returns the number of different ways to partition an n-element set into disjoint subsets of sizes k1, ..., km.

Parameters:sequence (tuple) – An un-normed tuple of form (k1, k2, k3).
Returns:nckm – The value of the multinomial coefficient for the series.
Return type:int
normed_seq(seq)

Normalizes the specified sequence using the powers of unexpanded terms in the multinomial.

Parameters:seq (list) – A list of exponents in an expanded term.
Returns:norm_seq – The normalized lust of exponents.
Return type:list
possible_powers = None

For each variable being considered, determines the possible powers based on the exponent in the multinomial.

powersum = None

Returns the integer value that all term exponents in the multinomial should sum to (or be less than).

class phenum.polyaburnside.Product(coefficient, targets)

Represents a product of multinomials for which only a single term is interesting.

coefficient

int – The scalar integer multiplying this product of multinomials.

targets

list – A list of exponents for the terms in the products.

multinoms

list – A list of multinomials.

coeff()

Returns the coefficient of the term with the target exponents if all the multinomials in the product were expanded and had their terms collected.

Returns:coefficient – The coefficient of th term with the target exponents.
Return type:int
class phenum.polyaburnside.Sequence(root, possibles, i, powersum, targets, parent=None)

Represents an exponent-limited sequence with a single root. Here, sequence represents a sequence of integer values k_1, k_2...k_j that are the exponents of a single term in a multinomial. The root of the sequence is one of the k_i; its children become sets of sequences including variables to the right of i.

used

int – The sum of the exponents that have already been used.

parent

Sequence – A Sequence instance for the variable previouse to this one.

kids

Sequence – The next possible Sequence instances.

varcount

int – The number of variables in the polynomial.

expand(depth=0)

Recursively generates a list of all relevant sequences for this multinomial term.

Parameters:depth (int, optional) – The current depth in the recursion.
Returns:sequences – The relevant sequences for the multinomial.
Return type:list
expand_noappend(sequences, start, varindex)

Implements an expansion that doesn’t use python’s append.

Parameters:
  • sequences (list) – A list of the relevant sequences for the multinomials.
  • start (int) – Starting index.
  • varindex (int) – The current variabl index.
Returns:

sequences – The relevant sequences for the multinomial.

Return type:

list

Raises:

ValueError – if self.kidcount is zero.

kidcount

Returns the number of children and grandchildren to the last generation.

Returns:_kidcount – The number of descendents of the parent.
Return type:int
phenum.polyaburnside.group(gen)

Generates an entire group using the specified generators by applying generators to each other recursively.

Parameters:gen (list) – A list of generators as integers.
Returns:groupi – A list of the symmetry group operators.
Return type:list
phenum.polyaburnside.polya(concentrations, group, arrowings=None, debug=False)

Uses a group and concentrations to find the number of unique arrangements as described by polya.

Parameters:
  • concentrations (list) – : A list of integers specifying how many of each coloring should be present in each of the enumerated lists.
  • group (list) – A list of group operations for permuting the colorings.
  • arrowings (int, optional) – The number of arrows present. Default is None.
  • debug (bool, optional) – True if the code is being debugged.
Returns:

polya – The polya number for the system.

Return type:

int

Raises:

ValueError – if the concentrations don’t sum to the size of the group operation.

This module contains methods used for visualizing the enuerated structures.

phenum.visualize.HNF_atoms(enum, lattice, show, testmode=False)

Plots the atomic positions of the atoms in the cells.

Parameters:
  • enum (str) – The enum.in style input file.
  • lattice (str) – The lattice.in style input file.
  • show (bool) – If true each HNF is plotted in an interactive window.
  • testmode (bool, optional) – True if unit tests are being run.
phenum.visualize.HNF_shapes(enum, lattice, show, testmode=False)

Plots the shape of each HNF.

Parameters:
  • enum (str) – The enum.in style input file.
  • lattice (str) – The lattice.in style input file.
  • show (bool) – If true each HNF is plotted in an interactive window. testmode (bool, optional): True if unittests are running.

Indices and tables