Winterbottom construction

class wulffpack.Winterbottom(surface_energies, interface_direction, interface_energy, primitive_structure=None, natoms=1000, symprec=1e-05, tol=1e-05)[source]

A Winterbottom object is a Winterbottom construction, i.e., the lowest energy shape adopted by a single crystalline particle in contact with an interface.

Parameters:
  • surface_energies (Dict[tuple, float]) – A dictionary with surface energies, where keys are Miller indices and values surface energies (per area) in a unit of choice, such as J/m^2.

  • interface_direction (tuple) – Miller indices for the interface facet.

  • interface_energy (float) – Energy per area for twin boundaries.

  • primtive_structure – primitive structure to implicitly define the point group as well as the atomic structure used if an atomic structure is requested. By default, an Au FCC structure is used.

  • natoms (int) – Together with primitive_structure, this parameter defines the volume of the particle. If an atomic structure is requested, the number of atoms will as closely as possible match this value.

  • symprec (float) – Numerical tolerance for symmetry analysis, forwarded to spglib.

  • tol (float) – Numerical tolerance parameter.

Example

The following example illustrates some possible uses of a Winterbottom object:

>>> from wulffpack import Winterbottom
>>> from ase.build import bulk
>>> from ase.io import write
>>> surface_energies = {(1, 1, 0): 1.0, (1, 0, 0): 1.08}
>>> prim = bulk('Fe', a=4.1, crystalstructure='bcc')
>>> particle = Winterbottom(surface_energies=surface_energies,
...                         interface_direction=(3, 2, 1),
...                         interface_energy=0.4,
...                         primitive_structure=prim)
>>> particle.view()
>>> write('winterbottom.xyz', particle.atoms)  # Writes atomic structure to file
property area: float

Returns total area of the surface of the particle (not including twin boundaries).

property atoms

Returns an ASE Atoms object

property average_surface_energy: float

Average surface energy for the Wulff construction, i.e., a weighted average over all the facets, where the weights are the area fraction of each facet.

property edge_length: float

Returns total edge length of the particle.

property facet_fractions: Dict[tuple, float]

Returns a dict specifying fraction of each form (not including twin boundaries).

property forms: List[Form]

List of inequivalent forms for the particle

get_continuous_color_scheme(base_colors=None, normalize=False)

Returns a dictionary with RGB colors for each form. The colors smoothly interpolate between three base colors, corresponding to (1, 1, 1), (1, 1, 0) and (1, 0, 0). Note that this is sensible primarily for cubic systems.

Parameters:
  • base_colors (Optional[dict]) – User chosen colors for one or several of (1, 1, 1), (1, 1, 0) and (1, 0, 0). To enforce, say, green (1, 1, 1), use base_colors={(1, 1, 1): 'g'}

  • normalize (bool) – If True, the norm of the RGB vectors will be 1. Note that this may affect the base_colors too.

Return type:

dict

make_plot(ax, alpha=0.85, linewidth=0.3, colors=None)

Plot a particle in an axis object. This function can be used to make customized plots of particles.

Parameters:
  • ax (matplotlib Axes3DSubplot) – An axis object with 3d projection

  • alpha (float) – Opacity of the faces

  • linewidth (float) – Thickness of lines between faces

  • colors (Optional[dict]) – Allows custom colors for facets of all or a subset of forms, example {(1, 1, 1): '#FF0000'}

Example

In the following example, three different particles are plotted in the same figure:

>>> from wulffpack import SingleCrystal, Decahedron, Icosahedron
>>> import matplotlib.pyplot as plt
>>> from mpl_toolkits.mplot3d import Axes3D
>>>
>>> surface_energies = {(1, 1, 1): 1.0,
...                     (1, 0, 0): 1.1,
...                     (1, 1, 0): 1.15,
...                     (3, 2, 1): 1.15}
>>> twin_energy = 0.05
>>>
>>> fig = plt.figure(figsize=(3*4.0, 4.0))
>>> ax = fig.add_subplot(131, projection='3d')
>>> particle = SingleCrystal(surface_energies)
>>> particle.make_plot(ax)
>>>
>>> ax = fig.add_subplot(132, projection='3d')
>>> particle = Decahedron(surface_energies,
...                       twin_energy=0.05)
>>> particle.make_plot(ax)
>>>
>>> ax = fig.add_subplot(133, projection='3d')
>>> particle = Icosahedron(surface_energies,
...                        twin_energy=0.05)
>>> particle.make_plot(ax)
>>>
>>> plt.subplots_adjust(top=1, bottom=0, left=0,
...                     right=1, wspace=0, hspace=0)
>>> plt.savefig('particles.png')
property natoms: List[int]

The approximate number of atoms in the particle (implicitly defining the volume).

property number_of_corners: float

Returns the number of corners (vertices) on the particle.

rotate_particle(rotation)

Rotate the particle.

Parameters:

rotation (ndarray) – Rotation matrix

property standardized_structure: Atoms

The standardized atomic structure that defines the geometry and thus the meaning of the Miller indices. Also forms the building blocks when particle.atoms is called.

property surface_energy: float

The total surface energy of the particle (including twin boundaries).

translate_particle(translation)

Translate the particle.

Parameters:

translation (list of 3 floats) – Translation vector

view(alpha=0.85, linewidth=0.3, colors=None, legend=True, save_as=None)

Use matplotlib to view a rendition of the particle.

Parameters:
  • alpha (float) – Opacity of the faces

  • linewidth (float) – Thickness of lines between faces

  • colors (Optional[dict]) – Allows custom colors for facets of all or a subset of forms, example {(1, 1, 1): '#FF0000'}

  • legend (bool) – Whether or not to show a legend with facet-color definitions

  • save_as (Optional[str]) – Filename to save figure as. If None, show the particle with the GUI instead.

property volume: float

Returns the volume of the particle

write(filename)

Write particle to file. The file format is derived from the filename. Currently supported fileformats are:

  • Wavefront .obj

Parameters:

filename (str) – Filename of file to write to