VASP is a density-functional theory code using pseudopotentials or the projector-augmented wave method and a plane wave basis set. This interface makes it possible to use VASP as a calculator in ASE, and also to use ASE as a post-processor for an already performed VASP calculation.
You need to write a script called run_vasp.py containing something like this:
import os
exitcode = os.system('vasp')
The environment variable VASP_SCRIPT must point to that file.
A directory containing the pseudopotential directories potpaw (LDA XC) potpaw_GGA (PW91 XC) and potpaw_PBE (PBE XC) is also needed, and it is to be put in the environment variable VASP_PP_PATH.
Set both environment variables in your shell configuration file:
$ export VASP_SCRIPT=$HOME/vasp/run_vasp.py
$ export VASP_PP_PATH=$HOME/vasp/mypps
The default setting used by the VASP interface is
Below follows a list with a selection of parameters
keyword | type | default value | description |
---|---|---|---|
restart | bool | None | Restart old calculation or use ASE for post-processing |
xc | str | 'PW91' | XC-functional |
setups | str | None | Additional setup option |
kpts | seq | \(\Gamma\)-point | k-point sampling |
gamma | bool | None | \(\Gamma\)-point centered k-point sampling |
reciprocal | bool | None | Use reciprocal units if k-points are specified explicitly |
prec | str | Accuracy of calculation | |
encut | float | Kinetic energy cutoff | |
ediff | float | Convergence break condition for SC-loop. | |
nbands | int | Number of bands | |
algo | str | Electronic minimization algorithm | |
ismear | int | Type of smearing | |
sigma | float | Width of smearing | |
nelm | int | Maximum number of SC-iterations |
seq: A sequence of three int‘s.
For parameters in the list without default value given, VASP will set the default value. Most of the parameters used in the VASP INCAR file are allowed keywords. See the official VASP manual for more details.
Note
Parameters can be changed after the calculator has been constructed by using the set() method:
>>> calc.set(prec='Accurate', ediff=1E-5)
This would set the precision to Accurate and the break condition for the electronic SC-loop to 1E-5 eV.
If the atoms object has non-zero magnetic moments, a spin-polarized calculation will be performed by default.
Here follows an example how to calculate the total magnetic moment of a sodium chloride molecule.
from ase import Atoms, Atom
from ase.calculators.vasp import Vasp
a = [6.5, 6.5, 7.7]
d = 2.3608
NaCl = Atoms([Atom('Na', [0, 0, 0], magmom=1.928),
Atom('Cl', [0, 0, d], magmom=0.75)],
cell=a)
calc = Vasp(prec='Accurate',
xc='PBE',
lreal=False)
NaCl.set_calculator(calc)
print(NaCl.get_magnetic_moment())
In this example the initial magnetic moments are assigned to the atoms when defining the Atoms object. The calculator will detect that at least one of the atoms has a non-zero magnetic moment and a spin-polarized calculation will automatically be performed. The ASE generated INCAR file will look like:
INCAR created by Atomic Simulation Environment
PREC = Accurate
LREAL = .FALSE.
ISPIN = 2
MAGMOM = 1*1.9280 1*0.7500
Note
It is also possible to manually tell the calculator to perform a spin-polarized calculation:
>>> calc.set(ispin=2)
This can be useful for continuation jobs, where the initial magnetic moment is read from the WAVECAR file.
To continue an old calculation which has been performed without the interface use the restart parameter when constructing the calculator
>>> calc = Vasp(restart=True)
Then the calculator will read atomic positions from the CONTCAR file, physical quantities from the OUTCAR file, k-points from the KPOINTS file and parameters from the INCAR file.
Note
Only Monkhorst-Pack and Gamma-centered k-point sampling are supported for restart at the moment. Some INCAR parameters may not be implemented for restart yet. Please report any problems to the ASE mailing list.
The restart parameter can be used , as the name suggest to continue a job from where a previous calculation finished. Furthermore, it can be used to extract data from an already performed calculation. For example, to get the total potential energy of the sodium chloride molecule in the previous section, without performing any additional calculations, in the directory of the previous calculation do:
>>> calc = Vasp(restart=True)
>>> atoms = calc.get_atoms()
>>> atoms.get_potential_energy()
-4.7386889999999999