Python implementation (2024), (c) Swinburne & Marinica 2024
C++ implementation (2023) available here
If using PAFI, please cite (bibtex). Details in Swinburne & Marinica PRL 2018.
- PAFI uses
mpi4py
,numpy
,scipy
,pandas
and LAMMPS-Python - If you can run the following python code
from lammps import lammps from mpi4py import MPI comm = MPI.COMM_WORLD lmp = lammps(comm=comm) lmp.close()
- You can try installing PAFI with
pip
(in an enviroment...) and run tests:cd /path/to/this/repo python -m pip install -e . # local install python unittests.py -v
- If this fails, see full installation instructions
-
Provide some initial pathway using e.g. LAMMPS NEB.
-
New
equal
style gives optimal spacing for force integration- e.g.fix neb all neb 1.0 parallel equal
-
Modify one of
examples/configuration_files/*_REAL.xml
to load in your pathway and potential -
Run with
mpirun
:
mpirun -np ${NUM_PROCS} python simple_run.py
where simple_run.py
:
from mpi4py import MPI
from pafi import PAFIManager
manager = PAFIManager(MPI.COMM_WORLD,"/path/to/config.xml")
manager.run()
manager.close()
- Alternatively, specify all custom options within python:
from mpi4py import MPI
from pafi import PAFIManager, PAFIParser
parameters = PAFIParser()
parameters.set_pathway("systems/EAM-SIA-Fe/image_*.dat") # NEB pathway
parameters.set_potential("systems/EAM-SIA-Fe/Fe.eam.fs","eam/fs",["Fe"]) # Potential
# typical sampling values
parameters.axes["Temperature"] = [100.*i for i in range(7)] # temperature range
parameters.set("CoresPerWorker",1)
parameters.set("SampleSteps",2000)
parameters.set("ThermSteps",1000)
parameters.set("ThermWindow",500)
manager = PAFIManager(MPI.COMM_WORLD,parameters=parameters)
manager.run()
manager.close()
PAFI uses mpi4py
, numpy
, scipy
, pandas
and LAMMPS-Python with at least MANYBODY
and ML-SNAP
If you have cmake and mpi installed:
export PREFIX=${HOME}/.local # example
export PYTHON=`which python` # to ensure same distribution
export MPICC=`which mpicc` # for mpi4py, your C++ MPI compiler (e.g. mpicc / mpiicc for intel)
# extract typical install location PLEASE CHECK THIS ON YOUR MACHINE!
# (see below for why this hack can be useful)
PYTHON_VERSION=`python --version | cut -f2 -d" " | cut -f2 -d"."`
export INSTALL_LOCATION=${PREFIX}/lib/python3.${PYTHON_VERSION}/site-packages
# get sources
git clone https://github.com/lammps/lammps.git
git clone https://github.com/tomswinburne/pafi.git
# install python packages
${PYTHON} -m pip install mpi4py numpy pandas
# LAMMPS build
cd /path/to/lammps
mkdir build
cd build
cmake -C ../../pafi/lammps_options.cmake ../cmake
make -j
# LAMMPS python install:
# whilst official command is 'make install python', can have env clashes
# instead, we do it "by hand":
cd ../python # within LAMMPS repository
${PYTHON} -m pip install -U .
# manually provide binary for LAMMPS package
cp ../build/liblammps.so ${INSTALL_LOCATION}/lammps
# Install and test PAFI
cd /path/to/pafi
pip install -e .
python unittests.py
-
See the tutorial for information on the
pafi-path-test
routine -
In general, we want a reference pathway with dense discretisation where energy gradients are large
-
The current non-smoothed spline implementation can oscillate between very similar image configurations, as a result, there should be non-negligible displacement between images
-
If your path isn't loading, try setting
LogLammps=1
inconfig.xml
to check for bugs inlog.lammps
-
If
SampleSteps
is too large workers will make thermally activated "jumps" to nearby paths in the hyperplane. This will return a warning messageReference path too unstable for sampling.
and increase error. If this happens, decreaseSampleSteps
and increasenRepeats
-
When running on
NPROCS
cores, we requireNPROCS%CoresPerWorker==0
, so we have an integer number of workers -
The total number of force calls per worker is
nPlanes * (ThermSteps+SampleSteps) * nRepeats
, spatially parallelised by LAMMPS acrossCoresPerWorker
cores for each worker. -
Each PAFI worker runs at the same speed as LAMMPS. Increasing
CoresPerWorker
will typically decrease execution time but also reducenWorkers
and increase error, as we have less samples. -
If you are core-limited, the
nRepeats
option forces workers to perform multiple independent sampling runs on each plane. For example, with all other parameters fixed, running on 32 cores withnRepeats=3
is equivalent to running on 3*32=96 cores withnRepeats=1
, but the latter will finish in a third of the time.
For more details please see our paper, citation:
@article{PhysRevLett.120.135503,
title = {Unsupervised Calculation of Free Energy Barriers in Large Crystalline Systems},
author = {Swinburne, Thomas D. and Marinica, Mihai-Cosmin},
journal = {Phys. Rev. Lett.},
volume = {120},
issue = {13},
pages = {135503},
numpages = {6},
year = {2018},
month = {Mar},
publisher = {American Physical Society},
doi = {10.1103/PhysRevLett.120.135503},
url = {https://link.aps.org/doi/10.1103/PhysRevLett.120.135503}
}
- Restart files from pathway deviations
- Smoothed spline interpolation for more general reference pathways