pyscal Reference

pyscal.core module

Main module of pyscal. This module contains definitions of the two major classes in pyscal - the System and Atom. Atom is a pure pybind11 class whereas System is a hybrid class with additional python definitions. For the ease of use, Atom class should be imported from the core module. The original pybind11 definitions of Atom and System can be found in catom and csystem respectively.

class pyscal.core.System[source]

Bases: pyscal.csystem.System

A python/pybind11 hybrid class for holding the properties of a system.

box

A list containing the dimensions of the simulation box in the format [[x_low, x_high], [y_low, y_high], [z_low, z_high]]

Type:list of list of floats
atoms
Type:list of Atom objects

Notes

A System consists of two major components - the simulation box and the atoms. All the associated variables are then calculated using this class.

Note

atoms can be accessed or set as atoms. However, due to technical reasons individual atoms should be accessed using the get_atom() method. An atom can be assigned to the atom using the set_atom() method.

Examples

>>> sys = System()
>>> sys.read_inputfile('atoms.dat')
calculate_angularcriteria()[source]

Calculate the angular criteria for each atom :param None:

Returns:
Return type:None

Notes

Calculates the angular criteria for each atom as defined in [1]_. Angular criteria is useful for identification of diamond cubic structures. Angular criteria is defined by,

\[A = \sum_{i=1}^6 (\cos(\theta_i) + \frac{1}{3})^2\]

where cos(theta) is the angle size suspended by each pair of neighbors of the central atom. A will have a value close to 0 for structures if the angles are close to 109 degrees. The calculated A parameter for each atom is stored in angular.

References

[1]Uttormark, MJ, Thompson, MO, Clancy, P, Phys. Rev. B 47, 1993
calculate_chiparams(angles=False)[source]

Calculate the chi param vector for each atom

Parameters:angles (bool, optional) – If True, return the list of cosines of all neighbor pairs
Returns:angles – list of all cosine values, returned only if angles is True.
Return type:array of floats

Notes

This method tries to distinguish between crystal structures by finding the cosines of angles formed by an atom with its neighbors. These cosines are then historgrammed with bins [-1.0, -0.945, -0.915, -0.755, -0.705, -0.195, 0.195, 0.245, 0.795, 1.0] to find a vector for each atom that is indicative of its local coordination. Compared to chi parameters from chi_0 to chi_7 in the associated publication, the vector here is from chi_0 to chi_8. This is due to an additional chi parameter which measures the number of neighbors between cosines -0.705 to -0.195.

Parameter nlimit specifies the number of nearest neighbors to be included in the analysis to find the cutoff. If parameter angles is true, an array of all cosine values is returned. The publication further provides combinations of chi parameters for structural identification which is not implemented here. The calculated chi params can be accessed using chiparams.

References

[1]Ackland, Jones, Phys. Rev. B 73, 2006
calculate_disorder(averaged=False, q=6)[source]

Calculate the disorder criteria for each atom

Parameters:
  • averaged (bool, optional) – If True, calculate the averaged disorder. Default False.
  • q (int, optional) – The Steinhardt parameter value over which the bonds have to be calculated. Default 6.
Returns:

Return type:

None

Notes

Calculate the disorder criteria as introduced in [1]. The disorder criteria value for each atom is defined by,

\[D_j = \frac{1}{N_b^j} \sum_{i=1}^{N_b} [ S_{jj} + S_{kk} -2S_{jk}]\]

where .. math:: S_{ij} = sum_{m=-6}^6 q_{6m}(i) q_{6m}^*(i)

The keyword averaged is True, the disorder value is averaged over the atom and its neighbors. The disorder value can be accessed using disorder and the averaged version can be accessed using avg_disorder. For ordered systems, the value of disorder would be zero which would increase and reach one for disordered systems.

References

[1]Kawasaki, T, Onuki, A, J. Chem. Phys. 135, 2011
calculate_nucsize(frenkelnums, threshold, avgthreshold)[source]

Calculate the size of the largest cluster in the given system.

Parameters:None
Returns:cluster size – size of the largest solid cluster in liquid (number of atoms)
Return type:int

Notes

Calculation of the the size of the largest solid cluster needs various prerequisites that can be set by the functions set_nucsize_parameters.

Warning

This function is deprecated and will be removed in a future release. Please use System.find_solids() instead.

calculate_q(q, averaged=False)[source]

Find the Steinhardt parameter q_l for all atoms.

Parameters:
  • q_l (int or list of ints) – A list of all Steinhardt parameters to be found from 2-12.
  • averaged (bool, optional) – If True, return the averaged q values, default False
Returns:

Return type:

None

Notes

Enables calculation of the Steinhardt parameters [1] q from 2-12. The type of q values depend on the method used to calculate neighbors. See the description find_neighbors() for more details. If the keyword average is set to True, the averaged versions of the bond order parameter [2] is returned.

References

[1]Steinhardt, PJ, Nelson, DR, Ronchetti, M. Phys Rev B 28, 1983
[2]Lechner, W, Dellago, C, J Chem Phys, 2013
calculate_rdf(histobins=100, histomin=0.0, histomax=None)[source]

Calculate the radial distribution function.

Parameters:
  • histobins (int) – number of bins in the histogram
  • histomin (float, optional) – minimum value of the distance histogram. Default 0.0.
  • histomax (float, optional) – maximum value of the distance histogram. Default, the maximum value in all pair distances is used.
Returns:

  • rdf (array of ints) – Radial distribution function
  • r (array of floats) – radius in distance units

calculate_solidneighbors()[source]

Find Solid neighbors of all atoms in the system.

Parameters:None
Returns:
Return type:None

Notes

A solid bond is considered between two atoms if the connection between them is greater than 0.6.

calculate_sro(reference_type=1, average=True)[source]

Calculate short range order

Parameters:
  • reference_type (int, optional) – type of the atom to be used a reference. default 1
  • average (bool, optional) – if True, average over all atoms of the reference type in the system. default True.
Returns:

vec – The short range order averaged over the whole system for atom of the reference type. Only returned if average is True. First value is SRO of the first neighbor shell and the second value corresponds to the second nearest neighbor shell.

Return type:

list of float

Notes

Calculates the short range order for an AB alloy using the approach by Cowley [1]. Short range order is calculated as,

\[\alpha_i = 1 - \frac{n_i}{m_A c_i}\]

where n_i is the number of atoms of the non reference type among the c_i atoms in the ith shell. m_A is the concentration of the non reference atom. Please note that the value is calculated for shells 1 and 2 by default. In order for this to be possible, neighbors have to be found first using the find_neighbors() method. The selected neighbor method should include the second shell as well. For this purpose method=cutoff can be chosen with a cutoff long enough to include the second shell. In order to estimate this cutoff, one can use the calculate_rdf() method.

References

[1]Cowley J. M., PR 77(5), 1950.
calculate_vorovector(edge_cutoff=0.05, area_cutoff=0.01, edge_length=False)[source]

get the voronoi structure identification vector.

Parameters:edge_cutoff (float, optional) – cutoff for edge length. Default 0.05.
area_cutoff : float, optional
cutoff for face area. Default 0.01.
edge_length : bool, optional
if True, a list of unrefined edge lengths are returned. Default false.
Returns:vorovector – array of the form (n3, n4, n5, n6)
Return type:array like, int

Notes

Returns a vector of the form (n3, n4, n5, n6), where n3 is the number of faces with 3 vertices, n4 is the number of faces with 4 vertices and so on. This can be used to identify structures [1] [2].

The keywords edge_cutoff and area_cutoff can be used to tune the values to minimise the effect of thermal distortions. Edges are only considered in the analysis if the edge_length/sum(edge_lengths) is at least edge_cutoff. Similarly, faces are only considered in the analysis if the face_area/sum(face_areas) is at least face_cutoff.

References

[1]Finney, JL, Proc. Royal Soc. Lond. A 319, 1970
[2]Tanemura, M, Hiwatari, Y, Matsuda, H,Ogawa, T, Ogita, N, Ueda, A. Prog. Theor. Phys. 58, 1977
cluster_atoms(condition, largest=True)[source]

Cluster atoms based on a property

Parameters:
  • condition (callable) – function which should take an Atom object, and give a True/False output
  • largest (bool, optional) – If True returns the size of the largest cluster. Default False.
Returns:

lc – Size of the largest cluster. Returned only if largest is True.

Return type:

int

Notes

This function helps to cluster atoms based on a defined property. This property is defined by the user through the function condition which is passed as a parameter. For each atom in the system, the condition should give a True/False values.

When clustering, the code loops over each atom and its neighbors. If the condition is true for both host atom and the neighbor, they are assigned to the same cluster. For example, a condition to cluster solid atoms would be,

def condition(atom):
    #if both atom is solid
    return (atom1.solid)

Check examples for more details.

find_clusters(recursive=True, largest=True)[source]

Find the clusters of all atoms in the system.

Parameters:
  • recursive (Bool, optional) – If True, use a recursive clustering algorithm, otherwise use an id based clustering. The difference in values between two methods can be upto 3 particles. Default True.
  • largest (Bool, optional) – If True, return the number of particles in the largest cluster. Default True.
Returns:

cluster – The size of the largest cluster in the system. Only returned if largest is set to True.

Return type:

int

Notes

Go through all the atoms in the system and cluster them together based on the issolid parameter of the atom. To cluster based on any user defined criteria, you can use set_solid method of Atom to explicitely set the issolid value.

Warning

This function is deprecated and will be removed in a future release. Please use cluster_atoms() instead.

find_largestcluster()[source]

Find the largest solid cluster of atoms in the system from all the clusters.

Parameters:None
Returns:cluster – the size of the largest cluster
Return type:int

Notes

pyscal.core.System.find_clusters() has to be used before using this function.

find_neighbors(method='cutoff', cutoff=None, threshold=2, filter=None, voroexp=1, padding=1.2, nlimit=6, cells=False)[source]

Find neighbors of all atoms in the System.

Parameters:
  • method ({'cutoff', 'voronoi'}) – cutoff method finds neighbors of an atom within a specified or adaptive cutoff distance from the atom. voronoi method finds atoms that share a Voronoi polyhedra face with the atom. Default, cutoff
  • cutoff ({ float, 'sann', 'adaptive'}) – the cutoff distance to be used for the cutoff based neighbor calculation method described above. If the value is specified as 0 or adaptive, adaptive method is used. If the value is specified as sann, sann algorithm is used.
  • threshold (float, optional) – only used if cutoff=adaptive. A threshold which is used as safe limit for calculation of cutoff.
  • filter ({'None', 'type'}, optional) – apply a filter to nearest neighbor calculation. If the filter keyword is set to type, only atoms of the same type would be included in the neighbor calculations. Default None.
  • voroexp (int, optional) – only used if method=voronoi. Power of the neighbor weight used to weight the contribution of each atom towards Steinhardt parameter values. Default 1.
  • padding (double, optional) – only used if cutoff=adaptive. A safe padding value used after an adaptive cutoff is found. Default 1.2.
  • nlimit (int, optional) – only used if cutoff=adaptive. The number of particles to be considered for the calculation of adaptive cutoff. Default 6.
Returns:

Return type:

None

Raises:
  • RuntimeWarning – raised when threshold value is too low. A low threshold value will lead to ‘sann’ algorithm not converging when finding a neighbor. This function will try to automatically increase threshold and check again.
  • RuntimeError – raised when neighbor search was unsuccessful. This is due to a low threshold value.

Notes

This function calculates the neighbors of each particle. There are several ways to do this. A complete description of the methods can be found here.

Method cutoff and specifying a cutoff radius uses the traditional approach being the one in which the neighbors of an atom are the ones that lie in the cutoff distance around it.

In order to reduce time during the distance sorting during the adaptive methods, pyscal sets an initial guess for a cutoff distance. This is calculated as,

\[r_{initial} = threshold * (simulation~box~volume/ number~of~particles)^{(1/3)}\]

threshold is a safe multiplier used for the guess value and can be set using the threshold keyword.

In Method cutoff, if cutoff='adaptive', an adaptive cutoff is found during runtime for each atom [1]. Setting the cutoff radius to 0 also uses this algorithm. The cutoff for an atom i is found using,

\[r_c(i) = padding * ((1/nlimit) * \sum_{j=1}^{nlimit}(r_{ij}))\]

padding is a safe multiplier to the cutoff distance that can be set through the keyword padding. nlimit keyword sets the limit for the top nlimit atoms to be taken into account to calculate the cutoff radius.

In Method cutoff, if cutoff='sann', sann algorithm is used [2]. There are no parameters to tune sann algorithm.

The second approach is using Voronoi polyhedra which also assigns a weight to each neighbor in the ratio of the face area between the two atoms. Higher powers of this weight can also be used [3]. The keyword voroexp can be used to set this weight.

Warning

Adaptive cutoff uses a padding over the intial guessed “neighbor distance”. By default it is 2. In case of a warning that threshold is inadequate, this parameter should be further increased. High/low value of this parameter will correspond to the time taken for finding neighbors.

References

[1]Stukowski, A, Model Simul Mater SC 20, 2012
[2]van Meel, JA, Filion, L, Valeriani, C, Frenkel, D, J Chem Phys 234107, 2012
[3]Haeberle, J, Sperl, M, Born, P, arxiv 2019
find_solids(bonds=0.5, threshold=0.5, avgthreshold=0.6, cluster=True, q=6)[source]

Distinguish solid and liquid atoms in the system.

Parameters:
  • bonds (int or float, optional) – Minimum number of solid bonds for an atom to be identified as a solid if the value is an integer. Minimum fraction of neighbors of an atom that should be solid for an atom to be solid if the value is float between 0-1. Default 0.5.
  • threshold (double, optional) – Solid bond cutoff value. Default 0.5.
  • avgthreshold (double, optional) – Value required for Averaged solid bond cutoff for an atom to be identified as solid. Default 0.6.
  • cluster (bool, optional) – If True, cluster the solid atoms and return the number of atoms in the largest cluster.
  • q (int, optional) – The Steinhardt parameter value over which the bonds have to be calculated. Default 6.
Returns:

solid – Size of the largest solid cluster. Returned only if cluster=True.

Return type:

int

Notes

The neighbors should be calculated before running this function. Check find_neighbors() method.

bonds define the number of solid bonds of an atom to be identified as solid. Two particles are said to be ‘bonded’ if [1],

\[s_{ij} = \sum_{m=-6}^6 q_{6m}(i) q_{6m}^*(i) \geq threshold\]

where threshold values is also an optional parameter.

If the value of bonds is a fraction between 0 and 1, at least that much of an atom’s neighbors should be solid for the atom to be solid.

An additional parameter avgthreshold is an additional parameter to improve solid-liquid distinction. In addition to having a the specified number of bonds,

\[\langle s_{ij} \rangle > avgthreshold\]

also needs to be satisfied. In case another q value has to be used for calculation of S_ij, it can be set used the q attribute.

References

[1]Auer, S, Frenkel, D. Adv Polym Sci 173, 2005
from_pickle(file)[source]

Read the contents of System object from a pickle file.

Parameters:file (string) – name of input file
Returns:
Return type:None

Notes

This function can be used to set up a system from a file. A System object needs to be created first.

Examples

>>> sys = System()
>>> sys.from_pickle(filename)
get_atom(index)[source]

Get the Atom object at the queried position in the list of all atoms in the System.

Parameters:index (int) – index of required atom in the list of all atoms.
Returns:atom – atom object at the queried position.
Return type:Atom object
get_distance(atom1, atom2, vector=False)[source]

Get the distance between two atoms.

Parameters:
  • atom1 (Atom object) – first atom
  • atom2 (Atom object) – second atom
  • vector (bool, optional) – If True, the displacement vector connecting the atoms is also returned. default false.
Returns:

distance – distance between the first and second atom.

Return type:

double

Notes

Periodic boundary conditions are assumed by default.

get_qvals(q, averaged=False)[source]

Get the required q_l (Steinhardt parameter) values of all atoms.

Parameters:
  • q_l (int or list of ints) – required q_l value with l from 2-12
  • averaged (bool, optional) – If True, return the averaged q values, default False
Returns:

qvals – list of q_l of all atoms.

Return type:

list of floats

Notes

The function returns a list of q_l values in the same order as the list of the atoms in the system.

prepare_pickle()[source]

Prepare the system for pickling and create a picklable system

Parameters:None
Returns:psys
Return type:picklable system object

Notes

This function prepares the system object for pickling. From a user perspective, the to_pickle() method should be used directly.

See also

to_file()

read_inputfile(filename, format='lammps-dump', frame=-1, compressed=False, customkeys=None)[source]

Read input file that contains the information of system configuration.

Parameters:
  • filename (string) – name of the input file.
  • format ({'lammps-dump', 'poscar'}) – format of the input file
  • compressed (bool, optional) – If True, force to read a gz compressed format, default False.
  • frame (int) –

    If the trajectory contains more than one time step, the slice can be specified using the frame option.

    Note

    works only with lammps-dump format.

  • customkeys (list) – A list containing names of headers of extra data that needs to be read in from the input file.
Returns:

Return type:

None

Notes

format keyword specifies the format of the input file. Currently only a lammps-dump and poscar files are supported. This function uses the traj_process() module to process a file which is then assigned to system.

compressed keyword is not required if a file ends with .gz extension, it is automatically treated as a compressed file.

frame keyword allows to read in a particular slice from a long trajectory. If all slices need to analysed, this is a very inefficient way. For handling multiple time slices, the traj_process() module offers a better set of tools.

Triclinic simulation boxes can also be read in for lammps-dump. No special keyword is necessary.

If custom_keys are provided, this extra information is read in from input files if available. This information is not passed to the C++ instance of atom, and is stored as a dictionary. It can be accessed directly as atom.custom[‘customval’]

reset_neighbors()[source]

Reset the neighbors of all atoms in the system.

Parameters:None
Returns:
Return type:None

Notes

It is used automatically when neighbors are recalculated.

set_atom(atom)[source]

Return the atom to its original location after modification.

Parameters:atom (Atom) – atom to be replaced
Returns:
Return type:None

Notes

For example, an Atom at location i in the list of all atoms in System can be queried by, atom = System.get_atom(i), then any kind of modification, for example, the position of the Atom can done by, atom.pos = [2.3, 4.5, 4.5]. After modification, the Atom can be set back to its position in System by set_atom().

Although the complete list of atoms can be accessed or set using atoms = sys.atoms, get_atom and set_atom functions should be used for accessing individual atoms. If an atom already exists at that index in the list, it will be overwritten and will lead to loss of information.

to_file(outfile, format='lammps-dump', custom=None, compressed=False)[source]

Save the system instance to a trajectory file.

Parameters:
  • outfile (string) – name of the output file
  • format (string, optional) – format of the output file, default lammps-dump Currently only lammps-dump format is supported.
  • custom (list of strings, optional) – a list of extra atom wise values to be written in the output file.
  • compressed (bool, optional) – If true, the output is written as a compressed file.
Returns:

Return type:

None

Notes

to_pickle(file)[source]

Save a system to pickle file

Parameters:file (string) – name of output file
Returns:
Return type:None

Notes

This function can be used to save a System object directly to file. This retains all the calculated quantities of the system, including the atoms and their properties. This can be useful to restart the calculation. The function uses numpy.save method to save the information. Hence pickling between different versions of python could lead to issues.

Warning

Pickling between different versions of numpy or python could be incompatible. Pickling is not secure. You should only unpickle objects that you trust.

pyscal.core.test()[source]

A simple function to test if the module works

Parameters:None
Returns:works – True if the module works and could create a System and Atom object False otherwise.
Return type:bool
class pyscal.catom.Atom

Bases: pybind11_builtins.pybind11_object

Class to store atom details.

Parameters:
  • pos (list of floats of length 3) – position of the Atom, default [0,0,0]
  • id (int) – id of the Atom, default 0
  • type (int) – type of the Atom, default 1

Notes

A pybind11 class for holding the properties of a single atom. Various properties of the atom can be accessed through the attributes and member functions which are described below in detail. Atoms can be created individually or directly by reading a file. Check the examples for more details on how atoms are created. For creating atoms directly from an input file check the documentation of System class.

Although an Atom object can be created independently, Atom should be thought of inherently as members of the System class. All the properties that define an atom are relative to the parent class. System has a list of all atoms. All the properties of an atom, hence should be calculated through System.

Examples

>>> #method 1 - individually
>>> atom = Atom()
>>> #now set positions of the atoms
>>> atom.pos = [23.0, 45.2, 34.2]
>>> #now set id
>>> atom.id = 23
>>> #now set type
>>> atom.type = 1
>>> #Setting through constructor
>>> atom = Atom([23.0, 45.2, 34.2], 23, 1)

References

Creation of atoms.

allaq

list of floats. list of all averaged q values of the atom.

allq

list of floats. list of all q values of the atom.

angular

Float. The value of angular parameter A of an atom. The angular parameter measures the tetrahedral coordination of an atom. Meaningful values are only returned if the property is calculated using calculate_angularcriteria().

avg_angular

Float. The average angular parameter value. Not used currently.

avg_connection

float. Value of averaged s_ij which is used for identification of solid atoms. s_ij is defined by

\[s_{ij} = \sum_{m=-6}^6 q_{6m}(i) q_{6m}^*(i)\]
avg_disorder

Float. The value of averaged disorder parameter.

avg_volume

float. Averaged version of the Voronoi volume which is calculated as an average over itself and its neighbors. Only calculated when the find_neighbors() using the method=’voronoi’ option is used.

bonds

Int. The number of solid bonds of an atom.

chiparams

Float. The value of chiparameter of an atom. The return value is a vector of length 8. Meaningful values are only returned if chi params are calculated using calculate_chiparams().

cluster

int. identification number of the cluster that the atom belongs to.

condition

int. condition that specifies if an atom is included in the clustering algorithm or not. Only atoms with the value of condition=1 will be used for clustering in cluster_atoms().

coordination

int. coordination number of the atom. Coordination will only be updated after neighbors are calculated using find_neighbors().

custom

dict. dictionary specfying custom values for an atom. The module only stores the id, type and position of the atom. If any extra values need to be stored, they can be stored in custom using atom.custom = {“velocity”:12}. read_inputfile() can also read in extra atom information. By default, custom values are treated as string.

disorder

Float. The value of disorder parameter.

edge_lengths

list of floats. For each face, this vector contains the lengths of edges that make up the Voronoi polyhedra of the atom. Only calculated when the find_neighbors() using the method=’voronoi’ option is used.

face_perimeters

list of floats. List consisting of the perimeters of each Voronoi face of an atom. Only calculated when the find_neighbors() using the method=’voronoi’ option is used.

face_vertices

list of floats. A list of the number of vertices shared between an atom and its neighbors. Only calculated when the find_neighbors() using the method=’voronoi’ option is used.

get_q()

Calculate the steinhardt parameter q_l value.

Parameters:
  • q (int or list of ints) – number of the required q_l - from 2-12
  • averaged (bool, optional) – If True, return the averaged q values, If False, return the non averaged ones default False
Returns:

q_l – the value(s) of the queried Steinhardt parameter(s).

Return type:

float or list of floats

Notes

Please check this link for more details about Steinhardts parameters and the averaged versions.

Meaningful values are only returned if calculate_q() is used.

id

int. Id of the atom.

largest_cluster

bool. True if the atom belongs to the largest cluster, False otherwise. Largest cluster is only identified after using the cluster_atoms() function.

loc

int. indicates the position of the atom in the list of all atoms.

neighbor_weights

List of floats. Used to weight the contribution of each neighbor atom towards the value of Steinhardt’s parameters. By default, each atom has a weight of 1 each. However, if find_neighbors() is used with method=’voronoi’, each neighbor gets a weight proportional to the area shared between the neighboring atom and host atom.

neighbors

List of ints. List of neighbors of the atom. The list contains indices of neighbor atoms which indicate their position in the list of all atoms.

pos

List of floats of the type [x, y, z], default [0, 0, 0]. Position of the atom.

set_q()

Set the value of steinhardt parameter q_l.

Parameters:
  • q (int or list of ints) – number of the required q_l - from 2-12
  • val (float or list of floats) – value(s) of Steinhardt parameter(s).
  • averaged (bool, optional) – If True, return the averaged q values, If False, return the non averaged ones default False
Returns:

Return type:

None

solid

bool. True if the atom is solid, False otherwise. Solid atoms are only identified after using the find_solids() function.

sro

Float. The value of short range order parameter.

structure

int. Indicates the structure of atom. Not used currently.

surface

bool. True if the atom has at least one liquid neighbor, False otherwise. Surface atoms are only identified after using the find_solids() function.

type

int. int specifying type of the atom.

vertex_numbers

list of floats. For each Voronoi face of the atom, this values includes a List of vertices that constitute the face. Only calculated when the find_neighbors() using the method=’voronoi’ option is used.

vertex_vectors

list of floats. A list of positions of each vertex of the Voronoi polyhedra of the atom. Only calculated when the find_neighbors() using the method=’voronoi’ option is used.

volume

float. Voronoi volume of the atom. The Voronoi volume is only calculated if neighbors are found using the find_neighbors() using the method=’voronoi’ option.

vorovector

list of ints. A vector of the form (n3, n4, n5, n6) where n3 is the number of faces with 3 vertices, n4 is the number of faces with 4 vertices and so on. This can be used to identify structures [1][2]. Vorovector is calculated if the calculate_vorovector() method is used.

References

[1]Finney, JL, Proc. Royal Soc. Lond. A 319, 1970
[2]Tanemura, M, Hiwatari, Y, Matsuda, H,Ogawa, T, Ogita, N, Ueda, A. Prog. Theor. Phys. 58, 1977

pyscal.crystal_structures module

pyscal module for creating crystal structures.

pyscal.crystal_structures.make_crystal(structure, lattice_constant=1.0, repetitions=None, ca_ratio=1.633, noise=0)[source]

Create a basic crystal structure and return it as a list of Atom objects and box dimensions.

Parameters:
  • structure ({'bcc', 'fcc', 'hcp', 'diamond' or 'l12'}) – type of the crystal structure
  • lattice_constant (float, optional) – lattice constant of the crystal structure, default 1
  • repetitions (list of ints of len 3, optional) – of type [nx, ny, nz], repetions of the unit cell in x, y and z directions. default [1, 1, 1].
  • ca_ratio (float, optional) – ratio of c/a for hcp structures, default 1.633
  • noise (float, optional) – If provided add normally distributed noise with standard deviation noise to the atomic positions.
Returns:

  • atoms (list of Atom objects) – list of all atoms as created by user input
  • box (list of list of floats) – list of the type [[xlow, xhigh], [ylow, yhigh], [zlow, zhigh]] where each of them are the lower and upper limits of the simulation box in x, y and z directions respectively.

Examples

>>> atoms, box = make_crystal('bcc', lattice_constant=3.48, repetitions=[2,2,2])
>>> sys = System()
>>> sys.assign_atoms(atoms, box)

pyscal.traj_process module

pyscal module containing methods for processing of a trajectory. Methods for reading of input files formats, writing of output files etc are provided in this module.

pyscal.traj_process.read_lammps_dump(infile, compressed=False, check_triclinic=False, box_vectors=False, customkeys=None)[source]

Function to read a lammps dump file format - single time slice.

Parameters:
  • infile (string) – name of the input file
  • compressed (bool, optional) – force to read a gz zipped file. If the filename ends with .gz, use of this keyword is not necessary. Default True.
  • check_triclinic (bool, optional) – If true check if the sim box is triclinic. Default False.
  • box_vectors (bool, optional) – If true, return the full box vectors along with boxdims which gives upper and lower bounds. default False.
  • customkeys (list of strings, optional) – A list of extra keywords to read from trajectory file.
Returns:

  • atoms (list of Atom objects) – list of all atoms as created by user input
  • boxdims (list of list of floats) – The dimensions of the box. This is of the form [[xlo, xhi],[ylo, yhi],[zlo, zhi]] where lo and hi are the upper and lower bounds of the simulation box along each axes. For triclinic boxes, this is scaled to [0, scalar length of the vector].
  • box (list of list of floats) – list of the type [[x1, x2, x3], [y1, y2, y3], [zz1, z2, z3]] which are the box vectors. Only returned if box_vectors is set to True.
  • triclinic (bool) – True if the box is triclinic. Only returned if check_triclinic is set to True
  • .. note:: – Values are always returned in the order atoms, boxdims, box, triclinic if all return keywords are selected. For example, ff check_triclinic is not selected, the return values would still preserve the order and fall back to atoms, boxdims, box.

Notes

Read a lammps-dump style snapshot that can have variable headers, reads in type and so on. Zipped files which end with a .gz can also be read automatically. However, if the file does not end with a .gz extension, keyword compressed = True can also be used.

Examples

>>> atoms, box = read_lammps_dump('conf.dump')
>>> atoms, box = read_lammps_dump('conf.dump.gz')
>>> atoms, box = read_lammps_dump('conf.d', compressed=True)
pyscal.traj_process.read_poscar(infile, compressed=False)[source]

Function to read a POSCAR format.

Parameters:
  • infile (string) – name of the input file
  • compressed (bool, optional) – force to read a gz zipped file. If the filename ends with .gz, use of this keyword is not necessary, Default False
Returns:

  • atoms (list of Atom objects) – list of all atoms as created by user input
  • box (list of list of floats) – list of the type [[xlow, xhigh], [ylow, yhigh], [zlow, zhigh]] where each of them are the lower and upper limits of the simulation box in x, y and z directions respectively.

Examples

>>> atoms, box = read_poscar('POSCAR')
>>> atoms, box = read_poscar('POSCAR.gz')
>>> atoms, box = read_poscar('POSCAR.dat', compressed=True)
pyscal.traj_process.split_traj_lammps_dump(infile, compressed=False)[source]

Read in a LAMMPS dump trajectory file and convert it to individual time slices.

Parameters:
  • filename (string) – name of input file
  • compressed (bool, optional) – force to read a gz zipped file. If the filename ends with .gz, use of this keyword is not necessary, Default False.
Returns:

snaps – a list of filenames which contain individual frames from the main trajectory.

Return type:

list of strings

pyscal.traj_process.split_trajectory(infile, format='lammps-dump', compressed=False)[source]

Read in a trajectory file and convert it to individual time slices.

Parameters:
  • filename (string) – name of input file
  • format (format of the input file) – only lammps-dump is supported now.
  • compressed (bool, optional) – force to read a gz zipped file. If the filename ends with .gz, use of this keyword is not necessary.
Returns:

snaps – a list of filenames which contain individual frames from the main trajectory.

Return type:

list of strings

Notes

This is a wrapper function around split_traj_lammps_dump function.

pyscal.traj_process.write_structure(sys, outfile, format='lammps-dump', compressed=False, customkey=None, customvals=None)[source]

Write the state of the system to a trajectory file.

Parameters:
  • sys (System object) – the system object to be written out
  • outfile (string) – name of the output file
  • format (string, optional) – the format of the output file, as of now only lammps-dump format is supported.
  • compressed (bool, default false) – write a .gz format
  • customkey (string or list of strings, optional) – If specified, it adds this custom column to the dump file. Default None.
  • customvals (list or list of lists, optional) – If customkey is specified, customvals take an array of the same length as number of atoms, which contains the values to be written out.
Returns:

Return type:

None

pyscal.pickle module