soprano.utils#
utils.py
Contains package-wide useful routines that don’t fall under any specific category. Many of these handle common operations involving periodicity, conversions between different representations etc.
Functions
|
Transforms an axes and angles representation of lattice parameters into a Cartesian one |
|
Find all the periodic equivalent vectors for a list of vectors and a given lattice falling within a given length. |
|
Average a list of quaternions. |
|
Transforms a Cartesian representation of lattice parameters into an axes and angles one |
|
Clebsch-Gordan cohefficients for given quantum numbers: j, m, j1, m1, j2, m2 |
|
Given a symmetric structure, compute a distance matrix for the given fractional coordinate points which contains only the distances between their closest symmetry equivalent sites in the asymmetric unit cell. |
|
|
|
Distance in bonds, returns -1 if the two are not connected |
|
|
|
Cluster points using given method if found in sklearn.clusters module. |
|
|
|
Check if the atoms object has CIF labels Does this simply by comparing the labels with the element symbols. |
|
Generate a matrix that turns hkl indices into inverse crystal plane distances for a given lattice in ABC form |
|
|
|
|
|
Calculate inverse planar distance for a given set of Miller indices h, k, l. |
|
Checks whether s is a string, with Python 2 and 3 compatibility |
|
Return an integer distance between two lists (number of differing elements) |
|
Maximum distance between two points achievable inside a single unit cell. |
Merge a list of arrays by concatenating them |
|
|
Merge a list of arrays by taking the first |
|
Merge a list of arrays by taking the mean |
|
Merge sites in a structure. |
|
Merge a list of arrays by summing them |
|
Find the shortest periodic equivalent vector for a list of vectors and a given lattice. |
|
Generate the bounds for a supercell containing a sphere of given radius, knowing the unit cell. |
|
Periodic version of the Bridson algorithm for generation of Poisson sphere distributions of points. |
|
A textual progress bar for the command line |
|
Creates a string for a given atom by traversing the molecular network starting with it. |
|
Repulsion algorithm, begins with a series of vectors v, finds a new one that is as far as possible, angle-wise, from all of them. |
|
Replace the folder of the given path with a new one |
|
Executes a Popen.communicate and returns output in a way that is compatible with Python 2 & 3 keeping input and output as strings (since Python 3 requires bytes objects otherwise) |
|
Ask for user input with Python 2 and 3 compatibility |
|
Get the filename (with no extension) from a full path |
|
Useful stdout/err silencer |
|
Generate a full linearized grid for a supercell with r_bounds and a base unit cell in Cartesian form. |
|
Perform a Swing*Twist decomposition of a Quaternion. |
|
Wigner 3j symbols for given quantum numbers: |
- soprano.utils.abc2cart(abc)[source]#
Transforms an axes and angles representation of lattice parameters into a Cartesian one
- soprano.utils.all_periodic(v, latt_cart, max_r, pbc=[True, True, True])[source]#
Find all the periodic equivalent vectors for a list of vectors and a given lattice falling within a given length.
Args:v (np.ndarray): list of 3-vectors representing points or vectors toproduce periodic versions oflatt_cart (np.ndarray): unit cell in cartesian formmax_r (float): maximum length of periodic copies of vectorspbc ([bool, bool, bool]): periodic boundary conditions - ifa boundary is not periodic therange returned will always be zeroin that dimensionReturns:v_period (np.ndarray): array with the same shape as v, containing thevectors in periodic reduced formv_index (np.ndarray): indices (referring to the original array v) ofthe array of which the corresponding element ofv_period is a copyv_cells (np.ndarray): array of triples of ints, corresponding to thecells from which the various periodic copies ofthe vectors were taken. For an unchanged vectorwill be all [0,0,0]
- soprano.utils.average_quaternions(quats)[source]#
Average a list of quaternions. Following: christophhagen/averaging-quaternions
Returns a normalized quaternion.
- Parameters:
quats (List[Quaternion])
- Return type:
Quaternion
- soprano.utils.cart2abc(cart)[source]#
Transforms a Cartesian representation of lattice parameters into an axes and angles one
- soprano.utils.clebsch_gordan(j, m, j1, m1, j2, m2)[source]#
Clebsch-Gordan cohefficients for given quantum numbers: j, m, j1, m1, j2, m2
The numbers passed can be arrays (just make sure they’re the same size).
- soprano.utils.compute_asymmetric_distmat(struct, points, linearized=True, return_images=False, images_centre=0, symprec=1e-05, spg=None)[source]#
Given a symmetric structure, compute a distance matrix for the given fractional coordinate points which contains only the distances between their closest symmetry equivalent sites in the asymmetric unit cell.
Distances are computed in fractional coordinates space and are not real distances. They do not necessarily map to real space distances either. Their purpose is merely to group points together by removing the effects of symmetry operations.
Parameters:struct (ase.Atoms): the structure from which to extract the symmetryoperationspoints (np.ndarray): list of points, in fractional coordinates, tocompute the distance matrix forlinearized (bool): if True, returns only the linearised uppertriangle of the distance matrix, in the formataccepted by scipy.cluster.hierarchy.linkage.Default is Truereturn_images (bool): if True, return also a list of the closestimages for each point. The closestsymmetric image to images_centre is chosen.images_centre (int or np.ndarray): centre to use to compute closestimages. If an integer, uses thegiven point of that index. If avector, uses the given coordinates.Should not be a high symmetry point.Returns:distmat (np.ndarray): distance matrix for pointsimages (np.ndarray): closest point images (only if return_images isTrue)
- soprano.utils.get_bonding_distance(bond_graph, i, j, nx=None)[source]#
Distance in bonds, returns -1 if the two are not connected
- soprano.utils.get_sklearn_clusters(points, method, params, sk=None)[source]#
Cluster points using given method if found in sklearn.clusters module. Use params as a dictionary of parameters passed to the method.
- soprano.utils.has_cif_labels(atoms)[source]#
Check if the atoms object has CIF labels Does this simply by comparing the labels with the element symbols. If they’re all the same, then so special labels are present.
- soprano.utils.hkl2d2_matgen(abc)[source]#
Generate a matrix that turns hkl indices into inverse crystal plane distances for a given lattice in ABC form
- soprano.utils.inv_plane_dist(hkl, hkl2d2)[source]#
Calculate inverse planar distance for a given set of Miller indices h, k, l.
- soprano.utils.list_distance(l1, l2)[source]#
Return an integer distance between two lists (number of differing elements)
- soprano.utils.max_distance_in_cell(cell)[source]#
Maximum distance between two points achievable inside a single unit cell. Relies on the fact that the isosurface for equal distance in cartesian space is an ellipsoid, meaning that the extreme has to be one of the 8 corners of the cube of side 1 in fractional space.
- soprano.utils.merge_sites(atoms, indices, merging_strategies={}, keep_all=False)[source]#
Merge sites in a structure.
Parameters:atoms (ase.Atoms): structure to merge sites inindices (list): list of lists of indices of sites to mergemerging_strategies (dict): dictionary of merging strategies forproperties, e.g. {‘positions’: lambda x: x.mean(axis=0)}if not specified, the default strategies will be used(see DEFAULT_MERGING_STRATEGIES)keep_all (bool): whether to keep all sites in the structure, default is False. If True, themerged sites will have the same properties (determined by merging_strategies)but will have a different index, so that they can be identified. If False, only one of the mergedsites will be kept, and the others will be removed from the Atoms object.Returns:atoms (ase.Atoms): structure with sites merged- Parameters:
atoms (Atoms)
- soprano.utils.minimum_periodic(v, latt_cart, exclude_self=False, pbc=[True, True, True])[source]#
Find the shortest periodic equivalent vector for a list of vectors and a given lattice.
Args:v (np.ndarray): list of 3-vectors representing points or vectors toreduce to their closest periodic versionlatt_cart (np.ndarray): unit cell in cartesian formexclude_self (bool): if True, any vector that is equal to zero will beexcluded, and its closest non-zero periodicversion will be considered instead. Default isFalsepbc ([bool, bool, bool]): periodic boundary conditions - ifa boundary is not periodic therange returned will always be zeroin that dimensionReturns:v_period (np.ndarray): array with the same shape as v, containing thevectors in periodic reduced formv_cells (np.ndarray): array of triples of ints, corresponding to thecells from which the various periodic copies ofthe vectors were taken. For an unchanged vectorwill be all [0,0,0]
- soprano.utils.minimum_supcell(max_r, latt_cart=None, r_matrix=None, pbc=[True, True, True])[source]#
Generate the bounds for a supercell containing a sphere of given radius, knowing the unit cell.
Args:max_r (float): radius of the sphere contained in the supercelllatt_cart (np.ndarray): unit cell in cartesian formr_matrix (np.ndarray): matrix for the quadratic form returningr^2 for this supercell.Alternative to latt_cart, for a directspace cell would be equal tonp.dot(latt_cart, latt_cart.T)pbc ([bool, bool, bool]): periodic boundary conditions - ifa boundary is not periodic therange returned will always be zeroin that dimensionReturns:shape (tuple[int]): shape of the supercell to be built.Raises:ValueError: if some of the arguments are invalid
- soprano.utils.periodic_bridson(cell, rmin, max_attempts=30, prepoints=None, prepoints_cuts=None, maxblock=1000)[source]#
Periodic version of the Bridson algorithm for generation of Poisson sphere distributions of points. This returns a generator.
Args:cell (np.ndarray): periodic cell in which to create the points.rmin (float): minimum distance between each generated point.max_attempts (int): maximum number of candidate neighbours generatedfor each point.prepoints (np.ndarray or list): pre-existing points to avoid duringgeneration. These must be infractional coordinates.prepoints_cuts (np.ndarray or list): custom cutoffs for each prepoint.If not included defaults to rmin.maxblock (int): maximum size of a grid block to be processed in asingle function (used to avoid excess memory use)Returns:bridsonGen (generator): an iterator producing Poisson-sphere likedistributed points in cell until space runsout or enough attempts fail.
- soprano.utils.progbar(i, i_max, bar_len=20, spinner=True, spin_rate=3.0)[source]#
A textual progress bar for the command line
Args:i (int): current progress indexmax_i (int): final progress indexbar_len (Optional[int]): length in characters of the bar (no brackets)spinner (Optional[bool]): show a spinner at the endspin_rate (Optional[float]): spinner rotation speed (turns per fullprogress)Returns:bar (str): a progress bar formatted as requested
- soprano.utils.recursive_mol_label(site_i, mol_indices, bonds, elems)[source]#
Creates a string for a given atom by traversing the molecular network starting with it.
Parameters:site_i (int): index of the atom for which the string must becalculatedmol_indices ([int]): indices of all atoms belonging to the moleculebonds (dict): dictionary of bonds within the molecule, containingthe indices of atoms as keys and lists of all theindices of atoms they’re bonded to as values.elems ([str]): list of element symbols for the entire system (mustbe returned by ase.Atoms’ get_chemical_symbols, NOTfollow the order the atoms appear in mol_indices)Returns:recursive_label (str): a string representing the molecular networktraversed starting from site i
- soprano.utils.rep_alg(v, iters=1000, attempts=10, step=0.1, simtol=1e-05)[source]#
Repulsion algorithm, begins with a series of vectors v, finds a new one that is as far as possible, angle-wise, from all of them. The process is that of treating the tips of the vectors as particles on a sphere which repel each other. Of course, it’s possible that multiple equilibrium positions exist, which is why multiple attempts can be made, and are judged identical or not based on a tolerance parameter
Parameters:v (np.ndarray): list of fixed vectors to avoiditers (int): number of iterations, default is 1000attempts (int): number of independent attempts with randomised starts,default is 10step (float): step by which the ‘particle’ describing the vector willbe displaced at each iteration during the search.Default is 1e-1.simtol (float): tolerance based on which two attempts are considered tobe effectively the same, checking if 1-v1.v2 < simtol.Default is 1e-5.Returns:out_v ([np.ndarray]): list of unit vectors of maximal distance fromthe ones in v
- soprano.utils.replace_folder(path, new_folder)[source]#
Replace the folder of the given path with a new one
- soprano.utils.safe_communicate(subproc, stdin='')[source]#
Executes a Popen.communicate and returns output in a way that is compatible with Python 2 & 3 keeping input and output as strings (since Python 3 requires bytes objects otherwise)
- soprano.utils.silence_stdio(silence_stdout=True, silence_stderr=True)[source]#
Useful stdout/err silencer
- soprano.utils.supcell_gridgen(latt_cart, shape)[source]#
Generate a full linearized grid for a supercell with r_bounds and a base unit cell in Cartesian form.
Args:latt_cart (np.ndarray): unit cell in cartesian formshape (tuple[int]): shape of the supercell to be built,as returned by minimum_supcell.Returns:neigh_i_grid (np.ndarray): supercell grid in fractional coordinatesneigh_grid (np.ndarray): supercell grid in cartesian coordinatesRaises:ValueError: if some of the arguments are invalid
- soprano.utils.swing_twist_decomp(quat, axis)[source]#
Perform a Swing*Twist decomposition of a Quaternion. This splits the quaternion in two: one containing the rotation around axis (Twist), the other containing the rotation around a vector parallel to axis (Swing).
Returns two quaternions: Swing, Twist.