matscipy.gamma_surface
Classes
|
A class for generating gamma surface/generalised stacking fault images & plots |
|
Class for stacking fault-specific generation & plotting |
- class matscipy.gamma_surface.GammaSurface(a, surface_direction, glide_direction=None, crystalstructure=None, symbol='C')
Bases:
object
A class for generating gamma surface/generalised stacking fault images & plots
Methods
generate_images
(nx, ny[, z_reps, z_offset, ...])Generate gamma surface images on an (nx, ny) grid
get_energy_densities
([calculator, relax])Get (self.nx, self.ny) grid of gamma surface energies from self.images
plot_energy_densities
([Es, ax, si, ...])Produce a matplotlib plot of the gamma surface energy, from the data gathered in self.generate_images and self.get_surface_energies
relax_images
([calculator, ftol, optimiser, ...])Utility function to relax gamma surface images using calculator
show
([CNA_color, plot_energies, si])Overload of GammaSurface.show() Plots an animation of the stacking fault structure, and optionally the associated energy densities (requires self.get_energy_densitities() to have been called)
- __init__(a, surface_direction, glide_direction=None, crystalstructure=None, symbol='C')
Initialise by cutting and rotating the input structure.
- Parameters:
a (float or ase.Atoms) – Lattice Constant or Starting structure to generate gamma surface from (Operates similarly to CubicCrystalDislocation) If lattice constant is provided, crystalstructure must also be set
surface_direction (np.array of int or subclass of) – matscipy.dislocation.CubicCrystalDissociatedDislocation Vector direction of gamma surface, in miller index notation EG: np.array([0, 0, 1]), np.array([-1, 1, 0]), np.array([1, 1, 1]) A subclass of matscipy.dislocation.CubicCrystalDissociatedDislocation (EG: DiamondGlideScrew or FCCEdge110Dislocation)
glide_direction (np.array of int or None) – Basis vector (in miller indices) to form the glide direction of the stacking fault, which is oriented along the y axis of generated images Should be orthogonal to surface_direction If None, a suitable glide_direction will be found automatically
crystalstructure (str) – Crystal Structure to use in building a base cubic cell Required when a is a lattice constant Current accepted values: “fcc”, “bcc”, “diamond”
symbol (str) – Chemical symbol to feed to ase.lattice.cubic cell builders Required when a is a lattice constant
- self.cut_at
Cut Atoms object used as base for gamma surface image generation
- Type:
ase.atoms.Atoms object
- self.surf_directions
Dict giving miller indices for “x”, “y” and “z” directions of gamma surface plot - self.surf_directions[“z”] = surface_direction - self.surf_directions[“y”] = glide_direction, if glide_direction was specified
- Type:
dict
- self.nx, self.ny
Dimensions of the (nx, ny) gamma surface grid
- Type:
int
- self.images
Generated gamma surface images (populated after self.generate_images is called)
- Type:
list of ase.atoms.Atoms objects
- generate_images(nx, ny, z_reps=1, z_offset=0.0, cell_strain=0.0, vacuum=0.0, path_xlims=[0, 1], path_ylims=None, cell_move=True)
Generate gamma surface images on an (nx, ny) grid
- Parameters:
nx (int) – Number of points in the x direction
ny (int) – Number of points in the y direction
z_reps (int) – Number of supercell copies in z (increases separation between periodic surfaces)
z_offset (float) – Offset in the z direction (in A) to apply to all atoms Used to select different stacking fault planes sharing the same normal direction (e.g. glide and shuffle planes in Diamond)
cell_strain (float) – Fractional strain to apply in z direction to cell only (0.1 = +10% strain; atoms aren’t moved) Helpful for issues when atoms are extremely close in unrelaxed images.
vacuum (float) – Additional vacuum layer (in A) to add between periodic images of the gamma surface
vacuum_offset (float) – Offset (in A) applied to the position of the vacuum layer in the cell The position of the vacuum layer is given by: vac_pos = self.cut_at.cell[2, 2] * (1 + cell_strain) / 2 + vacuum_offset
path_xlims (list or array of floats) – Limits (in fractional coordinates) of the stacking fault path in the x direction
path_ylims (list or array of floats) – Limits (in fractional coordinates) of the stacking fault path in the x direction If not supplied, will be set to either [0, 1], or will be set based on the glide distance of the supplied dislocation
cell_move (bool) – Toggles using the cell move method (True) or atom move method (False)
- Returns:
images – nx*ny length list of gamma surface images
- Return type:
list of ase.atoms.Atoms
- relax_images(calculator=None, ftol=0.001, optimiser=<class 'ase.optimize.bfgslinesearch.BFGSLineSearch'>, constrain_atoms=True, cell_relax=True, logfile=None, steps=200)
Utility function to relax gamma surface images using calculator
- Parameters:
calculator (ase calculator instance) – Calculator to give forces for relaxation
ftol (float) – force tolerance for optimisation
optimiser (ase optimiser) – optimiser method to use
constrain_atoms (bool) – Toggle constraint that all atoms only relax in z direction (normal to the gamma surface)
cell_relax (bool) – Allow z component of cell (cell[2, 2]) to relax (normal to the gamma surface)
steps (int) – Maximum number of iterations allowed for each image relaxation Fed directly to optimiser.run()
- get_energy_densities(calculator=None, relax=False, **relax_kwargs)
Get (self.nx, self.ny) grid of gamma surface energies from self.images
- Parameters:
calculator (ase calculator) – Calculator to use for finding surface energies (and optionally in the relaxation)
relax (bool) – Whether to additionally relax the images (through a call to self.relax_images) before finding energies
**relax_kwargs (keyword args) – Extra kwargs to be passed to self.relax_images if relax=True
- Returns:
Es (np.array)
(nx, ny) array of energy densities (energy per unit surface area)
for the gamma surface, in eV/A**2
- plot_energy_densities(Es=None, ax=None, si=True, interpolation='bicubic')
Produce a matplotlib plot of the gamma surface energy, from the data gathered in self.generate_images and self.get_surface_energies
Returns a matplotlib fig and ax object Es: np.array
(nx, ny) array of energy densities. If None, uses self.Es (if populated from self.get_energy_densities())
- ax: matplotlib axis object
Axis to draw plot
- si: bool
Use SI units (J/m^2) if True, else ASE units (eV/A^2)
- interpolation: str
arg to pass to matplotlib imshow interpolation
- show(CNA_color=True, plot_energies=False, si=False, **kwargs)
Overload of GammaSurface.show() Plots an animation of the stacking fault structure, and optionally the associated energy densities (requires self.get_energy_densitities() to have been called)
- CNA_color: bool
Toggle atom colours based on Common Neighbour Analysis (Structure identification)
- plot_energies:
Add additional plot showing energy density profile, alongside the structure
- si: bool
Plot energy densities in SI units (J/m^2), or “ASE” units eV/A^2
- **kwargs
extra kwargs passed to ase.visualize.plot.plot_atoms
- class matscipy.gamma_surface.StackingFault(a, surface_direction, glide_direction=None, crystalstructure=None, symbol='C')
Bases:
GammaSurface
Class for stacking fault-specific generation & plotting
Methods
generate_images
(n, *args, **kwargs)Generate gamma surface images on an (n) line n: int Number of images
get_energy_densities
([calculator, relax])Get (self.nx, self.ny) grid of gamma surface energies from self.images
plot_energy_densities
([Es, ax, si])Produce a matplotlib plot of the stacking fault energy, from the data gathered in self.generate_images and self.get_surface_energy
relax_images
([calculator, ftol, optimiser, ...])Utility function to relax gamma surface images using calculator
show
([CNA_color, plot_energies, si])Overload of GammaSurface.show() Plots an animation of the stacking fault structure, and optionally the associated energy densities (requires self.get_energy_densitities() to have been called)
- __init__(a, surface_direction, glide_direction=None, crystalstructure=None, symbol='C')
Initialise by cutting and rotating the input structure.
- Parameters:
a (float or ase.Atoms) – Lattice Constant or Starting structure to generate gamma surface from (Operates similarly to CubicCrystalDislocation) If lattice constant is provided, crystalstructure must also be set
surface_direction (np.array of int or subclass of) – matscipy.dislocation.CubicCrystalDissociatedDislocation Vector direction of gamma surface, in miller index notation EG: np.array([0, 0, 1]), np.array([-1, 1, 0]), np.array([1, 1, 1]) A subclass of matscipy.dislocation.CubicCrystalDissociatedDislocation (EG: DiamondGlideScrew or FCCEdge110Dislocation)
glide_direction (np.array of int or None) – Basis vector (in miller indices) to form the glide direction of the stacking fault, which is oriented along the y axis of generated images Should be orthogonal to surface_direction If None, a suitable glide_direction will be found automatically
crystalstructure (str) – Crystal Structure to use in building a base cubic cell Required when a is a lattice constant Current accepted values: “fcc”, “bcc”, “diamond”
symbol (str) – Chemical symbol to feed to ase.lattice.cubic cell builders Required when a is a lattice constant
- self.cut_at
Cut Atoms object used as base for gamma surface image generation
- Type:
ase.atoms.Atoms object
- self.surf_directions
Dict giving miller indices for “x”, “y” and “z” directions of gamma surface plot - self.surf_directions[“z”] = surface_direction - self.surf_directions[“y”] = glide_direction, if glide_direction was specified
- Type:
dict
- self.nx, self.ny
Dimensions of the (nx, ny) gamma surface grid
- Type:
int
- self.images
Generated gamma surface images (populated after self.generate_images is called)
- Type:
list of ase.atoms.Atoms objects
- generate_images(n, *args, **kwargs)
Generate gamma surface images on an (n) line n: int
Number of images
- get_energy_densities(calculator=None, relax=False, **relax_kwargs)
Get (self.nx, self.ny) grid of gamma surface energies from self.images
- Parameters:
calculator (ase calculator) – Calculator to use for finding surface energies (and optionally in the relaxation)
relax (bool) – Whether to additionally relax the images (through a call to self.relax_images) before finding energies
**relax_kwargs (keyword args) – Extra kwargs to be passed to self.relax_images if relax=True
- Returns:
Es (np.array)
(nx, ny) array of energy densities (energy per unit surface area)
for the gamma surface, in eV/A**2
- relax_images(calculator=None, ftol=0.001, optimiser=<class 'ase.optimize.bfgslinesearch.BFGSLineSearch'>, constrain_atoms=True, cell_relax=True, logfile=None, steps=200)
Utility function to relax gamma surface images using calculator
- Parameters:
calculator (ase calculator instance) – Calculator to give forces for relaxation
ftol (float) – force tolerance for optimisation
optimiser (ase optimiser) – optimiser method to use
constrain_atoms (bool) – Toggle constraint that all atoms only relax in z direction (normal to the gamma surface)
cell_relax (bool) – Allow z component of cell (cell[2, 2]) to relax (normal to the gamma surface)
steps (int) – Maximum number of iterations allowed for each image relaxation Fed directly to optimiser.run()
- plot_energy_densities(Es=None, ax=None, si=False)
Produce a matplotlib plot of the stacking fault energy, from the data gathered in self.generate_images and self.get_surface_energy
Returns a matplotlib fig and ax object Es: np.array
(nx, ny) array of energy densities. If None, uses self.Es (if populated from self.get_energy_densities())
- ax: matplotlib axis object
Axis to draw plot
- si: bool
Use SI units (J/m^2) if True, else ASE units (eV/A^2)
- show(CNA_color=True, plot_energies=False, si=False, **kwargs)
Overload of GammaSurface.show() Plots an animation of the stacking fault structure, and optionally the associated energy densities (requires self.get_energy_densitities() to have been called)
- CNA_color: bool
Toggle atom colours based on Common Neighbour Analysis (Structure identification)
- plot_energies:
Add additional plot showing energy density profile, alongside the structure
- si: bool
Plot energy densities in SI units (J/m^2), or “ASE” units eV/A^2
- rotation: str
rotation to apply to the structures, passed directly to ase.visualize.plot.plot_atoms
- **kwargs
extra kwargs passed to ase.visualize.plot.plot_atoms