Dislocation Quadrupoles

Dislocation Quadrupoles are structures designed with the goal of enabling the study of dislocation systems in periodic cells. In order for the system to be periodic, the cell must have a net Burgers vector of zero. To achieve this, we place two dislocations with Burgers vectors of equal magnitude and opposite sign (i.e. \(+\textbf{b}\) and \(-\textbf{b}\)) in the structure.

The code is designed to allow the generation of these structures in as small a structure is required for a given separation of the dislocations along the glide direction.

To start, lets look at the \(90^\circ\) partial dislocation in Diamond, using a dislocation cylinder approach documented in earlier docs and the same potential by D. Holland and M. Marder:

from matscipy.dislocation import DiamondGlide90degreePartial, get_elastic_constants
from matscipy.calculators.manybody.explicit_forms.stillinger_weber import StillingerWeber,\
                                                                Holland_Marder_PRL_80_746_Si
from matscipy.calculators.manybody import Manybody
calc = Manybody(**StillingerWeber(Holland_Marder_PRL_80_746_Si))

alat, C11, C12, C44 = get_elastic_constants(calculator=calc, symbol="Si", verbose=False)


Si_disloc = DiamondGlide90degreePartial(alat, C11, C12, C44, symbol="Si")

bulk, disloc = Si_disloc.build_cylinder(radius=20)

Si_disloc.view_cyl(disloc, scale=0.3, add_bonds=True)
/usr/lib/python3/dist-packages/scipy/__init__.py:146: UserWarning: A NumPy version >=1.17.3 and <1.25.0 is required for this version of SciPy (detected version 1.26.4
  warnings.warn(f"A NumPy version >={np_minversion} and <{np_maxversion}"

To build dislocation quadrupoles, we can use the Quadrupole class. The __init__() method of the class takes a dislocation class (e.g. DiamondGlide90degreePartial), followed by args and/or kwargs used to initialise that dislocation class (i.e, C11, C12, C44, and either alat or a bulk structure). From here, there are many class methods which generate different kinds of structures. Let’s first look at the basic quadrupole structure, using build_quadrupole

from matscipy.dislocation import Quadrupole

quad = Quadrupole(DiamondGlide90degreePartial, alat, C11, C12, C44, symbol="Si")

quad_bulk, quad_disloc = quad.build_quadrupole(glide_separation=4)

view = quad.view_quad(quad_disloc, scale=0.3, add_bonds=True)
view.control.zoom(0.2)
view
Periodic displacement iteration 0 -> max(|dr|) = 1.0963
Periodic displacement iteration 1 -> max(|dr|) = 0.1552
Periodic displacement iteration 2 -> max(|dr|) = 0.0059
Periodic displacement iteration 3 -> max(|dr|) = 0.0015
Periodic displacement iteration 4 -> max(|dr|) = 0.0006
Periodic displacements converged to r_tol=0.001 after 4 iterations

glide_separation is the most important argument to build_quadrupole, as it controls how far apart (in terms of glide distances) the quadrupole cores are, and thus controls the size of the quadrupole cell.

By default, the cell is constructed such that each dislocation core is surrounded by four cores of opposite sign, and the periodicity of the system tries to make the distances between the central core and it’s opposites as equal as possible (i.e. that the cores form two sublattices with square packing).

Dislocation Glide in quadrupoles

As quadrupoles are extremely useful for modelling dislocations using plane-wave DFT, it can be convenient to be able to set up initial guesses for complex processes such as dislocation glide. In this scenario, we assume that the full infinite dislocation line glides in unison, ignoring the true “kink” process.

We can use the function build_glide_quadrupoles to construct a set of images for this system, which can optionally model the glide of either the “left” (\(+\mathbf{b}\)) or “right” (\(-\mathbf{b}\)) dislocation cores, or both at once.

num_images = 11

glide_quads = quad.build_glide_quadrupoles(nims=num_images, 
                                            glide_left=True, # Allow left dislocation glide
                                            glide_right=True, # Allow right dislocation glide
                                            glide_separation=6,
                                            verbose=False)

To visualise the glide structures, we will combine ase’s plot_atoms to convert a structure to a matplotlib plot, and then use FuncAnimation to animate the glide structures:

def animate_glide(images, diamond=True):
    from ase.visualize.plot import plot_atoms
    import matplotlib.pyplot as plt
    from matplotlib.animation import FuncAnimation
    from matscipy.utils import get_structure_types
    from visualisation import show_HTML

    fig, ax = plt.subplots(figsize=(10, 15))
    ax.axis("off")

    # Add extra reps of start and end points for clarity
    anim_images = [images[0]] * 3 + images + [images[-1]] * 3

    def plot_frame(framedata):
        ax.clear()
        # Plot an individual frame of the animation 
        framenum, atoms = framedata

        # get CNA colours to enhance plot
        atom_labels, struct_names, colors = get_structure_types(atoms, 
                                                                diamond_structure=diamond)
        atom_colors = [colors[atom_label] for atom_label in atom_labels]

        plot_atoms(atoms, ax=ax, colors=atom_colors)


    animation = FuncAnimation(fig, plot_frame, frames=enumerate(anim_images),
                                save_count=len(anim_images),
                                init_func=lambda: None,
                                interval=200)
    
    # Need to convert animation to HTML in order for it to be visible on the docs
    return show_HTML(animation)

animate_glide(glide_quads)

Quadrupole Kink

Dislocation kink is commonly a more energetically accessible way for a dislocation to move in the glide direction. Rather than the full dislocation line gliding as one, the kink mechanism involves a small section of the dislocation line nucleating out by a glide vector, forming a pair of kinks. This pair of kinks can then diffuse along the direction of the dislocation line, advancing the line locally by a single glide vector.

In the quadrupole cells, we model an infinite chain of single kinks, where the dislocation line advances by one glide direction every periodic repetition of the structure in the Z direction. This is distinct from the kink-pair mechanism, which would need two kinks (which migrate the dislocation line in opposite directions), however the single kink structures can be useful in making kink behaviour accessible to the small scales required by expensive methods such as DFT.

To actually construct kink structures from a quadrupole object, we use the build_kink_quadrupole() method of Quadrupole. The main two important parameters to define are invert_direction, which changes the directionality of the kinks (some systems have distinct “left” kinks and “right” kinks), and z_reps, which both changes the size of the final structure, as well as the initial guess for the length of the kink.

Note

Because the kink advances the system by a single glide vector, and cells are typically multiple glide vectors long, the structure will not initially be periodic. To solve this, some atomic layers are removed (as controlled by layer_tol). In multispecies systems, this may disrupt the atomic composition of the end structure.

kink_quad = quad.build_kink_quadrupole(z_reps=5,
                                       glide_separation=4,
                                       verbose=False, layer_tol=1e-2)

quad.view_quad(kink_quad, scale=0.6, add_bonds=False)

Quadrupoles in other systems

The code has many other features to the general purpose methods described above, which are intended to support some common choices and desires which arise in other dislocation systems. Below are a few curated examples of dislocation systems, and how the output can be modified by additional parameters.

BCC Screw in W

Taking the perfect screw dislocation example in BCC W (and using the same Embedded Atom Potential from Marinica et. al. 2013 paper (version EAM4) as was used in previous dislocation documentation), we can see that the default behaviour of the code is to produce what can be called an “easy-hard” quadrupole, which is to say that one dislocation is in the “easy” core (with 3 atoms in white), and the other dislocation is in the “hard” core (with 6 atoms in white).

from matscipy.calculators.eam import EAM
from matscipy.dislocation import BCCScrew111Dislocation
eam_calc = EAM("../../tests/w_eam4.fs")
alat, C11, C12, C44 = get_elastic_constants(calculator=eam_calc, symbol="W", verbose=False)

quad = Quadrupole(BCCScrew111Dislocation, alat, C11, C12, C44, symbol="W")

W_bulk, W_quad = quad.build_quadrupole(glide_separation=8, verbose=False)

quad.view_quad(W_quad)

Whilst this kind of structure contains a lot of information about both kinds of cores, typically the “easy” core is much lower in energy, and thus it may be desirable to model an “easy-easy” quadrupole. To do this, we can just add an offset to the left core (half the glide distance in this case), using the left_offset parameter.

import numpy as np
W_bulk, W_quad = quad.build_quadrupole(glide_separation=8, verbose=False,
                                       left_offset=np.array([quad.glide_distance/2, 0, 0]))

quad.view_quad(W_quad)

Dissociation of Screw dislocations in FCC Ni Quadrupoles

Dissociation of dislocation cores is another phenomenon that me be desirable in a quadrupole structure, especially as quadrupoles of just the partials can’t give much insight into the physics controlling when the partial cores are very close.

We can build quadrupole structures containing these dissociated cores by passing the partial_distance argument to build_quadrupole, similarly to how we would when building a cylinder for that dissociated system.

from matscipy.dislocation import FCCScrew110Dislocation

eam_calc = EAM("../../tests/FeCuNi.eam.alloy")

# the function accepts any ASE type of calculator
alat, C11, C12, C44 = get_elastic_constants(calculator=eam_calc, symbol="Ni", verbose=False)

quad = Quadrupole(FCCScrew110Dislocation, alat, C11, C12, C44, symbol="Ni")

bulk, Ni_quad = quad.build_quadrupole(glide_separation=14, partial_distance=4, 
                                      disp_tol=1e-2) # Lower the tolerance for the periodic displacement convergence

quad.view_quad(Ni_quad)
Periodic displacement iteration 0 -> max(|dr|) = 1.2211
Periodic displacement iteration 1 -> max(|dr|) = 0.1928
Periodic displacement iteration 2 -> max(|dr|) = 0.0053
Periodic displacements converged to r_tol=0.01 after 2 iterations