Materials

Entry

class qmpy.Entry(*args, **kwargs)[source]

Base class for a database entry.

The core model for typical database entries. An Entry model represents an input structure to the database, and can be created from any input file. The Entry also ties together all of the associated qmpy.Structure, qmpy.Calculation, qmpy.Reference, qmpy.FormationEnergies, and other associated databas entries.

Relationships:
Calculation via calculation_set
DOS via dos_set
Entry via duplicate_of
Entry via duplicates
Element via element_set
FormationEnergy via formationenergy_set
Job via job_set
MetaData via meta_data
Project via project_set
Prototype via prototype
Species via species_set
Structure via structure_set
Task via task_set
Reference via reference
Composition via composition
Attributes:
id: Primary key (auto-incrementing int)
label: An identifying name for the structure. e.g. icsd-1001 or A3
exception DoesNotExist
exception MultipleObjectsReturned
property c

Dictionary of label:Calculation pairs.

property calculations

Dictionary of label:Calculation pairs.

property chg

Attempts to load the charge density of the final calculation, if it is done. If not, returns False.

static create(source, keywords=[], projects=[], prototype=None, **kwargs)[source]

Attempts to create an Entry object from a provided input file.

Processed in the following way:

  1. If an Entry exists at the specified path, returns that Entry.

  2. Create an Entry, and assign all fundamental attributes. (natoms, ntypes, input, path, elements, keywords, projects).

  3. If the input file is a CIF, and because CIF files have additional composition and reference information, if that file format is found, an additional test is performed to check that the reported composition matches the composition of the resulting structure. The reference for the work is also created and assigned to the entry.

  4. Attempt to identify another entry that this is either exactly equivalent to, or a defect cell of.

Keywords:

keywords: list of keywords to associate with the entry. projects: list of project names to associate with the entry.

do(module, *args, **kwargs)[source]

Looks for a computing script matching the first argument, and attempts to run it with itself as the first argument. Sends args and kwargs to the script. Should return a Calculation object, or list of Calculation objects.

Examples:

>>> e = Entry.objects.get(id=123)
>>> e.do('relaxation')
<Calculation: 523 @ relaxation settings>
property elements

List of Elements

property energy

If the structure has been relaxed, returns the formation energy of the final relaxed structure. Otherwise, returns None.

property errors

List of errors encountered in all calculations.

property hold_objects

Return list of holds (MetaData objects of type hold)

property holds

A note indicating a reason the entry should not be calculated

property html

HTML formatted name

property keyword_objects

Return list of keywords (MetaData objects of type keyword)

property keywords

Descriptive keyword for looking up entries

property latex

LaTeX formatted name

property mass

Return the mass of the entry, normalized to per atom.

move(path)[source]

Moves all calculation files to the specified path.

property name

Unformatted name

property projects

List of Projects

property red_comp

Composition dictionary, in reduced form.

reset()[source]

Deletes all calculations, removes all associated structures - returns the entry to a pristine state.

save(*args, **kwargs)[source]

Saves the Entry, as well as all associated objects.

property space

Return the set of elements in the input structure.

Examples:

>>> e = Entry.create("fe2o3/POSCAR") # an input containing Fe2O3
>>> e.space
set(["Fe", "O"])
property spec_comp

Composition dictionary, using species (element + oxidation state) instead of just the elements.

property species

List of Species

property total_energy

If the structure has been relaxed, returns the formation energy of the final relaxed structure. Otherwise, returns None.

property unit_comp

Composition dictionary, normalized to 1 atom.

visualize(structure='source')[source]

Attempts to open the input structure for visualization using VESTA

property volume

If the entry has gone through relaxation, returns the relaxed volume. Otherwise, returns the input volume.

Structure

class qmpy.Structure(*args, **kwargs)[source]

Structure model. Principal attributes are a lattice and basis set.

Relationships:
Entry via entry
Atom via atom_set: Atoms in the structure. More commonly
handled by the managed attributed atoms.
Calculation via calculated. Calculation objects
that the structure is an output from.
Calculation via calculation_set. Returns calculation
objects that the structure is an input to.
Composition via composition.
Element via element_set
Spacegroup via spacegroup
Species via species_set
Prototype via prototype. If the structure belongs to a
prototypical structure, it is referred to here.
Reference Original literature reference.
MetaData via meta_data
Attributes:
Identification
id
label: key in the Entry.structures dictionary.
natoms: Number of atoms.
nsites: Number of sites.
ntypes: Number of elements.
measured: Experimentally measured structure?
source: Name for source.

Lattice
x1, x2, x3
y1, y2, y3
z1, z2, z3: Lattice vectors of the cell. Accessed via cell.
volume
volume_pa

Calculated properties
delta_e: Formation energy (eV/atom)
meta_stability: Distance from the convex hull (eV/atom)
energy: Total DFT energy (eV/FU)
energy_pa: Total DFT energy (eV/atom)
magmom: Total magnetic moment (&Mu;<sub>b</sub>)
magmom_pa: Magnetic moment per atom.
sxx, sxy, syy
syx, szx, szz: Stresses on the cell. Accessed via stresses.

Examples:

>>> s = io.read(INSTALL_PATH+'/io/files/POSCAR_FCC')
>>> s.atoms
>>> s.cell
>>> s.magmoms
>>> s.forces
>>> s.stresses
exception DoesNotExist
exception MultipleObjectsReturned
add_atom(atom, tol=0.01)[source]

Adds atom to the structure if it isn’t already contained.

property atom_types

List of atomic symbols, length equal to number of atoms.

property atomic_numbers

List of atomic numbers, length equal to number of atoms.

property atoms

List of Atoms in the structure.

property cartesian_coords

Return atomic positions in cartesian coordinates.

property cell

Lattice vectors, 3x3 numpy.ndarray.

property comment_objects

Return list of comments (MetaData objects of type comment)

property comp

Composition dictionary.

compare(other, tol=0.01, atom_tol=10, volume=False, allow_distortions=False, check_spacegroup=False, wildcard=None)[source]

Credit to K. Michel for the algorithm.

  1. Check that elements are the same in both structures

  2. Convert both structures to primitive form

  3. Check that the total number of atoms in primitive cells are the same

4. Check that the number of atoms of each element are the same in primitive cells

4b. Check that the spacegroup is the same.

  1. If needed check that the primitive cell volumes are the same

6. Convert both primitive cells to reduced form There is one issue here - the reduce cell could be type I (all angles acute) or type II (all angles obtuse) and a slight difference in the initial cells could cause two structures to reduce to different types. So at this step, if the angles are not correct, the second cell is transformed as [[-1, 0, 0], [0, -1, 0], [0, 0, 1]].

7. Check that the cell internal angles are the same in both reduced cells.

8. Check that the ratios of reduced cell basis lengths are the same. ie a1/b1 = a2/b2, a1/c1 = a2/c2, and b1/c1 = b2/c2 where a1, b1, c1, are the lengths of basis vectors in cell 1 and a2, b2, c2 are the lengths of cell 2.

9. Get the lattice symmetry of the reduced cell 2 (this is a list of all rotations that leave the lattice itself unchanged). In turn, apply all rotations to reduced cell 2 and for each search for a vector that overlaps rotated cell positions with positions in reduced cell 1. If a rotation + translation overlaps reduced cells, then they are the same.

MODIFICATIONS: Only need one “base” atom from the first structure Get the distance from the origin for every atom first

Arguments:

other: Another Structure.

Keyword Arguments:

tol: Percent deviation in lattice parameters and angles.

Not Implemented Yet: wildcard: Elements of the specified type can match with any atom.

property coords

numpy.ndarray of atom coordinates.

copy()[source]

Create a complete copy of the structure, with any primary keys removed, so it is not associated with the original.

static create(cell, atoms=[], **kwargs)[source]

Creates a new Structure.

Arguments:

cell: 3x3 lattice vector array

Keyword Arguments:

atoms: List of Atom creation arguments. Can be a list of [element, coord], or a list of [element, coord, kwargs].

Examples:

>>> a = 3.54
>>> cell = [[a,0,0],[0,a,0],[0,0,a]]
>>> atoms = [('Cu', [0,0,0]),
             ('Cu', [0.5, 0.5, 0.5])]
>>> s = Structure.create(cell, atoms)
>>> atoms = [('Cu', [0,0,0], {'magmom':3}),
        ('Cu', [0.5, 0.5, 0.5], {'magmom':-3})]
>>> s2 = Structure.create(cell, atoms)
create_vacuum(direction, amount, in_place=True)[source]

Add vacuum along a lattice direction.

Arguments:

direction: direction to add the vacuum along. (0=x, 1=y, 2=z) amount: amount of vacuum in Angstroms.

Keyword Arguments:

in_place: apply change to current structure, or return a new one.

Examples:

>>> s = io.read(INSTALL_PATH+'/io/files/POSCAR_FCC')
>>> s.create_vacuum(2, 5)
property elements

List of Elements

find_lattice_points_by_transform(transform, tol=1e-06)[source]

Find the lattice points contained within the transformation.

find_lattice_points_within_distance(distance, tol=1e-06)[source]

Find the lattice points contained within radius distance from the origin.

find_nearest_neighbors(method='closest', tol=0.05, limit=5.0, **kwargs)[source]

Determine the nearest neighbors for all Atoms in Structure.

Calls get_nearest_neighbors() and assigns the nearest neighbor dictionary to Structure._neighbor_dict. Each atom is also given a list, nearest_neighbors that contains the nearest neighbor atoms. For atoms which have the “same” atom as a nearest neighbor across different periodic boundaries, a single atom may appear multiple times on the list.

Keyword Arguments:

limit: How far to look from each atom for nearest neighbors. Default=5.0.

tol: A tolerance which determines how much further than the closest atom a second atom can be and still be a part of the nearest neighbor shell.

Returns: None

property forces

numpy.ndarray of forces on atoms.

get_distance(atom1, atom2, limit=None, wrap_self=True)[source]

Calculate the shortest distance between two atoms.

Note

This is not as trivial a problem as it sounds. It is easy to demonstrate that for any non-cubic cell, the normal method of calculating the distance by wrapping the vector in fractional coordinates to the range (-0.5, 0.5) fails for cases near (0.5,0.5) in Type I cells and near (0.5, -0.5) for Type II.

To get the correct distance, the vector must be wrapped into the Wigner-Seitz cell.

Arguments:

atom1, atom2: (Atom, Site, int).

Keyword Arguments:
limit:

If a limit is provided, returns None if the distance is greater than the limit.

wrap_self:

If True, the distance from an atom to itself is 0, otherwise it is the distance to the shortest periodic image of itself.

get_kpoint_mesh_by_increment(kppra)[source]

DEPRECATED: Sometimes results in k-point meshes incommensurate with lattice symmetry. Use either get_TM_kpoint_mesh() or (if you have sympy installed) get_kpoint_mesh_with_sympy(kppra) instead.

get_kpoint_mesh_with_sympy(kppra)[source]

Generate the k-point mesh for a given KPPRA; requires sympy to be installed

get_lattice_network(elements=None, supercell=None, **kwargs)[source]

Constructs a networkx.Graph of lattice sites.

Keyword Arguments:
elements:

If elements is supplied, get_lattice_network will return the lattice of those elements only.

supercell:

Accepts any valid input to Structure.transform to construct a supercell, and return its lattice. Useful for finding AFM orderings for structures which have a smaller periodicity than their magnetic structure.

Returns:

A LatticeNetwork, which is a container for a lattice graph, which contains nodes which are atoms and edges which indicate nearest neighbors.

Examples:

>>> s = io.read(INSTALL_PATH+'/io/files/fe3o4.cif')
>>> sl = s.get_lattice_network(elements=['Fe'])
>>> sl.set_fraction(0.33333)
>>> sl.fraction
0.3333333333333333
>>> sl.run_MC()
get_sites(tol=0.1)[source]

From self.atoms, creates a list of Sites. Atoms which are closer than tol from one another are considered on the same site.

get_volume()[source]

Calculates the volume from the triple product of self.cell

group_atoms_by_symmetry()[source]

Sort self.atoms according to the site they occupy.

property inv

Precalculates the inverse of the lattice, for faster conversion between cartesian and direct coordinates.

is_buerger_cell(tol=1e-05)[source]

Tests whether or not the structure is a Buerger cell.

is_niggli_cell(tol=0.001)[source]

Tests whether or not the structure is a Niggli cell.

joggle_atoms(distance=0.001, in_place=True)[source]

Randomly displace all atoms in each direction by a distance up to +/- the distance keyword argument (in Angstroms).

Optional keyword arguments:
distanceRange within all internal coordinates are

displaced. Default=1e-3

in_placeIf True, returns an ndarray of the applied

translations. If False, returns (Structure, translations).

Examples:

>>> s = io.read('POSCAR')
>>> s2, trans = s.joggle_atoms(in_place=False)
>>> trans = s.joggle_atoms(1e-1)
>>> trans = s.joggle_atoms(distance=1e-1)
property keyword_objects

Return list of keywords (MetaData objects of type keyword)

property lat_param_dict

Dictionary of lattice parameters.

lat_param_string(format='screen')[source]

Generates a human friendly representation of the lattice parameters of a structure.

Keyword Arguments:

format: (‘screen’|’html’|’mathtype’)

property lat_params

Tuple of lattice parameters (a, b, c, alpha, beta, gamma).

property lp

Tuple of lattice parameters (a, b, c, alpha, beta, gamma).

property magmoms

numpy.ndarray of magnetic moments of shape (natoms,).

make_conventional(in_place=True, tol=0.001)[source]

Uses spglib to convert to the conventional cell.

Keyword Arguments:
in_place:

If True, changes the current structure. If false returns a new one

tol:

Symmetry precision for symmetry analysis

Examples:

>>> s = io.read(INSTALL_PATH+'io/files/POSCAR_FCC_prim')
>>> print len(s)
1
>>> s.make_conventional()
>>> print len(s)
4
make_perfect(in_place=True, tol=0.1)[source]

Constructs options for a ‘perfect’ lattice from the structure.

If a site is not fully occupied, but has only one atom type on it, it will be filled the rest of the way. If a site has two or more atom types on it, the higher fraction element will fill the site.

Keyword Arguments:

in_place: If False returns a new Structure, otherwise returns None

tol: maximum defect concentration.

Examples:

>>> s = io.read(INSTALL_PATH+'/io/files/partial_vac.cif')
>>> s
<Structure: Mn3.356Si4O16>
>>> s.make_perfect()
>>> s
<Structure: MnSiO4>
>>> s = io.read(INSTALL_PATH+'/io/files/partial_mix.cif')
>>> s
<Structure: Mn4.264Co3.736Si4O16>
>>> s2 = s.make_perfect(in_place=False)
>>> s2
<Structure: MnCoSiO4>
make_primitive(in_place=True, tol=0.001)[source]

Uses spglib to convert to the primitive cell.

Keyword Arguments:
in_place:

If True, changes the current structure. If false returns a new one

tol:

Symmetry precision for symmetry analysis

Examples:

>>> s = io.read(INSTALL_PATH+'io/files/POSCAR_FCC')
>>> print len(s)
4
>>> s.make_primitive()
>>> print len(s)
1
property metrical_matrix

np.dot(self.cell.T, self.cell)

property name

Unformatted name.

property nearest_neighbor_dict

Dict of Atom:[list of Atom] pairs.

pdf_compare(other, tol=0.01)[source]

Compute the PDF for each structure and evaluate the overlap integral for all pairs of species.

recenter(atom, in_place=True, middle=False)[source]

Translate the internal coordinates to center the specified atom. Atom can be an actual Atom from the Structure.atoms list, or can be identified by index.

Keyword Arguments:
in_place:

If False, return a new Structure with the transformation applied. defaults to True.

middle:

If False, “centers” the cell by putting the atom at the origin. If True, “centers” the cell by putting the atom at (0.5,0.5,0.5). defaults to False.

Examples:

>>> s = io.read('POSCAR')
>>> s.recenter(s[2])
>>> s2 = s.recenter(s[0], in_place=False)
>>> s2.recenter(2)
>>> s == s2
True
property reciprocal_lattice

Reciprocal lattice of the structure.

reduce(tol=0.001, limit=1000, in_place=True)[source]

Get the transformation matrix from unit to reduced cell Acta. Cryst. (1976) A32, 297 Acta. Cryst. (2003) A60, 1

Optional keyword arguments:
tol :

eps_rel in Acta. Cryst. 2003 above. Similar to tolerance for floating point comparisons. Defaults to 1e-5.

limit :

maximum number of loops through the algorithm. Defaults to 1000.

in_place :

Change the Structure or return a new one. If True, the transformation matrix is returned. If False, a tuple of (Structure, transformation_matrix) is returned.

Examples:

>>> s = io.read('POSCAR')
>>> s.reduce()
>>> s.reduce(in_place=False, get_transform=False)
refine(tol=0.001, recurse=True)[source]

Identify atoms that are close to high symmetry sites (within tol and shift them onto them.

Note:

“symprec” doesn’t appear to do anything with spglib, so I am unable to get “loose” symmetry operations. Without which, this doesn’t work.

Examples:

>>> s = io.read('POSCAR')
>>> s.symmetrize()
>>> print s.spacegroup
225L
>>> s.refine()
>>> print s.spacegroup
1L
replace(replace, rescale=True, rescale_method='relative', in_place=False, **kwargs)

Replace atoms, as specified in a dict of pairs.

Keyword Arguments:
rescale:

rescale the volume of the final structure based on the per atom volume of the new composition.

rescale_method:

How to rescale the

in_place:

change the species of the current Structure or return a new one.

Examples:

>>> s = io.read('POSCAR-Fe2O3')
>>> s2 = s.substitute({'Fe':'Ni', 'O':'F'} rescale=True)
>>> s2.substitute({'Ni':'Co'}, in_place=True, rescale=False)
save(*args, **kwargs)[source]

Save the current instance. Override this in a subclass if you want to control the saving process.

The ‘force_insert’ and ‘force_update’ parameters can be used to insist that the “save” must be an SQL insert or update (or equivalent for non-SQL backends), respectively. Normally, they should not be set.

set_magnetism(order, elements=None, scheme='primitive')[source]

Assigns magnetic moments to all atoms in accordance with the specified magnetism scheme.

Schemes:

Keyword

Description

None

all magnetic moments = None

“ferro”

atoms with partially filled d and f shells are assigned a magnetic moment of 5 mu_b and 7 mu_b respectively

“anti”

finds a highly ordererd arrangement arrangement of up and down spins. If only 1 magnetic atom is found a ferromagnetic arrangment is used. raises NotImplementedError

set_natoms(n)[source]

Set self.atoms to n blank Atoms.

set_nsites(n)[source]

Sets self.sites to n blank Sites.

set_volume(value)[source]

Rescales the unit cell to the specified volume, keeping the direction and relative magnitudes of all lattice vectors the same.

property site_comp_indices

List of site compositions, length equal to number of sites, each unique site composition identified by an integer.

property site_coords

numpy.ndarray of site coordinates.

property sites

List of Sites in the structure.

property spec_comp

Species composition dictionary.

property species

List of species

property species_id_types

List of species, length equal to number of atoms, each unique species identified by an integer.

property species_types

List of species, length equal to number of atoms.

property stresses

Calculated stresses, a numpy.ndarray of shape (6,)

sub(replace, rescale=True, rescale_method='relative', in_place=False, **kwargs)

Replace atoms, as specified in a dict of pairs.

Keyword Arguments:
rescale:

rescale the volume of the final structure based on the per atom volume of the new composition.

rescale_method:

How to rescale the

in_place:

change the species of the current Structure or return a new one.

Examples:

>>> s = io.read('POSCAR-Fe2O3')
>>> s2 = s.substitute({'Fe':'Ni', 'O':'F'} rescale=True)
>>> s2.substitute({'Ni':'Co'}, in_place=True, rescale=False)
substitute(replace, rescale=True, rescale_method='relative', in_place=False, **kwargs)[source]

Replace atoms, as specified in a dict of pairs.

Keyword Arguments:
rescale:

rescale the volume of the final structure based on the per atom volume of the new composition.

rescale_method:

How to rescale the

in_place:

change the species of the current Structure or return a new one.

Examples:

>>> s = io.read('POSCAR-Fe2O3')
>>> s2 = s.substitute({'Fe':'Ni', 'O':'F'} rescale=True)
>>> s2.substitute({'Ni':'Co'}, in_place=True, rescale=False)
symmetrize(tol=0.001, angle_tol=- 1)[source]

Analyze the symmetry of the structure. Uses spglib to find the symmetry.

symmetrize sets:
  • spacegroup -> Spacegroup

  • uniq_sites -> list of unique Sites

  • orbits -> lists of equivalent Atoms

  • rotations -> List of rotation operations

  • translations -> List of translation operations

  • operatiosn -> List of (rotation,translation) pairs

  • for each atom: atom.wyckoff -> WyckoffSite

  • for each site: site.multiplicity -> int

t(transform, in_place=True, tol=1e-05)

Apply lattice transform to the structure. Accepts transformations of shape (3,) and (3,3).

Optional keyword arguments:
in_placeIf False, return a new Structure with the

transformation applied.

Examples:

>>> s = io.read('POSCAR')
>>> s.transform([2,2,2]) # 2x2x2 supercell
>>> s.transform([[0,1,0],[1,0,0],[0,0,1]]) # swap axis 1 for 2
>>> s2 = s.transform([2,2,2], in_place=False)
transform(transform, in_place=True, tol=1e-05)[source]

Apply lattice transform to the structure. Accepts transformations of shape (3,) and (3,3).

Optional keyword arguments:
in_placeIf False, return a new Structure with the

transformation applied.

Examples:

>>> s = io.read('POSCAR')
>>> s.transform([2,2,2]) # 2x2x2 supercell
>>> s.transform([[0,1,0],[1,0,0],[0,0,1]]) # swap axis 1 for 2
>>> s2 = s.transform([2,2,2], in_place=False)
translate(cv, cartesian=True, in_place=True)[source]

Shifts the contents of the structure by a vector.

Optional keyword arguments:
cartesianIf True, translation vector is taken to be

cartesian coordinates. If False, translation vector is taken to be in fractional coordinates. Default=True

in_placeIf False, return a new Structure with the

transformation applied.

Examples:

>>> s = io.read('POSCAR')
>>> s.translate([1,2,3])
>>> s.translate([0.5,0.5, 0.5], cartesian=False)
>>> s2 = s.translate([-1,2,1], in_place=False)
property unit_comp

Composition dict, where sum(self.unit_comp.values()) == 1

class qmpy.Prototype(*args, **kwargs)[source]

Base class for a prototype structure.

Relationships:
Composition via composition_set
Structure via structure_set
Entry via entry_set
Attributes:
name: Prototype name.
exception DoesNotExist
exception MultipleObjectsReturned
classmethod get(name)[source]

Retrieves a Prototype named name if it exists. If not, creates a new one.

Examples:

>>> proto = Prototype.get('Corundum')

Atom

class qmpy.Atom(*args, **kwargs)[source]

Model for an Atom.

Relationships:
Structure via structure
Element via element
Site via site
WyckoffSite via wyckoff
Attributes:
id
x, y, z: Coordinate of the atom
fx, fy, fz: Forces on the atom
magmom: Magnetic moment on the atom (in &Mu;<sub>b</sub>)
occupancy: Occupation fraction (0-1).
ox: Oxidation state of the atom (can be different from charge)
charge: Charge on the atom
volume: Volume occupied by the atom
exception DoesNotExist
exception MultipleObjectsReturned
property cart_coord

Cartesian coordinates of the Atom.

property coord

[x,y,z] coordinates.

copy()[source]

Creates an exact copy of the atom, only without the matching primary key.

Examples:

>>> a = Atom.get('Fe', [0,0,0])
>>> a.save()
>>> a.id
1
>>> a.copy()
>>> a
<Atom: Fe - 0.000, 0.000, 0.000>
>>> a.id
None
classmethod create(element, coord, **kwargs)[source]

Creates a new Atom object.

Arguments:

element (str or Element): Specifies the element of the Atom. coord (iterable of floats): Specifies the coordinate of the Atom.

Keyword Arguments:
forces:

Specifies the forces on the atom.

magmom:

The magnitude of the magnetic moment on the atom.

charge:

The charge on the Atom.

volume:

The atomic volume of the atom (Angstroms^3).

Examples:

>>> Atom.create('Fe', [0,0,0])
>>> Atom.create('Ni', [0.5, 0.5, 0.5], ox=2, magmom=5,
>>>                                 forces=[0.2, 0.2, 0.2],
>>>                                 volume=101, charge=1.8,
>>>                                 occupancy=1)
property forces

Forces on the Atom in [x, y, z] directions.

property index

None if not in a Structure, otherwise the index of the atom in the structure.

is_on(site, tol=0.001)[source]

Tests whether or not the Atom is on the specified Site.

Examples:

>>> a = Atom.create('Fe', [0,0,0])
>>> s = a.get_site()
>>> a2 = Atom.create('Ni', [0,0,0])
>>> a2.is_on(s)
True
property species

Formatted Species string. e.g. Fe3+, O2-

Site

class qmpy.Site(*args, **kwargs)[source]

A lattice site.

A site can be occupied by one Atom, many Atoms or no Atoms.

Relationships:
Structure via structure
Atom via atom_set
WyckoffSite via wyckoff
Attributes:
id
x, y, z: Coordinate of the Site
exception DoesNotExist
exception MultipleObjectsReturned
add_atom(atom, tol=None)[source]

Adds Atom to Site.atoms.

Notes:

If the Site being assigned to doens’t have a coordinate, it is assigned the coordinate of atom.

Arguments:

atom (Atom): Atom to add to the structure.

Keyword Arguments:
tol (float): Distance between atom and the Site for the Atom to be

assigned to the Site. Raises a SiteError if the distance is greater than tol.

Raises:

SiteError: If atom is more than tol from the Site.

Examples:

>>> s = Site.create([0,0,0])
>>> a = Atom.create('Fe', [0,0,0])
>>> s.add_atom(a)
>>> s2 = Site()
>>> s2.add_atom(a)
property atoms

List of Atoms on the Site.

property cart_coord

Cartesian coordinates of the Atom.

property comp

Composition dictionary of the Site.

Returns:

dict: of (element, occupancy) pairs.

Examples:

>>> a1 = Atom('Fe', [0,0,0], occupancy=0.2)
>>> a2 = Atom('Ni', [0,0,0], occupancy=0.8)
>>> s = Site.from_atoms([a1,a2])
>>> s.comp
{'Fe':0.2, 'Ni':0.8}
property coord

[Site.x, Site.y, Site.z]

static create(coord, comp=None)[source]

Constructs a Site from a coordinate.

Note:

The Site is created without any Atoms occupying it.

Arguments:
coord (length 3 iterable): Assigns the x, y, and z coordinates of

the Site.

Keyword Arguments:
comp (dict, string, or qmpy.Element): Composition dictionary.

Flexible about input forms. Options include: <Element: Fe>, ‘Fe’, {“Fe”:0.5, “Co”:0.5}, and {<Element: Ni>:0.5, <Element: Co>:0.5}.

Raises:

TypeError: if comp isn’t a string, Atom, Element.

Examples:

>>> s = Site.create([0.5,0.5,0.5])
classmethod from_atoms(atoms, tol=0.0001)[source]

Constructs a Site from an iterable of Atoms.

Notes:

Site.coord is set as the average coord of all assigned Atoms.

Checks that the Atoms are close together. If the Atoms are further apart than tol, raises a SiteError

Arguments:

atoms (iterable of Atom): List of Atoms to occupy the Site.

Keyword Arguments:

tol (float): Atoms must be within tol of each other to be assigned to the same Site. Defaults to 1e-4.

Examples:

>>> a1 = Atom.create('Fe', [0,0,0])
>>> a2 = Atom.create('Ni', [1e-5, -1e-5, 0])
>>> s = Site.from_atoms([a1,a1])
property label

Assigns a human friendly label for the Site, based on its atomic composition. If singly occupied, returns the symbol of the atom on the site. If multiply occupied, returns a comma seperated string

Examples:

>>> a1 = Atom.create('Fe', [0,0,0], occupancy=0.2)
>>> a2 = Atom.create('Ni', [0,0,0], occupancy=0.8)
>>> s = Site.from_atoms([a1,a2])
property magmom

Calculates the composition weighted average magnetic moment of the atoms on the Site.

Returns:

float or None

property occupancy

Calculates the total occupancy of the site.

Returns:

float or None

property ox

Calculates the composition weighted average oxidation state of the atoms on the Site.

Returns:

float or None

save(*args, **kwargs)[source]

Save the current instance. Override this in a subclass if you want to control the saving process.

The ‘force_insert’ and ‘force_update’ parameters can be used to insist that the “save” must be an SQL insert or update (or equivalent for non-SQL backends), respectively. Normally, they should not be set.

property spec_comp

Composition dictionary of the Site.

Returns:

dict: of (species, occupancy) pairs.

Examples:

>>> a1 = Atom('Fe', [0,0,0], occupancy=0.2)
>>> a2 = Atom('Ni', [0,0,0], occupancy=0.8)
>>> s = Site.from_atoms([a1,a2])
>>> s.comp
{'Fe':0.2, 'Ni':0.8}

Element

class qmpy.Element(*args, **kwargs)[source]

Core model for an element.

Relationships:
Atom via atom_set
Species via species_set
Structure via structure_set
Entry via entry_set
Composition via composition_set
Calculation via calculation_set
Potential via potential_set
Hubbard via hubbards
HubbardCorrection via hubbardcorrection_set
ReferenceEnergy via referenceenergy_set
Attributes:
Identification
z: atomic number
name: full atomic name
symbol: atomic symbol
group: group in the periodic table
period: period in the periodic table

Physical properties
mass: Atomic mass, in AMU (float)
density: Density at STP, in g/cm^3 (float)
volume: Atomic volume at STP, in A^3/atom (float)
atomic_radii: in A (float)
van_der_waals radii: in A (float)
covalent_radii: in A (float)
scattering_factors: A dictionary of scattering factor coeffs.

Thermodynamic properties
melt: melting point in K
boil: boiling point in K
specific_heat: C_p in J/K

Electronic properties
electronegativity: Pauling electronegativity
ion_energy: First ionization energy. (eV)
s_elec: # of s electrons
p_elec: # of p electrons
d_elec: # of d electrons
f_elec: # of f electrons

Additional information
production: Annual tons of element produced.
abundance: Amount in earths crust (ppm)
radioactive: Are all isotopes unstable?
HHI_P: Herfindahl-Hirschman Index for production.
HHI_R: Herfindahl-Hirschman Index for reserve
Note:

HHI values from Gaultois, M. et al. Chem. Mater. 25, 2911-2920 (2013).

exception DoesNotExist
exception MultipleObjectsReturned
classmethod get(value)[source]

Return an element object. Accepts symbols and atomic numbers, or a list of symbols/atomic numbers.

Examples:

>>> Element.get('Fe')
>>> Element.get(26)
>>> Element.get(['Fe', 'O'])
class qmpy.Species(*args, **kwargs)[source]

Base model for an atomic species. (Element + charge state).

Relationships:
Element via element
Entry via entry_set
Structure via structure_set
Attributes:
name: Species name. e.g. Fe3+, O2-
ox: Oxidation state (float)
exception DoesNotExist
exception MultipleObjectsReturned
classmethod get(value)[source]

Gets or creates the specified species.

Arguments:
value:

Accepts multiple input types. Can be a string, e.g. Fe3+ or a tuple of (symbol, oxidation state) pairs, e.g. (Fe, 3).

Return:

A Species or list of Species.

Examples:

>>> Species.get('Fe3+')
>>> Species.get('Fe3')
>>> Species.get(('Fe', 3))
>>> Species.get([ 'Fe3+', 'O2-', 'Li1+'])

Composition

class qmpy.Composition(*args, **kwargs)[source]

Base class for a composition.

Relationships:
Calculation via calculation_set
Element via element_set
Entry via entry_set
ExptFormationEnergy via exptformationenergy_set
FormationEnergy via formationenergy_set
MetaData via meta_data
Structure via structure_set
Prototype via prototype_set
Attributes:
formula: Electronegativity sorted and normalized composition string.
e.g. Fe2O3, LiFeO2
generic: Genericized composition string. e.g. A2B3, ABC2.
mass: Mass per atom in AMUs
meidema: Meidema model energy for the composition
ntypes: Number of elements.
exception DoesNotExist
exception MultipleObjectsReturned
property comp

Return an element:amount composition dictionary.

property delta_e

Return the lowest formation energy.

property experiment

Return the lowest experimantally measured formation energy at the compositoin.

classmethod get(composition)[source]

Classmethod for getting Composition objects - if the Composition existsin the database, it is returned. If not, a new Composition is created.

Examples:

>>> Composition.get('Fe2O3')
<Composition: Fe2O3>
classmethod get_list(bounds, calculated=False, uncalculated=False)[source]

Classmethod for finding all compositions within the space bounded by a sequence of compositions.

Examples:

>>> from pprint import pprint
>>> comps = Composition.get_list(['Fe','O'], calculated=True)
>>> pprint(list(comps))
[<Composition: Fe>,
 <Composition: FeO>,
 <Composition: FeO3>,
 <Composition: Fe2O3>,
 <Composition: Fe3O4>,
 <Composition: Fe4O5>,
 <Composition: O>]
property ground_state

Return the most stable entry at the composition.

property icsd_delta_e

Return the lowest formation energy calculated from experimentally measured structures - i.e. excluding prototypes.

property ndistinct

Return the number of distinct entries.

property space

Return the set of element symbols

property unit_comp

Return an element:amoutn composition dictionary normalized to a unit composition.