Python API of the DFT-D3 project
Project description
Python interface for the D3 dispersion model. This Python project is targeted at developers who want to interface their project via Python with s-dftd3.
This interface provides access to the C-API of s-dftd3 via the CFFI module. The low-level CFFI interface is available in the dftd3.libdftd3 module and only required for implementing other interfaces. A more pythonic interface is provided in the dftd3.interface module which can be used to build more specific interfaces.
from dftd3.interface import RationalDampingParam, DispersionModel
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],
])
model = DispersionModel(numbers, positions)
res = model.get_dispersion(RationalDampingParam(method="pbe0"), grad=False)
print(res.get("energy")) # Results in atomic units
# => -0.029489232932494884
QCSchema Integration
This Python API natively understands QCSchema and the QCArchive infrastructure. If the QCElemental package is installed the dftd3.qcschema module becomes importable and provides the run_qcschema function.
from dftd3.qcschema import run_qcschema
import qcelemental as qcel
atomic_input = qcel.models.AtomicInput(
molecule = qcel.models.Molecule(
symbols = ["O", "H", "H"],
geometry = [
0.00000000000000, 0.00000000000000, -0.73578586109551,
1.44183152868459, 0.00000000000000, 0.36789293054775,
-1.44183152868459, 0.00000000000000, 0.36789293054775
],
),
driver = "energy",
model = {
"method": "tpss",
},
keywords = {
"level_hint": "d3bj",
},
)
atomic_result = run_qcschema(atomic_input)
print(atomic_result.return_result)
# => -0.0004204244108151285
Installing
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 dftd3-python
Now you are ready to use dftd3, check if you can import it with
>>> import dftd3
>>> from dftd3.libdftd3 import get_api_version
>>> get_api_version()
'0.6.0'
Building the extension module
To perform an out-of-tree build some version of s-dftd3 has to be available on your system and preferably findable by pkg-config. Try to find a s-dftd3 installation you build against first with
pkg-config --modversion s-dftd3
Adjust the PKG_CONFIG_PATH environment variable to include the correct directories to find the installation if necessary.
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 .
If you already have a s-dftd3 installation, e.g. from conda-forge, you can build the Python extension module directly without cloning this repository
pip install "https://github.com/awvwgk/simple-dftd3/archive/refs/heads/main.zip#egg=dftd3-python&subdirectory=python"
Using meson
This directory contains a separate meson build file to allow the out-of-tree build of the CFFI extension module. The out-of-tree build requires
C compiler to build the C-API and compile the extension module
meson version 0.53 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 -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/dftd3/_libdftd3.*.so. You can install as usual with
meson configure _build --prefix=/path/to/install
meson install -C _build
Project details
Release history Release notifications | RSS feed
Download files
Download the file for your platform. If you're not sure which to choose, learn more about installing packages.
Source Distributions
Built Distributions
Hashes for dftd3-0.6.0-cp310-cp310-manylinux_2_12_x86_64.manylinux2010_x86_64.whl
Algorithm | Hash digest | |
---|---|---|
SHA256 | 13a29035c1a1ce761d687bac5141ab457dd5e02caf482a34a81b13924c81242b |
|
MD5 | f4b6a0ea844f680de3fa24f705b6a092 |
|
BLAKE2b-256 | b0cac745e17878a219fbfe7a2cf3d4e6a4cfd779c1a93d826814cc58d7cd056d |
Hashes for dftd3-0.6.0-cp39-cp39-manylinux_2_12_x86_64.manylinux2010_x86_64.whl
Algorithm | Hash digest | |
---|---|---|
SHA256 | b76107d85dc57cf8d38e7ba8e9f91b75225a25b029f73c2b8d613ce3b4142527 |
|
MD5 | e68003ccad60f7bfcac94b11683d8baf |
|
BLAKE2b-256 | 7246ee26ce9e07456aff3a63bf3ad6dde5bba2269509cd5c7b1858b535e7622d |
Hashes for dftd3-0.6.0-cp38-cp38-manylinux_2_12_x86_64.manylinux2010_x86_64.whl
Algorithm | Hash digest | |
---|---|---|
SHA256 | 9c2b0576845b111773619788a5765bf8c8178a68a3938f07e0a42e24cc87f501 |
|
MD5 | 7eeef4517c4625f75e35ef1610d2e474 |
|
BLAKE2b-256 | b5f93add6c6138dd43be53059ed5bd496b620809ac7eb1973d7988439cf1bfab |
Hashes for dftd3-0.6.0-cp37-cp37m-manylinux_2_12_x86_64.manylinux2010_x86_64.whl
Algorithm | Hash digest | |
---|---|---|
SHA256 | f59a34d007ddd0f4c27872c8da6f67b10ed7488ccc336d7ece1d4dba4de3d976 |
|
MD5 | 2b25e4ed62dd7435c08eb10de4e0acd5 |
|
BLAKE2b-256 | 8d6ca46933906019fc3a8c395f38e79eaa58caab24c023335d96c8ab8f63571d |
Hashes for dftd3-0.6.0-cp36-cp36m-manylinux2010_x86_64.whl
Algorithm | Hash digest | |
---|---|---|
SHA256 | 295218c3fcbe10e49ce590bc25e75c4413b1ae364ed14d63b2ed0f5dc272c0f3 |
|
MD5 | ddcbc1fdc8a326c613af782be0e9b654 |
|
BLAKE2b-256 | 510dbda4cfab8af7f27b29f79dcd10b0b0dd62bda576dfb2e75d8b2120d6e5aa |