4. Get started

4.1. Python API

phonoLAMMPS has been designed to have a very similar python interface to phonopy. For this reason the python objects generated by phonoLAMMPS can be used with phonopy functions. The procedure to obtain the force constants from LAMMPS forces is the following:

  1. Loading phonoLAMMPS module

    from phonolammps import Phonolammps
    

Creating an instance of the main phonoLAMMPS class. In this call you have to introduce a lammps input that contains the definition of the crystal unit cell the unitcell and empirical potential/s. For this you have two options:

2a. Use a LAMMPS input file (in.lammps)

phlammps = Phonolammps('in.lammps',
                       supercell_matrix=[[2, 0, 0],
                                         [0, 2, 0],
                                         [0, 0, 2]],
                       primitive_matrix=[[0.0, 0.5 ,0.5],
                                         [0.5, 0.0, 0.5],
                                         [0.5, 0.5, 0.0])

Where supercell_matrix defines the supercell expansion used to calculate the force constants using the finite displacements method, and primitive_matrix (optional) defines the primitive cell axis matrix. If primitive_matrix is not defined, identity matrix is used (primitive cell = unit cell).

2b. Use a python list containing LAMMPS input commands (one command per line)

list_of_commands = open('in.lammps').read().split('\n')
phlammps = Phonolammps(list_of_commands,
                       supercell_matrix=[[2, 0, 0],
                                         [0, 2, 0],
                                         [0, 0, 2]])
                       primitive_matrix=[[0.0, 0.5 ,0.5],
                                         [0.5, 0.0, 0.5],
                                         [0.5, 0.5, 0.0])

This second option may be handy if you want to generate/modify LAMMPS commands in a python script.

  1. Get the data needed for phonopy
unitcell = phlammps.get_unitcell()
force_constants = phlammps.get_force_constants()
supercell_matrix = phlammps.get_supercell_matrix()
  1. From this you have all the information you need for phonopy calculations
from phonopy import Phonopy
phonon = Phonopy(unitcell,
                 supercell_matrix)

phonon.set_force_constants(force_constants)
phonon.set_mesh([20, 20, 20])

phonon.set_total_DOS()
phonon.plot_total_DOS().show()

phonon.set_thermal_properties()
phonon.plot_thermal_properties().show()

4.2. Command line interface

To use phonoLAMMPS from command line interface you need a LAMMPS input file containing the definition of the unit cell and potentials.

example:

units           metal

boundary        p p p

box tilt large

atom_style      atomic

read_data       data.si

pair_style      tersoff
pair_coeff      * * SiCGe.tersoff  Si(C)

neighbor    0.3 bin

Notice that run command, as well as other MD related commands (thermostat, velocities, etc..) should not be included in the input file.

Note

In this example read data command is used to define the atoms coordinates in a different file (refer to LAMMPS manual for further information).

Phonolammps script uses argparse to provide a clean command line interface using flags All options available are displayed by using -h

$ phonolammps -h
usage: phonolammps [-h] [-o file] [--dim N N N] [-p] [-c file]
                   [-pa F F F F F F F F F] [-t F] [--optimize] [--force_tol F]
                   [--amplitude F] [--total_time F] [--relaxation_time F]
                   [--timestep F] [--logshow] [--no_symmetrize] [--use_NAC]
                   [--write_force_sets]
                   lammps_file

phonoLAMMPS options

positional arguments:
  lammps_file           lammps input file

optional arguments:
  -h, --help            show this help message and exit
  -o file               force constants output file [default: FORCE_CONSTANTS]
  --dim N N N           dimensions of the supercell
  -p                    plot phonon band structure
  -c file, --cell file  generates a POSCAR type file containing the unit cell
  -pa F F F F F F F F F, --primitive_axis F F F F F F F F F
                        primitive axis
  -t F                  temperature in K
  --optimize            optimize atoms position of unitcell
  --force_tol F         forces tolerance for optimization
  --amplitude F         displacement distance [default: 0.01 angstrom]
  --total_time F        total MD time in picoseconds [default: 20 ps]
  --relaxation_time F   MD relaxation time in picoseconds [default: 5 ps]
  --timestep F          MD time step in picoseconds [default: 0.001 ps]
  --logshow             show LAMMPS & dynaphopy log on screen
  --no_symmetrize       deactivate force constant symmetrization
  --use_NAC             include non analytical corrections (Requires BORN file
                        in work directory)
  --write_force_sets    write FORCE_SETS file
  --version             print version

A simple example for crystalline silicon using a 2x2x2 supercell would be

phonolammps in.lammps --dim 2 2 2 -pa 0.0 0.5 0.5 0.5 0.0 0.5 0.5 0.5 0.0 -c POSCAR_unitcell -p

where in.lammps is a LAMMPS input containing the unit cell, –dim defines the supercell, –pa are the primitive axis in matrix format written in one line (phonopy style), -c FILENAME (optional) requests to write the unitcell (the same written in LAMMPS input) in VASP format on the disk to be used in phonopy calculations, and -p requests to show a preview of the phonon band structure in a matplotlib plot. The output of this script is a file named FORCE_CONSTANTS that contains the interatomic 2nd order force constants in phonopy format.

To use the obtained FORCE_CONSTANTS in more advanced calculations, or to have more control on the displayed data it is recommended to use PHONOPY by reading the FORCE_CONSTANTS file using –readfc option

phonopy --dim="2 2 2" --pa="0.0 0.5 0.5 0.5 0.0 0.5 0.5 0.5 0.0" --readfc -c POSCAR_unitcell band.conf