dislocation_monopole - Methodology and code

Python imports

[1]:
# Standard library imports
from pathlib import Path
import shutil
import datetime
from copy import deepcopy
from math import floor
from typing import Optional, Tuple, Union
import random

# http://www.numpy.org/
import numpy as np

# https://ipython.org/
from IPython.display import display, Markdown

# https://github.com/usnistgov/atomman
import atomman as am
import atomman.lammps as lmp
import atomman.unitconvert as uc
from atomman.tools import filltemplate

# https://github.com/usnistgov/iprPy
import iprPy
from iprPy.tools import read_calc_file

print('Notebook last executed on', datetime.date.today(), 'using iprPy version', iprPy.__version__)
Notebook last executed on 2023-07-31 using iprPy version 0.11.6

1. Load calculation and view description

1.1. Load the calculation

[2]:
# Load the calculation being demoed
calculation = iprPy.load_calculation('dislocation_monopole')

1.2. Display calculation description and theory

[3]:
# Display main docs and theory
display(Markdown(calculation.maindoc))
display(Markdown(calculation.theorydoc))

dislocation_monopole calculation style

Lucas M. Hale, lucas.hale@nist.gov, Materials Science and Engineering Division, NIST.

Introduction

The dislocation_monopole calculation style calculation inserts a dislocation monopole into an atomic system using the anisotropic elasticity solution for a perfectly straight dislocation. The system is divided into two regions: a boundary region at the system’s edges perpendicular to the dislocation line, and an active region around the dislocation. After inserting the dislocation, the atoms within the active region of the system are relaxed while the positions of the atoms in the boundary region are held fixed at the elasticity solution. The relaxed dislocation system and corresponding dislocation-free base systems are retained in the calculation’s archived record. Various properties associated with the dislocation’s elasticity solution are recorded in the calculation’s results record.

Version notes
  • 2018-09-25: Notebook added

  • 2019-07-30: Method and Notebook updated for iprPy version 0.9.

  • 2020-09-22: Notebook updated to reflect that calculation method has changed to now use atomman.defect.Dislocation.

  • 2022-03-11: Notebook updated to reflect version 0.11.

Additional dependencies
Disclaimers
  • NIST disclaimers

  • This calculation method holds the boundary atoms fixed during the relaxation process which results in a mismatch stress at the border between the active and fixed regions that interacts with the dislocation. Increasing the distance between the dislocation and the boundary region, i.e. increasing the system size, will reduce the influence of the mismatch stresses.

  • The boundary atoms are fixed at the elasticity solution, which assumes the dislocation to be compact (not spread out, dissociated into partials) and to remain at the center of the system. An alternate simulation design or boundary region method should be used if you want accurate simulations of dislocations with wide cores and/or in which the dislocation moves long distances.

  • The calculation allows for the system to be relaxed either using only static energy/force minimizations or with molecular dynamic steps followed by a minimization. Only performing a static relaxation is considerably faster than performing a dynamic relaxation, but the dynamic relaxation is more likely to result in a better final dislocation structure.

Method and Theory

Stroh theory

A detailed description of the Stroh method to solve the Eshelby solution for an anisotropic straight dislocation can be found in the atomman documentation.

Orientation

One of the most important considerations in defining an atomistic system containing a dislocation monopole system is the system’s orientation. In particular, care is needed to properly align the system’s Cartesian axes, \(x, y, z\), the system’s box vectors, \(a, b, c\), and the Stroh solution’s axes, \(u, m, n\).

  • In order for the dislocation to remain perfectly straight in the atomic system, the dislocation line, \(u\), must correspond to one of the system’s box vectors. The resulting dislocation monopole system will be periodic along the box direction corresponding to \(u\), and non-periodic in the other two box directions.

  • The Stroh solution is invariant along the dislocation line direction, \(u\), thereby the solution is 2 dimensional. \(m\) and \(n\) are the unit vectors for the 2D axis system used by the Stroh solution. As such, \(u, m\) and \(n\) are all normal to each other. The plane defined by the \(um\) vectors is the dislocation’s slip plane, i.e. \(n\) is normal to the slip plane.

  • For LAMMPS simulations, the system’s box vectors are limited such that \(a\) is parallel to the \(x\)-axis, and \(b\) is in the \(xy\)-plane (i.e. has no \(z\) component). Based on this and the previous two points, the most convenient, and therefore the default, orientation for a generic dislocation is to

    • Make \(u\) and \(a\) parallel, which also places \(u\) parallel to the \(x\)-axis.

    • Select \(b\) such that it is within the slip plane. As \(u\) and \(a\) must also be in the slip plane, the plane itself is defined by \(a \times b\).

    • Given that neither \(a\) nor \(b\) have \(z\) components, the normal to the slip plane will be perpendicular to \(z\). From this, it naturally follows that \(m\) can be taken as parallel to the \(y\)-axis, and \(n\) parallel to the \(z\)-axis.

Calculation methodology
  1. An initial system is generated based on the loaded system and uvw, atomshift, and sizemults input parameters. This initial system is saved as base.dump.

  2. The atomman.defect.Stroh class is used to obtain the dislocation solution based on the defined \(m\) and \(n\) vectors, \(C_{ij}\), and the dislocation’s Burgers vector, \(b_i\).

  3. The dislocation is inserted into the system by displacing all atoms according to the Stroh solution’s displacements. The dislocation line is placed parallel to the specified dislocation_lineboxvector and includes the Cartesian point (0, 0, 0).

  4. The box dimension parallel to the dislocation line is kept periodic, and the other two box directions are made non-periodic. A boundary region is defined that is at least bwidth thick at the edges of the two non-periodic directions, such that the shape of the active region corresponds to the bshape parameter. Atoms in the boundary region are identified by altering their integer atomic types.

  5. The dislocation is relaxed by performing a LAMMPS simulation. The atoms in the active region are allowed to relax either with nvt integration followed by an energy/force minimization, or with just an energy/force minimization. The atoms in the boundary region are held fixed at the elastic solution. The resulting relaxed system is saved as disl.dump.

  6. Parameters associated with the Stroh solution are saved to the results record.

2. Define calculation functions and generate files

This section defines the calculation functions and associated resource files exactly as they exist inside the iprPy package. This allows for the code used to be directly visible and modifiable by anyone looking to see how it works.

2.1. dislocation_monopole()

This is the primary function for the calculation. The version of this function built in iprPy can be accessed by calling the calc() method of an object of the associated calculation class.

[4]:
def dislocation_monopole(lammps_command: str,
                         ucell: am.System,
                         potential: lmp.Potential,
                         C: am.ElasticConstants,
                         burgers: Union[list, np.ndarray],
                         ξ_uvw: Union[list, np.ndarray],
                         slip_hkl: Union[list, np.ndarray],
                         mpi_command: Optional[str] = None,
                         m: Union[list, np.ndarray] = [0,1,0],
                         n: Union[list, np.ndarray] = [0,0,1],
                         sizemults=None,
                         amin: float = None,
                         bmin: float = None,
                         cmin: float = None,
                         shift: Union[list, np.ndarray, None] = None,
                         shiftscale: bool = False,
                         shiftindex: int = None,
                         tol: float = 1e-8,
                         etol: float = 0.0,
                         ftol: float = 0.0,
                         maxiter: int = 10000,
                         maxeval: int = 100000,
                         dmax: float = uc.set_in_units(0.01, 'angstrom'),
                         annealtemp: float = 0.0,
                         annealsteps: Optional[int] = None,
                         randomseed: Optional[int] = None,
                         boundaryshape: str = 'cylinder',
                         boundarywidth: float = 0.0,
                         boundaryscale: bool = False) -> dict:
    """
    Creates and relaxes a dislocation monopole system.

    Parameters
    ----------
    lammps_command :str
        Command for running LAMMPS.
    ucell : atomman.System
        The unit cell to use as the seed for generating the dislocation
        monopole system.
    potential : atomman.lammps.Potential
        The LAMMPS implemented potential to use.
    C : atomman.ElasticConstants
        The elastic constants associated with the bulk crystal structure
        for ucell.
    burgers : array-like object
        The dislocation's Burgers vector given as a Miller or
        Miller-Bravais vector relative to ucell.
    ξ_uvw : array-like object
        The dislocation's line direction given as a Miller or
        Miller-Bravais vector relative to ucell.
    slip_hkl : array-like object
        The dislocation's slip plane given as a Miller or Miller-Bravais
        plane relative to ucell.
    mpi_command : str or None, optional
        The MPI command for running LAMMPS in parallel.  If not given, LAMMPS
        will run serially.
    m : array-like object, optional
        The m unit vector for the dislocation solution.  m, n, and ξ
        (dislocation line) should be right-hand orthogonal.  Default value
        is [0,1,0] (y-axis).
    n : array-like object, optional
        The n unit vector for the dislocation solution.  m, n, and ξ
        (dislocation line) should be right-hand orthogonal.  Default value
        is [0,0,1] (z-axis). n is normal to the dislocation slip plane.
    sizemults : tuple, optional
        The size multipliers to use when generating the system.  Values are
        limited to being positive integers.  The multipliers for the two
        non-periodic directions must be even.  If not given, the default
        multipliers will be 2 for the non-periodic directions and 1 for the
        periodic direction.
    amin : float, optional
        A minimum thickness to use for the a box vector direction of the
        final system.  Default value is 0.0.  For the non-periodic
        directions, the resulting vector multiplier will be even.  If both
        amin and sizemults is given, then the larger multiplier for the two
        will be used.
    bmin : float, optional
        A minimum thickness to use for the b box vector direction of the
        final system.  Default value is 0.0.  For the non-periodic
        directions, the resulting vector multiplier will be even.  If both
        bmin and sizemults is given, then the larger multiplier for the two
        will be used.
    cmin : float, optional
        A minimum thickness to use for the c box vector direction of the
        final system.  Default value is 0.0.  For the non-periodic
        directions, the resulting vector multiplier will be even.  If both
        cmin and sizemults is given, then the larger multiplier for the two
        will be used.
    shift : float, optional
        A rigid body shift to apply to the rotated cell prior to inserting
        the dislocation.  Should be selected such that the ideal slip plane
        does not correspond to any atomic planes.  Is taken as absolute if
        shiftscale is False, or relative to the rotated cell's box vectors
        if shiftscale is True.  Cannot be given with shiftindex.  If
        neither shift nor shiftindex is given then shiftindex = 0 is used.
    shiftindex : float, optional
        The index of the identified optimum shifts based on the rotated
        cell to use.  Different values allow for the selection of different
        atomic planes neighboring the slip plane.  Note that shiftindex
        values only apply shifts normal to the slip plane; best shifts for
        non-planar dislocations (like bcc screw) may also need a shift in
        the slip plane.  Cannot be given with shiftindex.  If neither shift
        nor shiftindex is given then shiftindex = 0 is used.
    shiftscale : bool, optional
        If False (default), a given shift value will be taken as absolute
        Cartesian.  If True, a given shift will be taken relative to the
        rotated cell's box vectors.
    tol : float
        A cutoff tolerance used with obtaining the dislocation solution.
        Only needs to be changed if there are issues with obtaining a
        solution.
    etol : float, optional
        The energy tolerance for the structure minimization. This value is
        unitless. Default is 0.0.
    ftol : float, optional
        The force tolerance for the structure minimization. This value is in
        units of force. Default is 0.0.
    maxiter : int, optional
        The maximum number of minimization iterations to use. Default is
        10000.
    maxeval : int, optional
        The maximum number of minimization evaluations to use. Default is
        100000.
    dmax : float, optional
        The maximum distance in length units that any atom is allowed to relax
        in any direction during a single minimization iteration. Default is
        0.01 Angstroms.
    annealtemp : float, optional
        The temperature to perform a dynamic relaxation at. Default is 0.0,
        which will skip the dynamic relaxation.
    annealsteps : int, optional
        The number of time steps to run the dynamic relaxation for.  Default
        is None, which will run for 10000 steps if annealtemp is not 0.0.
    randomseed : int or None, optional
        Random number seed used by LAMMPS in creating velocities and with
        the Langevin thermostat.  Default is None which will select a
        random int between 1 and 900000000.
    boundaryshape : str, optional
        Indicates the shape of the boundary region to use.  Options are
        'cylinder' (default) and 'box'.  For 'cylinder', the non-boundary
        region is defined by a cylinder with axis along the dislocation
        line and a radius that ensures the boundary is at least
        boundarywidth thick.  For 'box', the boundary region will be
        exactly boundarywidth thick all around.
    boundarywidth : float, optional
        The width of the boundary region to apply.  Default value is 0.0,
        i.e. no boundary region.  All atoms in the boundary region will
        have their atype values changed.
    boundaryscale : bool, optional
        If False (Default), the boundarywidth will be taken as absolute.
        If True, the boundarywidth will be taken relative to the magnitude
        of the unit cell's a box vector.

    Returns
    -------
    dict
        Dictionary of results consisting of keys:

        - **'dumpfile_base'** (*str*) - The filename of the LAMMPS dump file
          for the relaxed base system.
        - **'symbols_base'** (*list of str*) - The list of element-model
          symbols for the Potential that correspond to the base system's
          atypes.
        - **'dumpfile_disl'** (*str*) - The filename of the LAMMPS dump file
          for the relaxed dislocation monopole system.
        - **'symbols_disl'** (*list of str*) - The list of element-model
          symbols for the Potential that correspond to the dislocation
          monopole system's atypes.
        - **'dislocation'** (*atomman.defect.Dislocation*) - The Dislocation
          object used to generate the monopole system.
        - **'E_total_disl'** (*float*) - The total potential energy of the
          dislocation monopole system.
    """
    # Construct dislocation configuration generator
    dislocation = am.defect.Dislocation(ucell, C, burgers, ξ_uvw, slip_hkl,
                                        m=m, n=n, shift=shift, shiftindex=shiftindex,
                                        shiftscale=shiftscale, tol=tol)

    # Generate the base and dislocation systems
    base_system, disl_system = dislocation.monopole(sizemults=sizemults,
                                                    amin=amin, bmin=bmin, cmin=cmin,
                                                    shift=shift,
                                                    shiftindex=shiftindex,
                                                    shiftscale=shiftscale,
                                                    boundaryshape=boundaryshape,
                                                    boundarywidth=boundarywidth,
                                                    boundaryscale=boundaryscale,
                                                    return_base_system=True)

    # Initialize results dict
    results_dict = {}

    # Save initial perfect system
    base_system.dump('atom_dump', f='base.dump')
    results_dict['dumpfile_base'] = 'base.dump'
    results_dict['symbols_base'] = base_system.symbols

    # Save dislocation generator
    results_dict['dislocation'] = dislocation

    # Relax system
    relaxed = disl_relax(lammps_command, disl_system, potential,
                         mpi_command=mpi_command, annealtemp=annealtemp,
                         annealsteps=annealsteps, randomseed=randomseed,
                         etol=etol, ftol=ftol, maxiter=maxiter,
                         maxeval=maxeval, dmax=dmax)

    # Save relaxed dislocation system with original box vects
    system_disl = am.load('atom_dump', relaxed['dumpfile'], symbols=disl_system.symbols)

    system_disl.box_set(vects=disl_system.box.vects, origin=disl_system.box.origin)
    system_disl.dump('atom_dump', f='disl.dump')
    results_dict['dumpfile_disl'] = 'disl.dump'
    results_dict['symbols_disl'] = system_disl.symbols

    results_dict['E_total_disl'] = relaxed['E_total']

    # Cleanup files
    Path('0.dump').unlink()
    Path(relaxed['dumpfile']).unlink()
    for dumpjsonfile in Path('.').glob('*.dump.json'):
        dumpjsonfile.unlink()

    return results_dict

2.2. disl_relax()

[5]:
def disl_relax(lammps_command: str,
               system: am.System,
               potential: lmp.Potential,
               mpi_command: Optional[str] = None,
               annealtemp: float = 0.0,
               annealsteps: Optional[int] = None,
               randomseed: Optional[int] = None,
               etol: float = 0.0,
               ftol: float = 1e-6,
               maxiter: int = 10000,
               maxeval: int = 100000,
               dmax: float = uc.set_in_units(0.01, 'angstrom')) -> dict:
    """
    Sets up and runs the disl_relax.in LAMMPS script for relaxing a
    dislocation monopole system.

    Parameters
    ----------
    lammps_command :str
        Command for running LAMMPS.
    system : atomman.System
        The system to perform the calculation on.
    potential : atomman.lammps.Potential
        The LAMMPS implemented potential to use.
    mpi_command : str, optional
        The MPI command for running LAMMPS in parallel.  If not given, LAMMPS
        will run serially.
    annealtemp : float, optional
        The temperature to perform a dynamic relaxation at. Default is 0.0,
        which will skip the dynamic relaxation.
    annealsteps : int, optional
        The number of time steps to run the dynamic relaxation for.  Default
        is None, which will run for 10000 steps if annealtemp is not 0.0.
    randomseed : int or None, optional
        Random number seed used by LAMMPS in creating velocities and with
        the Langevin thermostat.  Default is None which will select a
        random int between 1 and 900000000.
    etol : float, optional
        The energy tolerance for the structure minimization. This value is
        unitless. Default is 0.0.
    ftol : float, optional
        The force tolerance for the structure minimization. This value is in
        units of force. Default is 0.0.
    maxiter : int, optional
        The maximum number of minimization iterations to use default is
        10000.
    maxeval : int, optional
        The maximum number of minimization evaluations to use default is
        100000.
    dmax : float, optional
        The maximum distance in length units that any atom is allowed to relax
        in any direction during a single minimization iteration default is
        0.01 Angstroms.

    Returns
    -------
    dict
        Dictionary of results consisting of keys:

        - **'logfile'** (*str*) - The name of the LAMMPS log file.
        - **'dumpfile'** (*str*) - The name of the LAMMPS dump file
          for the relaxed system.
        - **'E_total'** (*float*) - The total potential energy for the
          relaxed system.
    """

    # Get lammps units
    lammps_units = lmp.style.unit(potential.units)

    #Get lammps version date
    lammps_date = lmp.checkversion(lammps_command)['date']

    # Define lammps variables
    lammps_variables = {}
    system_info = system.dump('atom_data', f='system.dat',
                              potential=potential)
    lammps_variables['atomman_system_pair_info'] = system_info
    lammps_variables['anneal_info'] = anneal_info(annealtemp, annealsteps,
                                                  randomseed, potential.units)
    lammps_variables['etol'] = etol
    lammps_variables['ftol'] = uc.get_in_units(ftol, lammps_units['force'])
    lammps_variables['maxiter'] = maxiter
    lammps_variables['maxeval'] = maxeval
    lammps_variables['dmax'] = dmax
    lammps_variables['group_move'] = ' '.join(np.array(range(1, system.natypes // 2 + 1), dtype=str))

    # Set dump_modify format based on dump_modify_version
    if lammps_date < datetime.date(2016, 8, 3):
        lammps_variables['dump_modify_format'] = '"%d %d %.13e %.13e %.13e %.13e"'
    else:
        lammps_variables['dump_modify_format'] = 'float %.13e'

    # Write lammps input script
    lammps_script = 'disl_relax.in'
    template = read_calc_file('iprPy.calculation.dislocation_monopole',
                              'disl_relax.template')
    with open(lammps_script, 'w') as f:
        f.write(filltemplate(template, lammps_variables, '<', '>'))

    # Run LAMMPS
    output = lmp.run(lammps_command, script_name=lammps_script,
                     mpi_command=mpi_command)
    thermo = output.simulations[-1]['thermo']

    # Extract output values
    results = {}
    results['logfile'] = 'log.lammps'
    results['dumpfile'] = '%i.dump' % thermo.Step.values[-1]
    results['E_total'] = uc.set_in_units(thermo.PotEng.values[-1],
                                         lammps_units['energy'])

    return results

2.3. anneal_info()

[6]:
def anneal_info(temperature: float = 0.0,
                runsteps: Optional[int] = None,
                randomseed: Optional[int] = None,
                units: str = 'metal') -> str:
    """
    Generates LAMMPS commands for thermo anneal.

    Parameters
    ----------
    temperature : float, optional
        The temperature to relax at (default is 0.0).
    randomseed : int or None, optional
        Random number seed used by LAMMPS in creating velocities and with
        the Langevin thermostat.  (Default is None which will select a
        random int between 1 and 900000000.)
    units : str, optional
        The LAMMPS units style to use (default is 'metal').

    Returns
    -------
    str
        The generated LAMMPS input lines for performing a dynamic relax.
        Will be '' if temperature==0.0.
    """
    # Return nothing if temperature is 0.0 (don't do thermo anneal)
    if temperature == 0.0:
        return ''

    # Generate velocity, fix nvt, and run LAMMPS command lines
    else:
        if randomseed is None:
            randomseed = random.randint(1, 900000000)
        if runsteps is None:
            runsteps = 10000

        start_temp = 2 * temperature
        tdamp = 100 * lmp.style.timestep(units)
        timestep = lmp.style.timestep(units)
        info = '\n'.join([
            'velocity move create %f %i mom yes rot yes dist gaussian' % (start_temp, randomseed),
            'fix nvt all nvt temp %f %f %f' % (temperature, temperature,
                                               tdamp),
            'timestep %f' % (timestep),
            'thermo %i' % (runsteps),
            'run %i' % (runsteps),
            ])

    return info

2.4. disl_relax.template file

[7]:
with open('disl_relax.template', 'w') as f:
    f.write("""#LAMMPS input script for relaxing a dislocation

<atomman_system_pair_info>

group move type <group_move>
group hold subtract all move

compute peatom all pe/atom

dump first all custom <maxeval> *.dump id type x y z c_peatom
dump_modify first format <dump_modify_format>
thermo_style custom step pe

fix nomove hold setforce 0.0 0.0 0.0

<anneal_info>

min_modify dmax <dmax>

minimize <etol> <ftol> <maxiter> <maxeval>""")

3. Specify input parameters

3.1. System-specific paths

  • lammps_command is the LAMMPS command to use (required).

  • mpi_command MPI command for running LAMMPS in parallel. A value of None will run simulations serially.

[8]:
lammps_command = 'lmp'
mpi_command = None

# Optional: check that LAMMPS works and show its version
print(f'LAMMPS version = {am.lammps.checkversion(lammps_command)["version"]}')
LAMMPS version = 15 Sep 2022

3.2. Interatomic potential

  • potential_name gives the name of the potential_LAMMPS reference record in the iprPy library to use for the calculation.

  • potential is an atomman.lammps.Potential object (required).

[10]:
potential_name = '1999--Mishin-Y--Ni--LAMMPS--ipr1'

# Retrieve potential and parameter file(s) using atomman
potential = am.load_lammps_potential(id=potential_name, getfiles=True)

3.3. Initial unit cell system

  • ucell is an atomman.System representing a fundamental unit cell of the system (required). Here, this is generated using the load parameters and symbols.

[11]:
# Create ucell by loading prototype record
ucell = am.load('crystal', potential=potential, family='A1--Cu--fcc')

print(ucell)
avect =  [ 3.520,  0.000,  0.000]
bvect =  [ 0.000,  3.520,  0.000]
cvect =  [ 0.000,  0.000,  3.520]
origin = [ 0.000,  0.000,  0.000]
natoms = 4
natypes = 1
symbols = ('Ni',)
pbc = [ True  True  True]
per-atom properties = ['atype', 'pos']
     id   atype  pos[0]  pos[1]  pos[2]
      0       1   0.000   0.000   0.000
      1       1   0.000   1.760   1.760
      2       1   1.760   0.000   1.760
      3       1   1.760   1.760   0.000

3.4. Elastic constants

  • C_dict is a dictionary containing the unique elastic constants for the potential and crystal structure defined above.

  • C is an atomman.ElasticConstants object built from C_dict.

[12]:
C_dict = {}
C_dict['C11'] = uc.set_in_units(247.86, 'GPa')
C_dict['C12'] = uc.set_in_units(147.83, 'GPa')
C_dict['C44'] = uc.set_in_units(124.84, 'GPa')

# -------------- Derived parameters -------------- #
# Build ElasticConstants object from C_dict terms
C = am.ElasticConstants(**C_dict)

3.5. Defect parameters

  • burgers is the crystallographic Miller Burgers vector for the dislocation.

  • ξ_uvw is the Miller [uvw] line vector direction for the dislocation. The angle between burgers and ξ_uvw determines the dislocation’s character

  • slip_hkl is the Miller (hkl) slip plane for the dislocation.

  • m is the Cartesian vector of the final system that the dislocation solution’s m vector (in-plane, perpendicular to ξ) should align with. Limited to being parallel to one of the three Cartesian axes.

  • n is the Cartesian vector of the final system that the dislocation solution’s n vector (slip plane normal) should align with. Limited to being parallel to one of the three Cartesian axes.

  • shift is a rigid body shift to apply to the atoms in the system. This controls how the atomic positions align with the ideal position of the dislocation core, which is at coordinates (0,0) for the two Cartesian axes aligned with m and n.

  • shiftscale allows for shift to be defined relative to the cell created by rotating ucell to coincide with the dislocation solution orientation. This is useful as it allows for shift values to be defined relative to the defect type and crystal prototype rather than on a per-crystal basis.

  • shiftindex alternate to specifying shift values, the shiftindex allows for one of the identified suggested shift values to be used that will position the slip plane halfway between two planes of atoms. Note that shiftindex values only shift atoms in the slip plane normal direction and may not be the ideal positions for some dislocation cores.

[13]:
# fcc a/2 <110>{111} dislocations
burgers = 0.5 * np.array([ 1., -1., 0.])
slip_hkl = np.array([1, 1, 1])

# Line direction determines dislocation character
ξ_uvw = np.array([ 1, 1,-2]) # 90 degree edge
#ξ_uvw = np.array([ 1, 0,-1]) # 60 degree mixed
#ξ_uvw = np.array([ 1,-2, 1]) # 30 degree mixed
#ξ_uvw = np.array([ 1,-1, 0]) # 0 degree screw

# Best choice for m + n as it works for non-cubic systems
m = [0,1,0]
n = [0,0,1]

# Specify shift or shiftindex
shift = None
shiftscale = True
shiftindex = 0

3.6. Calculation-specific parameters

  • boundarywidth sets the minimum width of the fixed-atom boundary region.

  • boundaryscale flag indicating if boundarywidth is absolute (False) or scaled relative to ucell’s a lattice parameter (True).

  • boundaryshape specifies what shape to make the fixed-atom boundary region.

    • cylinder will create a cylindrical active region.

    • box will create a rectangular active region.

  • annealtemperature is the temperature at which to relax (anneal) the dislocation system. If annealtemperature is 0.0, then only a static relaxation will be performed. Default value is 0.0.

  • annealsteps is the number of nvt iteration steps to perform at the given temperature. Default value is 0 if annealtemperature is zero, 10000 otherwise.

  • randomseed allows for the random seed used in generating initial atomic velocities for a dynamic relaxation to be specified. This is an integer between 1 and 900000000. Default value is None, which will randomly pick a number in that range.

  • energytolerance is the energy tolerance to use during the minimizations. This is unitless.

  • forcetolerance is the force tolerance to use during the minimizations. This is in energy/length units.

  • maxiterations is the maximum number of minimization iterations to use.

  • maxevaluations is the maximum number of minimization evaluations to use.

  • maxatommotion is the largest distance that an atom is allowed to move during a minimization iteration. This is in length units.

[14]:
# Specify the boundary shape to use
boundarywidth = 3
boundaryscale = True
boundaryshape = 'cylinder'

# Specify MD anneal parameters
annealtemperature = 50.0
annealsteps = 10000
randomseed = None

# Specify minimization parameters
energytolerance = 0.0
forcetolerance = uc.set_in_units(1e-6, 'eV/angstrom')
maxiterations = 10000
maxevaluations = 100000
maxatommotion = uc.set_in_units(0.01, 'angstrom')

3.7. System modifications

  • sizemults list of three sets of two integers specifying how many times the ucell vectors of \(a\), \(b\) and \(c\) are replicated in positive and negative directions when creating system.

  • amin, bmin, cmin specifies minimum widths in the three directions. The corresponding sizemults will be increased to ensure these widths if needed.

[15]:
sizemults = [1, 40, 40]
amin = uc.set_in_units(0.0, 'angstrom')
bmin = uc.set_in_units(0.0, 'angstrom')
cmin = uc.set_in_units(0.0, 'angstrom')

4. Run calculation and view results

4.1. Run calculation

All primary calculation method functions take a series of inputs and return a dictionary of outputs.

[16]:
results_dict = dislocation_monopole(lammps_command, ucell, potential, C,
                                    burgers, ξ_uvw, slip_hkl,
                                    mpi_command = mpi_command,
                                    m = m,
                                    n = n,
                                    shift = shift,
                                    shiftscale = shiftscale,
                                    shiftindex = shiftindex,
                                    sizemults = sizemults,
                                    amin = amin,
                                    bmin = bmin,
                                    cmin = cmin,
                                    etol = energytolerance,
                                    ftol = forcetolerance,
                                    maxiter = maxiterations,
                                    maxeval = maxevaluations,
                                    dmax = maxatommotion,
                                    annealtemp = annealtemperature,
                                    annealsteps = annealsteps,
                                    randomseed = randomseed,
                                    boundaryshape = boundaryshape,
                                    boundarywidth = boundarywidth,
                                    boundaryscale = boundaryscale)
print(results_dict.keys())
dict_keys(['dumpfile_base', 'symbols_base', 'dislocation', 'dumpfile_disl', 'symbols_disl', 'E_total_disl'])

4.2. Report results

Values returned in the results_dict:

  • ‘dumpfile_base’ (str) - The filename of the LAMMPS dump file for the relaxed base system.

  • ‘symbols_base’ (list of str) - The list of element-model symbols for the Potential that correspond to the base system’s atypes.

  • ‘dumpfile_disl’ (str) - The filename of the LAMMPS dump file for the relaxed dislocation monopole system.

  • ‘symbols_disl’ (list of str) - The list of element-model symbols for the Potential that correspond to the dislocation monopole system’s atypes.

  • ‘dislocation’ (atomman.defect.Dislocation) - The Dislocation object used to generate the monopole system.

  • ‘E_total_disl’ (float) - The total potential energy of the dislocation monopole system.

[17]:
length_unit = 'angstrom'
energy_unit = 'eV'
pressure_unit = 'GPa'
energy_per_length_unit = energy_unit+'/'+length_unit
[18]:
print('pre-ln factor (alpha) =', uc.get_in_units(results_dict['dislocation'].dislsol.preln, energy_per_length_unit), energy_per_length_unit)
print('K_tensor =')
print(uc.get_in_units(results_dict['dislocation'].dislsol.K_tensor, pressure_unit), pressure_unit)
pre-ln factor (alpha) = 0.38115202169620926 eV/angstrom
K_tensor =
[[ 81.07544348   0.         -12.76027037]
 [  0.         123.86918865   0.        ]
 [-12.76027037   0.         127.31323026]] GPa
[19]:
print('Perfect system is saved as', results_dict['dumpfile_base'])
print('with symbols', results_dict['symbols_base'])
print('Defect system is saved as',  results_dict['dumpfile_disl'])
print('with symbols', results_dict['symbols_disl'])
Perfect system is saved as base.dump
with symbols ('Ni',)
Defect system is saved as disl.dump
with symbols ('Ni', 'Ni')