Welcome to aimsgb’s documentation!
aimsgb is an efficient and open-source Python library for generating atomic coordinates in periodic grain boundary models. It is designed to construct various grain boundary structures from cubic and non-cubic initial configurations. A convenient command line tool has also been provided to enable easy and fast construction of tilt and twist boundaries by assigining the degree of fit (Σ), rotation axis, grain boundary plane and initial crystal structure.
We also provide a web-based GUI to access aimsgb framework: aimsgb.org
Getting started
Installation
Method 1: Use Pip
The easiest way to install aimsgb is to simply run a one-liner in pip:
pip install aimsgb
Method 2: Use Git to install
Clone the latest version from github:
git clone https://github.com/ksyang2013/aimsgb.git
Navigate to aimsgb folder:
cd aimsgb
Type in the root of the repo:
pip install .
or to install the package in development mode:
pip install -e .
How to cite aimsgb
If you use aimsgb in your research, please consider citing the following work:
Jianli Cheng, Jian Luo, Kesong Yang. Aimsgb: An algorithm and open-source python library to generate periodic grain boundary structures. Computational Materials Science, 2018, 155, 92-103. doi:10.1016/j.commatsci.2018.08.029
GUI example
BCC Fe
Let’s take Fe Σ5[001]/(120) tilt grain boundary as an example. The initial structure is a BCC lattice structure. Enter 001
for the Rotation axis,
5
for the Σ, and 1 2 0
for the GB plane, and upload an initial structure
.
Remember to add the space between the integers in 1 2 0
. The default output formate is poscar
for VASP calculations.

Then click Build, you will see the following page. You can either copy & past the output or directly Download the file.

FCC Ni
Let’s take Ni Σ7[111]/(1 -3 2) tilt grain boundary as an example. The initial structure is an FCC lattice structure. Enter 111
for the Rotation axis,
7
for the Σ, and 1 -3 0
for the GB plane, and upload an initial structure
.
Remember to add the space between the integers in 1 -3 0
. The default output formate is poscar
for VASP calculations.

Then click Build, you will see the following page. You can either copy & past the output or directly Download the file.

SrTiO3
Let’s take SrTiO3 Σ3[110]/(1-12) tilt grain boundary as an example. The initial structure is a cubic structure. Enter 110
for the Rotation axis,
3
for the Σ, and 1 -1 2
for the GB plane, and upload an initial structure
.
Remember to add the space between the integers in 1 -1 2
. The default output formate is poscar
for VASP calculations.

Then click Build, you will see the following page. You can either copy & past the output or directly Download the file.

In VASP calculations of interface models, we often want to fix some layers in the bulk. The web app provides a method to fix some layers in each grain.
For example, if we want to fix 3 center layers in each grain, we can enter 3
for the Fix center layers.

Then click Build, you will see the following page.

Please note that the default tolerance factor to determnine if sites are in the same layer is 0.25 Angstrom.
More to come.
Python examples
BCC Fe
Let’s take Fe Σ5[001]/(120) tilt grain boundary as an example. The initial structure is a BCC lattice structure:
from aimsgb import GrainBoundary, Grain
s_input = Grain.from_file("POSCAR_Fe")
gb = GrainBoundary([0, 0, 1], 5, [1, 2, 0], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

s_input
can also be created using from_mp_id
method by giving an mp_id
from Materials Project.
The mp_id
of a BCC Fe is mp-13:
s_input = Grain.from_mp_id("mp-13")
To construct a bigger supercell for each grain, say 2 unit cells:
gb = GrainBoundary([0, 0, 1], 5, [1, 2, 0], s_input, uc_a=2, uc_b=2)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

To construct a Fe Σ5[001]/(001) twist grain boundary:
gb = GrainBoundary([0, 0, 1], 5, [0, 0, 1], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

It’s possible to shift one of the grain along the grain boundary plane. For example, let’s say we want to have one of the grain shifted for 0.2 Å along a-axis:
gb = GrainBoundary([0, 0, 1], 5, [1, 2, 0], s_input)
gb.grain_a.translate_sites(range(len(gb.grain_a)), [0.2, 0, 0])
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

FCC Ni
Let’s take Ni Σ7[111]/(1 -3 2) tilt grain boundary as an example. The initial structure is an FCC lattice structure:
from aimsgb import GrainBoundary, Grain
s_input = Grain.from_file("POSCAR_Ni")
gb = GrainBoundary([1, 1, 1], 7, [1, -3, 2], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

s_input
can also be created using from_mp_id
method by giving an mp_id
from Materials Project.
The mp_id
of an FCC Ni is mp-23:
s_input = Grain.from_mp_id("mp-23")
MoS2
Let’s take MoS2 Σ5[001]/(130) tilt grain boundary as an example. The initial structure can be downloaded from this Materials Project link:
from aimsgb import GrainBoundary, Grain
s_input = Grain.from_file("POSCAR_MoS2")
gb = GrainBoundary([0, 0, 1], 5, [1, 3, 0], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

s_input
can also be created using from_mp_id
method by giving an mp_id
from Materials Project.
The mp_id
of MoS2 is mp-2815:
s_input = Grain.from_mp_id("mp-2815")
Notice that some atoms at the grain boundary sit too close to each other. Let’s add 1 Å distance between the grains:
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, gap=1, direction=gb.direction)

Graphene
Let’s take graphene Σ5[001]/(120) tilt grain boundary as an example. The initial structure can be found from Materials Project:
from aimsgb import GrainBoundary, Grain
s_input = Grain.from_file("POSCAR_C")
gb = GrainBoundary([0, 0, 1], 5, [1, 2, 0], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)


s_input
can also be created using from_mp_id
method by giving an mp_id
from Materials Project.
The mp_id
of C is mp-48:
s_input = Grain.from_mp_id("mp-48")
The structure above has two layers of graphenes. To remove the top layer along the c-axis:
structure.delete_bt_layer('t', direction=2)

Here, t
means the top layer. direction
is the direction to sort the sites by layers. 0: a, 1: b, 2: c
SrTiO3
Let’s take SrTiO3 Σ3[110]/(1-12) tilt grain boundary as an example. The initial structure is a cubic structure:
from aimsgb import GrainBoundary, Grain
s_input = Grain.from_file("POSCAR_SrTiO3")
gb = GrainBoundary([1, 1, 0], 3, [1, -1, 2], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

s_input
can also be created using from_mp_id
method by giving an mp_id
from Materials Project.
The mp_id
of a cubic SrTiO3 is mp-5229:
s_input = Grain.from_mp_id("mp-5229")
In VASP calculations of interface models, we often want to fix some layers in the bulk. aimsgb
provides a method to fix some layers in each grain.
First, let’s check the number of layers in each grain:
layers = gb.grain_a.sort_sites_in_layers(direction=gb.direction)
print(len(layers))
We can see that there are 12 layers in each grain. Please note that the default tolerance factor to determnine if sites are in the same layer is 0.25 Angstrom.
You can change the tolerance factor by giving a tol
argument:
layers = gb.grain_a.sort_sites_in_layers(direction=gb.direction, tol=1)
print(len(layers))
We can see that there are 6 layers in each grain.

Next, let’s fix layer 2 and 3 in each grain:
gb.grain_a.fix_sites_in_layers([2, 3], direction=gb.direction, tol=1)
gb.grain_b.fix_sites_in_layers([2, 3], direction=gb.direction, tol=1)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)
structure.to(filename="POSCAR")
We can see that Selective dynamics
has been added in the POSCAR
file.
LiNiO2
Let’s take LiNiO2 Σ13[001]/(230) tilt grain boundary as an example. The initial structure can be downloaded from this Materials Project link:
from aimsgb import GrainBoundary, Grain
s_input = Grain.from_file("POSCAR_LiNiO2")
gb = GrainBoundary([0, 0, 1], 13, [2, 3, 0], s_input)
structure = Grain.stack_grains(gb.grain_a, gb.grain_b, direction=gb.direction)

s_input
can also be created using from_mp_id
method by giving an mp_id
from Materials Project.
The mp_id
of LiNiO2 is mp-25411:
s_input = Grain.from_mp_id("mp-25411")
Find Grain Boundaries
In this example, I will show you how to print out available Σ, rotation angles, grain boundary (GB) planes, and coincident site lattice (CSL) matrices from a
given rotation axis. The purpose of this method is to provide you the grain boundary information. For example, if you are not sure what the GB plane is for a
GB you want to generate, you can find this information in GBInformation
, which is a dict
object. Let’s say you want to create a Σ3[110] GB, but you are not sure what the available
GB planes are. To find this information
from aimsgb import GBInformation
print(GBInformation([1, 1, 0], 30).__str__())
Grain boundary information for rotation axis: 110
Show the sigma values up to 30 (Note: * means twist GB, Theta is the rotation angle)
| Sigma | Theta | GB plane | CSL matrix |
|---------+---------+------------+--------------|
| 3 | 70.53 | (-1 1 1) | -1 1 1 |
| | | (1 -1 2) | 1 -1 1 |
| | | (1 1 0)* | 1 2 0 |
| 9 | 38.94 | (-1 1 -4) | -1 -2 1 |
| | | (-2 2 1) | 1 2 1 |
| | | (1 1 0)* | -4 1 0 |
| 11 | 50.48 | (-3 3 2) | -3 1 1 |
| | | (1 -1 3) | 3 -1 1 |
| | | (1 1 0)* | 2 3 0 |
| 17 | 86.63 | (3 -3 -4) | 3 -2 1 |
| | | (-2 2 -3) | -3 2 1 |
| | | (1 1 0)* | -4 -3 0 |
| 19 | 26.53 | (-3 3 1) | -3 1 1 |
| | | (1 -1 6) | 3 -1 1 |
| | | (1 1 0)* | 1 6 0 |
| 27 | 31.59 | (-5 5 2) | -5 1 1 |
| | | (1 -1 5) | 5 -1 1 |
| | | (1 1 0)* | 2 5 0 |
Here, I print out the GB information with Σ up to 30. As you can see, the available GB planes for a Σ3[110] GB are (-1 1 1)
, (1 -1 2)
, and (1 1 0)
.
(-1 1 1)
and (1 -1 2)
are for tilt GB and (1 1 0)
is for twist GB.
More to come.