Light-weight tight-binding framework
Contents
Light-weight tight-binding framework#
This pages describe the usage and functionality of tblite tight-binding framework. The tblite project aims to provide an efficent and uniform interface to the (extended) tight-binding Hamiltonians.
Important
This project does not provide a complete atomistic simulation environment, but the tools to build one. If you are looking for a simulation environment that provides the functionality available in tblite, checkout the available integrations.
Warning
The current state of this project should be considered as highly experimental.
Installation#
This project is currently in a highly experimental stage, however we can already offer binary distributions via some channels.
Installing from conda#
This project is packaged for the conda package manager and available on the conda-forge channel. To install the conda package manager we recommend the miniforge installer. If the conda-forge channel is not yet enabled, add it to your channels with
conda config --add channels conda-forge
Once the conda-forge channel has been enabled, this project can be installed with:
conda install tblite
It is possible to list all of the versions available on your platform with:
conda search tblite --channel conda-forge
Now you are ready to use the tblite
executable, find the tblite.h
header and link against the tblite
library.
FreeBSD Port#
A port for FreeBSD is available
pkg install science/tblite
In case no package is available build the port using
cd /usr/ports/science/tblite
make install clean
For more information see the tblite port details.
Building from source#
This library depends on several Fortran modules to provide the desired functionality
mctc-lib: Modular computation tool chain library
dftd4: Reference implementation of the generally applicable charge-dependent London-dispersion correction, DFT-D4
s-dftd3: Reimplementation of the DFT-D3 dispersion correction
mstore: Molecular structure store (testing only)
toml-f: Library for processing and emitting TOML data
Meson based build#
The primary build system of this project is meson. For the full feature-complete build it is highly recommended to perform the build and development with meson. To setup a build the following software is required
A Fortran 2008 compliant compiler, like GCC Fortran and Intel Fortran classic
meson, version 0.57.2 or newer
ninja, version 1.10 or newer
a linear algebra backend, like MKL or OpenBLAS
Optional dependencies are
asciidoctor, to build the manual pages
A C compiler to test the C API and compile the Python extension module
Python 3.6 or newer with the CFFI package installed to build the Python API
To setup a new build run
meson setup _build --prefix=$HOME/.local
The Fortran and C compiler can be selected with the FC
and CC
environment variable, respectively.
The installation location is selected using the --prefix
option.
The required Fortran modules will be fetched automatically from the upstream repositories and checked out in the subprojects directory.
Note
For Intel Fortran oneAPI (2021 or newer) builds with MKL backend the -Dfortran_link_args=-qopenmp
option has to be added.
Tip
To produce statically linked binaries set --default-library=static
and add -Dfortran_link_args=-static
as well.
To compile the project run
meson compile -C _build
Verify that the resulting projects is working correctly by running the testsuite with
meson test -C _build --print-errorlogs
In case meson is spawning too many concurrent test jobs limit the number of processes with the --num-processes
option when running the test command.
By default the whole library and its subprojects are tested, to limit the testing to the project itself add --suite tblite
as option.
To verify the included parametrizations are working correctly run the extra testsuite by passing the --benchmark
argument
meson test -C _build --print-errorlogs --benchmark
Finally, you can make the project available by installation with
meson install -C _build
CMake based build#
This project also provides support for CMake to give projects using it as build system an easier way to interface. The CMake build files usually do not provide a feature-complete build, but contributions are more than welcome. To setup a build the following software is required
A Fortran 2008 compliant compiler, like GCC Fortran and Intel Fortran classic
cmake, version 3.14 or newer
ninja, version 1.10 or newer
a linear algebra backend, like MKL or OpenBLAS
Configure a new build with
cmake -B _build -G Ninja -DCMAKE_INSTALL_PREFIX=$HOME/.local
You can set the Fortran compiler in the FC
environment variable.
The installation location can be selected with the CMAKE_INSTALL_PREFIX
, GNU install directories are supported by default.
CMake will automatically fetch the required Fortran modules, you can provide specific version in the subprojects directory which will be used instead.
To run a build use
cmake --build _build
After a successful build make sure the testsuite passes
pushd _build && ctest --output-on-failure && popd
To make the project available install it with
cmake --install _build
Fpm based build#
This projects supports building with the Fortran package manager (fpm). Create a new build by running
fpm build
You can adjust the Fortran compiler with the --compiler
option and select the compilation profile with --profile release
.
To test the resulting build run the testsuite with
fpm test
The command line driver can be directly used from fpm wih
fpm run --profile release -- --help
To make the installation accessible install the project with
fpm install --profile release --prefix $HOME/.local
Python extension module#
The Python API is available as Python extension module.
The easiest way to setup is to add -Dpython=true
to a meson tree build and follow the meson installation instructions.
The extension module will become available once the project is installed.
Important
When building with Intel compilers make sure to use the real-time version of the MKL.
Add -Dlapack=mkl-rt
when configuring the build.
Otherwise, when using the normal MKL libraries dynamically loading the tblite library from Python will fail.
This section describes alternative ways to build the Python API
Using pip#
This project support installation with pip as an easy way to build the Python API.
C compiler to build the C-API and compile the extension module (the compiler name should be exported in the
CC
environment variable)Python 3.6 or newer
The following Python packages are required additionally
Make sure to have your C compiler set to the CC
environment variable
export CC=gcc
Install the project with pip
pip install .
To install extra dependencies as well use
pip install '.[ase]'
Using meson#
The Python extension module can be built on-top of an existing installation, either provided by meson or CMake.
Building requires against an existing tblite installation requires
C compiler to build the C-API and compile the extension module
meson version 0.55 or newer
a build-system backend, i.e. ninja version 1.7 or newer
Python 3.6 or newer with the CFFI package installed
Setup a build with
meson setup _build_python python -Dpython_version=$(which python3)
The Python version can be used to select a different Python version, it defaults to 'python3'
.
Python 2 is not supported with this project, the Python version key is meant to select between several local Python 3 versions.
Compile the project with
meson compile -C _build
The extension module is now available in _build_python/tblite/_libtblite.*.so
.
You can install as usual with
meson configure _build --prefix=/path/to/install
meson install -C _build
Supported Compilers#
This is a non-comprehensive list of tested compilers for tblite. Compilers with the label latest are tested with continuous integration for each commit.
Compiler |
Version |
Platform |
Architecture |
tblite |
---|---|---|---|---|
GCC |
11.1, 10.3 |
Ubuntu 20.04 |
x86_64 |
0.2.0, 0.2.1, 0.3.0, latest |
GCC |
10.3, 9.4 |
MacOS 11.6.5 |
x86_64 |
0.2.0, 0.2.1, 0.3.0, latest |
GCC |
9.4 |
MacOS 10.15.7 |
x86_64 |
0.2.0, 0.2.1 |
GCC |
11.0 |
MacOS 11.0 |
arm64 |
0.2.0, 0.2.1 |
GCC |
10.3 |
CentOS 7 |
aarch64, ppc64le |
0.2.0, 0.2.1 |
GCC/MinGW |
11.2 |
Windows 2022 |
x86_64 |
0.2.1 0.3.0, latest |
Intel |
2021.2 |
Ubuntu 20.04 |
x86_64 |
0.2.0, 0.2.1, 0.3.0, latest |
NAG |
7.1 |
AlmaLinux 8.5 |
x86_64 |
0.2.0, 0.2.1 |
Compiler known to fail are documented here, together with the last commit where this behaviour was encountered. If available an issue in on the projects issue tracker or the issue tracker of the dependencies is linked. Usually, it safe to assume that older versions of the same compiler will fail to compile as well and this failure is consistent over platforms and/or architectures.
Compiler |
Version |
Platform |
Architecture |
Reference |
---|---|---|---|---|
GCC |
6.4.0 |
MacOS 10.15.7 |
x86_64 |
|
Intel |
19.0.5 |
AlmaLinux 8.5 |
x86_64 |
|
Intel |
17.0.1 |
OpenSuse 42.1 |
x86_64 |
|
Intel |
16.0.3 |
CentOS 7.3 |
x86_64 |
|
Flang |
20190329 |
Ubuntu 20.04 |
x86_64 |
|
NVHPC |
20.9 |
Manjaro Linux |
x86_64 |
Tutorials#
This section contains tutorials with step by step instructions for using tblite and is meant for beginners not yet familiar with the command line interace.
Note
This section is currently under constructions more tutorials will be added in the future.
Single point calculations#
Evaluating the tight-binding methods is possible using the tblite-run command.
Calculating with the standalone#
The tblite standalone allows to perform single point calculations with the built-in methods and parameter files.
Energy evaluation with different methods#
To run a first calculation we only need a geometry input in one of the supported formats.
As a simple example we choose the caffeine molecule and save it as struc.xyz
:
24
# caffeine
C 1.07317 0.04885 -0.07573
N 2.51365 0.01256 -0.07580
C 3.35199 1.09592 -0.07533
N 4.61898 0.73028 -0.07549
C 4.57907 -0.63144 -0.07531
C 3.30131 -1.10256 -0.07524
C 2.98068 -2.48687 -0.07377
O 1.82530 -2.90038 -0.07577
N 4.11440 -3.30433 -0.06936
C 5.45174 -2.85618 -0.07235
O 6.38934 -3.65965 -0.07232
N 5.66240 -1.47682 -0.07487
C 7.00947 -0.93648 -0.07524
C 3.92063 -4.74093 -0.06158
H 0.73398 1.08786 -0.07503
H 0.71239 -0.45698 0.82335
H 0.71240 -0.45580 -0.97549
H 2.99301 2.11762 -0.07478
H 7.76531 -1.72634 -0.07591
H 7.14864 -0.32182 0.81969
H 7.14802 -0.32076 -0.96953
H 2.86501 -5.02316 -0.05833
H 4.40233 -5.15920 0.82837
H 4.40017 -5.16929 -0.94780
To perform a single point calculation we use the tblite-run subcommand:
tblite run struc.xyz
The output is written to the standard output, no further files are produced at this point. By default the GFN2-xTB method is used.
repulsion energy 4.9222938485600E-01 Eh
dispersion energy 9.6582008183032E-05 Eh
number of electrons 7.4000000000000E+01 e
integral cutoff 1.7347787504999E+01 bohr
------------------------------------------------------------
cycle total energy energy error density error
------------------------------------------------------------
1 -41.75569301832 -4.2248019E+01 1.9479952E-01
2 -42.10609780052 -3.5040478E-01 7.5972187E-02
3 -42.14069903418 -3.4601234E-02 4.6343369E-02
4 -42.14507107704 -4.3720429E-03 1.2550655E-02
5 -42.14724061029 -2.1695333E-03 6.1305110E-03
6 -42.14770623967 -4.6562938E-04 2.4358040E-03
7 -42.14738299681 3.2324286E-04 1.0726483E-03
8 -42.14747017326 -8.7176447E-05 3.6117531E-04
9 -42.14746157709 8.5961684E-06 1.6698164E-04
10 -42.14746263709 -1.0599990E-06 5.4777834E-05
11 -42.14746337487 -7.3778077E-07 2.5306995E-05
12 -42.14746403618 -6.6130801E-07 8.2333411E-06
------------------------------------------------------------
electronic energy -4.2639790003039E+01 Eh
total energy -4.2147464036175E+01 Eh
------------------------------------------------------------------
# Occupation Energy/Eh Energy/eV
------------------------------------------------------------------
1 2.0000 -0.7201238 -19.5956
... ... ... ...
29 2.0000 -0.4792385 -13.0407
30 2.0000 -0.4764138 -12.9639
31 2.0000 -0.4564304 -12.4201
32 2.0000 -0.4379226 -11.9165
33 2.0000 -0.4264752 -11.6050
34 2.0000 -0.4193573 -11.4113
35 2.0000 -0.4148063 -11.2875
36 2.0000 -0.4035253 -10.9805
37 2.0000 -0.3892737 -10.5927 (HOMO)
38 -0.2657013 -7.2301 (LUMO)
39 -0.2265846 -6.1657
40 -0.2098114 -5.7093
41 -0.1225595 -3.3350
42 -0.0559053 -1.5213
43 -0.0080016 -0.2177
44 0.0133961 0.3645
... ... ...
66 0.6592542 17.9392
------------------------------------------------------------------
HL-Gap 0.1235725 Eh 3.3626 eV
------------------------------------------------------------------
We can select a method by adding the --method
argument to the command line and rerun
tblite run --method gfn1 struc.xyz
This will perform the calculation again using the GFN1-xTB method.
Calculating derivatives#
So far only energies are calculated, to enable the evaluation of molecular gradients we add the --grad
keyword to the command line for our GFN1-xTB calculation
tblite run --grad --method gfn1 struc.xyz
Note
The gradient option takes one optional argument, make sure to separate the geometry input file by the --
separator in case you do not provide the output file for the derivatives:
tblite run --grad -- struc.xyz
The output will show that the gradient output has been written to a file
[Info] Tight-binding results written to 'tblite.txt'
Tip
Generally, tblite will never write to a file without explicitly stating it in the standard output.
We can inspect the result using a pager like less or any text editor of your choice, the output should look similar to this:
energy :real:0:
-4.4509702550094509E+01
gradient :real:2:3,24
3.8751207233845646E-03 1.0589576505981196E-03 4.0101833008853537E-06
-6.8197644662543388E-03 -2.5469121433516095E-02 -4.7470927704034663E-05
1.0648728891662693E-02 8.8366335898543858E-03 4.7996917287539728E-05
3.1902722527341862E-04 7.4242248622642722E-03 -2.8496918953763304E-05
-2.6208973446986007E-02 -9.0131378921360604E-03 -3.8669871255270171E-05
-1.7509318247989529E-03 3.6705299785918994E-03 3.4118097773531792E-05
1.9484148290568090E-02 4.7415207908131207E-03 4.0471496678204742E-05
-1.1722356880465069E-02 7.4828408528079754E-03 -1.1612153192208481E-04
-1.2241607563977187E-03 -1.9132992220688570E-02 1.2461589545076379E-04
-1.6766717498318165E-02 1.2103660828186829E-02 -3.4573645420344302E-05
2.4844203069064843E-02 -1.9287321306216326E-02 -4.9197479201307158E-05
1.1084471631753511E-02 1.7004159137354469E-02 -3.2214753172841759E-06
-6.4592191454479423E-03 -4.2860776915391351E-03 -4.0858618678683595E-06
-6.5899183170761598E-03 1.3221070732562403E-02 -9.0165446293088262E-05
-1.3268732296168747E-03 2.1245633439957766E-03 4.0929253563542303E-06
1.7227511511327858E-03 1.1717306153688867E-03 2.6844186750092148E-03
1.7254183949399506E-03 1.1710614866087398E-03 -2.6879103256391306E-03
-1.5499406234334675E-03 3.9313477434577580E-03 2.5224714308249690E-05
6.3570530436816057E-03 2.2899202031480908E-03 3.4807215781027190E-06
5.1739877369726952E-04 -5.0616912429584586E-04 1.9423398812888478E-03
5.4061192657724807E-04 -4.8581923304342556E-04 -1.9455409758072379E-03
-1.2675930676672822E-03 -7.6587507428546708E-03 8.4916914284767576E-05
2.7792577845120209E-04 -1.8293010741325709E-04 1.4003154522332375E-03
2.8959035627490622E-04 -2.0990206390925764E-04 -1.3505474151682963E-03
virial :real:2:3,3
7.8713512660909479E-02 -3.0635463233997452E-02 -1.6134089053178151E-04
-3.0635463233997452E-02 8.4807095398895124E-02 -1.1909491256157945E-04
-1.6134089053176757E-04 -1.1909491256159344E-04 2.0323518845924864E-02
The ASCII output format is adopted from the DFTB+ project and called tagged output (see tblite-tag). It has the advantage of being easily human-readable and a custom parser is straight-forward to write.
Instead of relying on a custom parser, a more standardized format like JSON might be preferable for automating calculations with tblite.
We can request to create to use JSON to store the results with the --json
argument.
tblite run --grad --json --method gfn1 struc.xyz
In the output we will now find
[Info] JSON dump of results written to 'tblite.json'
Note
The tagged output file will not be written by default when requesting a JSON output, to enable both provide the name of the output files in the command line
tblite run --grad caffeine.txt --json caffeine.json struc.xyz
The output will show that both output are now available
[Info] Tight-binding results written to 'caffeine.txt'
[Info] JSON dump of results written to 'caffeine.json'
We find that tblite also included meta data like the version number of the program in the JSON dump, which can be useful when automatically processing the data later.
Parameter optimization#
Optimization of parameters is available with the tblite-fit subcommand.
Setting up a fit#
In this example we try to create a parametrization for titanium oxides based on GFN2-xTB
Important
We will assume a data set to fit against is already existing and prepared here.
First, we create our base parameter file using the tblite-param command.
tblite param --method gfn2 --output gfn2-xtb.toml
Now, we can inspect the parameter file, we will find six main sections named hamiltonian, dispersion, repulsion, charge, thirdorder, and multipole as well as the element sections.
In the next step we create our input file for the parameter optimization, the input file is written in TOML:
script = "exec ./run.sh"
[mask]
hamiltonian = {}
dispersion = {}
repulsion = {}
charge = {}
thirdorder = {}
multipole = {}
[mask.element]
O = {lgam=[false, true]}
Ti = {lgam=[true, false, true]}
The most important part is the mask section as it selects the parameters to optimize. Here we include all sections from above to enable them in the element records automatically. To actually select the elements we add them in mask.element, at this place we can also overwrite the automatic defaults. We will only disable the shell hardnesses for the s shells here, since those should use the atomic hardnesses unmodified.
With our exported parameter file and the fit input we start a test run to check if everything is working correctly
tblite fit --dry-run gfn2-xtb.toml tbfit.toml
If tblite returns an error, check the input file again for typos.
The next step is to create the runner, usually a shell script is sufficient for this purpose.
The script will get its input via environment variables.
The fit driver will set the name of the tblite executable in TBLITE_EXE
, we can use this variable to use the same version of tblite to read the modified parameter file again and perform single points.
Also, we get the name of the created parameter file in TBLITE_PAR
and the name of the data output we are supposed to generate in the command in TBLITE_OUT
.
We use the script shown below
#!/usr/bin/env bash
# Uncomment for debugging
#set -ex
# Input from fit runner with defaults for standalone use
fitpar=$(realpath ${TBLITE_PAR:-./fitpar.toml})
output=$(realpath ${TBLITE_OUT:-.data})
tblite=$(which ${TBLITE_EXE:-tblite})
# Balance to get nproc == njob * nthreads
nthreads=1
njob=$(nproc | awk "{print int(\$1/$nthreads)}")
# Temporary data file
data=.data
# Temporary wrapper
wrapper=./wrapped_runner
# Arguments for tblite runner
tblite_args="run --param \"$fitpar\" --grad results.tag coord"
# Ad hoc error in case the Hamiltonian does not work
# (SCC does not converge or similar)
penalty="1.0e3"
# Create our wrapper script
cat > "$wrapper" <<-EOF
#!/usr/bin/env bash
if [ -d \$1 ]; then
pushd "\$1" > /dev/null 2>&1
test -f "$data" && rm "$data"
OMP_NUM_THREADS=1 "$tblite" $tblite_args > tblite.out 2> tblite.err \
|| echo "0.0 $penalty # run: \$1" > "$data"
"$tblite" tagdiff --fit results.tag reference.tag >> "$data" \
|| echo "0.0 $penalty # diff: \$1" >> "$data"
fi
EOF
chmod +x "$wrapper"
# Create the actual multiprocessing queue
printf "%s\0" data/*/ | xargs -n 1 -P $njob -0 "$wrapper"
# Collect the data
cat data/*/$data > "$output"
# Cleanup
rm data/*/$data
rm "$wrapper"
In short this script will process multiple structures in parallel and create the final data output.
The structures we want to evaluate with this script are expected to be in the subdirectories of the data
directory stored as Turbomole coord
files.
Always run the command you want to use before start the fit driver in production mode, for this purpose also uncomment the debugging line in the script and run it with
TBLITE_PAR=./gfn2-xtb.toml TBLITE_OUT=data.txt time ./run.sh
The most important step is to check the output generated in data.txt
.
Failed runs will be marked with the directory name, inspect the output of tblite to check whether there is a setup error that has to be fixed first before starting into production.
We also use the time
command here to determine the approximate runtime of our script.
A good target is a maximum runtime of one second.
Note
While long evaluation times are possible, it makes the fit more difficult in practice. A good rule of thumb is that everything above 10 seconds will become problematic.
If necessary the data set has to be cut down by a bit or smaller input have to be used. Alternatively, a more powerful machine should be used if available. Especially for large numbers of data more cores can help to reduce the wall time.
Finally, if you are happy with the setup start the actual fit in verbose mode.
tblite fit -v gfn2-xtb.toml tbfit.toml
Using the verbose printout will show the objective function in every step as well as the relative change in all parameters. Initially, the verbose printout is useful to track the stability of the fit.
Tip
To create long-running fits detach the fit driver from the shell using
nohup tblite fit -v gfn2-xtb.toml tbfit.toml > fitlog.txt 2>&1 &
You can check the output in fitlog.txt
while the fit driver is running in the background.
Once the fit finishes the optimized parameters are found in fitpar.toml
.
You should update the meta section to specify how those parameters were fitted.
[meta]
version = 1
name = "...-xTB"
reference = "..."
Important
The method we are fitting is an extended tight-binding (xTB) Hamiltonian. A parametrization like the second geometry, frequency, non-covalent interaction parametrization (GFN2) results in the final method name GFN2-xTB. Parametrizations are always given as prefix to xTB, it is highly discouraged to add something to the xTB part itself, i.e. there is no xTB2.
Available integrations#
The following packages are currently using or in the process of integrating tblite:
Usage in xtb#
The tblite project originated from the xtb program package.1
Note
Support in xtb
is planned for the version 6.5.0 release.
Literature#
- 1
Christoph Bannwarth, Eike Caldeweyher, Sebastian Ehlert, Andreas Hansen, Philipp Pracht, Jakob Seibert, Spicher Spicher, and Stefan Grimme. Extended tight‐binding quantum chemistry methods. WIREs Comput. Mol. Sci., 11:e01493, 2020. doi:10.1002/wcms.1493.
Usage in DFTB+#
The tblite project started as an effort to make the xTB Hamiltonian available in the DFTB+ program package.1
DFTB+ supports a communication bridge to the tblite library when setting -DWITH_TBLITE=true
at the CMake configuration step.
Note
Support in DFTB+ is available since version 21.2, released in fall 2021.
Installing DFTB+ and the tblite library#
A complete DFTB+ package is available via conda-forge. To install the conda package manager we recommend the miniforge installer. If the conda-forge channel is not yet enabled, add it to your channels with
conda config --add channels conda-forge
Once the conda-forge channel has been enabled, DFTB+ can be installed with:
conda install 'dftbplus=*=nompi_*'
Or to install the MPI enabled version, you can add mpich or openmpi as MPI provider or just let conda choose:
conda install 'dftbplus=*=mpi_*'
It is possible to list all of the versions available on your platform with:
conda search dftbplus --channel conda-forge
Now you are ready to use the dftb+
executable, find the tblite library is installed as well as part of the DFTB+ dependencies.
Input structure#
Input to DFTB+ is provided by creating a file named dftb_in.hsd
in the custom human-friendly structured data (HSD) format, which is inspired by XML.
The geometry input is provided in the Geometry group, either as xyzFormat for molecular geometries or as vaspFormat for periodic geometries.
To enable the xTB methods the Hamiltonian group is set to xTB.
In the xTB group the Method keyword is provided to select between GFN2-xTB, GFN1-xTB and IPEA1-xTB.
Geometry = xyzFormat {
<<< "struc.xyz"
}
Hamiltonian = xTB {
Method = "GFN2-xTB"
}
To perform periodic calculations with the xTB Hamiltonian only the k-point sampling has to be added to the xTB group with kPointsAndWeights. Using the SuperCellFolding provides Monkhorst–Pack k-point sampling.
Geometry = vaspFormat {
<<< "POSCAR"
}
Hamiltonian = xTB {
Method = "GFN1-xTB"
kPointsAndWeights = SuperCellFolding {
2 0 0
0 2 0
0 0 2
0.5 0.5 0.5
}
}
Instead of providing a Method the xTB method can be initialized from a parameter file by providing its path in ParameterFile.
# ...
Hamiltonian = xTB {
ParameterFile = "gfn2-xtb.toml"
# ...
}
Finally, to perform more than just single point calculations, the Driver group has to be provided.
Best is to use the new GeometryOptimization driver, which defaults to a rational function minimizer as present in the xtb
program package.
For periodic structures the lattice optimization can be enabled by setting LatticeOpt to Yes.
Driver = GeometryOptimization {
LatticeOpt = Yes
}
For the full capabilities for DFTB+ check the reference manual.
Literature#
- 1
B. Hourahine, B. Aradi, V. Blum, F. Bonafé, A. Buccheri, C. Camacho, C. Cevallos, M. Y. Deshaye, T. Dumitrică, A. Dominguez, S. Ehlert, M. Elstner, T. van der Heide, J. Hermann, S. Irle, J. J. Kranz, C. Köhler, T. Kowalczyk, T. Kubař, I. S. Lee, V. Lutsker, R. J. Maurer, S. K. Min, I. Mitchell, C. Negre, T. A. Niehaus, A. M. N. Niklasson, A. J. Page, A. Pecchia, G. Penazzi, M. P. Persson, J. Řezáč, C. G. Sánchez, M. Sternberg, M. Stöhr, F. Stuckenberg, A. Tkatchenko, V. W.-z. Yu, and T. Frauenheim. DFTB+, a software package for efficient approximate density functional theory based atomistic simulations. J. Chem. Phys., 152(12):124101, 2020. doi:10.1063/1.5143190.
Usage in QCxMS#
The QCxMS project for the quantum chemical simulation of mass spectra12 uses tblite as backend to implement the GFN1-xTB, GFN2-xTB and IPEA1-xTB Hamiltonians.
Note
Support in QCxMS is available with version 5.1.0 and later. For more details see the release notes of version 5.1.0.
Literature#
- 1
Vilhjálmur Ásgeirsson, Christoph A Bauer, and Stefan Grimme. Quantum chemical calculation of electron ionization mass spectra for general organic and inorganic molecules. Chem. Sci., 8(7):4879–4895, 2017. doi:10.1039/c7sc00601b.
- 2
Jeroen Koopman and Stefan Grimme. Calculation of electron ionization mass spectra with semiempirical gfnn-xtb methods. ACS Omega, 4(12):15120–15133, 2019. doi:10.1021/acsomega.9b02011.
Usage in ASE#
The Python API of tblite natively support integration with the atomic simulation environment (ASE). By constructing a calculator most functionality of ASE is readily available. For details on building the Python API checkout the installation guide.
ASE calculator#
- class tblite.ase.TBLite(*args: Any, **kwargs: Any)[source]#
ASE calculator for using xTB Hamiltonians from the tblite library. Supported properties by this calculator are:
energy (free_energy)
forces
stress
dipole
charges
Supported keywords are
Keyword
Default
Description
method
“GFN2-xTB”
Underlying method for energy and forces
accuracy
1.0
Numerical accuracy of the calculation
electronic_temperature
300.0
Electronic temperatur in Kelvin
max_iterations
250
Iterations for self-consistent evaluation
cache_api
True
Reuse generate API objects (recommended)
verbosity
1
Set verbosity of printout
Example
An ASE calculator can be constructed by using the TBLite class provided by the tblite.ase module. For example to perform a single point calculation for a CO2 crystal use
>>> from tblite.ase import TBLite >>> from ase.atoms import Atoms >>> import numpy as np >>> atoms = Atoms( ... symbols="C4O8", ... positions=np.array( ... [ ... [0.9441259872, 0.9437851680, 0.9543505632], ... [3.7179966528, 0.9556570368, 3.7316862240], ... [3.7159517376, 3.7149292800, 0.9692330016], ... [0.9529872864, 3.7220864832, 3.7296981120], ... [1.6213905408, 1.6190616096, 1.6313879040], ... [0.2656685664, 0.2694175776, 0.2776540416], ... [4.3914553920, 1.6346256864, 3.0545920000], ... [3.0440834880, 0.2764611744, 4.4080419264], ... [4.3910577696, 3.0416409504, 0.2881058304], ... [3.0399936576, 4.3879335936, 1.6497353376], ... [0.2741322432, 4.4003734944, 3.0573754368], ... [1.6312174944, 3.0434586528, 4.4023048032], ... ] ... ), ... cell=np.array([5.68032, 5.68032, 5.68032]), ... pbc=np.array([True, True, True]), ... ) >>> atoms.calc = TBLite(method="GFN1-xTB") >>> atoms.get_potential_energy() # result in eV -1257.0943962462964
The resulting calculator can be used like most ASE calculator, e.g. for optimizing geometries.
- calculate(atoms: Optional[ase.atoms.Atoms] = None, properties: Optional[List[str]] = None, system_changes: List[str] = ase.calculators.calculator.all_changes) None [source]#
Perform actual calculation with by calling the TBLite API
Example
>>> from ase.build import molecule >>> from tblite.ase import TBLite >>> calc = TBLite(method="GFN2-xTB") >>> calc.calculate(molecule("H2O")) >>> calc.get_potential_energy() -137.96777625229421 >>> calc.calculate(molecule("CH4")) >>> calc.get_potential_energy() -113.60956621093894
- Raises
ase.calculators.calculator.InputError – on invalid input passed to the interface module
ase.calculators.calculator.CalculationFailed – in case of an RuntimeError in the library
- reset() None [source]#
Clear all information from old calculation. This will only remove the cached API objects in case the cache_api is set to False.
- set(**kwargs) dict [source]#
Set new parameters to TBLite. Will automatically reconstruct the underlying model in case critical parameters change.
Example
>>> from ase.build import molecule >>> from tblite.ase import TBLite >>> atoms = molecule("H2O") >>> atoms.calc = TBLite(method="GFN2-xTB") >>> atoms.get_potential_energy() -137.96777625229421 >>> atoms.calc.set(method="GFN1-xTB") {'method': 'GFN1-xTB'} >>> atoms.get_potential_energy() -156.9675057724589
Parameter Specification#
This section contains the detailed specification of all contributions entering the tight-binding Hamiltonian. Additionally, the format to construct the Hamiltonian is described for writing parameter files. A list of built-in methods is available as well.
Built-in methods#
The library supports GFN1-xTB1, GFN2-xTB2, and IPEA1-xTB3 as built-in parametrizations.
The parameters can be loaded in the standalone using the --method
keyword and can be loaded directly from the API without requiring to create a parametrization first.
method |
keyword |
---|---|
GFN1-xTB |
gfn1 |
GFN2-xTB |
gfn2 |
IPEA1-xTB |
ipea1 |
The built-in parameters can be exported using the tblite-param command.
tblite param --method gfn2 --output gfn2-xtb.toml
For more information on the format of the parameter file see the parameter specification.
Hamiltonian Specification#
Contents
Basis set definition#
The Hamiltonian is defined in a minimal basis set using the linear combination of atomic orbitals (LCAO) approach. Atomic orbitals are represented by contracted Gaussian basis functions of the STO-NG type basis sets. This implementation support overlap integrals up to g-functions and one moment lower, f-functions, for derivatives.
Spherical harmonics ordering#
Standard spherical harmonic sorting is used for the basis functions.
angular momentum |
order |
---|---|
0 (s) |
0 |
1 (p) |
-1, 0, 1 |
2 (d) |
-2, -1, 0, 1, 2 |
3 (f) |
-3, -2, -1, 0, 1, 2, 3 |
4 (g) |
-4, -3, -2, -1, 0, 1, 2, 3, 4 |
Potential definition#
The potential shift vectors are defined with respect the partial charge and cumulative atomic moments. This is the negative of the potential with respect to the populations. The sign is applied when adding the potential shifts to the Hamiltonian.
Note
This convention is in line with the sign convention for partial charges, but is unintuitive from a population / density perspective. The sign convention might change to the latter case as the latter convention is more in line with usual density functional theory.
Overlap integrals#
The overlap integral is defined as
Normalisation of the basis function in practise is accurate to 10:sup:-10, due to the precision of the tabulated STO-NG contraction coefficients.1
Atomic partial charges#
Orbital populations are obtained by Mulliken population analysis
and summed up to respective shell-resolved or atomic partial charges.
Dipole moment integrals#
The dipole moment integral is defined origin independent as
The dipole operator is always placed on the aufpunkt of the basis function in ket.
Cumulative dipole moments#
The dipole moments are evaluated from the dipole moment integrals using Mulliken population analysis
and summed up to respective cumulative shell-resolved or atomic dipole moments.
Quadrupole moment integrals#
The quadrupole moment integral is defined origin independent as
The quadrupole operator is always placed on the aufpunkt of the basis function in ket. Only the lower triangle of the traceless quadrupole integrals are saved packed (ordering: xx, xy, yy, xz, yz, zz).
Cumulative quadrupole moments#
The quadrupole moments are evaluated from the quadrupole moment integrals using Mulliken population analysis
and summed up to respective cumulative shell-resolved or atomic quadrupole moments.
Literature#
- 1
Robert F. Stewart. Small gaussian expansions of Slater-type orbitals. J. Chem. Phys., 52:431–438, 1970. doi:10.1063/1.1672702.
Repulsive contributions#
Effective coulombic repulsion#
The xTB Hamiltonian is using an effective coulomb-like repulsion term
where RAB is the interatomic distance and ZAB, αAB, rAB and kAB are pairwise parameters for each species pair. The pairwise parameters are usually formed from element-specific parameters to limit the actual amount of free parameters. The effective nuclear charges are formed by
Repulsion exponents are obtained as geometric mean from element specific exponents by
Finally, the exponents, rAB and kAB, of the distance dependent contributions are usually set globally for all element-pairs.
Electrostatic interactions#
Contents
Second order#
Isotropic electrostatics#
The isotropic electrostatic in a shell-resolved formulation is given by the parametrized Coulomb interaction between shellwise partial charges
The interaction potential is parametrized by a Klopman–Ohno type potential in the xTB Hamiltonian or the γ-functional as used in the DFTB Hamiltonian.
Klopman–Ohno kernel#
The interaction kernel for the Klopman–Ohno electrostatic is given by
where η:sub:A/B are the chemical hardness parameters of the respective shells and g is the exponent to manipulate the potential shape.
γ-functional kernel#
The interaction kernel for the DFTB γ-functional is derived from the integral of two exponential densities
where τ:sub:A/B are scaled Hubbard parameters of the respective shells and R is the distance between the atomic sides.
Anisotropic electrostatics#
The anisotropic electrostatic in an atom-resolved formulation is given by the multipole interactions between the different moments:
Third order#
The isotropic third-order contributions are included as the trace of the on-site shell-resolved Hubbard derivatives.
London dispersion interactions#
D3 dispersion model#
The D3 dispersion model is provided by the s-dftd3 library. The DFT-D3 dispersion correction1 is used together with the rational damping function.2
D4 dispersion model#
The D4 dispersion model is provided by the dftd4 project and extended in this project as a self-consistent dispersion model.3
The self-consistent D4 dispersion model can be classified as infinite order charge-dependent contribution in the density expansion, due to the non-linear dependency on the atomic partial charges in the charge scaling function.
Literature#
- 1
S. Grimme, J. Antony, S. Ehrlich, and H. Krieg. A consistent and accurate ab initio parametrization of density functional dispersion correction (DFT-D) for the 94 elements H-Pu. J. Chem. Phys., 132:154104, 2010.
- 2
S. Grimme, S. Ehrlich, and L. Goerigk. Effect of the damping function in dispersion corrected density functional theory. J. Comput. Chem., 32:1456–1465, 2011.
- 3
Eike Caldeweyher, Sebastian Ehlert, Andreas Hansen, Hagen Neugebauer, Sebastian Spicher, Christoph Bannwarth, and Stefan Grimme. A generally applicable atomic-charge dependent London dispersion correction. J. Chem. Phys., 150(15):154122, 2019. doi:10.1063/1.5090222.
Parameter file#
The tblite parameter files are written in TOML. Each major contribution to the energy is specified in its own table. The presence of the table enables the contribution in the method. The parameters are splitted in global and element specific parameters, global parameters are specified in the respective method table, while element specific parameters are collected together in the element records. The exception are the optional element pairwise parameters for the scaling of the Hamiltonian elements which are specified in the Hamiltonian section.
TOML quickstart#
Note
This guide will only cover the aspects of TOML relevant for representing the parameter file in this project. For an extensive guide on TOML visit its homepage.
We choose TOML to represent our parameter files since it produces both human and machine readable files and has good support over all major progamming languages. TOML is a configuration file format, which allows to represent data in tables, arrays and values. Values can have data types of integer, boolean, real, and string.
Note
Integers do not cast implicitly to real numbers in TOML, make sure to actually specify a real number in the parametrization file if a real is required. Some TOML parsers offer an implicit conversion but this behaviour is not generally available in standard-compliant TOML parser.
Tables are specified by headers enclosed in brackets, nesting is represented by dots in the header name.
[hamiltonian]
# ...
[hamiltonian.xtb]
# ...
[hamiltonian.xtb.kpair]
# ...
The toplevel header is optional as long as the table only contains other tables, the hamiltonian header can be left away without changing the data structure. Small tables can be represented by inline tables with curly braces:
[dispersion]
d4 = {sc=true, s8=2.70, a1=0.52, a2=5.00, s9=5.00}
The above structure is equivalent to
[dispersion.d4]
s6 = 1.00
s8 = 2.70
a1 = 0.52
a2 = 5.00
s9 = 5.00
sc = true
Arrays can be created by using brackets after a key
[element.C]
shells = [ "2s", "2p" ]
levels = [ -13.970922, -10.063292 ]
slater = [ 2.096432, 1.80 ]
ngauss = [ 4, 4 ]
refocc = [ 1.0, 3.0 ]
# ...
Note
In TOML 0.5.0 arrays were constrained to equivalent data only, this restriction was lifted in TOML 1.0.0. In this format no mixed data arrays are required, which allows the usage of TOML 0.5.0 compliant parsers as well, if no updated library is available yet.
For libraries to handle TOML documents visit the TOML wiki.
Meta section#
The meta section allows to specify data describing the parametrization, like the name of the method, the version of the parametrization and the format used to specify it as well as relevant publications for this parametrizations. Only the data type of the entries is checked and whether the format version is supported by the library. The current format version is 1, specifying a higher version will result in an error.
[meta]
format = 1
name = "GFN2-xTB"
version = 1
reference = "DOI: 10.1021/acs.jctc.8b01176"
Allowed entries:
Keyword |
Description |
Type |
---|---|---|
format |
version of the parameter file format |
integer |
name |
name of the parametrization |
string |
version |
version of the parametrization data |
integer |
reference |
relevant publications for this parameters |
string |
Hamiltonian section#
The Hamiltonian section is used to declare the details of the used Hamiltonian. Currently, only xTB type Hamiltonians can be declared in the xtb subtable.
[hamiltonian.xtb]
wexp = 0.5
enscale = 2.0e-2
cn = "gfn"
shell = {ss=1.85, pp=2.23, dd=2.23, sd=2.0, pd=2.0}
Pair parameters can be specified for all elements of the element records. Specifying an X-Y entry implies the Y-X entry and only one will be read and used since the pair parameters cannot be asymetric. The default value is one for each pair.
[hamiltonian.xtb.kpair]
H-H = 0.96
H-B = 0.95
H-N = 1.04
N-Si = 1.01
B-P = 0.97
Sc-Sc = 1.10
Sc-Ti = 1.10
# ...
Note
The kpair section can be used to tune the Hamiltonian to better describe certain bonding situations. The suitable range parameters is close to one and large deviations from unity will create an instable Hamiltonian.
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
wexp |
scaling from slater exponents |
real |
dimensionless |
enscale |
scaling from EN difference |
real |
dimensionless |
cn |
CN type for selfenergy shift |
string |
|
shell |
shell-specific scaling |
table of reals |
dimensionless |
kpair |
pairwise scaling for elements |
table of reals |
dimensionless |
Dispersion section#
The dispersion section supports the d4 subtable for DFT-D4 type dispersion corrections and the d3 subtable for DFT-D3 type dispersion corrections. The rational (Becke–Johnson) damping scheme is always used.
[dispersion]
d4 = {sc=true, s8=2.70, a1=0.52, a2=5.00, s9=5.00}
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
s6 |
scaling for C6 dispersion terms |
real |
dimensionless |
s8 |
scaling for C8 dispersion terms |
real |
dimensionless |
a1 |
scaling of critical radius |
real |
dimensionless |
a2 |
offset for critical radius |
real |
Bohr |
s9 |
scaling for triple-dipole terms |
real |
dimensionless |
sc |
use self-consistent dispersion |
logical |
Repulsion section#
The xTB repulsion term can be specified in the effective subtable.
[repulsion]
effective = {kexp=1.5}
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
kexp |
exponent for repulsion |
real |
dimensionless |
klight |
exponent for light atom pairs |
real |
dimensionless |
Halogen section#
The GFN1-xTB specific halogen bonding correction can be specified in the classical subtable.
[halogen]
classical = {rscale=1.3, damping=0.44}
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
rscale |
scaling parameter for radii |
real |
dimensionless |
damping |
damping parameter |
real |
dimensionless |
Charge section#
The Klopman–Ohno parametrized electrostatic model is available with the effective subtable, while the DFTB γ-functional electrostatic can be enabled with the gamma subtable. Only one electrostatic model can be active at a time.
Klopman–Ohno electrostatic#
The gam parameter in the element records is used as atomic Hubbard parameters and scaled with the lgam parameter to obtain shell-resolved values. The exponent of the kernel can be modified as well as the averaging scheme for the Hubbard parameters. Available averaging schemes are arithmetic (GFN2-xTB), harmonic (GFN1-xTB) and geometric. The electrostatic is always constructed shell-resolved.
[charge]
effective = {gexp=2.0, average="arithmetic"}
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
gexp |
exponent of kernel |
real |
dimensionless |
average |
averaging scheme |
string |
DFTB γ-functional electrostatic#
The gam parameter in the element records is used as Hubbard parameter and scaled with the lgam parameter to obtain shell-resolved values. The electrostatic is always constructed shell-resolved.
[charge]
gamma = {}
Thirdorder section#
An on-site thirdorder charge is supported, to use atomic Hubbard derivatives the shell keyword can be set to false, while for shell-resolved Hubbard derivatives the scaling parameters for the respective shells have to specified. The highest specified angular momentum is implicitly used for all but absent higher momenta.
[thirdorder]
shell = {s=1.00, p=0.50, d=0.25}
Important
While the following setup uses the atomic Hubbard derivative for all shells
[thirdorder]
shell.s = 1.0
it is fundamentally different from using an atom-resolved third-order model.
Multipole section#
The anisotropic electrostatic of GFN2-xTB can be enabled using the damped subtable. It requires five parameters to setup the damping function to reduce the short-range contributions from the multipole electrostatics.
[multipole]
damped = {dmp3=3.0, dmp5=4.0, kexp=4.0, shift=1.2, rmax=5.0}
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
dmp3 |
damping for quadratic terms |
real |
dimensionless |
dmp5 |
damping for cubic terms |
real |
dimensionless |
kexp |
exponent for multipole radii |
real |
dimensionless |
shift |
shift for valence CN value |
real |
dimensionless |
rmax |
maximum multipole radius |
real |
dimensionless |
Element records#
The main body of the parameter file contains of element records. The parameters here are used to initialize contributions from the tables other tables, but are collected in the element records for easy usage. Most keywords require entries, even if the respective contribution is not used in the method.
Note
Each record is identified by its symbol, which allows to have multiple parameter sets for the same element. Input elements which do not match any symbol, will use the parametrization of the first element record with the same atomic number. To ensure that the right element is used as fallback an ordered dictionary is recommended to represent the element records.
[element.H]
shells = [ "1s" ]
levels = [ -10.707211 ]
slater = [ 1.23 ]
ngauss = [ 3 ]
refocc = [ 1.0 ]
kcn = [ -5.0e-2 ]
gam = 0.405771
lgam = [ 1.0 ]
gam3 = 0.08
zeff = 1.105388
arep = 2.213717
en = 2.20
dkernel = 5.563889e-2
qkernel = 2.7431e-4
mprad = 1.4
mpvcn = 1.0
[element.C]
shells = [ "2s", "2p" ]
levels = [ -13.970922, -10.063292 ]
slater = [ 2.096432, 1.80 ]
ngauss = [ 4, 4 ]
refocc = [ 1.0, 3.0 ]
kcn = [ -1.02144e-2, 1.61657e-2 ]
gam = 5.38015e-1
lgam = [ 1.0, 1.1056358 ]
gam3 = 1.50e-1
zeff = 4.231078
arep = 1.247655
en = 2.55
dkernel = -4.11674e-3
qkernel = 2.13583e-3
mprad = 3.0
mpvcn = 3.0
Allowed entries:
Keyword |
Description |
Type |
Unit |
---|---|---|---|
shells |
included valence shells |
array of strings |
dimensionless |
levels |
atomic self-energies |
array of reals |
eV |
slater |
exponents of basis functions |
array of reals |
1/Bohr² |
ngauss |
number of STO-NG primitives |
array of integers |
dimensionless |
refocc |
reference occupation of atom |
array of reals |
Unitcharge |
kcn |
CN dependent self-energy shift |
array of reals |
eV |
shpoly |
polynomial enhancement factor |
array of reals |
dimensionless |
gam |
atomic Hubbard parameter |
real |
Hartree/Unitcharge² |
lgam |
relative shell hardness |
array of reals |
dimensionless |
gam3 |
atomic Hubbard derivative |
real |
Hartree/Unitcharge³ |
zeff |
effective nuclear charge |
real |
Unitcharge |
arep |
repulsion exponent |
real |
dimensionless |
dkernel |
on-site dipole kernel |
real |
Hartree |
qkernel |
on-site quadrupole kernel |
real |
Hartree |
mprad |
critical multipole radius |
real |
Bohr |
mpvcn |
multipole valence CN |
real |
dimensionless |
xbond |
halogen bonding strength |
real |
Hartree |
en |
atomic electronegativity |
real |
dimensionless |
Public API documentation#
Fortran API#
The tblite library seamlessly integrates with other Fortran projects via module interfaces,
Note
Generally, all quantities used in the library are stored in atomic units.
Handling of geometries and structure#
The basic infrastructure to handle molecular and periodic structures is provided by the modular computation tool chain library. The library provides a structure type which is used to represent all geometry related informations in tblite. A structure type can be constructed from arrays or read from a file.
The constructor is provided with the generic interface new
and takes an array of atomic numbers (integer
) or element symbols (character(len=*)
) as well as the cartesian coordinates in Bohr.
Additionally, the molecular charge and the number of unpaired electrons can be provided the charge
and uhf
keyword, respectively.
To create a periodic structure the lattice parameters can be passed as 3 by 3 matrix with the lattice
keyword.
An example for using the constructor is given here
subroutine example
use mctc_env, only : wp
use mctc_io, only : structure_type, new
implicit none
type(structure_type) :: mol
real(wp), allocatable :: xyz(:, :)
integer, allocatable :: num(:)
num = [6, 1, 1, 1, 1]
xyz = reshape([ &
& 0.00000000000000_wp, -0.00000000000000_wp, 0.00000000000000_wp, &
& -1.19220800552211_wp, 1.19220800552211_wp, 1.19220800552211_wp, &
& 1.19220800552211_wp, -1.19220800552211_wp, 1.19220800552211_wp, &
& -1.19220800552211_wp, -1.19220800552211_wp, -1.19220800552211_wp, &
& 1.19220800552211_wp, 1.19220800552211_wp, -1.19220800552211_wp],&
& [3, size(num)])
call new(mol, num, xyz, charge=0.0_wp, uhf=0)
! ...
end subroutine example
To interact with common input file formats for structures the read_structure
procedure is available.
The file type is inferred from the name of the file automatically or if a file type hint is provided directly from the enumerator of available file types.
The read_structure
routine can also use an already opened unit, but in this case the file type hint is mandatory to select the correct format to read from.
subroutine example
use mctc_env, only : error_type
use mctc_io, only : structure_type, read_structure, file_type
implicit none
type(structure_type) :: mol
type(error_type), allocatable :: error
character(len=:), allocatable :: input
input = "struc.xyz"
call read_structure(mol, input, error, file_type%xyz)
if (allocated(error)) then
print '(a)', error%message
stop 1
end if
! ...
end subroutine example
The structure type as well as the error type are using only allocatable members and can therefore be used without requiring explicit deconstruction.
Certain members of the structure type should be considered immutable, like the number of atoms (nat
), the identifiers for unique atoms (id
) and the boundary conditions (periodic
).
To change those specific structure parameters the structure type and all dependent objects should be reconstructed to ensure a consistent setup.
Other properties, like the geometry (xyz
), molecular charge (charge
), number of unpaired electrons (uhf
) and lattice parameters (lattice
) can be changed without requiring to reconstruct dependent objects like calculators or restart data.
Error handling#
The basic error handler is an allocatable derived type, available from mctc_env
as error_type
, which signals an error by its allocation status.
use mctc_env, only : error_type, fatal_error
implicit none
type(error_type), allocatable :: error
call always_ok(error)
if (allocated(error)) then
print '(a)', "Unexpected failure:", error%message
end if
call always_failed(error)
if (allocated(error)) then
print '(a)', "Error:", error%message
end if
contains
subroutine always_ok(error)
type(error_type), allocatable, intent(out) :: error
end subroutine always_ok
subroutine always_failed(error)
type(error_type), allocatable, intent(out) :: error
call fatal_error(error, "Message associated with this error")
end subroutine always_failed
end
An unhandled error might get dropped by the next procedure call.
Calculation context#
The calculation context is available with the context_type
from the tblite_context
module.
The context stores error messages generated while running which can be queried using the type bound function failed
.
To access the actual errors the messages can be removed using the type bound subroutine get_error
.
An output verbosity is available in the context as the member verbosity, all procedures with access to the context will default to the verbosity of the context unless the verbosity level is overwritten by an argument.
To cutomize the output the context_logger
abstract base class is available.
It must implement a type bound message
procedure, which is used by the context to create output.
This type can be used to create callbacks for customizing or redirecting the output of the library.
High-level interface#
The high-level interface is defined by the calculation context, the calculator instance and its restart data.
The calculation context is defined with the context_type
, which stores general settings regarding the overall method independent setup of the calculation.
The actual parametrisation data is stored in the xtb_calculator
type.
An instance of the calculator can be used in a thread-safe way to perform calculations for a specific structure (defined by its number of atoms, unique elements and boundary conditions).
Changing the specific structure parameters requires to reconstruct the calculator.
Finally the specific persient data for a geometry is stored in a wavefunction_type
, which allows to restart calculations based on previous results.
C API#
The C API bindings are provided by using the iso_c_binding
intrinsic module.
Generally, objects are exported as opaque pointers and can only be manipulated within the library.
The API user is required delete all objects created in the library by using the provided deconstructor functions to avoid memory leaks.
Overall four classes of objects are provided by the library
error handlers (
tblite_error
andtblite_context
), used to communicate exceptional conditions and errors from the library to the user, also usable to communicate general setting to the library in case of a context handlestructure containers (
tblite_structure
), used to represent the system specific information and geometry data, only the latter are mutable for the usercalculator objects (
tblite_calculator
), used to store actual method parametrisationresult and restart data (
tblite_result
) used to hold calculated properties and persistent data for restarting between calculations
Note
Generally, all quantities provided to the library are assumed to be in atomic units.
Main header#
#include "tblite.h"
Provides convenience functionality for working with the C API bindings for tblite by including all relevant headers and defining general type generic macros for working with the object handles.
Defines
-
tblite_delete(ptr)#
Generic macro to free an object handle.
- Parameters
ptr – Opaque object handle
Version query#
#include "tblite/version.h"
Provides access to the version, compatibility and features exported by this API.
Functions
-
TBLITE_API_ENTRY int TBLITE_API_CALL tblite_get_version(void)#
Retrieve version of library used
- Returns
Compact version number in the format 10000 * major + 100 * minor + patch
Macro definitions#
#include "tblite/macros.h"
General macro definitions for the tblite C API bindings.
Error handler#
#include "tblite/error.h"
Provides a light-weight error handler for communicating with the library.
The library provides two different kinds handlers, a light error handle type tblite_error is defined here.
The error handle is used in the context of simple tasks and requires only small overhead to construct and use. It is mainly used in the context of retrieving data or building structures.
Typedefs
-
typedef struct _tblite_error *tblite_error#
Error instance.
Functions
-
TBLITE_API_ENTRY tblite_error TBLITE_API_CALL tblite_new_error(void)#
Create new error handle object.
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_error(tblite_error*)#
Delete an error handle object.
-
TBLITE_API_ENTRY int TBLITE_API_CALL tblite_check_error(tblite_error)#
Check error handle status.
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_clear_error(tblite_error)#
Clear error handle status.
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_error(tblite_error, char*, const int*)#
Get error message from error handle.
Context handler#
#include "tblite/context.h"
Provides a persistent configuration object to modify the behaviour of a calculation. Acts as an error handler.
The environment context tblite_context can be considered a persistent setup for all calculations performed with the library, it is usually used together with calculator objects tblite_calculator. While the error handle can only contain a single error, multiple errors can be accumulated in a context object, which allows storing more complex error information like they can occur in an actual calculation.
Typedefs
-
typedef struct _tblite_context *tblite_context#
Context manager for the library usage.
-
typedef void (*tblite_logger_callback)(char*, int, void*)#
Define callback function for use in custom logger.
Functions
-
TBLITE_API_ENTRY tblite_context TBLITE_API_CALL tblite_new_context(void)#
Create new calculation environment object
- Returns
New context handle
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_context(tblite_context *ctx)#
Delete a calculation environment object
- Parameters
ctx – Context handle
-
TBLITE_API_ENTRY int TBLITE_API_CALL tblite_check_context(tblite_context ctx)#
Check calculation environment status
- Parameters
ctx – Context handle
- Returns
Current status, 0: okay, 1: error
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_context_error(tblite_context ctx, char *buffer, const int *buffersize)#
Get error message from calculation environment
- Parameters
ctx – Context handle
buffer – Character buffer for writing error message to
buffersize – Maximum length of buffer (optional)
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_context_logger(tblite_context ctx, tblite_logger_callback callback, void *userdata)#
Set custom logger function
- Parameters
ctx – Context handle
callback – Procedure pointer implementing logger
userdata – Passthrough data pointer
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_context_color(tblite_context ctx, int color)#
Enable colorful output
- Parameters
ctx – Context handle
color – Set color support, 0: disabled, 1: enabled
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_context_verbosity(tblite_context ctx, int verbosity)#
Set verbosity level of printout
- Parameters
ctx – Context handle
verbosity – Printout verbosity
Molecular structure data#
#include "tblite/structure.h"
The structure data tblite_structure is used to represent the system of interest in the library.
It contains immutable system specific information like the number of atoms, the unique atom groups and the boundary conditions as well as mutable geometry data like cartesian coordinates and lattice parameters.
To change immutable parameters of the tblite_structure data the object should be reinstantiated as well as all dependent objects, like tblite_calculator or tblite_result instances.
Typedefs
-
typedef struct _tblite_structure *tblite_structure#
Molecular structure data class.
The structure data is used to represent the system of interest in the library. It contains immutable system specific information like the number of atoms, the unique atom groups and the boundary conditions as well as mutable geometry data like cartesian coordinates and lattice parameters.
Functions
-
TBLITE_API_ENTRY tblite_structure TBLITE_API_CALL tblite_new_structure(tblite_error error, const int natoms, const int *numbers, const double *positions, const double *charge, const int *uhf, const double *lattice, const bool *periodic)#
Create new molecular structure data
- Parameters
error – Handle for error messages
natoms – Number of atoms
numbers – Atomic numbers for each atom, shape [natoms]
positions – Cartesian coordinates in Bohr for each atom, shape [natoms][3]
charge – Total charge of the system, (optional)
uhf – Number of unpaired electrons, (optional)
lattice – Lattice parameters in Bohr, shape [3][3], (optional)
periodic – Periodic dimensions, shape [3], (optional)
- Returns
New molecular structure data
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_structure(tblite_structure *mol)#
Delete molecular structure data
- Parameters
mol – Molecular structure data
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_update_structure_geometry(tblite_error error, tblite_structure mol, const double *positions, const double *lattice)#
Update coordinates and lattice parameters (quantities in Bohr)
- Parameters
error – Handle for error messages
mol – Molecular structure data
positions – Cartesian coordinates in Bohr for each atom, shape [natoms][3]
lattice – Lattice parameters in Bohr, shape [3][3], (optional)
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_update_structure_charge(tblite_error error, tblite_structure mol, const double *charge)#
Update total charge in structure object
- Parameters
error – Handle for error messages
mol – Molecular structure data
charge – Total charge of the system
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_update_structure_uhf(tblite_error error, tblite_structure mol, const int *uhf)#
Update number of unpaired electrons in structure object
- Parameters
error – Handle for error messages
mol – Molecular structure data
uhf – Number of unpaired electrons
Calculator instance#
#include "tblite/calculator.h"
Provides a single point calculator for performing extended tight-binding computations.
Provides a parametrized type tblite_calculator storing the method parameters required for the evaluation of single point calculations. The calculator is parametrized for a method as well as for the molecular structure data tblite_structure it was created with. Changes in the tblite_structure required to reinstantiate the structure data also require to reinstantiate the tblite_calculator object.
Typedefs
-
typedef struct _tblite_calculator *tblite_calculator#
Single point calculator.
Enums
Functions
-
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_gfn2_calculator(tblite_context ctx, tblite_structure mol)#
Construct calculator with GFN2-xTB parametrisation loaded
- Parameters
ctx – Context handle
mol – Molecular structure data
- Returns
New calculator instance
-
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_gfn1_calculator(tblite_context ctx, tblite_structure mol)#
Construct calculator with GFN1-xTB parametrisation loaded
- Parameters
ctx – Context handle
mol – Molecular structure data
- Returns
New calculator instance
-
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_ipea1_calculator(tblite_context ctx, tblite_structure mol)#
Construct calculator with IPEA1-xTB parametrisation loaded
- Parameters
ctx – Context handle
mol – Molecular structure data
- Returns
New calculator instance
-
TBLITE_API_ENTRY tblite_calculator TBLITE_API_CALL tblite_new_xtb_calculator(tblite_context ctx, tblite_structure mol, tblite_param param)#
Construct calculator from parametrization records
- Parameters
ctx – Context handle
mol – Molecular structure data
param – Parametrization records
- Returns
New calculator instance
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_calculator(tblite_calculator *calc)#
Delete calculator
- Parameters
calc – Calculator instance
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_accuracy(tblite_context ctx, tblite_calculator calc, double acc)#
Set calculation accuracy for the calculator object
- Parameters
ctx – Context handle
calc – Calculator instance
acc – Accuracy value for numerical thresholds
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_max_iter(tblite_context ctx, tblite_calculator calc, int max_iter)#
Set maximum number of allowed iterations in calculator object
- Parameters
ctx – Context handle
calc – Calculator instance
max_iter – Maximum number of allowed iterations
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_mixer_damping(tblite_context ctx, tblite_calculator calc, double damping)#
Set damping parameter for mixer in calculator object
- Parameters
ctx – Context handle
calc – Calculator instance
damping – Value for mixer damping
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_guess(tblite_context ctx, tblite_calculator calc, tblite_guess guess)#
Set initial guess for creating new wavefunction objects
- Parameters
ctx – Context handle
calc – Calculator instance
guess – Guess for initializing the wavefunction
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_temperature(tblite_context ctx, tblite_calculator calc, double etemp)#
Set electronic temperature for the calculator object (in Hartree)
- Parameters
ctx – Context handle
calc – Calculator instance
etemp – Electronic temperature in Hartree
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_set_calculator_save_integrals(tblite_context ctx, tblite_calculator calc, int save_integrals)#
Set the flag in the calculator to retain the integral matrices
- Parameters
ctx – Context handle
calc – Calculator instance
save_integrals – Flag to enable storing of integrals
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_calculator_push_back(tblite_context ctx, tblite_calculator calc, tblite_container *cont)#
Add container to calculator object.
Note: Ownership is transferred and container handle is destroyed after function call
- Parameters
ctx – Context handle
calc – Calculator instance
cont – Interaction container
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_shell_count(tblite_context ctx, tblite_calculator calc, int *nsh)#
Query calculator for the number of shells
- Parameters
ctx – Context handle
calc – Calculator instance
nsh – Number of shells
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_shell_map(tblite_context ctx, tblite_calculator calc, int *sh2at)#
Query calculator for index mapping from shells to atomic centers
- Parameters
ctx – Context handle
calc – Calculator instance
sh2at – Index mapping from shells to atomic centers
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_angular_momenta(tblite_context ctx, tblite_calculator calc, int *am)#
Query calculator for angular momenta of shells
- Parameters
ctx – Context handle
calc – Calculator instance
am – Angular momenta of each shell
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_orbital_count(tblite_context ctx, tblite_calculator calc, int *nao)#
Query calculator for the number of orbitals
- Parameters
ctx – Context handle
calc – Calculator instance
nao – Number of atomic orbitals
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_calculator_orbital_map(tblite_context ctx, tblite_calculator calc, int *ao2sh)#
Query calculator for index mapping from atomic orbitals to shells
- Parameters
ctx – Context handle
calc – Calculator instance
ao2sh – Index mapping from atomic orbitals to shells
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_singlepoint(tblite_context ctx, tblite_structure mol, tblite_calculator calc, tblite_result res)#
Perform single point calculation
- Parameters
ctx – Context handle
mol – Molecular structure data
calc – Calculator instance
res – Result container
Interaction container#
#include "tblite/container.h"
Provides an interaction container which can be added to a tblite_calculator.
Typedefs
-
typedef struct _tblite_container *tblite_container#
Interaction container.
Functions
-
TBLITE_API_ENTRY tblite_container TBLITE_API_CALL tblite_new_electric_field(double *efield)#
Create new electric field container
- Parameters
efield – Electric field in atomic units (Hartree/(Bohr*e)), shape: [3]
- Returns
New interaction container
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_container(tblite_container *cont)#
Delete container handle
- Parameters
cont – Container handle
Result container#
#include "tblite/result.h"
Provides a storage container tblite_result for all persistent quantities produced in a calculation by tblite_get_singlepoint. The data stored in tblite_result can be used as restart data for subsequent calculations.
The individual entries of the tblite_result container can be queried and retrieved. For some quantities, like the count of the shells and orbitals duplicated queries are implemented to allow the usage of the tblite_result container without requiring access to the tblite_calculator which produced it.
Typedefs
-
typedef struct _tblite_result *tblite_result#
Container to for storing and handling calculation results.
Functions
-
TBLITE_API_ENTRY tblite_result TBLITE_API_CALL tblite_new_result(void)#
Create new result container
- Returns
New result container
-
TBLITE_API_ENTRY tblite_result TBLITE_API_CALL tblite_copy_result(tblite_result res)#
Create new result container from existing container
- Parameters
res – Existing result container
- Returns
Copy of the results in a new container
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_result(tblite_result *res)#
Delete a calculation environment object
- Parameters
res – Result container
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_number_of_atoms(tblite_error error, tblite_result res, int *nat)#
Retrieve number of atoms from result container
- Parameters
error – Handle for error messages
res – Result container
nat – Number of atoms
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_number_of_shells(tblite_error error, tblite_result res, int *nsh)#
Retrieve number of shells from result container
- Parameters
error – Handle for error messages
res – Result container
nsh – Number of shells
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_number_of_orbitals(tblite_error error, tblite_result res, int *nao)#
Retrieve number of orbitals from result container
- Parameters
error – Handle for error messages
res – Result container
nao – Number of orbitals
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_energy(tblite_error error, tblite_result res, double *energy)#
Retrieve energy from result container
- Parameters
error – Handle for error messages
res – Result container
energy – Total energy
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_energies(tblite_error error, tblite_result res, double *energies)#
Retrieve atom-resolved energies from result container
- Parameters
error – Handle for error messages
res – Result container
energies – Atom-resolved energies, shape [nat]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_gradient(tblite_error error, tblite_result res, double *gradient)#
Retrieve gradient from result container
- Parameters
error – Handle for error messages
res – Result container
gradient – Cartesian gradient, shape [nat][3]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_virial(tblite_error error, tblite_result res, double *sigma)#
Retrieve virial from result container
- Parameters
error – Handle for error messages
res – Result container
sigma – Strain derivatives, shape [3][3]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_charges(tblite_error error, tblite_result res, double *charges)#
Retrieve atomic charges from result container
- Parameters
error – Handle for error messages
res – Result container
charges – Atomic partial charges, shape [nat]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_bond_orders(tblite_error error, tblite_result res, double *mbo)#
Retrieve bond orders from result container
- Parameters
error – Handle for error messages
res – Result container
mbo – Bond orders, shape [nat][nat]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_dipole(tblite_error error, tblite_result res, double dipole[3])#
Retrieve dipole moment from result container (order x, y, z)
- Parameters
error – Handle for error messages
res – Result container
dipole – Total dipole moment, shape [3]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_quadrupole(tblite_error error, tblite_result res, double quadrupole[6])#
Retrieve traceless quadrupole moment from result container (packed xx, xy, yy, xz, yz, zz)
- Parameters
error – Handle for error messages
res – Result container
quadrupole – Total traceless quadrupole moment, shape [6]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_orbital_energies(tblite_error error, tblite_result res, double *emo)#
Retrieve orbital energies from result container
- Parameters
error – Handle for error messages
res – Result container
emo – Eigenvalues for each orbital
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_orbital_occupations(tblite_error error, tblite_result res, double *occ)#
Retrieve orbital occupations from result container
- Parameters
error – Handle for error messages
res – Result container
occ – Occupation numbers for each orbital
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_orbital_coefficients(tblite_error error, tblite_result res, double *cmat)#
Retrieve orbital coefficients from result container
- Parameters
error – Handle for error messages
res – Result container
cmat – Orbital coefficient matrix, shape [nspin][nao][nao]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_density_matrix(tblite_error error, tblite_result res, double *pmat)#
Retrieve density matrix from result container
- Parameters
error – Handle for error messages
res – Result container
pmat – Density matrix, shape [nspin][nao][nao]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_overlap_matrix(tblite_error error, tblite_result res, double *smat)#
Retrieve overlap matrix from result container
- Parameters
error – Handle for error messages
res – Result container
smat – Overlap matrix, shape [nao][nao]
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_get_result_hamiltonian_matrix(tblite_error error, tblite_result res, double *hmat)#
Retrieve Hamiltonian matrix from result container
- Parameters
error – Handle for error messages
res – Result container
hmat – Hamiltonian matrix, shape [nao][nao]
Parametrization records#
#include "tblite/param.h"
Provides a representation of a parametrization of an xTB Hamiltonian which can be used to instantiate a tblite_calculator object.
The parametrization data itself can be represented as a tblite_table data structure, which provides the possibility to customize the parametrization programmatically.
Typedefs
-
typedef struct _tblite_param *tblite_param#
Parametrization records.
Functions
-
TBLITE_API_ENTRY tblite_param TBLITE_API_CALL tblite_new_param(void)#
Create new parametrization records object
- Returns
Parametrization records
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_param(tblite_param *param)#
Delete a parametrization records object
- Parameters
param – Parametrization records
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_load_param(tblite_error error, tblite_param param, tblite_table table)#
Load parametrization records from data table
- Parameters
error – Handle for error messages
param – Parametrization records
table – Table data structure
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_dump_param(tblite_error error, tblite_param param, tblite_table table)#
Dump parametrization records to data table
- Parameters
error – Handle for error messages
param – Parametrization records
table – Table data structure
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_export_gfn2_param(tblite_error error, tblite_param param)#
Export GFN2-xTB parametrization records
- Parameters
error – Handle for error messages
param – Parametrization records
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_export_gfn1_param(tblite_error error, tblite_param param)#
Export GFN1-xTB parametrization records
- Parameters
error – Handle for error messages
param – Parametrization records
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_export_ipea1_param(tblite_error error, tblite_param param)#
Export IPEA1-xTB parametrization records
- Parameters
error – Handle for error messages
param – Parametrization records
Table data structure#
#include "tblite/table.h"
Provides a representation of a generic table data structure.
Used to mirror the data available in the tblite_param object. It aims to provide a programmatic accessible representation of the parametrization records.
Defines
-
tblite_table_set_value(error, table, key, value, ...)#
Generic setter based on the type of the value
- Parameters
error – Error handle
table – Table data structure
key – Key to set value at
value – Value to set at key
Typedefs
-
typedef struct _tblite_table *tblite_table#
Handle for holding a table data structure
The table can either own its data or reference another table. Table references can be created from existing table data structures or by adding new table entries into an existing table, which returns a reference to the newly created table.
Functions
-
TBLITE_API_ENTRY tblite_table TBLITE_API_CALL tblite_new_table(tblite_table *table)#
Create new data table object
- Parameters
table – Table object to reference in new table (optional)
- Returns
: New table data structure
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_delete_table(tblite_table *table)#
Delete a data table object
- Parameters
table – Table object to be deleted
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_double(tblite_error error, tblite_table table, char key[], double *value, int n)#
Set floating point number to data table
- Parameters
error – Error handle
table – Table data structure
key – Key to set value at
value – Double value or array to add to table
n – Number of entries to add, 0 for adding scalars
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_int64_t(tblite_error error, tblite_table table, char key[], int64_t *value, int n)#
Set integer number to data table (use int64_t rather than long)
- Parameters
error – Error handle
table – Table data structure
key – Key to set value at
value – Integer value or array to add to table
n – Number of entries to add, 0 for adding scalars
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_bool(tblite_error error, tblite_table table, char key[], bool *value, int n)#
Set boolean value to data table
- Parameters
error – Error handle
table – Table data structure
key – Key to set value at
value – Boolean value or array to add to table
n – Number of entries to add, 0 for adding scalars
-
TBLITE_API_ENTRY void TBLITE_API_CALL tblite_table_set_char(tblite_error error, tblite_table table, char key[], char (*value)[], int n)#
Set character string to data table
- Parameters
error – Error handle
table – Table data structure
key – Key to set value at
value – Character value or array to add to table
n – Number of entries to add, 0 for adding scalars
-
TBLITE_API_ENTRY tblite_table TBLITE_API_CALL tblite_table_add_table(tblite_error error, tblite_table table, char key[])#
Create new subtable in existing data table
- Parameters
error – Error handle
table – Table data structure
key – Key to add new subtable at
- Returns
New table data structure
Python API#
Python API of the tblite project
Note
The documentation of the ASE integration can be found in Usage in ASE.
Library interface#
Definition of the basic interface to library for most further integration in other Python frameworks. The classes defined here allow a more Pythonic usage of the library in actual workflows than the low-level access provided in the CFFI generated wrappers.
Structure#
- class tblite.interface.Structure(numbers: numpy.ndarray, positions: numpy.ndarray, charge: Optional[float] = None, uhf: Optional[int] = None, lattice: Optional[numpy.ndarray] = None, periodic: Optional[numpy.ndarray] = None)[source]#
Represents a wrapped structure object in
tblite
. The molecular structure data object has a fixed number of atoms and immutable atomic identifiers.Example
>>> from tblite.interface import Structure >>> import numpy as np >>> mol = Structure( ... positions=np.array([ ... [+0.00000000000000, +0.00000000000000, -0.73578586109551], ... [+1.44183152868459, +0.00000000000000, +0.36789293054775], ... [-1.44183152868459, +0.00000000000000, +0.36789293054775], ... ]), ... numbers = np.array([8, 1, 1]), ... ) ... >>> len(mol) 3
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
- update(positions: numpy.ndarray, lattice: Optional[numpy.ndarray] = None) None [source]#
Update coordinates and lattice parameters, both provided in atomic units (Bohr). The lattice update is optional also for periodic structures.
Generally, only the cartesian coordinates and the lattice parameters can be updated, every other modification, regarding total charge, total spin, boundary condition, atomic types or number of atoms requires the complete reconstruction of the object.
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
Calculator#
- class tblite.interface.Calculator(method: str, numbers: numpy.ndarray, positions: numpy.ndarray, charge: Optional[float] = None, uhf: Optional[int] = None, lattice: Optional[numpy.ndarray] = None, periodic: Optional[numpy.ndarray] = None)[source]#
Represents a wrapped calculator object in
tblite
and the associated structure data. The calculator is instantiated for the respective structure and is immutable once created. The cartesian coordinates and the lattice parameters of the structure can be updated, while changing the boundary conditions, atomic types or number of atoms require the complete reconstruction of the calculator instance.Available methods to parametrization of the calculator are:
GFN2-xTB:
Self-consistent extended tight binding Hamiltonian with anisotropic second order electrostatic contributions, third order on-site contributions and self-consistent D4 dispersion.
Geometry, frequency and non-covalent interactions parametrisation for elements up to Z=86.
Cite as:
C. Bannwarth, S. Ehlert and S. Grimme., J. Chem. Theory Comput. (2019), 15, 1652-1671. DOI: 10.1021/acs.jctc.8b01176
GFN1-xTB:
Self-consistent extended tight binding Hamiltonian with isotropic second order electrostatic contributions and third order on-site contributions.
Geometry, frequency and non-covalent interactions parametrisation for elements up to Z=86.
Cite as:
S. Grimme, C. Bannwarth, P. Shushkov, J. Chem. Theory Comput. (2017), 13, 1989-2009. DOI: 10.1021/acs.jctc.7b00118
IPEA1-xTB:
Special parametrisation for the GFN1-xTB Hamiltonian to improve the description of vertical ionisation potentials and electron affinities. Uses additional diffuse s-functions on light main group elements. Parametrised up to Z=86.
Cite as:
V. Asgeirsson, C. Bauer and S. Grimme, Chem. Sci. (2017), 8, 4879. DOI: 10.1039/c7sc00601b
Example
>>> from tblite.interface import Calculator >>> import numpy as np >>> numbers = np.array([1, 1, 6, 5, 1, 15, 8, 17, 13, 15, 5, 1, 9, 15, 1, 15]) >>> positions = np.array([ # Coordinates in Bohr ... [+2.79274810283778, +3.82998228828316, -2.79287054959216], ... [-1.43447454186833, +0.43418729987882, +5.53854345129809], ... [-3.26268343665218, -2.50644032426151, -1.56631149351046], ... [+2.14548759959147, -0.88798018953965, -2.24592534506187], ... [-4.30233097423181, -3.93631518670031, -0.48930754109119], ... [+0.06107643564880, -3.82467931731366, -2.22333344469482], ... [+0.41168550401858, +0.58105573172764, +5.56854609916143], ... [+4.41363836635653, +3.92515871809283, +2.57961724984000], ... [+1.33707758998700, +1.40194471661647, +1.97530004949523], ... [+3.08342709834868, +1.72520024666801, -4.42666116106828], ... [-3.02346932078505, +0.04438199934191, -0.27636197425010], ... [+1.11508390868455, -0.97617412809198, +6.25462847718180], ... [+0.61938955433011, +2.17903547389232, -6.21279842416963], ... [-2.67491681346835, +3.00175899761859, +1.05038813614845], ... [-4.13181080289514, -2.34226739863660, -3.44356159392859], ... [+2.85007173009739, -2.64884892757600, +0.71010806424206], ... ]) >>> calc = Calculator("GFN2-xTB", numbers, positions) >>> res = calc.singlepoint() >>> res.get("energy") # Results in atomic units -31.716159156026254
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
- get(attribute: str) Any [source]#
Set an attribute from the calculator instance. Supported attributes are
name
description
angular-momenta
Angular momenta of shells
orbital-map
Mapping from orbitals to shells
shell-map
Mapping from shells to atoms
- Raises
ValueError – on invalid attributes
- set(attribute: str, value) None [source]#
Set an attribute in the calculator instance. Supported attributes are
name
description
default
accuracy
Numerical thresholds for SCC
1.0
guess
Initial guess for wavefunction
0 (SAD)
max-iter
Maximum number of SCC iterations
250
mixer-damping
Parameter for the SCC mixer
0.4
save-integrals
Keep integral matrices in results
0 (False)
temperature
Electronic temperature for filling
9.500e-4
verbosity
Set verbosity of printout
1
Note
The electronic temperature is given in Hartree, rather than Kelvin. The conversion factor from Kelvin to Hartree is the Boltzmann constant in Hartree/Kelvin (3.166808578545117e-6).
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
- singlepoint(res: Optional[tblite.interface.Result] = None, copy: bool = False) tblite.interface.Result [source]#
Perform actual single point calculation in the library backend. The output of the library will be forwarded to the standard output.
The routine returns an object containing the results, which can be used to restart the calculation. Unless specified the restart object will be updated inplace rather than copied.
- Raises
RuntimeError – in case the calculation fails
Result#
- class tblite.interface.Result(other=None)[source]#
Container for calculation results, can be passed to a single point calculation as restart data. Allows to retrieve individual results or export all results as dict.
Example
>>> from tblite.interface import Calculator >>> import numpy as np >>> calc = Calculator( ... method="GFN2-xTB", ... numbers=np.array([14, 1, 1, 1, 1]), ... positions=np.array([ ... [ 0.000000000000, 0.000000000000, 0.000000000000], ... [ 1.617683897558, 1.617683897558,-1.617683897558], ... [-1.617683897558,-1.617683897558,-1.617683897558], ... [ 1.617683897558,-1.617683897558, 1.617683897558], ... [-1.617683897558, 1.617683897558, 1.617683897558], ... ]), ... ) >>> res = calc.singlepoint() ------------------------------------------------------------ cycle total energy energy error density error ------------------------------------------------------------ 1 -3.710873811182 -3.7424104E+00 1.5535536E-01 2 -3.763115011046 -5.2241200E-02 4.8803250E-02 3 -3.763205560205 -9.0549159E-05 2.0456319E-02 4 -3.763213200846 -7.6406409E-06 2.7461272E-03 5 -3.763250843634 -3.7642788E-05 2.4071582E-03 6 -3.763236758520 1.4085114E-05 2.3635321E-04 7 -3.763237287151 -5.2863152E-07 3.5812281E-04 8 -3.763233829618 3.4575334E-06 8.5040576E-06 9 -3.763233750684 7.8934355E-08 1.1095285E-07 ------------------------------------------------------------ >>> res.get("energy") -3.7632337506836944 >>> res = calc.singlepoint(res) ------------------------------------------------------------ cycle total energy energy error density error ------------------------------------------------------------ 1 -3.763233752583 -3.7947704E+00 1.3823013E-07 2 -3.763233751557 1.0256169E-09 1.4215249E-08 ------------------------------------------------------------
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
RuntimeError – if a requested quantity is not available in the container
- dict() dict [source]#
Return all quantities inside the result container as dict. In case no results are present an empty dict is returned.
- get(attribute: str)[source]#
Get a quantity stored instade the result container. The following quantities are available
property
dimension
unit
energy
scalar
Hartree
energies
nat
Hartree
gradient
nat, 3
Hartree/Bohr
virial
3, 3
Hartree
charges
n
e
dipole
3
e·Bohr
quadrupole
6
e·Bohr²
orbital-energies
norb
Hartree
orbital-occupations
norb
e
orbital-coefficients
norb
unitless
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
RuntimeError – if a requested quantity is not available in the container
- set(attribute: str, value)[source]#
Get a quantity stored instade the result container. Currently, no quantities can be set in the result container.
- Raises
ValueError – on invalid input, like incorrect shape / type of the passed arrays
RuntimeError – if a requested quantity cannot be set in the container