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)

You can toggle the following to investigate different models:

[1]:
model_params = {"arch": "mace_mp", "model_path": "medium-0b3"}
# model_params = {"arch": "mace_mp", "model_path": "medium-mpa-0"}
# model_params = {"arch": "mace_mp", "model_path": "medium-omat-0"}
# model_params = {"arch": "chgnet"}
# model_params = {"arch": "sevennet"}

Set up environment (optional)

These steps are required for Google Colab, but may work on other systems too:

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

# ! pip uninstall torch torchaudio torchvision numpy -y
# ! uv pip install janus-core[all] data-tutorials torch==2.5.1 --system
# get_ipython().kernel.do_shutdown(restart=True)
[3]:
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:

[4]:
get_data(
    # url="https://raw.githubusercontent.com/stfc/janus-core/main/tests/data/",
    url="https://raw.githubusercontent.com/stfc/janus-tutorials/main/neb/data/",
    filename="LiFePO4_supercell.cif",
    folder="data",
)
try to download LiFePO4_supercell.cif from https://raw.githubusercontent.com/stfc/janus-tutorials/main/neb/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):

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

First, we will relax the supercell:

[6]:
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.
/home/runner/work/janus-core/janus-core/.venv/lib/python3.12/site-packages/mace/calculators/mace.py:139: 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 17:12:40     -762.842666        0.654854
LBFGS:    1 17:12:43     -763.038742        0.437508
LBFGS:    2 17:12:46     -763.090096        0.392867
LBFGS:    3 17:12:48     -763.119488        0.374220
LBFGS:    4 17:12:50     -763.172830        0.346620
LBFGS:    5 17:12:52     -763.213488        0.334927
LBFGS:    6 17:12:55     -763.256842        0.327520
LBFGS:    7 17:12:57     -763.297792        0.328964
LBFGS:    8 17:12:59     -763.342419        0.332009
LBFGS:    9 17:13:01     -763.385108        0.328513
LBFGS:   10 17:13:03     -763.427187        0.308144
LBFGS:   11 17:13:05     -763.474267        0.274918
LBFGS:   12 17:13:07     -763.528185        0.234548
LBFGS:   13 17:13:09     -763.579927        0.246328
LBFGS:   14 17:13:11     -763.627605        0.245408
LBFGS:   15 17:13:13     -763.684752        0.231582
LBFGS:   16 17:13:15     -763.766360        0.283626
LBFGS:   17 17:13:17     -763.855346        0.296599
LBFGS:   18 17:13:19     -763.931149        0.188692
LBFGS:   19 17:13:21     -763.965778        0.154808
LBFGS:   20 17:13:23     -763.996599        0.213159
LBFGS:   21 17:13:26     -764.039520        0.255107
LBFGS:   22 17:13:27     -764.108589        0.255570
LBFGS:   23 17:13:29     -764.175608        0.200210
LBFGS:   24 17:13:31     -764.208816        0.096659
[6]:

Next, we will create the start and end structures:

[7]:
# 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.

[8]:
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.
/home/runner/work/janus-core/janus-core/.venv/lib/python3.12/site-packages/mace/calculators/mace.py:139: 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)
[9]:
results = neb_b.run()
       Step     Time          Energy          fmax
LBFGS:    0 17:13:33     -758.975432        1.930990
LBFGS:    1 17:13:35     -759.140669        0.711957
LBFGS:    2 17:13:37     -759.209571        0.514600
LBFGS:    3 17:13:38     -759.298948        0.404890
LBFGS:    4 17:13:40     -759.316304        0.265548
LBFGS:    5 17:13:42     -759.336518        0.252446
LBFGS:    6 17:13:44     -759.351338        0.330410
LBFGS:    7 17:13:45     -759.368546        0.309163
LBFGS:    8 17:13:47     -759.378014        0.229741
LBFGS:    9 17:13:49     -759.386280        0.221552
LBFGS:   10 17:13:50     -759.395290        0.278156
LBFGS:   11 17:13:52     -759.404981        0.291370
LBFGS:   12 17:13:54     -759.411852        0.192426
LBFGS:   13 17:13:56     -759.415355        0.075502
       Step     Time          Energy          fmax
LBFGS:    0 17:13:57     -758.975436        1.930961
LBFGS:    1 17:13:59     -759.140671        0.711930
LBFGS:    2 17:14:01     -759.209572        0.514582
LBFGS:    3 17:14:03     -759.298949        0.404887
LBFGS:    4 17:14:04     -759.316305        0.265544
LBFGS:    5 17:14:06     -759.336519        0.252446
LBFGS:    6 17:14:08     -759.351338        0.330395
LBFGS:    7 17:14:10     -759.368546        0.309158
LBFGS:    8 17:14:11     -759.378015        0.229745
LBFGS:    9 17:14:13     -759.386281        0.221556
LBFGS:   10 17:14:15     -759.395291        0.278147
LBFGS:   11 17:14:16     -759.404982        0.291359
LBFGS:   12 17:14:18     -759.411852        0.192423
LBFGS:   13 17:14:20     -759.415355        0.075502
                   Step     Time         fmax
NEBOptimizer[ode]:    0 17:14:56       1.5090
NEBOptimizer[ode]:    1 17:15:07       1.1032
NEBOptimizer[ode]:    2 17:15:19       0.9983
NEBOptimizer[ode]:    3 17:15:31       0.8704
NEBOptimizer[ode]:    4 17:15:43       0.8013
NEBOptimizer[ode]:    5 17:15:55       0.7712
NEBOptimizer[ode]:    6 17:16:07       0.7441
NEBOptimizer[ode]:    7 17:16:19       0.6738
NEBOptimizer[ode]:    8 17:16:31       0.4346
NEBOptimizer[ode]:    9 17:16:55       0.4299
NEBOptimizer[ode]:   10 17:17:07       0.4226
NEBOptimizer[ode]:   11 17:17:19       0.3946
NEBOptimizer[ode]:   12 17:17:31       0.2843
NEBOptimizer[ode]:   13 17:17:55       0.2726
NEBOptimizer[ode]:   14 17:18:07       0.2652
NEBOptimizer[ode]:   15 17:18:19       0.2612
NEBOptimizer[ode]:   16 17:18:31       0.2552
NEBOptimizer[ode]:   17 17:18:42       0.2323
NEBOptimizer[ode]:   18 17:19:06       0.2219
NEBOptimizer[ode]:   19 17:19:18       0.2141
NEBOptimizer[ode]:   20 17:19:30       0.2087
NEBOptimizer[ode]:   21 17:19:42       0.2055
NEBOptimizer[ode]:   22 17:19:54       0.2005
NEBOptimizer[ode]:   23 17:20:06       0.1809
NEBOptimizer[ode]:   24 17:20:30       0.1729
NEBOptimizer[ode]:   25 17:20:42       0.1664
NEBOptimizer[ode]:   26 17:20:53       0.1624
NEBOptimizer[ode]:   27 17:21:05       0.1598
NEBOptimizer[ode]:   28 17:21:17       0.1559
NEBOptimizer[ode]:   29 17:21:29       0.1409
NEBOptimizer[ode]:   30 17:21:41       0.2612
NEBOptimizer[ode]:   31 17:22:05       0.0826

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

[10]:
print(results)
{'barrier': 0.28471736968049277, 'delta_E': -1.6499620869581122e-07, 'max_force': 0.09009922955812238}

We can also plot the band:

[11]:
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

[11]:
../../_images/tutorials_python_neb_24_1.png

Path c

We can calculate the barrier height along path c similarly

[12]:
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,
)
[13]:
results = neb_c.run()
       Step     Time          Energy          fmax
LBFGS:    0 17:22:09     -759.415355        0.075502
       Step     Time          Energy          fmax
LBFGS:    0 17:22:11     -758.975431        1.930990
LBFGS:    1 17:22:12     -759.140669        0.711958
LBFGS:    2 17:22:14     -759.209571        0.514601
LBFGS:    3 17:22:16     -759.298948        0.404890
LBFGS:    4 17:22:17     -759.316304        0.265548
LBFGS:    5 17:22:19     -759.336518        0.252446
LBFGS:    6 17:22:21     -759.351338        0.330410
LBFGS:    7 17:22:23     -759.368546        0.309163
LBFGS:    8 17:22:24     -759.378014        0.229741
LBFGS:    9 17:22:26     -759.386280        0.221553
LBFGS:   10 17:22:28     -759.395290        0.278156
LBFGS:   11 17:22:29     -759.404981        0.291370
LBFGS:   12 17:22:31     -759.411852        0.192426
LBFGS:   13 17:22:33     -759.415355        0.075502
                   Step     Time         fmax
NEBOptimizer[ode]:    0 17:23:08       2.3204
NEBOptimizer[ode]:    1 17:23:20       1.6499
NEBOptimizer[ode]:    2 17:23:32       0.9977
NEBOptimizer[ode]:    3 17:23:43       0.6537
NEBOptimizer[ode]:    4 17:24:07       0.4186
NEBOptimizer[ode]:    5 17:24:19       0.3969
NEBOptimizer[ode]:    6 17:24:31       0.3176
NEBOptimizer[ode]:    7 17:24:54       0.2998
NEBOptimizer[ode]:    8 17:25:18       0.2690
NEBOptimizer[ode]:    9 17:25:30       0.2608
NEBOptimizer[ode]:   10 17:25:41       0.2293
NEBOptimizer[ode]:   11 17:26:05       0.2186
NEBOptimizer[ode]:   12 17:26:17       0.2622
NEBOptimizer[ode]:   13 17:26:29       0.2009
NEBOptimizer[ode]:   14 17:26:40       0.1946
NEBOptimizer[ode]:   15 17:26:52       0.1883
NEBOptimizer[ode]:   16 17:27:04       0.1653
NEBOptimizer[ode]:   17 17:27:27       0.1600
NEBOptimizer[ode]:   18 17:27:39       0.1538
NEBOptimizer[ode]:   19 17:27:51       0.1475
NEBOptimizer[ode]:   20 17:28:03       0.1379
NEBOptimizer[ode]:   21 17:28:15       0.1346
NEBOptimizer[ode]:   22 17:28:26       0.1288
NEBOptimizer[ode]:   23 17:28:38       0.1263
NEBOptimizer[ode]:   24 17:28:50       0.1169
NEBOptimizer[ode]:   25 17:29:02       0.1806
NEBOptimizer[ode]:   26 17:29:26       0.0821
[14]:
print(results)
{'barrier': 1.7875754406807047, 'delta_E': 4.249613994034007e-09, 'max_force': 0.09762768246406749}
[15]:
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
[15]:
../../_images/tutorials_python_neb_30_1.png