dscribe.core package

Submodules

dscribe.core.lattice module

Copyright 2019 DScribe developers

Licensed under the Apache License, Version 2.0 (the “License”); you may not use this file except in compliance with the License. You may obtain a copy of the License at

Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an “AS IS” BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License.

class dscribe.core.lattice.Lattice(matrix)[source]

Bases: object

A lattice object. Essentially a matrix with conversion matrices. In general, it is assumed that length units are in Angstroms and angles are in degrees unless otherwise stated.

Create a lattice from any sequence of 9 numbers. Note that the sequence is assumed to be read one row at a time. Each row represents one lattice vector.

Parameters:

matrix – Sequence of numbers in any form. Examples of acceptable input. i) An actual numpy array. ii) [[1, 0, 0], [0, 1, 0], [0, 0, 1]] iii) [1, 0, 0 , 0, 1, 0, 0, 0, 1] iv) (1, 0, 0, 0, 1, 0, 0, 0, 1) Each row should correspond to a lattice vector. E.g., [[10, 0, 0], [20, 10, 0], [0, 0, 30]] specifies a lattice with lattice vectors [10, 0, 0], [20, 10, 0] and [0, 0, 30].

property abc

Lengths of the lattice vectors, i.e. (a, b, c)

get_cartesian_coords(fractional_coords)[source]

Returns the cartesian coordinates given fractional coordinates.

Parameters:

fractional_coords (3x1 array) – Fractional coords.

Returns:

Cartesian coordinates

get_fractional_coords(cart_coords)[source]

Returns the fractional coordinates given cartesian coordinates.

Parameters:

cart_coords (3x1 array) – Cartesian coords.

Returns:

Fractional coordinates.

get_points_in_sphere(frac_points, center, r, zip_results=True)[source]

Find all points within a sphere from the point taking into account periodic boundary conditions. This includes sites in other periodic images.

Algorithm:

  1. place sphere of radius r in crystal and determine minimum supercell (parallelpiped) which would contain a sphere of radius r. for this we need the projection of a_1 on a unit vector perpendicular to a_2 & a_3 (i.e. the unit vector in the direction b_1) to determine how many a_1”s it will take to contain the sphere.

    Nxmax = r * length_of_b_1 / (2 Pi)

  2. keep points falling within r.

Parameters:
  • frac_points – All points in the lattice in fractional coordinates.

  • center – Cartesian coordinates of center of sphere.

  • r – radius of sphere.

  • zip_results (bool) – Whether to zip the results together to group by point, or return the raw fcoord, dist, index arrays

Returns:

[(fcoord, dist, index) …] since most of the time, subsequent

processing requires the distance.

else:

fcoords, dists, inds

Return type:

if zip_results

property inv_matrix

Inverse of lattice matrix.

property lengths
property matrix

Copy of matrix representing the Lattice

property reciprocal_lattice

Return the reciprocal lattice. Note that this is the standard reciprocal lattice used for solid state physics with a factor of 2 * pi. If you are looking for the crystallographic reciprocal lattice, use the reciprocal_lattice_crystallographic property. The property is lazily generated for efficiency.

property reciprocal_lattice_crystallographic

Returns the crystallographic reciprocal lattice, i.e., no factor of 2 * pi.

dscribe.core.system module

Copyright 2019 DScribe developers

Licensed under the Apache License, Version 2.0 (the “License”); you may not use this file except in compliance with the License. You may obtain a copy of the License at

Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an “AS IS” BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License.

class dscribe.core.system.System(symbols=None, positions=None, numbers=None, tags=None, momenta=None, masses=None, magmoms=None, charges=None, scaled_positions=None, cell=None, pbc=None, celldisp=None, constraint=None, calculator=None, info=None, wyckoff_positions=None, equivalent_atoms=None)[source]

Bases: Atoms

Represents atomic systems that are used internally by the package. Inherits from the ase.Atoms class, but adds the possibility to cache various time-consuming quantities that can be shared when creating multiple descriptors.

static from_atoms(atoms)[source]

Creates a System object from ASE.Atoms object.

get_cell_inverse()[source]

Get the matrix inverse of the lattice matrix.

get_displacement_tensor()[source]

A matrix where the entry A[i, j, :] is the vector self.cartesian_pos[j] - self.cartesian_pos[i].

For periodic systems the distance of an atom from itself is the smallest displacement of an atom from one of it’s periodic copies, and the distance of two different atoms is the distance of two closest copies.

Returns:

3D matrix containing the pairwise distance vectors.

Return type:

np.array

get_distance_matrix()[source]

Calculates the distance matrix A defined as:

\[A_{ij} = \lvert \mathbf{r}_i - \mathbf{r}_j \rvert\]

For periodic systems the distance of an atom from itself is the smallest displacement of an atom from one of it’s periodic copies, and the distance of two different atoms is the distance of two closest copies.

Returns:

Symmetric 2D matrix containing the pairwise distances.

Return type:

np.array

get_distance_matrix_within_radius(radius, pos=None, output_type='coo_matrix')[source]

Calculates a sparse distance matrix by only considering distances within a certain cutoff. Uses a k-d tree to reach O(n log(N)) time complexity.

Parameters:
  • radius (float) – The cutoff radius within which distances are calculated. Distances outside this radius are not included.

  • output_type (str) – Which container to use for output data. Options: “dok_matrix”, “coo_matrix”, “dict”, or “ndarray”. Default: “dok_matrix”.

Returns:

Symmetric sparse 2D matrix containing the pairwise distances.

Return type:

dok_matrix | np.array | coo_matrix | dict

get_inverse_distance_matrix()[source]

Calculates the inverse distance matrix A defined as:

\[A_{ij} = \frac{1}{\lvert \mathbf{r}_i - \mathbf{r}_j \rvert }\]

For periodic systems the distance of an atom from itself is the smallest displacement of an atom from one of it’s periodic copies, and the distance of two different atoms is the distance of two closest copies.

Returns:

Symmetric 2D matrix containing the pairwise inverse distances.

Return type:

np.array

set_cell(cell, scale_atoms=False)[source]

Set unit cell vectors.

Parameters:

cell: 3x3 matrix or length 3 or 6 vector

Unit cell. A 3x3 matrix (the three unit cell vectors) or just three numbers for an orthorhombic cell. Another option is 6 numbers, which describes unit cell with lengths of unit cell vectors and with angles between them (in degrees), in following order: [len(a), len(b), len(c), angle(b,c), angle(a,c), angle(a,b)]. First vector will lie in x-direction, second in xy-plane, and the third one in z-positive subspace.

scale_atoms: bool

Fix atomic positions or move atoms with the unit cell? Default behavior is to not move the atoms (scale_atoms=False).

apply_constraint: bool

Whether to apply constraints to the given cell.

Examples:

Two equivalent ways to define an orthorhombic cell:

>>> atoms = Atoms('He')
>>> a, b, c = 7, 7.5, 8
>>> atoms.set_cell([a, b, c])
>>> atoms.set_cell([(a, 0, 0), (0, b, 0), (0, 0, c)])

FCC unit cell:

>>> atoms.set_cell([(0, b, b), (b, 0, b), (b, b, 0)])

Hexagonal unit cell:

>>> atoms.set_cell([a, a, c, 90, 90, 120])

Rhombohedral unit cell:

>>> alpha = 77
>>> atoms.set_cell([a, a, a, alpha, alpha, alpha])
set_pbc(pbc)[source]

Set periodic boundary condition flags.

set_positions(newpositions, apply_constraint=True)[source]

Set positions, honoring any constraints. To ignore constraints, use apply_constraint=False.

set_scaled_positions(scaled)[source]

Set positions relative to unit cell.

to_cartesian(scaled_positions, wrap=False)[source]

Used to transform a set of relative positions to the cartesian basis defined by the cell of this system.

Parameters:
  • positions (numpy.ndarray) – The positions to scale

  • wrap (numpy.ndarray) – Whether the positions should be wrapped inside the cell.

Returns:

The cartesian positions

Return type:

numpy.ndarray

to_scaled(positions, wrap=False)[source]

Used to transform a set of positions to the basis defined by the cell of this system.

Parameters:
  • positions (numpy.ndarray) – The positions to scale

  • wrap (numpy.ndarray) – Whether the positions should be wrapped inside the cell.

Returns:

The scaled positions

Return type:

numpy.ndarray

Module contents

Copyright 2019 DScribe developers

Licensed under the Apache License, Version 2.0 (the “License”); you may not use this file except in compliance with the License. You may obtain a copy of the License at

Unless required by applicable law or agreed to in writing, software distributed under the License is distributed on an “AS IS” BASIS, WITHOUT WARRANTIES OR CONDITIONS OF ANY KIND, either express or implied. See the License for the specific language governing permissions and limitations under the License.