Open In Colab

Nudged Elastic Band

In this tutorial, we will determine the activation energies of Li diffusion along the [010] and [001] directions (referred to as paths b and c respectively) in lithium iron phosphate (LiFePO_4), a cathode material for lithium ion batteries.

DFT references energies are:

  • Barrier heights:

    • path b = 0.27 eV

    • path c = 2.5 eV

(see table 1 in https://doi.org/10.1039/C5TA05062F)

Set up environment (optional)

These steps are required to run this tutorial with Google Colab. To do so, uncomment and run the cell below.

This will replace pre-installed versions of numpy and torch in Colab with versions that are known to be compatible with janus-core.

It may be possible to skip the steps that uninstall and reinstall torch, which will save a considerable amount of time.

These instructions but may work for other systems too, but it is typically preferable to prepare a virtual environment separately before running this notebook if possible.

[1]:
# import locale
# locale.getpreferredencoding = lambda: "UTF-8"

# ! pip uninstall numpy -y # Uninstall pre-installed numpy

# ! pip uninstall torch torchaudio torchvision transformers -y # Uninstall pre-installed torch
# ! uv pip install torch==2.5.1 # Install pinned version of torch

# ! uv pip install janus-core[mace,orb,chgnet,visualise] data-tutorials --system # Install janus-core with MACE, Orb, CHGNet, and WeasWidget, and data-tutorials

# get_ipython().kernel.do_shutdown(restart=True) # Restart kernel to update libraries. This may warn that your session has crashed.

To ensure you have the latest version of janus-core installed, compare the output of the following cell to the latest version available at https://pypi.org/project/janus-core/

[2]:
from janus_core import __version__

print(__version__)
0.8.6

Prepare data, modules, and model parameters

You can toggle the following to investigate different models:

[3]:
model_params = {"arch": "mace_mp", "model": "medium-0b3"}
# model_params = {"arch": "mace_mp", "model": "medium-mpa-0"}
# model_params = {"arch": "mace_mp", "model": "medium-omat-0"}
# model_params = {"arch": "chgnet"}
# model_params = {"arch": "orb"}
[4]:
from weas_widget import WeasWidget
from ase.io import read
from data_tutorials.data import get_data

from janus_core.calculations.geom_opt import GeomOpt
from janus_core.calculations.neb import NEB

Use data_tutorials to get the data required for this tutorial:

[5]:
get_data(
    url="https://raw.githubusercontent.com/stfc/janus-core/main/docs/source/tutorials/data/",
    filename="LiFePO4_supercell.cif",
    folder="data",
)
try to download LiFePO4_supercell.cif from https://raw.githubusercontent.com/stfc/janus-core/main/docs/source/tutorials/data/ and save it in data/LiFePO4_supercell.cif
saved in data/LiFePO4_supercell.cif

Preparing end structures

The initial structure can be downloaded from the Materials Project (mp-19017):

[6]:
LFPO = read("data/LiFePO4_supercell.cif")
v=WeasWidget()
v.from_ase(LFPO)
v.avr.model_style = 1
v.avr.show_hydrogen_bonds = True
v
[6]:

First, we will relax the supercell:

[7]:
GeomOpt(struct=LFPO, **model_params).run()

v1=WeasWidget()
v1.from_ase(LFPO)
v1.avr.model_style = 1
v1.avr.show_hydrogen_bonds = True
v1
/home/runner/work/janus-core/janus-core/.venv/lib/python3.12/site-packages/e3nn/o3/_wigner.py:10: UserWarning: Environment variable TORCH_FORCE_NO_WEIGHTS_ONLY_LOAD detected, since the`weights_only` argument was not explicitly passed to `torch.load`, forcing weights_only=False.
  _Jd, _W3j_flat, _W3j_indices = torch.load(os.path.join(os.path.dirname(__file__), 'constants.pt'))
cuequivariance or cuequivariance_torch is not available. Cuequivariance acceleration will be disabled.
Downloading MACE model from 'https://github.com/ACEsuit/mace-mp/releases/download/mace_mp_0b3/mace-mp-0b3-medium.model'
Cached MACE model to /home/runner/.cache/mace/macemp0b3mediummodel
Using Materials Project MACE for MACECalculator with /home/runner/.cache/mace/macemp0b3mediummodel
Using float64 for MACECalculator, which is slower but more accurate. Recommended for geometry optimization.
Using head default out of ['default']
/home/runner/work/janus-core/janus-core/.venv/lib/python3.12/site-packages/mace/calculators/mace.py:197: UserWarning: Environment variable TORCH_FORCE_NO_WEIGHTS_ONLY_LOAD detected, since the`weights_only` argument was not explicitly passed to `torch.load`, forcing weights_only=False.
  torch.load(f=model_path, map_location=device)
       Step     Time          Energy          fmax
LBFGS:    0 16:58:52     -762.842666        0.654854
LBFGS:    1 16:58:54     -763.038742        0.437508
LBFGS:    2 16:58:56     -763.090096        0.392867
LBFGS:    3 16:58:58     -763.119488        0.374220
LBFGS:    4 16:59:00     -763.172830        0.346620
LBFGS:    5 16:59:02     -763.213488        0.334927
LBFGS:    6 16:59:04     -763.256842        0.327520
LBFGS:    7 16:59:06     -763.297792        0.328964
LBFGS:    8 16:59:08     -763.342419        0.332009
LBFGS:    9 16:59:10     -763.385108        0.328513
LBFGS:   10 16:59:12     -763.427187        0.308144
LBFGS:   11 16:59:14     -763.474267        0.274918
LBFGS:   12 16:59:16     -763.528185        0.234548
LBFGS:   13 16:59:17     -763.579927        0.246328
LBFGS:   14 16:59:19     -763.627605        0.245408
LBFGS:   15 16:59:21     -763.684752        0.231582
LBFGS:   16 16:59:23     -763.766360        0.283626
LBFGS:   17 16:59:25     -763.855346        0.296599
LBFGS:   18 16:59:27     -763.931149        0.188692
LBFGS:   19 16:59:29     -763.965778        0.154808
LBFGS:   20 16:59:30     -763.996599        0.213159
LBFGS:   21 16:59:32     -764.039520        0.255107
LBFGS:   22 16:59:34     -764.108589        0.255570
LBFGS:   23 16:59:36     -764.175608        0.200210
LBFGS:   24 16:59:38     -764.208816        0.096659
[7]:

Next, we will create the start and end structures:

[8]:
# NEB path along b and c directions have the same starting image.
# For start bc remove site 5
LFPO_start_bc = LFPO.copy()
del LFPO_start_bc[5]

# For end b remove site 11
LFPO_end_b = LFPO.copy()
del LFPO_end_b[11]

# For end c remove site 4
LFPO_end_c = LFPO.copy()
del LFPO_end_c[4]

Path b

We can now calculate the barrier height along path b.

This also includes running geometry optimization on the end points of this path.

[9]:
n_images = 7
interpolator="pymatgen" # ASE interpolation performs poorly in this case

neb_b = NEB(
    init_struct=LFPO_start_bc,
    final_struct=LFPO_end_b,
    n_images=n_images,
    interpolator=interpolator,
    minimize=True,
    fmax=0.1,
    **model_params,
)
Using Materials Project MACE for MACECalculator with /home/runner/.cache/mace/macemp0b3mediummodel
Using float64 for MACECalculator, which is slower but more accurate. Recommended for geometry optimization.
Using head default out of ['default']
/home/runner/work/janus-core/janus-core/.venv/lib/python3.12/site-packages/mace/calculators/mace.py:197: UserWarning: Environment variable TORCH_FORCE_NO_WEIGHTS_ONLY_LOAD detected, since the`weights_only` argument was not explicitly passed to `torch.load`, forcing weights_only=False.
  torch.load(f=model_path, map_location=device)
[10]:
results = neb_b.run()
       Step     Time          Energy          fmax
LBFGS:    0 16:59:40     -758.975432        1.930990
LBFGS:    1 16:59:42     -759.140669        0.711957
LBFGS:    2 16:59:44     -759.209571        0.514600
LBFGS:    3 16:59:46     -759.298948        0.404890
LBFGS:    4 16:59:47     -759.316304        0.265548
LBFGS:    5 16:59:49     -759.336518        0.252446
LBFGS:    6 16:59:51     -759.351338        0.330410
LBFGS:    7 16:59:53     -759.368546        0.309163
LBFGS:    8 16:59:55     -759.378014        0.229741
LBFGS:    9 16:59:56     -759.386280        0.221552
LBFGS:   10 16:59:58     -759.395290        0.278156
LBFGS:   11 17:00:00     -759.404981        0.291370
LBFGS:   12 17:00:02     -759.411852        0.192426
LBFGS:   13 17:00:03     -759.415355        0.075502
       Step     Time          Energy          fmax
LBFGS:    0 17:00:05     -758.975436        1.930961
LBFGS:    1 17:00:07     -759.140671        0.711930
LBFGS:    2 17:00:09     -759.209572        0.514582
LBFGS:    3 17:00:10     -759.298949        0.404887
LBFGS:    4 17:00:12     -759.316305        0.265544
LBFGS:    5 17:00:14     -759.336519        0.252446
LBFGS:    6 17:00:16     -759.351338        0.330395
LBFGS:    7 17:00:17     -759.368546        0.309158
LBFGS:    8 17:00:19     -759.378015        0.229745
LBFGS:    9 17:00:21     -759.386281        0.221556
LBFGS:   10 17:00:23     -759.395291        0.278147
LBFGS:   11 17:00:25     -759.404982        0.291359
LBFGS:   12 17:00:26     -759.411852        0.192423
LBFGS:   13 17:00:28     -759.415355        0.075502
                   Step     Time         fmax
NEBOptimizer[ode]:    0 17:01:06       1.5090
NEBOptimizer[ode]:    1 17:01:18       1.1032
NEBOptimizer[ode]:    2 17:01:30       0.9983
NEBOptimizer[ode]:    3 17:01:42       0.8704
NEBOptimizer[ode]:    4 17:01:54       0.8013
NEBOptimizer[ode]:    5 17:02:06       0.7712
NEBOptimizer[ode]:    6 17:02:18       0.7441
NEBOptimizer[ode]:    7 17:02:30       0.6738
NEBOptimizer[ode]:    8 17:02:43       0.4346
NEBOptimizer[ode]:    9 17:03:07       0.4299
NEBOptimizer[ode]:   10 17:03:20       0.4226
NEBOptimizer[ode]:   11 17:03:32       0.3946
NEBOptimizer[ode]:   12 17:03:45       0.2843
NEBOptimizer[ode]:   13 17:04:10       0.2726
NEBOptimizer[ode]:   14 17:04:22       0.2652
NEBOptimizer[ode]:   15 17:04:35       0.2612
NEBOptimizer[ode]:   16 17:04:47       0.2552
NEBOptimizer[ode]:   17 17:05:00       0.2323
NEBOptimizer[ode]:   18 17:05:24       0.2219
NEBOptimizer[ode]:   19 17:05:37       0.2141
NEBOptimizer[ode]:   20 17:05:50       0.2087
NEBOptimizer[ode]:   21 17:06:02       0.2055
NEBOptimizer[ode]:   22 17:06:14       0.2005
NEBOptimizer[ode]:   23 17:06:26       0.1809
NEBOptimizer[ode]:   24 17:06:51       0.1729
NEBOptimizer[ode]:   25 17:07:03       0.1664
NEBOptimizer[ode]:   26 17:07:15       0.1624
NEBOptimizer[ode]:   27 17:07:27       0.1598
NEBOptimizer[ode]:   28 17:07:40       0.1559
NEBOptimizer[ode]:   29 17:07:52       0.1409
NEBOptimizer[ode]:   30 17:08:05       0.2612
NEBOptimizer[ode]:   31 17:08:29       0.0826

The results include the barrier (without any interpolation between highest images) and maximum force at the point in the simulation

[11]:
print(results)
{'barrier': np.float64(0.2847173696816294), 'delta_E': np.float64(-1.6499575394846033e-07), 'max_force': np.float64(0.09009922957482004)}

We can also plot the band:

[12]:
fig = neb_b.nebtools.plot_band()
v1=WeasWidget()
v1.from_ase(neb_b.nebtools.images)
v1.avr.model_style = 1
v1.avr.show_hydrogen_bonds = True
v1

[12]:
../../_images/tutorials_python_neb_27_1.png

Path c

We can calculate the barrier height along path c similarly

[13]:
n_images = 7
interpolator="pymatgen"

neb_c = NEB(
    init_struct=LFPO_start_bc,
    final_struct=LFPO_end_c,
    n_images=n_images,
    interpolator=interpolator,
    minimize=True,
    fmax=0.1,
    **model_params,
)
Using Materials Project MACE for MACECalculator with /home/runner/.cache/mace/macemp0b3mediummodel
Using float64 for MACECalculator, which is slower but more accurate. Recommended for geometry optimization.
Using head default out of ['default']
/home/runner/work/janus-core/janus-core/.venv/lib/python3.12/site-packages/mace/calculators/mace.py:197: UserWarning: Environment variable TORCH_FORCE_NO_WEIGHTS_ONLY_LOAD detected, since the`weights_only` argument was not explicitly passed to `torch.load`, forcing weights_only=False.
  torch.load(f=model_path, map_location=device)
[14]:
results = neb_c.run()
       Step     Time          Energy          fmax
LBFGS:    0 17:08:35     -759.415355        0.075502
       Step     Time          Energy          fmax
LBFGS:    0 17:08:37     -758.975431        1.930990
LBFGS:    1 17:08:38     -759.140669        0.711958
LBFGS:    2 17:08:40     -759.209571        0.514601
LBFGS:    3 17:08:42     -759.298948        0.404890
LBFGS:    4 17:08:44     -759.316304        0.265548
LBFGS:    5 17:08:46     -759.336518        0.252446
LBFGS:    6 17:08:47     -759.351338        0.330410
LBFGS:    7 17:08:49     -759.368546        0.309163
LBFGS:    8 17:08:51     -759.378014        0.229741
LBFGS:    9 17:08:53     -759.386280        0.221553
LBFGS:   10 17:08:55     -759.395290        0.278156
LBFGS:   11 17:08:57     -759.404981        0.291370
LBFGS:   12 17:08:58     -759.411852        0.192426
LBFGS:   13 17:09:00     -759.415355        0.075502
                   Step     Time         fmax
NEBOptimizer[ode]:    0 17:09:37       2.3204
NEBOptimizer[ode]:    1 17:09:49       1.6499
NEBOptimizer[ode]:    2 17:10:01       0.9977
NEBOptimizer[ode]:    3 17:10:14       0.6537
NEBOptimizer[ode]:    4 17:10:38       0.4186
NEBOptimizer[ode]:    5 17:10:50       0.3969
NEBOptimizer[ode]:    6 17:11:02       0.3176
NEBOptimizer[ode]:    7 17:11:26       0.2998
NEBOptimizer[ode]:    8 17:11:51       0.2690
NEBOptimizer[ode]:    9 17:12:04       0.2608
NEBOptimizer[ode]:   10 17:12:16       0.2293
NEBOptimizer[ode]:   11 17:12:40       0.2186
NEBOptimizer[ode]:   12 17:12:52       0.2622
NEBOptimizer[ode]:   13 17:13:04       0.2009
NEBOptimizer[ode]:   14 17:13:17       0.1946
NEBOptimizer[ode]:   15 17:13:29       0.1883
NEBOptimizer[ode]:   16 17:13:41       0.1653
NEBOptimizer[ode]:   17 17:14:06       0.1600
NEBOptimizer[ode]:   18 17:14:19       0.1538
NEBOptimizer[ode]:   19 17:14:31       0.1475
NEBOptimizer[ode]:   20 17:14:44       0.1379
NEBOptimizer[ode]:   21 17:14:56       0.1346
NEBOptimizer[ode]:   22 17:15:09       0.1288
NEBOptimizer[ode]:   23 17:15:21       0.1263
NEBOptimizer[ode]:   24 17:15:33       0.1169
NEBOptimizer[ode]:   25 17:15:45       0.1806
NEBOptimizer[ode]:   26 17:16:09       0.0821
[15]:
print(results)
{'barrier': np.float64(1.787575440772347), 'delta_E': np.float64(4.249386620358564e-09), 'max_force': np.float64(0.09762768267248513)}
[16]:
fig = neb_c.nebtools.plot_band()
v2=WeasWidget()
v2.from_ase(neb_c.nebtools.images)
v2.avr.model_style = 1
v2.avr.show_hydrogen_bonds = True
v2
[16]:
../../_images/tutorials_python_neb_33_1.png