Dislocation
- matscipy.dislocation.make_screw_cyl(alat, C11, C12, C44, cylinder_r=10, cutoff=5.5, hard_core=False, center=[0.0, 0.0, 0.0], l_extend=[0.0, 0.0, 0.0], symbol='W')
Makes screw dislocation using atomman library
- Parameters:
alat (float) – Lattice constant of the material.
C11 (float) – C11 elastic constant of the material.
C12 (float) – C12 elastic constant of the material.
C44 (float) – C44 elastic constant of the material.
cylinder_r (float) – radius of cylinder of unconstrained atoms around the dislocation in angstrom
cutoff (float) – Potential cutoff for Marinica potentials for FS cutoff = 4.4
hard_core (bool) – Description of parameter hard_core.
center (type) –
- The position of the dislocation core and the center of the
cylinder with FixAtoms condition
l_extend (float) – extension of the box. used for creation of initial dislocation position with box equivalent to the final position
symbol (string) – Symbol of the element to pass to ase.lattice.cubic.SimpleCubicFactory default is “W” for tungsten
- Returns:
disloc (ase.Atoms object) – screw dislocation cylinder.
bulk (ase.Atoms object) – bulk disk used to generate dislocation
u (np.array) – displacement per atom.
- matscipy.dislocation.make_edge_cyl(alat, C11, C12, C44, cylinder_r=10, cutoff=5.5, symbol='W')
makes edge dislocation using atomman library
- cylinder_r - radius of cylinder of unconstrained atoms around the
dislocation in angstrom
cutoff - potential cutoff for Marinica potentials for FS cutoff = 4.4
- symbolstring
Symbol of the element to pass to ase.lattice.cubic.SimpleCubicFactory default is “W” for tungsten
- matscipy.dislocation.plot_vitek(dislo, bulk, alat=3.16, plot_axes=None, xyscale=10)
Plots vitek map from ase configurations.
- Parameters:
dislo (ase.Atoms) – Dislocation configuration.
bulk (ase.Atoms) – Corresponding bulk configuration for calculation of displacements.
alat (float) – Lattice parameter for calculation of neighbour list cutoff.
plot_axes (matplotlib.Axes.axes object) – Existing axes to plot on, allows to pass existing matplotlib axes have full control of the graph outside the function. Makes possible to plot multiple differential displacement maps using subplots. Default is None, then new graph is created by plt.subplots() Description of parameter plot_axes.
xyscale (float) – xyscale of the graph
- Return type:
None
- matscipy.dislocation.show_NEB_configurations(images, bulk, xyscale=7, show=True, core_positions=None)
- Plots Vitek differential displacement maps for the list of images
for example along the NEB path.
- Parameters:
images (list of ase.Atoms) – List of configurations with dislocations.
bulk (ase.Atoms) – Corresponding bulk configuration for calculation of displacements.
xyscale (float) – xyscale of the graph
show (bool) – Show the figure after plotting. Default is True.
core_positions (list) – [x, y] position of dislocation core to plot
- Returns:
If the show is False else returns None
- Return type:
figure
- matscipy.dislocation.show_configuration(disloc, bulk, u, fixed_mask=None)
shows the displacement fixed atoms.
- matscipy.dislocation.get_elastic_constants(pot_path=None, calculator=None, delta=0.01, symbol='W')
return lattice parameter, and cubic elastic constants: C11, C12, 44 using matscipy function pot_path - path to the potential
- symbolstring
Symbol of the element to pass to ase.lattice.cubic.SimpleCubicFactory default is “W” for tungsten
- matscipy.dislocation.make_barrier_configurations(elastic_param=None, pot_path=None, calculator=None, cylinder_r=10, hard_core=False, **kwargs)
- Creates the initial and final configurations for the NEB calculation
The positions in FixedAtoms constrained region are average between final and initial configurations
- Parameters:
pot_path (string) – Path to the potential file.
calculator (type) – Description of parameter calculator.
cylinder_r (float) –
- Radius of cylinder of unconstrained atoms around the
dislocation in angstrom.
hard_core (bool) – Type of the core hard or soft. If hard is chosen the displacement field is reversed.
**kwargs – Keyword arguments to pass to make_screw_cyl() function.
- Returns:
disloc_ini (ase.Atoms) – Initial dislocation configuration.
disloc_fin (ase.Atoms) – Final dislocation configuration.
bulk (ase.Atoms) – Perfect bulk configuration for Vitek displacement maps
- matscipy.dislocation.make_screw_cyl_kink(alat, C11, C12, C44, cylinder_r=40, kink_length=26, kind='double', **kwargs)
Function to create kink configuration based on make_screw_cyl() function. Double kink configuration is in agreement with quadrupoles in terms of formation energy. Single kink configurations provide correct and stable structure, but formation energy is not accessible?
- Parameters:
alat (float) – Lattice constant of the material.
C11 (float) – C11 elastic constant of the material.
C12 (float) – C12 elastic constant of the material.
C44 (float) – C44 elastic constant of the material.
cylinder_r (float) – radius of cylinder of unconstrained atoms around the dislocation in angstrom
kink_length (int) – Length of the cell per kink along b in unit of b, must be even.
kind (string) – kind of the kink: right, left or double
**kwargs – Keyword arguments to pass to make_screw_cyl() function.
- Returns:
kink (ase.atoms) – kink configuration
reference_straight_disloc (ase.atoms) – reference straight dislocation configuration
large_bulk (ase.atoms) – large bulk cell corresponding to the kink configuration
- matscipy.dislocation.slice_long_dislo(kink, kink_bulk, b)
- Function to slice a long dislocation configuration to perform
dislocation structure and core position analysis
- Parameters:
kink (ase.Atoms) – kink configuration to slice
kink_bulk (ase.Atoms) – corresponding bulk configuration to perform mapping for slicing
b (float) – burgers vector b should be along z direction
- Returns:
sliced_kink (list of [sliced_bulk, sliced_kink]) – sliced configurations 1 b length each
disloc_z_positions (float) – positions of each sliced configuration (center along z)
- matscipy.dislocation.compare_configurations(dislo, bulk, dislo_ref, bulk_ref, alat, cylinder_r=None, print_info=True, remap=True, bulk_neighbours=None, origin=(0.0, 0.0))
- Compares two dislocation configurations based on the gradient of
the displacements along the bonds.
- Parameters:
dislo (ase.Atoms) – Dislocation configuration.
bulk (ase.Atoms) – Corresponding bulk configuration for calculation of displacements.
dislo_ref (ase.Atoms) – Reference dislocation configuration.
bulk_ref (ase.Atoms) – Corresponding reference bulk configuration for calculation of displacements.
alat (float) – Lattice parameter for calculation of neghbour list cutoff.
cylinder_r (float or None) – Radius of region of comparison around the dislocation coreself. If None makes global comparison based on the radius of dislo configuration, else compares the regions with cylinder_r around the dislocation core position.
print_info (bool) – Flag to switch print statement about the type of the comparison
remap (bool) – Flag to swtich off remapping of atoms between deformed and reference configurations. Only set this to true if atom order is the same!
bulk_neighbours – Optionally pass in bulk neighbours as a tuple (bulk_i, bulk_j)
origin (tuple) – Optionally pass in coordinate origin (x0, y0)
- Returns:
The Du norm of the differences per atom.
- Return type:
float
- matscipy.dislocation.cost_function(pos, dislo, bulk, cylinder_r, elastic_param, hard_core=False, print_info=True, remap=True, bulk_neighbours=None, origin=(0, 0))
- Cost function for fitting analytical displacement field
and detecting dislocation core position. Uses compare_configurations function for the minimisation of the core position.
- Parameters:
pos (list of float) – Positions of the core to build the analytical solution [x, y].
dislo (ase.Atoms) – Dislocation configuration.
bulk (ase.Atoms) – Corresponding bulk configuration for calculation of displacements.
cylinder_r (float or None) – Radius of region of comparison around the dislocation coreself. If None makes global comparison based on the radius of dislo configuration, else compares the regions with cylinder_r around the dislocation core position.
elastic_param (list of float) – List containing alat, C11, C12, C44
hard_core (bool) – type of the core True for hard
print_info (bool) – Flag to switch print statement about the type of the comparison
bulk_neighbours (tuple or None) – Optionally pass in neighbour list for bulk reference config to save computing it each time.
origin (tuple) – Optionally pass in coordinate origin (x0, y0)
- Returns:
Error for optimisation (result from compare_configurations function)
- Return type:
float
- matscipy.dislocation.fit_core_position(dislo_image, bulk, elastic_param, hard_core=False, core_radius=10, current_pos=None, bulk_neighbours=None, origin=(0, 0))
Use cost_function() to fit atomic positions to Stroh solution with
scipy.optimize.minimize is used to perform the fit using Powell’s method.
- Parameters:
dislo_image (ase.atoms.Atoms) –
bulk (ase.atoms.Atoms) –
elastic_param (array-like) – [alat, C11, C12, C44]
hard_core (bool) –
core_radius (float) –
current_pos (array-like) – array [core_x, core_y] containing initial guess for core position
bulk_neighbours (tuple) – cache of bulk neigbbours to speed up calcualtion. Should be a tuple (bulk_I, bulk_J) as returned by matscipy.neigbbours.neighbour_list(‘ij’, bulk, alat).
origin (tuple) – Optionally pass in coordinate origin (x0, y0)
- Return type:
core_pos - array [core_x, core_y]
- matscipy.dislocation.fit_core_position_images(images, bulk, elastic_param, bulk_neighbours=None, origin=(0, 0))
Call fit_core_position() for a list of Atoms objects, e.g. NEB images
- Parameters:
images (list) – list of Atoms object for dislocation configurations
bulk (ase.atoms.Atoms) – bulk reference configuration
elastic_param (list) – as for fit_core_position().
bulk_neighbours – as for fit_core_position().
origin (tuple) – Optionally pass in coordinate origin (x0, y0)
- Returns:
core_positions
- Return type:
array of shape (len(images), 2)
- matscipy.dislocation.screw_cyl_tetrahedral(alat, C11, C12, C44, scan_r=15, symbol='W', imp_symbol='H', hard_core=False, center=(0.0, 0.0, 0.0))
- Generates a set of tetrahedral positions with scan_r radius.
Applies the screw dislocation displacement for creating an initial guess for the H positions at dislocation core.
- Parameters:
alat (float) – Lattice constant of the material.
C11 (float) – C11 elastic constant of the material.
C12 (float) – C12 elastic constant of the material.
C44 (float) – C44 elastic constant of the material.
scan_r (float) – Radius of the region to create tetrahedral positions.
symbol (string) – Symbol of the element to pass to ase.lattuce.cubic.SimpleCubicFactory default is “W” for tungsten
imp_symbol (string) – Symbol of the elemnt to pass creat Atoms object default is “H” for hydrogen
hard_core (float) – Type of the dislocatino core if True then -u (sign of displacement is flipped) is applied. Default is False i.e. soft core is created.
- centertuple of floats
Coordinates of dislocation core (center) (x, y, z). Default is (0., 0., 0.)
- Returns:
Atoms object with predicted tetrahedral positions around dislocation core.
- Return type:
ase.Atoms object
- matscipy.dislocation.screw_cyl_octahedral(alat, C11, C12, C44, scan_r=15, symbol='W', imp_symbol='H', hard_core=False, center=(0.0, 0.0, 0.0))
- Generates a set of octahedral positions with scan_r radius.
Applies the screw dislocation displacement for creating an initial guess for the H positions at dislocation core.
- Parameters:
alat (float) – Lattice constant of the material.
C11 (float) – C11 elastic constant of the material.
C12 (float) – C12 elastic constant of the material.
C44 (float) – C44 elastic constant of the material.
symbol (string) – Symbol of the element to pass to ase.lattuce.cubic.SimpleCubicFactory default is “W” for tungsten
imp_symbol (string) – Symbol of the elemnt to pass creat Atoms object default is “H” for hydrogen
symbol – Symbol of the elemnt to pass creat Atoms object
hard_core (float) – Type of the dislocatino core if True then -u (sign of displacement is flipped) is applied. Default is False i.e. soft core is created.
center (tuple of floats) – Coordinates of dislocation core (center) (x, y, z). Default is (0., 0., 0.)
- Returns:
Atoms object with predicted tetrahedral positions around dislocation core.
- Return type:
ase.Atoms object
- class matscipy.dislocation.BodyCenteredCubicTetrahedralFactory
Bases:
SimpleCubicFactory
A factory for creating tetrahedral lattices in bcc structure
- xtal_name = 'bcc_tetrahedral'
- bravais_basis: Sequence[Sequence[float]] | None = [[0.0, 0.5, 0.25], [0.0, 0.5, 0.75], [0.0, 0.25, 0.5], [0.0, 0.75, 0.5], [0.5, 0.0, 0.75], [0.25, 0.0, 0.5], [0.75, 0.0, 0.5], [0.5, 0.0, 0.25], [0.5, 0.25, 0.0], [0.5, 0.75, 0.0], [0.25, 0.5, 0.0], [0.75, 0.5, 0.0]]
- align()
Align the first axis along x-axis and the second in the x-y plane.
- atoms_in_unit_cell = 1
- basis_factor = 1.0
- calc_num_atoms()
- check_basis_volume()
Check the volume of the unit cell.
- chop_tolerance = 1e-10
- convert_to_natural_basis()
Convert directions and miller indices to the natural basis.
- element_basis: Sequence[int] | None = None
- find_directions(directions, miller)
Find missing directions and miller indices from the specified ones.
- find_ortho(idx)
Replace keyword ‘ortho’ or ‘orthogonal’ with a direction.
- get_lattice_constant()
Get the lattice constant of an element with cubic crystal structure.
- inside(point)
Is a point inside the unit cell?
- int_basis = array([[1, 0, 0], [0, 1, 0], [0, 0, 1]])
- inverse_basis = array([[1, 0, 0], [0, 1, 0], [0, 0, 1]])
- inverse_basis_factor = 1.0
- make_crystal_basis()
Make the basis matrix for the crystal unit cell and the system unit cell.
- make_list_of_atoms()
Repeat the unit cell.
- make_unit_cell()
Make the unit cell.
- other = {0: (1, 2), 1: (2, 0), 2: (0, 1)}
- print_directions_and_miller(txt='')
Print direction vectors and Miller indices.
- process_element(element)
Extract atomic number from element
- put_atom(point)
Place an atom given its integer coordinates.
- class matscipy.dislocation.BodyCenteredCubicOctahedralFactory
Bases:
SimpleCubicFactory
A factory for creating octahedral lattices in bcc structure
- xtal_name = 'bcc_octahedral'
- bravais_basis: Sequence[Sequence[float]] | None = [[0.5, 0.5, 0.0], [0.0, 0.0, 0.5], [0.5, 0.0, 0.0], [0.5, 0.0, 0.5], [0.0, 0.5, 0.0], [0.0, 0.5, 0.5]]
- align()
Align the first axis along x-axis and the second in the x-y plane.
- atoms_in_unit_cell = 1
- basis_factor = 1.0
- calc_num_atoms()
- check_basis_volume()
Check the volume of the unit cell.
- chop_tolerance = 1e-10
- convert_to_natural_basis()
Convert directions and miller indices to the natural basis.
- element_basis: Sequence[int] | None = None
- find_directions(directions, miller)
Find missing directions and miller indices from the specified ones.
- find_ortho(idx)
Replace keyword ‘ortho’ or ‘orthogonal’ with a direction.
- get_lattice_constant()
Get the lattice constant of an element with cubic crystal structure.
- inside(point)
Is a point inside the unit cell?
- int_basis = array([[1, 0, 0], [0, 1, 0], [0, 0, 1]])
- inverse_basis = array([[1, 0, 0], [0, 1, 0], [0, 0, 1]])
- inverse_basis_factor = 1.0
- make_crystal_basis()
Make the basis matrix for the crystal unit cell and the system unit cell.
- make_list_of_atoms()
Repeat the unit cell.
- make_unit_cell()
Make the unit cell.
- other = {0: (1, 2), 1: (2, 0), 2: (0, 1)}
- print_directions_and_miller(txt='')
Print direction vectors and Miller indices.
- process_element(element)
Extract atomic number from element
- put_atom(point)
Place an atom given its integer coordinates.
- matscipy.dislocation.dipole_displacement_angle(W_bulk, dislo_coord_left, dislo_coord_right, shift=0.0, mode=1.0)
Generates a simple displacement field for two dislocations in a dipole configuration uding simple Voltera solution as u = b/2 * angle
- matscipy.dislocation.get_u_img(W_bulk, dislo_coord_left, dislo_coord_right, n_img=10, n1_shift=0, n2_shift=0)
Function for getting displacemnt filed for images of quadrupole cells used by make_screw_quadrupole
- matscipy.dislocation.make_screw_quadrupole(alat, left_shift=0, right_shift=0, n1u=5, symbol='W')
- Generates a screw dislocation dipole configuration
for effective quadrupole arrangement. Works for BCC systems.
- Parameters:
alat (float) – Lattice parameter of the system in Angstrom.
left_shift (float, optional) – Shift of the left dislocation core in number of dsitances to next equivalent disocation core positions needed for creation for final configuration for NEB. Default is 0.
right_shift (float, optional) – shift of the right dislocation core in number of dsitances to next equivalent disocation core positions needed for creation for final configuration for NEB. Default is 0.
n1u (int, odd number) – odd number! length of the cell a doubled distance between core along x. Main parameter to calculate cell vectors
symbol (string) – Symbol of the element to pass to ase.lattuce.cubic.SimpleCubicFactory default is “W” for tungsten
- Returns:
disloc_quadrupole (ase.Atoms) – Resulting quadrupole configuration.
W_bulk (ase.Atoms) – Perfect system.
dislo_coord_left (list of float) – Coodrinates of left dislocation core [x, y]
dislo_coord_right (list of float) – Coodrinates of right dislocation core [x, y]
Notes
Calculation of cell vectors
From [1] we take:
Unit vectors for the cell are:
\[u = \frac{1}{3}[1 \bar{2} 1];\]\[v = \frac{1}{3}[2 \bar{1} \bar{1}];\]\[z = b = \frac{1}{2}[1 1 1];\]Cell vectors are:
\[C_1 = n^u_1 u + n^v_1 v + C^z_1 z;\]\[C_2 = n^u_2 u + n^v_2 v + C^z_2 z;\]\[C_3 = z\]For quadrupole arrangement n1u needs to be odd number, for 135 atoms cell we take n1u=5
To have quadrupole as as close as possible to a square one has to take:
\[2 n^u_2 + n^v_2 = n^u_1\]\[n^v_2 \approx \frac{n^u_1}{\sqrt{3}}\]for n1u = 5:
\[n^v_2 \approx \frac{n^u_1}{\sqrt{3}} = 2.89 \approx 3.0\]\[n^u_2 = \frac{1}{2} (n^u_1 - n^v_2) \approx \frac{1}{2} (5-3)=1\]Following [2] cell geometry is optimized by ading tilt compomemts
Cz1 and Cz2 for our case of n1u = 3n - 1:
Easy core
\[C^z_1 = + \frac{1}{3}\]\[C^z_2 = + \frac{1}{6}\]Hard core
\[C^z_1 = + \frac{1}{3}\]\[C^z_2 = + \frac{1}{6}\]may be typo in the paper check the original!
References:
- matscipy.dislocation.make_screw_quadrupole_kink(alat, kind='double', n1u=5, kink_length=20, symbol='W')
- Generates kink configuration using make_screw_quadrupole() function
works for BCC structure. The method is based on paper https://doi.org/10.1016/j.jnucmat.2008.12.053
- Parameters:
alat (float) – Lattice parameter of the system in Angstrom.
kind (string) – kind of the kink: right, left or double
n1u (int) – Number of lattice vectors for the quadrupole cell (make_screw_quadrupole() function)
kink_length (int) – Length of the cell per kink along b in unit of b, must be even.
symbol (string) – Symbol of the element to pass to ase.lattuce.cubic.SimpleCubicFactory default is “W” for tungsten
- Returns:
kink (ase.atoms) – kink configuration
reference_straight_disloc (ase.atoms) – reference straight dislocation configuration
large_bulk (ase.atoms) – large bulk cell corresponding to the kink configuration
- matscipy.dislocation.make_edge_cyl_001_100(a0, C11, C12, C44, cylinder_r, cutoff=5.5, tol=1e-06, symbol='W')
Function to produce consistent edge dislocation configuration.
- Parameters:
alat (float) – Lattice constant of the material.
C11 (float) – C11 elastic constant of the material.
C12 (float) – C12 elastic constant of the material.
C44 (float) – C44 elastic constant of the material.
cylinder_r (float) – Radius of cylinder of unconstrained atoms around the dislocation in angstrom.
cutoff (float) – Potential cutoff for determenition of size of fixed atoms region (2*cutoff)
tol (float) – Tolerance for generation of self consistent solution.
symbol (string) – Symbol of the element to pass to ase.lattuce.cubic.SimpleCubicFactory default is “W” for tungsten
- Returns:
bulk (ase.Atoms object) – Bulk configuration.
disloc (ase.Atoms object) – Dislocation configuration.
disp (np.array) – Corresponding displacement.
- matscipy.dislocation.read_dislo_QMMM(filename=None, image=None)
Reads extended xyz file with QMMM configuration Uses “region” for mapping of QM, MM and fixed atoms Sets ase.constraints.FixAtoms constraint on fixed atoms
- Parameters:
filename (path to xyz file) –
image (image with "region" array to set up constraint and extract qm_mask) –
- Returns:
dislo_QMMM (Output ase.Atoms object) – Includes “region” array and FixAtoms constraint
qm_mask (array mask for QM atoms mapping)
- matscipy.dislocation.plot_bulk(atoms, n_planes=3, ax=None, ms=200)
Plots x, y coordinates of atoms colored according to non-equivalent planes in z plane
- matscipy.dislocation.ovito_dxa_straight_dislo_info(disloc, structure='BCC', replicate_z=3)
A function to extract information from ovito dxa analysis. Current version works for 1b thick configurations containing straight dislocations.
- Parameters:
disloc (ase.Atoms) – Atoms object containing the atomic configuration to analyse
replicate_z (int) – Specifies number of times to replicate the configuration along the dislocation line. Ovito dxa analysis needs at least 3b thick cell to work.
- Returns:
Results
- Return type:
np.array(position, b, line, angle)
- matscipy.dislocation.get_centering_mask(atoms, radius, core_position=[0.0, 0.0, 0.0], extension=[0.0, 0.0, 0.0])
- matscipy.dislocation.check_duplicates(atoms, distance=0.1)
Returns a mask of atoms that have at least one other atom closer than distance
- class matscipy.dislocation.CubicCrystalDislocation(unit_cell, alat, C11, C12, C44, axes, burgers, unit_cell_core_position=None, parity=None, glide_distance=None, n_planes=None, self_consistent=None)
Bases:
object
- __init__(unit_cell, alat, C11, C12, C44, axes, burgers, unit_cell_core_position=None, parity=None, glide_distance=None, n_planes=None, self_consistent=None)
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- init_stroh()
- set_burgers(burgers)
- plot_unit_cell(ms=250, ax=None)
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- class matscipy.dislocation.BCCScrew111Dislocation(alat, C11, C12, C44, symbol='W')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='W')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.BCCEdge111Dislocation(alat, C11, C12, C44, symbol='W')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='W')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.BCCMixed111Dislocation(alat, C11, C12, C44, symbol='W')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='W')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.BCCEdge100Dislocation(alat, C11, C12, C44, symbol='W')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='W')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.BCCEdge100110Dislocation(alat, C11, C12, C44, symbol='W')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='W')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.DiamondGlide30degreePartial(alat, C11, C12, C44, symbol='C')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='C')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.DiamondGlide90degreePartial(alat, C11, C12, C44, symbol='C')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='C')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.CubicCrystalDissociatedDislocation(left_dislocation, right_dislocation, burgers)
Bases:
CubicCrystalDislocation
- __init__(left_dislocation, right_dislocation, burgers)
This class represents a dissociated dislocation in a cubic crystal with burgers vercor b = b_left + b_right.
- Parameters:
left_dislocation (CubicCrystalDislocation) – dislocation with b_left
right_dislocation (CubicCrystalDislocation) – dislocation with b_right
burgers (ndarray of 3 floats) – resulting burgers vector
- Raises:
ValueError – If resulting burgers vector burgers is not a sum of burgers vectors of left and right dislocations.
ValueError – If one of the properties of left and righ dislocations are not the same.
- build_cylinder(radius, partial_distance=0, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
Overloaded function to make dissociated dislocations. Partial distance is provided as an integer to define number of glide distances between two partials.
- Parameters:
radius (float) – radius of the cell
partial_distance (int) – distance between partials (SF length) in number of glide distances. Default is 0 -> non dissociated dislocation with b = b_left + b_right is produced
- displacements(bulk_positions, center, partial_distance=0, **kwargs)
Overloaded function to provide correct displacements for the dissociated dislocation. Partial distance is provided as an integer to define number of glide distances between two partials.
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.DiamondGlideScrew(alat, C11, C12, C44, symbol='C')
Bases:
CubicCrystalDissociatedDislocation
- __init__(alat, C11, C12, C44, symbol='C')
This class represents a dissociated dislocation in a cubic crystal with burgers vercor b = b_left + b_right.
- Parameters:
left_dislocation (CubicCrystalDislocation) – dislocation with b_left
right_dislocation (CubicCrystalDislocation) – dislocation with b_right
burgers (ndarray of 3 floats) – resulting burgers vector
- Raises:
ValueError – If resulting burgers vector burgers is not a sum of burgers vectors of left and right dislocations.
ValueError – If one of the properties of left and righ dislocations are not the same.
- build_cylinder(radius, partial_distance=0, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
Overloaded function to make dissociated dislocations. Partial distance is provided as an integer to define number of glide distances between two partials.
- Parameters:
radius (float) – radius of the cell
partial_distance (int) – distance between partials (SF length) in number of glide distances. Default is 0 -> non dissociated dislocation with b = b_left + b_right is produced
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, partial_distance=0, **kwargs)
Overloaded function to provide correct displacements for the dissociated dislocation. Partial distance is provided as an integer to define number of glide distances between two partials.
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.DiamondGlide60Degree(alat, C11, C12, C44, symbol='C')
Bases:
CubicCrystalDissociatedDislocation
- __init__(alat, C11, C12, C44, symbol='C')
This class represents a dissociated dislocation in a cubic crystal with burgers vercor b = b_left + b_right.
- Parameters:
left_dislocation (CubicCrystalDislocation) – dislocation with b_left
right_dislocation (CubicCrystalDislocation) – dislocation with b_right
burgers (ndarray of 3 floats) – resulting burgers vector
- Raises:
ValueError – If resulting burgers vector burgers is not a sum of burgers vectors of left and right dislocations.
ValueError – If one of the properties of left and righ dislocations are not the same.
- build_cylinder(radius, partial_distance=0, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
Overloaded function to make dissociated dislocations. Partial distance is provided as an integer to define number of glide distances between two partials.
- Parameters:
radius (float) – radius of the cell
partial_distance (int) – distance between partials (SF length) in number of glide distances. Default is 0 -> non dissociated dislocation with b = b_left + b_right is produced
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, partial_distance=0, **kwargs)
Overloaded function to provide correct displacements for the dissociated dislocation. Partial distance is provided as an integer to define number of glide distances between two partials.
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.FCCScrewShockleyPartial(alat, C11, C12, C44, symbol='Fe')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='Fe')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.FCCScrew110Dislocation(alat, C11, C12, C44, symbol='Fe')
Bases:
CubicCrystalDissociatedDislocation
- __init__(alat, C11, C12, C44, symbol='Fe')
This class represents a dissociated dislocation in a cubic crystal with burgers vercor b = b_left + b_right.
- Parameters:
left_dislocation (CubicCrystalDislocation) – dislocation with b_left
right_dislocation (CubicCrystalDislocation) – dislocation with b_right
burgers (ndarray of 3 floats) – resulting burgers vector
- Raises:
ValueError – If resulting burgers vector burgers is not a sum of burgers vectors of left and right dislocations.
ValueError – If one of the properties of left and righ dislocations are not the same.
- build_cylinder(radius, partial_distance=0, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
Overloaded function to make dissociated dislocations. Partial distance is provided as an integer to define number of glide distances between two partials.
- Parameters:
radius (float) – radius of the cell
partial_distance (int) – distance between partials (SF length) in number of glide distances. Default is 0 -> non dissociated dislocation with b = b_left + b_right is produced
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, partial_distance=0, **kwargs)
Overloaded function to provide correct displacements for the dissociated dislocation. Partial distance is provided as an integer to define number of glide distances between two partials.
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.FCCEdgeShockleyPartial(alat, C11, C12, C44, symbol='Fe')
Bases:
CubicCrystalDislocation
- __init__(alat, C11, C12, C44, symbol='Fe')
This class represents a dislocation in a cubic crystal
The dislocation is defined by the crystal unit cell, elastic constants C11, C12 and C44, crystal axes, burgers vector and optional shift and parity vectors.
- Parameters:
unit_cell (unit cell to build the dislocation configuration) –
alat (lattice constant) –
C11 (elastic constants) –
C12 –
C44 –
axes (cell axes (b is normally along z direction)) –
burgers (burgers vector of the dislocation) –
unit_cell_core_position (dislocation core position in the unit cell) – used to shift atomic positions to make the dislocation core the center of the cell
parity –
glide_distance (distance to the next equivalent) – core position in the glide direction
n_planes (int) – number of non equivalent planes in z direction
self_consistent (float) – default value for the displacement calculation
- build_cylinder(radius, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, self_consistent=True, tol=1e-06, max_iter=100, verbose=True)
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.FCCEdge110Dislocation(alat, C11, C12, C44, symbol='Fe')
Bases:
CubicCrystalDissociatedDislocation
- __init__(alat, C11, C12, C44, symbol='Fe')
This class represents a dissociated dislocation in a cubic crystal with burgers vercor b = b_left + b_right.
- Parameters:
left_dislocation (CubicCrystalDislocation) – dislocation with b_left
right_dislocation (CubicCrystalDislocation) – dislocation with b_right
burgers (ndarray of 3 floats) – resulting burgers vector
- Raises:
ValueError – If resulting burgers vector burgers is not a sum of burgers vectors of left and right dislocations.
ValueError – If one of the properties of left and righ dislocations are not the same.
- build_cylinder(radius, partial_distance=0, core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), fix_width=10.0, self_consistent=None)
Overloaded function to make dissociated dislocations. Partial distance is provided as an integer to define number of glide distances between two partials.
- Parameters:
radius (float) – radius of the cell
partial_distance (int) – distance between partials (SF length) in number of glide distances. Default is 0 -> non dissociated dislocation with b = b_left + b_right is produced
- build_glide_configurations(radius, average_positions=False, **kwargs)
- build_impurity_cylinder(disloc, impurity, radius, imp_symbol='H', core_position=array([0., 0., 0.]), extension=array([0., 0., 0.]), self_consistent=False, extra_bulk_at_core=False, core_radius=0.5, shift=array([0., 0., 0.]))
- displacements(bulk_positions, center, partial_distance=0, **kwargs)
Overloaded function to provide correct displacements for the dissociated dislocation. Partial distance is provided as an integer to define number of glide distances between two partials.
- init_stroh()
- plot_unit_cell(ms=250, ax=None)
- set_burgers(burgers)
- class matscipy.dislocation.FixedLineAtoms(a, direction)
Bases:
object
Constrain atoms to move along a given direction only.
- __init__(a, direction)
- adjust_positions(atoms, newpositions)
- adjust_forces(atoms, forces)
- matscipy.dislocation.gamma_line(unit_cell, calc=None, shift_dir=0, surface=2, size=[2, 2, 2], factor=15, n_dots=11, relax=True, fmax=0.01, return_images=False)
This function performs a calculation of a cross-sections in ‘shift_dir` of the generalized stacking fault (GSF) gamma surface with surface orientation. A gamma surface is defined as the energy variation when the crystal is cut along a particular plane and then one of the resulting parts is displaced along a particular direction. This quantity is related to the energy landscape of dislocations and provides data out of equilibrium, preserving the crystal state. For examples for the case of W and more details see section 4.2 and figure 2 in [J. Phys.: Condens. Matter 25 (2013) 395502 (15pp)] (http://iopscience.iop.org/0953-8984/25/39/395502)
- Parameters:
unit_cell (ase.Atoms) – Unit cell to construct gamma surface from. Should have a ase.calculator attached as calc in order to perform relaxation.
calc (ase.calculator) – if unit_cell.calc is None set unit_cell.calc to calc
shift_dir (int) – index of unit_cell axes to shift atoms
surface (int) – index of unit_cell axes to be the surface normal direction
size (list of ints) – start size of the cell
factor (int) – factor to increase the size of the cell along the surface normal direction
n_dots (int) – number of images along the gamma line
relax (bool) – flag to perform relaxation
fmax (float) – maximum force value for relaxation
return_images (bool) – flag to control if the atomic configurations are returned together with the results
- Returns:
deltas (np.array) – shift distance of every image in Angstroms
totens (np.array) – gamma surface energy in eV / Angstroms^2
images (list of ase.Atoms) – images along the gamma surface. Returned if return_images is True