Nanoparticles and clusters

There are modules for creating nanoparticles (clusters) with a given crystal structure by specifying either the number of layers in different directions, or by making a Wulff construction.

Examples

Layer specification

This example sets up a nanoparticle of copper in the FCC crystal structure, by specifying 6 layers in the (100) directions, 9 in the (110) directions and 5 in the (111) directions:

import ase
from ase.cluster.cubic import FaceCenteredCubic

surfaces = [(1, 0, 0), (1, 1, 0), (1, 1, 1)]
layers = [6, 9, 5]
lc = 3.61000
atoms = FaceCenteredCubic('Cu', surfaces, layers, latticeconstant=lc)

culayer

Wulff construction

To set up a Wulff construction, the surface energies should be specified, in units of energy per area (not energy per atom). The actual unit used does not matter, as only the ratio between surface energies is important. In addition, the approximate size of the nanoparticle should be given. As the Wulff construction is build from whole layers, it is not possible to hit the desired particles size exactly:

from ase.cluster import wulff_construction

surfaces = [(1, 0, 0), (1, 1, 0), (1, 1, 1)]
esurf = [1.0, 1.1, 0.9]   # Surface energies.
lc = 3.61000
size = 1000  # Number of atoms
atoms = wulff_construction('Cu', surfaces, esurf,
                           size, 'fcc',
                           rounding='above', latticeconstant=lc)

Note that the Wulff construction currently only work with cubic lattices.

Creating a nanoparticle

The ase.cluster module contains a number of sub-modules for defining clusters, one for each crystal structure. They are all called the same way, by specifying the element, the number of layers in different directions, and optionally the lattice constant.

The layer specification is the only part that may not be intuitive. It is given as two arrays, one specifying the Miller indices of the surfaces, and one specifying the number of layers from the center of the cluster to the respective surfaces.

The surface specification allows for one or more surfaces of a given family of surfaces to be different from the other surfaces. This can be used e.g. to create a cluster where one part has been truncated by a substrate. This is done by first specifying the number of layers for the family of surfaces, and later specifying the number of layers for a given surface. Consider the surface specification

surfaces = [(1, 0, 0), (1, 1, 1), (1, -1, 1)]
layers = [6, 5, -1]
atoms = FaceCenteredCubic('Cu', surfaces, layers)

Here we first ask for 6 layers for the {100} surface family, i.e. the directions (100), (010), (001), (-1,0,0), etc. Then we ask for 5 layers for the {111} family of surfaces. Finally, we change the number of layers for the (-1,1,1) surface. This is interpreted as a single surface, since it is part of a family that has already been specified. Asking for a negative number of layers is allowed, this cause the particle to be truncated before its center point. The result is seen below.

truncated

The functions for creating nanoparticles take the following arguments:

symbols: A string specifying the element (or a tuple of strings for compounds).

surfaces: A list of surfaces, as explained above.

layers: A corresponding list of the number of layers to be included.

vacuum=0.0: The amount of vacuum to include around the particle. Defaults to 0.

latticeconstant=None: The lattice constant of the lattice. If not specified, the experimental value from ase.data is used.

Possible crystal structures

You select the crystal structure by selecting the right function for creating the nanoparticle. Currently, these modules only work for the three cubic crystal structures: FaceCenteredCubic, BodyCenteredCubic, and SimpleCubic. Other structures are implemented, but do currently not work correctly.

Module for creating clusters.

class ase.cluster.Cluster(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, velocities=None)[source]
copy()[source]

Return a copy.

get_diameter(method='volume')[source]

Returns an estimate of the cluster diameter based on two different methods.

Parameters:

method ({'volume', 'shape'}) – ‘volume’ (default) returns the diameter of a sphere with the same volume as the atoms. ‘shape’ returns the averaged diameter calculated from the directions given by the defined surfaces.

get_layers()[source]

Return number of atomic layers in stored surfaces directions.

get_surfaces()[source]

Returns the miller indexs of the stored surfaces of the cluster.

ase.cluster.Decahedron(symbol, p, q, r, latticeconstant=None)[source]

Return a decahedral cluster.

Parameters:
  • symbol (Chemical symbol (or atomic number) of the element.) –

  • p (Number of atoms on the (100) facets perpendicular to the five) –

  • axis. (fold) –

  • q (Number of atoms on the (100) facets parallel to the five fold) –

  • facets. (axis. q = 1 corresponds to no visible (100)) –

  • r (Depth of the Marks re-entrence at the pentagon corners. r = 0) –

  • re-entrence. (corresponds to no) –

  • (optional) (latticeconstant) –

  • ase.data. (then it is extracted form) –

ase.cluster.Icosahedron(symbol, noshells, latticeconstant=None)[source]

Returns a cluster with the icosahedra symmetry.

Parameters:
  • symbol (str or int) – The chemical symbol (or atomic number) of the element.

  • noshells (int) – The number of shells (>= 1).

  • latticeconstant (float, optional) – The lattice constant. If not given, then it is extracted from \(ase.data\).

ase.cluster.Octahedron(symbol, length, cutoff=0, latticeconstant=None, alloy=False)[source]

Returns Face Centered Cubic clusters of the octahedral class depending on the choice of cutoff.

Type

Condition

Regular octahedron

cutoff = 0

Truncated octahedron

cutoff > 0

Regular truncated octahedron

length = 3 * cutoff + 1

Cuboctahedron

length = 2 * cutoff + 1

Parameters:
  • symbol (str or list) – The chemical symbol or atomic number of the element(s).

  • length (int) – Number of atoms on the square edges of the complete octahedron.

  • cutoff (int, default 0) – Number of layers cut at each vertex.

  • latticeconstant (float, optional) – The lattice constant. If not given, then it is extracted from \(ase.data\).

  • alloy (bool, default False) – If True the L1_2 structure is used.

ase.cluster.wulff_construction(symbol, surfaces, energies, size, structure, rounding='closest', latticeconstant=None, debug=False, maxiter=100)[source]

Create a cluster using the Wulff construction.

A cluster is created with approximately the number of atoms specified, following the Wulff construction, i.e. minimizing the surface energy of the cluster.

Parameters:
  • symbol (str or int) – The chemical symbol (or atomic number) of the desired element.

  • surfaces (list) – A list of surfaces. Each surface is an (h, k, l) tuple or list of integers.

  • energies (list) – A list of surface energies for the surfaces.

  • size (int) – The desired number of atoms.

  • structure ({'fcc', bcc', 'sc'}) – The desired crystal structure.

  • rounding ({'closest', 'above', 'below'}) – Specifies what should be done if no Wulff construction corresponds to exactly the requested number of atoms. ‘above’, ‘below’, and ‘closest’ mean that the nearest cluster above or below - or the closest one - is created instead.

  • latticeconstant (float (optional)) – The lattice constant. If not given, extracted from \(ase.data\).

  • debug (bool, default False) – If True, information about the iteration towards the right cluster size is printed.