AbinitInput object#

The creation of the Abinit input file is one of the most repetitive and error-prone operations we have to perform before running our calculations. To facilitate the creation of the input files, AbiPy provides the AbinitInput object, a dict-like object storing the Abinit variables and providing methods to automate the specification of multiple parameters.

This notebook discusses how to create an AbinitInput and how to define the parameters of the calculation. In the last part, we introduce the MultiDataset object that is mainly designed for the generation of multiple inputs sharing the same structure and the same list of pseudopotentials. In another notebook, we briefly discuss how to use factory functions to generate automatically input objects for typical calculations.

See e.g the abipy.abio.inputs.AbinitInput

Creating an AbinitInput object#

import os
import warnings
warnings.filterwarnings("ignore") # to get rid of deprecation warnings

import abipy.data as abidata
import abipy.abilab as abilab
abilab.enable_notebook() # This line tells AbiPy we are running inside a notebook
from abipy.abilab import AbinitInput

# This line configures matplotlib to show figures embedded in the notebook.
# Replace `inline` with `notebook` in classic notebook
%matplotlib inline

# Option available in jupyterlab. See https://github.com/matplotlib/jupyter-matplotlib
#%matplotlib widget

To create an Abinit input, we must specify the paths of the pseudopotential files. In this case, we have a pseudo named 14si.pspnc located in the abidata.pseudo_dir directory.

inp = AbinitInput(structure=abidata.cif_file("si.cif"),
                  pseudos="14si.pspnc", pseudo_dir=abidata.pseudo_dir)

print(inp) returns a string with our input. In this case, the input is almost empty since only the structure and the pseudos have been specified.

print(inp)
##############################################
####                SECTION: files
##############################################
 indata_prefix "indata/in"
 tmpdata_prefix "tmpdata/tmp"
 outdata_prefix "outdata/out"
 pseudos "/usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/abipy/data/pseudos/14si.pspnc"
##############################################
####                  STRUCTURE
##############################################
 natom 2
 ntypat 1
 typat 1 1
 znucl 14
 xred
    0.0000000000    0.0000000000    0.0000000000
    0.2500000000    0.2500000000    0.2500000000
 acell    1.0    1.0    1.0
 rprim
    6.3285005244    0.0000000000    3.6537614813
    2.1095001748    5.9665675141    3.6537614813
    0.0000000000    0.0000000000    7.3075229627


#<JSON>
#{
#    "pseudos": [
#        {
#            "@module": "pymatgen.io.abinit.pseudos",
#            "@class": "NcAbinitPseudo",
#            "basename": "14si.pspnc",
#            "type": "NcAbinitPseudo",
#            "symbol": "Si",
#            "Z": 14,
#            "Z_val": 4.0,
#            "l_max": 2,
#            "md5": "3916b143991b1cfa1542b130be320e5e",
#            "filepath": "/usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/abipy/data/pseudos/14si.pspnc"
#        }
#    ]
#}
#</JSON>

Inside the jupyter notebook, it is possible to visualize the input in HTML including the links to the official ABINIT documentation:

inp
##############################################
#### SECTION: files
##############################################
indata_prefix indata/in
tmpdata_prefix tmpdata/tmp
outdata_prefix outdata/out
pseudos /usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/abipy/data/pseudos/14si.pspnc
##############################################
#### STRUCTURE
##############################################
natom 2
ntypat 1
typat 1 1
znucl 14
xred
0.0000000000 0.0000000000 0.0000000000
0.2500000000 0.2500000000 0.2500000000
acell 1.0 1.0 1.0
rprim
6.3285005244 0.0000000000 3.6537614813
2.1095001748 5.9665675141 3.6537614813
0.0000000000 0.0000000000 7.3075229627

The input has a structure:

print(inp.structure)
Full Formula (Si2)
Reduced Formula: Si
abc   :   3.866975   3.866975   3.866975
angles:  60.000000  60.000000  60.000000
pbc   :       True       True       True
Sites (2)
  #  SP       a     b     c
---  ----  ----  ----  ----
  0  Si    0     0     0
  1  Si    0.25  0.25  0.25

and a list of Pseudo objects:

for pseudo in inp.pseudos:
    print(pseudo)
<NcAbinitPseudo: 14si.pspnc>
  summary: Troullier-Martins psp for element  Si        Thu Oct 27 17:31:21 EDT 1994
  number of valence electrons: 4.0
  maximum angular momentum: d
  angular momentum for local part: d
  XC correlation: LDA_XC_TETER93
  supports spin-orbit: False
  radius for non-linear core correction: 1.80626423934776
  hint for low accuracy: ecut: 0.0, pawecutdg: 0.0
  hint for normal accuracy: ecut: 0.0, pawecutdg: 0.0
  hint for high accuracy: ecut: 0.0, pawecutdg: 0.0

that have been constructed by parsing the pseudopotential files passed to AbinitInput.

Use set_vars to set the value of several variables with a single call:

inp.set_vars(ecut=8, paral_kgb=0)
{'ecut': 8, 'paral_kgb': 0}

AbinitInput is a dict-like object, hence one can test for the presence of a variable in the input:

"ecut" in inp
True

To list all the variables that have been defined, use:

list(inp.keys())
['ecut', 'paral_kgb']

To access the value of a particular variable use the syntax:

inp["ecut"]
8

To iterate over keywords and values:

for varname, varvalue in inp.items():
    print(varname, "-->", varvalue)
ecut --> 8
paral_kgb --> 0

Use lists, tuples or numpy arrays when Abinit expects arrays

inp.set_vars(kptopt=1,
             ngkpt=[2, 2, 2],
             nshiftk=2,
             shiftk=[0.0, 0.0, 0.0, 0.5, 0.5, 0.5]  # 2 shifts in one list
            )

# It is possible to use strings but use them only for special cases such as:
inp["istwfk"] = "*1"
inp
##############################################
#### SECTION: basic
##############################################
ecut 8
kptopt 1
ngkpt 2 2 2
nshiftk 2
shiftk
0.0 0.0 0.0
0.5 0.5 0.5
##############################################
#### SECTION: dev
##############################################
istwfk *1
##############################################
#### SECTION: files
##############################################
indata_prefix indata/in
tmpdata_prefix tmpdata/tmp
outdata_prefix outdata/out
pseudos /usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/abipy/data/pseudos/14si.pspnc
##############################################
#### SECTION: paral
##############################################
paral_kgb 0
##############################################
#### STRUCTURE
##############################################
natom 2
ntypat 1
typat 1 1
znucl 14
xred
0.0000000000 0.0000000000 0.0000000000
0.2500000000 0.2500000000 0.2500000000
acell 1.0 1.0 1.0
rprim
6.3285005244 0.0000000000 3.6537614813
2.1095001748 5.9665675141 3.6537614813
0.0000000000 0.0000000000 7.3075229627

If you mistype the name of the variable, AbinitInput raises an Exception:

try:
    inp.set_vars(perl=0)
except Exception as exc:
    print(exc)
Cannot find variable `perl` in the internal database. If you believe this is not a typo, use:

    input.set_spell_check(False)

to disable spell checking. Perhaps the internal database is not in synch
with the Abinit version you are using. Please contact the AbiPy developers.

Warning

The AbinitInput is a mutable object so changing it will aftect all the references to the object. See this page for further info.

a = {"foo": "bar"}
b = a
c = a.copy()
a["hello"] = "world"
print("a dict:", a)
print("b dict:", b)
print("c dict:", c)
a dict: {'foo': 'bar', 'hello': 'world'}
b dict: {'foo': 'bar', 'hello': 'world'}
c dict: {'foo': 'bar'}

The set_structure method sets the value of the ABINIT variables:

It is always a good idea to set the structure immediately after the creation of an AbinitInput because several methods use this information to facilitate the specification of other variables. For instance, the set_kpath method uses the structure to generate the high-symmetry \(k\)-path for band structure calculations.

Warning

typat must be consistent with the list of pseudopotentials passed to AbinitInput

Creating a structure from Abinit variables#

It is possible to create a structure in different ways.

The most explicit (and verbose) consists in passing a dictionary with ABINIT variables provided one uses python lists (or lists or lists) when ABINIT expects a 1D (or a multidimensional array):

si_struct = dict(
    ntypat=1,
    natom=2,
    typat=[1, 1],
    znucl=14,
    acell=3*[10.217],
    rprim=[[0.0,  0.5,  0.5],
           [0.5,  0.0,  0.5],
           [0.5,  0.5,  0.0]],
    xred=[[0.0 , 0.0 , 0.0],
          [0.25, 0.25, 0.25]]
)
print(si_struct)
{'ntypat': 1, 'natom': 2, 'typat': [1, 1], 'znucl': 14, 'acell': [10.217, 10.217, 10.217], 'rprim': [[0.0, 0.5, 0.5], [0.5, 0.0, 0.5], [0.5, 0.5, 0.0]], 'xred': [[0.0, 0.0, 0.0], [0.25, 0.25, 0.25]]}

If you already have a string with the Abinit variable, you can use the from_abistring class method:

lif_struct = abilab.Structure.from_abistring("""
acell      7.7030079150    7.7030079150    7.7030079150 Angstrom
rprim      0.0000000000    0.5000000000    0.5000000000
           0.5000000000    0.0000000000    0.5000000000
           0.5000000000    0.5000000000    0.0000000000
natom      2
ntypat     2
typat      1 2
znucl      3 9
xred       0.0000000000    0.0000000000    0.0000000000
           0.5000000000    0.5000000000    0.5000000000
""")

print(lif_struct)
Full Formula (Li1 F1)
Reduced Formula: LiF
abc   :   5.446849   5.446849   5.446849
angles:  60.000000  60.000000  60.000000
pbc   :       True       True       True
Sites (2)
  #  SP      a    b    c
---  ----  ---  ---  ---
  0  Li    0    0    0
  1  F     0.5  0.5  0.5

This approach requires less input yet we still need to specify ntypat, znucland {[typat}}. Fortunately, from_abistring supports another Abinit-specific format in which the fractional coordinates and the element symbol are specified via the xred_symbols variable. In this case ntypat, znucl and typat do not need to be specified as they are automatically computed from xred_symbols:

mgb2_struct = abilab.Structure.from_abistring("""

# MgB2 lattice structure.
natom   3
acell   2*3.086  3.523 Angstrom
rprim   0.866025403784439  0.5  0.0
       -0.866025403784439  0.5  0.0
        0.0                0.0  1.0

# Atomic positions
xred_symbols
 0.0  0.0  0.0 Mg
 1/3  2/3  0.5 B
 2/3  1/3  0.5 B
""")

print(mgb2_struct)
Full Formula (Mg1 B2)
Reduced Formula: MgB2
abc   :   3.086000   3.086000   3.523000
angles:  90.000000  90.000000 120.000000
pbc   :       True       True       True
Sites (3)
  #  SP           a         b    c
---  ----  --------  --------  ---
  0  Mg    0         0         0
  1  B     0.333333  0.666667  0.5
  2  B     0.666667  0.333333  0.5

Structure from file#

From a CIF file:

inp.set_structure(abidata.cif_file("si.cif"))
Structure Summary
Lattice
    abc : 3.86697462 3.86697462 3.86697462
 angles : 59.99999999999999 59.99999999999999 59.99999999999999
 volume : 40.88829179346891
      A : 3.3488982567096763 0.0 1.9334873100000005
      B : 1.1162994189032256 3.1573715557642927 1.9334873100000005
      C : 0.0 0.0 3.86697462
    pbc : True True True
PeriodicSite: Si1 (Si) (0.0, 0.0, 0.0) [0.0, 0.0, 0.0]
PeriodicSite: Si2 (Si) (1.116, 0.7893, 1.933) [0.25, 0.25, 0.25]

From a Netcdf file produced by ABINIT:

inp.set_structure(abidata.ref_file("si_scf_GSR.nc"))
Structure Summary
Lattice
    abc : 3.866974636909136 3.866974636913999 3.866974636891966
 angles : 60.000000000188486 60.00000000014689 60.00000000040088
 volume : 40.88829232994348
      A : 3.3488982713583737 0.0 1.933487318445983
      B : 1.1162994237684851 3.157371569587409 1.933487318445983
      C : 0.0 0.0 3.866974636891966
    pbc : True True True
PeriodicSite: Si (0.0, 0.0, 0.0) [0.0, 0.0, 0.0]
PeriodicSite: Si (1.116, 0.7893, 1.933) [0.25, 0.25, 0.25]

Supported formats include:

  • CIF

  • POSCAR/CONTCAR

  • CHGCAR

  • LOCPOT

  • vasprun.xml

  • CSSR

  • ABINIT netcdf files

  • pymatgen’s JSON serialized structures

From the Materials Project database:#

# https://www.materialsproject.org/materials/mp-149/
inp.set_structure(abilab.Structure.from_mpid("mp-149"))
---------------------------------------------------------------------------
ValueError                                Traceback (most recent call last)
Cell In[20], line 2
      1 # https://www.materialsproject.org/materials/mp-149/
----> 2 inp.set_structure(abilab.Structure.from_mpid("mp-149"))

File /usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/abipy/core/structure.py:387, in Structure.from_mpid(cls, material_id)
    385 # Get pytmatgen structure and convert it to abipy structure
    386 from abipy.core import restapi
--> 387 with restapi.get_mprester() as rest:
    388     new = rest.get_structure_by_material_id(material_id)
    389     return cls.as_structure(new)

File /usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/abipy/core/restapi.py:33, in get_mprester()
     29 def get_mprester():
     30     """
     31     Args:
     32     """
---> 33     rester = MPRester()
     34     #print(f"{type(rester)=}")
     35     return rester

File /usr/share/miniconda/envs/abipy/lib/python3.12/site-packages/pymatgen/ext/matproj.py:414, in MPRester.__new__(cls, *args, **kwargs)
    411     kwargs["api_key"] = api_key
    413 if not api_key:
--> 414     raise ValueError("Please supply an API key. See https://materialsproject.org/api for details.")
    416 if len(api_key) != 32:
    417     from pymatgen.ext.matproj_legacy import _MPResterLegacy

ValueError: Please supply an API key. See https://materialsproject.org/api for details.

Remember to set the PMG_MAPI_KEY in ~/.pmgrc.yaml as described here

Note that you can avoid the call to set_structure if the structure argument is passed to AbiniInput:

AbinitInput(structure=abidata.cif_file("si.cif"), pseudos=abidata.pseudos("14si.pspnc"))

Brillouin zone sampling#

There are two different types of sampling of the BZ: homogeneous and high-symmetry k-path. The later is mainly used for band structure calculations and requires the specification of:

  • kptopt

  • kptbounds

  • ndivsm

whereas the homogeneous sampling is needed for all the calculations in which we have to compute integrals in the Brillouin zone e.g. total energy calculations, DOS, etc. The \(k\)-mesh is usually specified via:

Explicit \(k\)-mesh#

inp = AbinitInput(structure=abidata.cif_file("si.cif"), pseudos=abidata.pseudos("14si.pspnc"))

# Set ngkpt, shiftk explicitly
inp.set_kmesh(ngkpt=(1, 2, 3), shiftk=[0.0, 0.0, 0.0, 0.5, 0.5, 0.5])

Automatic \(k\)-mesh#

# Define a homogeneous k-mesh.
# nksmall is the number of divisions to be used to sample the smallest lattice vector,
# shiftk is automatically selected from an internal database.

inp.set_autokmesh(nksmall=4)

High-symmetry \(k\)-path#

# Generate a high-symmetry k-path (taken from an internal database)
# Ten points are used to sample the smallest segment,
# the other segments are sampled so that proportions are preserved.
# A warning is issued by pymatgen about the structure not being standard.
# Be aware that this might possibly affect the automatic labelling of the boundary k-points on the k-path.
# So, check carefully the k-point labels on the figures that are produced in such case.

inp.set_kpath(ndivsm=10)

Utilities#

Once the structure has been defined, one can compute the number of valence electrons with:

print("The number of valence electrons is: ", inp.num_valence_electrons)

If we need to change a particular (scalar) variable to generate inputs for convergence studies:

# When using a non-integer step, such as 0.1, the results will often not
# be consistent.  It is better to use ``linspace`` for these cases.
# See also numpy.arange and numpy.linspace

ecut_inps = inp.arange("ecut", start=2, stop=5, step=2)

print([i["ecut"] for i in ecut_inps])
tsmear_inps = inp.linspace("tsmear", start=0.001, stop=0.003, num=3)
print([i["tsmear"] for i in tsmear_inps])

Invoking Abinit with AbinitInput#

Once you have an AbinitInput, you can call Abinit to get useful information or simply to validate the input file before running the calculation. All the method that invoke Abinit starts with the abi prefix followed by a verb e.g. abiget or abivalidate.

inp = AbinitInput(structure=abidata.cif_file("si.cif"), pseudos=abidata.pseudos("14si.pspnc"))

inp.set_vars(ecut=-2)
inp.set_autokmesh(nksmall=4)

v = inp.abivalidate()
if v.retcode != 0:
    # If there is a mistake in the input, one can acces the log file of the run with the log_file object
    print("".join(v.log_file.readlines()[-10:]))

Let’s fix the problem with the negative ecut and rerun abivalidate!

inp["ecut"] = 2
inp["toldfe"] = 1e-10

v = inp.abivalidate()

if v.retcode == 0:
    print("All ok")
else:
    print(v)

At this point, we have a valid input file and we can get the k-points in the irreducible zone with:

ibz = inp.abiget_ibz()
print("number of k-points:", len(ibz.points))
print("k-points:", ibz.points)
print("weights:", ibz.weights)
print("weights are normalized to:", ibz.weights.sum())

We can also call the Abinit spacegroup finder with:

abistruct = inp.abiget_spacegroup()
print("spacegroup found by Abinit:", abistruct.abi_spacegroup)

To get the list of possible parallel configurations for this input up to 5 max_ncpus

inp["paral_kgb"] = 1
pconfs = inp.abiget_autoparal_pconfs(max_ncpus=5)
print("best efficiency:\n", pconfs.sort_by_efficiency()[0])
print("best speedup:\n", pconfs.sort_by_speedup()[0])

To get the list of irreducible phonon perturbations at Gamma (Abinit notation)

inp.abiget_irred_phperts(qpt=(0, 0, 0))

Multiple datasets#

Multiple datasets are handy when you have to generate several input files sharing several common variables e.g. the crystalline structure, the value of ecut etc… In this case, one can use the MultiDataset object that is essentially a list of AbinitInput objects. Note however that Abipy workflows do not support input files with more than one dataset.

# A MultiDataset object with two datasets (a.k.a. AbinitInput)
multi = abilab.MultiDataset(structure=abidata.cif_file("si.cif"),
                            pseudos="14si.pspnc", pseudo_dir=abidata.pseudo_dir, ndtset=2)

# A MultiDataset is essentially a list of AbinitInput objects
# with handy methods to perform global modifications.
# i.e. changes that will affect all the inputs in the MultiDataset
# For example:
multi.set_vars(ecut=4)

# is equivalent to
#
#   for inp in multi: inp.set_vars(ecut=4)
#
# and indeed:

for inp in multi:
    print(inp["ecut"])
# To change the values in a particular dataset use:
multi[0].set_vars(ngkpt=[2, 2, 2], tsmear=0.004)
multi[1].set_vars(ngkpt=[4, 4, 4], tsmear=0.008)

To build a table with the values of ngkpt and tsmear:

multi.get_vars_dataframe("ngkpt", "tsmear")
multi

Warning

Remember that in python we start to count from zero hence the first dataset has index 0.

Calling set_structure on MultiDataset will set the structure of the inputs:

multi.set_structure(abidata.cif_file("si.cif"))

# The structure attribute of a MultiDataset returns a list of structures
# equivalent to [inp.structure for inp in multi]
print(multi.structure)

The function split_datasets return the list of AbinitInput stored in MultiDataset

inp0, inp1 = multi.split_datasets()
inp0

Note

You can use MultiDataset to build your input files but remember that Abipy workflows will never support input files with more than one dataset. As a consequence, you should always pass an AbinitInput to the AbiPy functions that are building Tasks, Works or Flows.

print("Number of datasets:", multi.ndtset)
# To create and append a new dataset (initialized from dataset number 1)
multi.addnew_from(1)
multi[-1].set_vars(ecut=42)
print("Now multi has", multi.ndtset, "datasets and the ecut in the last dataset is:",
      multi[-1]["ecut"])