Core components¶
Base particle¶
- class wulffpack.core.BaseParticle(forms, standardized_structure, natoms, ngrains=1, volume_scale=1.0, tol=1e-05)[source]¶
Base class for Wulff constructions, useful for systems of (almost) any symmetry.
The Wulff construction calculates the vertices of the dual, calculates the convex hull of those vertices using scipy’s ConvexHull, and then converts the dual back to the Wulff shape.
- Parameters:
forms (
List
[Form
]) – List of forms belonging to the system at handstandardized_structure (
Atoms
) – The structural unit implicitly carrying information about basis vectors etc.natoms (
int
) – Used as a target number of atoms if an atomistic structure is requested. Together withstandardized_structure
,natoms
also defines the volume of the particle.ngrains (
int
) – If only a part of a complete particle is to be constructed, for example a fifth of a decahedron, this should be reflected by this parameter.volume_scale (
int
) – Decahedral and icosahedral particle need to be strained to fill space. With this parameter, the volume is scaled such that the volume of the particle is correct after this straining.tolerance – Numerical tolerance
- property area: float¶
Returns total area of the surface of the particle (not including twin boundaries).
- 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).
- get_continuous_color_scheme(base_colors=None, normalize=False)[source]¶
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), usebase_colors={(1, 1, 1): 'g'}
normalize (
bool
) – If True, the norm of the RGB vectors will be 1. Note that this may affect thebase_colors
too.
- Return type:
dict
- make_plot(ax, alpha=0.85, linewidth=0.3, colors=None)[source]¶
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 faceslinewidth (
float
) – Thickness of lines between facescolors (
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)[source]¶
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)[source]¶
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)[source]¶
Use matplotlib to view a rendition of the particle.
- Parameters:
alpha (
float
) – Opacity of the faceslinewidth (
float
) – Thickness of lines between facescolors (
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 definitionssave_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
Form¶
- class wulffpack.core.Form(miller_indices, energy, cell, symmetries, parent_miller_indices)[source]¶
A
Form
object contains all facets that are equivalent by the symmetry of the system under consideration. For a cubic system, this means that the form with Miller indices {100} will contain the facets with normals [100], [-100], [010] etc.- Parameters:
miller_indices (
tuple
) – Miller indices for a representative member of the form.energy (
float
) – The surface energy per area for a facet in this form.cell (
ndarray
) – The cell, with the basis vectors as columns.symmetries (
List
[ndarray
]) – Symmetry elements of the system. If symmetry is broken by, e.g., the presence of an interface, this list should only contain the symmetries surviving.parent_miller_indices (
tuple
) – If symmetry is broken, it may still be of interest to know what its form would be had not the symmetry been broken. This attribute contains the Miller indices for a representative of such a form.
- property area: float¶
Returns the total area of the facets in this form.
- property edge_length: float¶
Returns the length of all facets in the form.
- property surface_energy: float¶
Returns the total surface energy of the facets in this form.
- property volume: float¶
Returns the total volume formed by pyramids formed by each facet and the origin.
Facet¶
- class wulffpack.core.Facet(normal, energy, symmetry, original_grain=True, tol=1e-05)[source]¶
A
Facet
carries information about a particular facet, i.e., only one facet out of potentially several that are equivalent by symmetry.- Parameters:
normal (
Union
[List
[float
],Tuple
[float
]]) – The normal to the surface in Cartesian coordintesenergy (
float
) – Surface energy (in units of energy per surface area)rotation – Symmetry operation (acting on scaled coordinates) that carried an arbtriraily defined “parent facet” to this facet
original_grain (
bool
) – If True, this facet is part of an “original grain”, i.e., in the case of decahedron/icosahedron, it is part of the first created tetrahedrontol (
float
) – Numerical tolerance
- add_vertex(vertex)[source]¶
Add a vertex to the list of vertices. This function checks that there is not already such a vertex, to ensure that no vertex is duplicated. It is also checked that the vertex lies in the expected plane (as defined by the first vertex and the surface normal).
- Parameters:
vertex (
Union
[List
,tuple
]) – The vertex in Cartesian coordinates
- property area: float¶
Returns the total area of the facet.
- property distance_from_origin: float¶
Returns the distance from the origin to the facet.
- property face_as_triangles: List[List[ndarray]]¶
Given N 3D coordinates, ABCD…N, split them into triangles that do not overlap. Two of the vertices of each triangle will be vertices of the face, the third vertex will be the origin.
- get_volume(origin=None)[source]¶
Calculate the volume formed by the pyramid having the face as its base and the specified origin as its tip.
- Parameters:
origin (list of 3 ints) – The origin
- Return type:
The total volume
- property ordered_vertices: List[ndarray]¶
Returns the vertices ordered such that they enclose a polygon. The first/last vertex occurs twice.
It is assumed (but not checked) that the coordinates lie in a plane and that they form a convex polygon.
- property perimeter_length: float¶
Returns the length of the perimeter
- remove_redundant_vertices()[source]¶
Remove any vertex which is midway between two other vertices. This would not work if the polygon were concave, but that is not allowed anyway.
- property surface_energy: float¶
Returns the total surface energy of the facet.