Coverage for wulffpack/winterbottom.py: 100%

25 statements  

« prev     ^ index     » next       coverage.py v7.9.1, created at 2025-06-29 11:40 +0000

1from typing import Dict 

2import numpy as np 

3from ase import Atoms 

4from .core import BaseParticle 

5from .core.form import setup_forms 

6from .core.geometry import (get_symmetries, 

7 break_symmetry, 

8 get_angle, 

9 get_rotation_matrix, 

10 get_standardized_structure) 

11 

12 

13class Winterbottom(BaseParticle): 

14 """ 

15 A `Winterbottom` object is a Winterbottom construction, i.e., 

16 the lowest energy shape adopted by a single crystalline particle 

17 in contact with an interface. 

18 

19 Parameters 

20 ---------- 

21 surface_energies 

22 A dictionary with surface energies, where keys are 

23 Miller indices and values surface energies (per area) 

24 in a unit of choice, such as J/m^2. 

25 interface_direction 

26 Miller indices for the interface facet. 

27 interface_energy 

28 Energy per area for interfaces. 

29 Caution: The interface energy here is equivalent to the generalized 

30 surface energy described in the original Winterbottom paper. 

31 Thus, this interface energy is the difference between the 

32 substrate-particle surface energy and the substrate-vapor surface energy. 

33 Note that this parameter can be negative. 

34 primtive_structure 

35 primitive structure to implicitly define the point group as 

36 well as the atomic structure used if an atomic structure 

37 is requested. By default, an Au FCC structure is used. 

38 natoms 

39 Together with `primitive_structure`, this parameter 

40 defines the volume of the particle. If an atomic structure 

41 is requested, the number of atoms will as closely as possible 

42 match this value. 

43 symprec 

44 Numerical tolerance for symmetry analysis, forwarded to spglib. 

45 tol 

46 Numerical tolerance parameter. 

47 

48 Example 

49 ------- 

50 The following example illustrates some possible uses of a 

51 `Winterbottom` object:: 

52 

53 >>> from wulffpack import Winterbottom 

54 >>> from ase.build import bulk 

55 >>> from ase.io import write 

56 >>> surface_energies = {(1, 1, 0): 1.0, (1, 0, 0): 1.08} 

57 >>> prim = bulk('Fe', a=4.1, crystalstructure='bcc') 

58 >>> particle = Winterbottom(surface_energies=surface_energies, 

59 ... interface_direction=(3, 2, 1), 

60 ... interface_energy=0.4, 

61 ... primitive_structure=prim) 

62 >>> particle.view() 

63 >>> write('winterbottom.xyz', particle.atoms) # Writes atomic structure to file 

64 

65 """ 

66 

67 def __init__(self, 

68 surface_energies: Dict[tuple, float], 

69 interface_direction: tuple, 

70 interface_energy: float, 

71 primitive_structure: Atoms = None, 

72 natoms: int = 1000, 

73 symprec: float = 1e-5, 

74 tol: float = 1e-5): 

75 standardized_structure = get_standardized_structure(primitive_structure, symprec=symprec) 

76 full_symmetries = get_symmetries(standardized_structure, symprec=symprec) 

77 broken_symmetries = break_symmetry(full_symmetries, 

78 [interface_direction]) 

79 

80 if abs(interface_energy) > min(surface_energies.values()): 

81 raise ValueError('The construction expects an absolute interface energy ' 

82 'that is smaller than than the smallest ' 

83 'surface energy.') 

84 surface_energies = surface_energies.copy() 

85 surface_energies['interface'] = interface_energy 

86 forms = setup_forms(surface_energies, 

87 standardized_structure.cell.T, 

88 broken_symmetries, 

89 full_symmetries, 

90 interface=interface_direction) 

91 

92 super().__init__(forms=forms, 

93 standardized_structure=standardized_structure, 

94 natoms=natoms, 

95 tol=tol) 

96 

97 # Rotate particle such that the interface aligns with the plane z=0 

98 target = (0, 0, -1) 

99 rotation_axis = np.cross(interface_direction, target) 

100 angle = get_angle(interface_direction, target) 

101 R = get_rotation_matrix(angle, rotation_axis.astype(float)) 

102 self.rotate_particle(R) 

103 

104 @property 

105 def atoms(self): 

106 """ 

107 Returns an ASE Atoms object 

108 """ 

109 return self._get_atoms()