Supercell
Supercell:
The supercell module defines the supercell
class for building supercells from crystal.Crystal
classes.
Supercell class
Class to store supercells of crystals. A supercell is a lattice model of a crystal, with periodically repeating unit cells. In that framework we can
add/remove/substitute atoms
find the transformation map between two different representations of the same supercell
output POSCAR format (possibly other formats?)
- class supercell.ClusterSupercell(crys, superlatt, spectator=())[source]
A class that defines a Supercell of a crystal for purposes of evaluating a cluster expansion. We intend to use this with Monte Carlo sampling.
Takes in a crystal, a supercell (3x3 integer matrix). We can identify sites as spectator sites (that is, they can have different occupancies, but we do not intend for those to change during a Monte Carlo simulation.
- Supercell_occ(sup, chemmapping=None)[source]
Takes in a Supercell object (that is assumed to be consistent with this supercell!) and produces the corresponding occupancy vectors for this supercell, using a specific chemical mapping (described below).
In a Supercell object, each site has a “native” chemistry; moreover, those sites may be occupied by everything from a vacancy (-1) to a different chemical element (>=0). We need to define how that happens, since ClusterSupercells only have occupancies of 0 or 1.
chemmapping
is a dictionary of dictionaries.chemmapping[csite][cocc]
= 0 or 1 to dictate what the occupancy for a site should be if chemistry of typecocc
occurs on a site with native chemistrycsite
.If the
chemmapping
is None, we use a default “defect” occupancy mapping; namely, ifcsite
!= Interstitial, then we use 0 whencsite==cocc
, 1 otherwise; and ifcsite
== Interstitial, we use 0 whencsite==-1
, 1 otherwise. SeeSupercell.defect_chemmapping()
- Parameters:
sup – Supercell object, with appropriate chemical occupancies
chemmapping – mapping of chemical identities to occupancy variables.
- Return mocc:
mobile occupancy vector
- Return socc:
spectator occupancy vector
- __init__(crys, superlatt, spectator=())[source]
Initialize our supercell to an empty supercell.
- Parameters:
crys – crystal object
superlatt – 3x3 integer matrix
spectator – list of indices of chemistries that will be considered “spectators”
- __weakref__
list of weak references to the object
- ciR(ind, mobile=True)[source]
Return the chem/index and lattice vector for a specific indexed position
- Parameters:
ind – index of site
mobile – True if mobile; false if spectator
- Return ci:
(c, i) index
- Return R:
lattice vector
- clusterevaluator(socc, clusters, values)[source]
Construct the information necessary for an (efficient) cluster evaluator, for a given spectator occupancy, set of clusters, and values for those clusters.
- Parameters:
socc – spectator occupancy vector (0 or 1 only)
clusters – list of lists (or sets) of Cluster objects
values – vector of values for the clusters; if it is longer than the list of clusters by one, the last values is assumed to be the constant value.
- Return siteinteract:
list of lists of interactions for each site
- Return interact:
list of interaction values
- evalTScluster(mocc, socc, TSclusters, initial, final, dx)[source]
Evaluate a TS cluster expansion for a given mobile occupancy and spectator occupancy. Indexing corresponds to
mobilepos
andspecpos
. The clusters are input as a list of lists of clusters (where it is assumed that all of the clusters in a given sublist have equal coefficients (i.e., grouped by symmetry). We return a vector of length Nclusters; each entry is the number of times each cluster appears. This can then be dotted into the vector of values to get the cluster expansion value. This is evaluated for the transition where the mobile species atinitial
jumps to the position atfinal
. Requires mocc[initial] == 1 and mocc[final] == 0- Parameters:
mocc – mobile occupancy vector (0 or 1 only)
socc – spectator occupancy vector (0 or 1 only)
TSclusters – list of lists (or sets) of (transition state) Cluster objects
initial – index of initial state
final – index of final state
dx – displacement vector (necessary to deal with PBC)
- Returns:
clustercount: count of how many of each cluster is in this supercell.
- evalcluster(mocc, socc, clusters)[source]
Evaluate a cluster expansion for a given mobile occupancy and spectator occupancy. Indexing corresponds to
mobilepos
andspecpos
. The clusters are input as a list of lists of clusters (where it is assumed that all of the clusters in a given sublist have equal coefficients (i.e., grouped by symmetry). We return a vector of length Nclusters + 1; each entry is the number of times each cluster appears, and the last entry is equal to the size of the supercell (which would be an “empty” cluster). This can then be dotted into the vector of values to get the cluster expansion value.- Parameters:
mocc – mobile occupancy vector (0 or 1 only)
socc – spectator occupancy vector (0 or 1 only)
clusters – list of lists (or sets) of Cluster objects
- Returns:
clustercount: count of how many of each cluster is in this supercell.
- expandcluster_matrices(socc, clusters)[source]
Expand a cluster expansion for a given spectator occupancy into matrices of indices. This is designed for rapid evaluation for a fixed spectator occupancy. The clusters are input as a list of lists of clusters (i.e., grouped by symmetry). We return a list of lists of integer matrices of indices. This can then be used to efficiently evaluate cluster counts for a given mobile occupancy. The given row of indices must be all 1 in order to increment the particular cluster count.
- Parameters:
socc – spectator occupancy vector (0 or 1 only)
clusters – list of lists (or sets) of Cluster objects
- Returns:
clustermatrices: list of lists of matrices of indices
- expiqx()[source]
Construct a Fourier transform matrix for our mobile species
- Returns:
exp( I q[i] . x[j]) as a matrix
- Returns:
gamma_index of the gamma point (0)
- index(R, ci)[source]
Return the index that corresponds to a position specified by a lattice vector R and a chem/index (c,i). We also need to specify if its in the spectator basis or mobile basis.
- Parameters:
R – lattice vector
ci – (c, i) index
- Return ind:
index of site in our position
- Return mobile:
boolean; True if mobile, False if spectator
- indexpos(pos, threshold=1.0, CARTESIAN=False)[source]
Return the index that corresponds to the position closest to pos in the supercell. Done in direct coordinates of the supercell, using periodic boundary conditions.
- Parameters:
pos – 3-vector
threshold – (optional) minimum squared “distance” in supercell for a match; default=1.
- Return index:
index of closest position
- Return mobile:
boolean; True if mobile, False if spectator
- jumpnetworkevaluator(socc, clusters, values, chem, jumpnetwork, KRAvalues=0, TSclusters=(), TSvalues=(), siteinteract=(), interact=())[source]
Build out an efficient jump network evaluator. Similar inputs to
clusterevaluator
, with the addition of a jumpnetwork and energies. The interactions can be appended onto existing interactions, if included. The information about all of the transitions is: initial state, final state, delta x.- Parameters:
socc – spectator occupancy vector (0 or 1 only)
clusters – list of lists (or sets) of Cluster objects
values – vector of values for the clusters; if it is longer than the list of clusters by one, the last values is assumed to be the constant value.
chem – index of species that transitions
jumpnetwork – list of lists of jumps; each is ((i, j), dx) where
i
andj
are unit cell indices for specieschem
KRAvalues – list of “KRA” values for barriers (relative to average energy of endpoints); if
TSclusters
are used, choosing 0 is more straightforward.TSclusters – (optional) list of transition state cluster expansion terms; this is always added on to KRAvalues (thus using 0 is recommended if TSclusters are also used)
TSvalues – (optional) values for TS cluster expansion entries
siteinteract – (optional) list of lists of interactions for each site, to append
interact – (optional) list of interaction values, to append
- Return siteinteract:
list of lists of interactions for each site
- Return interact:
list of interaction values
- Return jumps:
list of ((initial, final), dx)
- Return interactrange:
range of indices to count in interact for each jump; for the nth jump, sum over interactrange[n-1]:interactrange[n]; interactrange[-1] == range for energy
- jumpnetworkevaluator_vacancy(socc, clusters, values, chem, jumpnetwork, KRAvalues=0, TSclusters=(), TSvalues=(), siteinteract=(), interact=())[source]
Build out an efficient jump network evaluator for a vacancy. Similar inputs to
jumpnetworkevaluator
. This is designed for a “stationary” vacancy, where we’ll just look at its jumps in the supercell.- Parameters:
socc – spectator occupancy vector (0 or 1 only)
clusters – list of lists (or sets) of Cluster objects
values – vector of values for the clusters; if it is longer than the list of clusters by one, the last values is assumed to be the constant value.
chem – index of species that transitions
jumpnetwork – list of lists of jumps; each is ((i, j), dx) where
i
andj
are unit cell indices for specieschem
KRAvalues – list of “KRA” values for barriers (relative to average energy of endpoints); if
TSclusters
are used, choosing 0 is more straightforward.TSclusters – (optional) list of transition state cluster expansion terms; this is always added on to KRAvalues (thus using 0 is recommended if TSclusters are also used)
TSvalues – (optional) values for TS cluster expansion entries
siteinteract – (optional) list of lists of interactions for each site, to append
interact – (optional) list of interaction values, to append
- Return siteinteract:
list of lists of interactions for each site
- Return interact:
list of interaction values
- Return jumps:
list of ((initial, final), dx)
- Return interactrange:
range of indices to count in interact for each jump; for the nth jump, sum over interactrange[n-1]:interactrange[n]; interactrange[-1] == range for energy
- class supercell.Supercell(crys, superlatt, interstitial=(), Nsolute=0, empty=False, NOSYM=False)[source]
A class that defines a Supercell of a crystal.
Takes in a crystal, a supercell (3x3 integer matrix). We can identify sites as interstitial sites, and specify if we’ll have solutes.
- KrogerVink()[source]
Attempt to make a “simple” string based on the defectindices, using Kroger-Vink notation. That is, we identify: vacancies, antisites, and interstitial sites, and return a string. NOTE: there is no relative charges, so this is a pseudo-KV notation.
- Return KV:
string representation
- POSCAR(name=None, stoichiometry=True)[source]
Return a VASP-style POSCAR, returned as a string.
- Parameters:
name – (optional) name to use for first list
stoichiometry – (optional) if True, append stoichiometry to name
- Return POSCAR:
string
- POSCAR_occ(POSCAR_str, EMPTY_SUPER=True, disp_threshold=-1, latt_threshold=-1)[source]
Takes in a POSCAR_str, and sets the occupancy of the supercell accordingly. Note: if we want to read a POSCAR from a file instead, the proper usage is
with open(POSCAR_filename, "r") as f: sup.from_POSCAR(f.read())
Warning: there is only minimal validity checking; it makes a strong assumption that a reasonable POSCAR file is being given, and that the sites should correspond to the supercell object. Should that not be the case, the behavior is unspecified.
- Parameters:
POSCAR_str – string form of a POSCAR
EMPTY_SUPER – initialize supercell by emptying it first (default=True)
disp_threshold – threshold for difference in position to raise an error, in unit cell coordinates. For a negative value, we compute a value that is ~0.1A. If you do not want this check to be meaningful, choose a value >sqrt(3)=1.733
latt_threshold – threshold for supercell lattice check, in units of strain. For a negative values do not check.
- Returns:
name from the POSCAR
- __eq__(other)[source]
Return True if two supercells are equal; this means they should have the same occupancy. and the same ordering
- Parameters:
other – supercell for comparison
- Returns:
True if same crystal, supercell, occupancy, and ordering; False otherwise
- __getitem__(key)[source]
Index into supercell
- Parameters:
key – index (either an int, a slice, or a position)
- Returns:
chemical occupation at that point
- __hash__ = None
- __imul__(other)[source]
Multiply by a GroupOp, in place.
- Parameters:
other – must be a GroupOp (and should be a GroupOp of the supercell!)
- Returns:
self
- __init__(crys, superlatt, interstitial=(), Nsolute=0, empty=False, NOSYM=False)[source]
Initialize our supercell to an empty supercell.
- Parameters:
crys – crystal object
superlatt – 3x3 integer matrix
interstitial – (optional) list/tuple of indices that correspond to interstitial sites
Nsolute – (optional) number of substitutional solute elements to consider; default=0
empty – (optional) designed to allow “copy” to work–skips all derived info
NOSYM – (optional) does not do symmetry analysis (intended ONLY for testing purposes)
- __mul__(other)[source]
Multiply by a GroupOp; returns a new supercell (constructed via copy).
- Parameters:
other – must be a GroupOp (and should be a GroupOp of the supercell!)
- Returns:
rotated supercell
- __rmul__(other)[source]
Multiply by a GroupOp; returns a new supercell (constructed via copy).
- Parameters:
other – must be a GroupOp (and should be a GroupOp of the supercell!)
- Returns:
rotated supercell
- __setitem__(key, value)[source]
Set specific composition for site; uses same indexing as __getitem__
- Parameters:
key – index (either an int, a slice, or a position)
value – chemical occupation at that point
- __weakref__
list of weak references to the object
- copy()[source]
Make a copy of the supercell; initializes, then copies over
__copyattr__
and__eqattr__
.- Returns:
new supercell object, copy of the original
- defectindices()[source]
Return a dictionary that corresponds to the “defect” content of the supercell.
- Return defects:
dictionary, keyed by defect type, with a set of indices of corresponding defects
- definesolute(c, chemistry)[source]
Set the name of the chemistry of chemical index c. Only works for substitutional solutes.
- Parameters:
c – index
chemistry – string
- equivalencemap(other)[source]
Given the superlatt
other
we want to find a group operation that transformsself
into other. This is a GroupOp along with an index mapping of chemorder. The index mapping is to get the occposlist to match up:(g*self).occposlist()[c][mapping[c][i]] == other.occposlist()[c][i]
(We can write a similar expression using chemorder, since chemorder indexes into pos). We’re going to return both g and mapping.Remember:
g
does not change the presentation ordering;mapping
is necessary for full equivalence. If no such equivalence, returnNone,None
.- Parameters:
other – Supercell
- Return g:
GroupOp to transform sites from
self
toother
- Return mapping:
list of maps, such that (g*self).chemorder[c][mapping[c][i]] == other.chemorder[c][i]
- fillperiodic(ci, Wyckoff=True)[source]
Occupies all of the (Wyckoff) sites corresponding to chemical index with the appropriate chemistry.
- Parameters:
ci – tuple of (chem, index) in crystal
Wyckoff – (optional) if False, only occupy the specific tuple, but still periodically
- Return self:
- gengroup()[source]
Generate the group operations internal to the supercell
- Return G:
set of GroupOps
- index(pos, threshold=1.0)[source]
Return the index that corresponds to the position closest to pos in the supercell. Done in direct coordinates of the supercell, using periodic boundary conditions.
- Parameters:
pos – 3-vector
threshold – (optional) minimum squared “distance” in supercell for a match; default=1.
- Return index:
index of closest position
- makesites()[source]
Generate the array corresponding to the sites; the indexing is based on the translations and the atomindices in crys. These may not all be filled when the supercell is finished.
- Return pos:
array [N*size, 3] of supercell positions in direct coordinates
- static maketrans(superlatt)[source]
Takes in a supercell matrix, and returns a list of all translations of the unit cell that remain inside the supercell
- Parameters:
superlatt – 3x3 integer matrix
- Return size:
integer, corresponding to number of unit cells
- Return invsuper:
integer matrix inverse of supercell (needs to be divided by size)
- Return translist:
list of integer vectors (to be divided by
size
) corresponding to unit cell positions- Return transdict:
dictionary of tuples and their corresponding index (inverse of trans)
- occposlist()[source]
Returns a list of lists of occupied positions, in (chem)order.
- Return occposlist:
list of lists of supercell coord. positions
- reorder(mapping)[source]
Reorder (in place) the occupied sites. Does not change the occupancies, only the ordering for “presentation”.
- Parameters:
mapping – list of maps; will make newchemorder[c][i] = chemorder[c][mapping[c][i]]
- Return self:
If mapping is not a proper permutation, raises ValueError.